2015-06-29 10:35:29 +02:00
|
|
|
subroutine get_s2(key_i,key_j,s2,Nint)
|
2015-04-20 16:45:06 +02:00
|
|
|
implicit none
|
|
|
|
use bitmasks
|
|
|
|
BEGIN_DOC
|
|
|
|
! Returns <S^2>
|
|
|
|
END_DOC
|
|
|
|
integer, intent(in) :: Nint
|
|
|
|
integer(bit_kind), intent(in) :: key_i(Nint,2)
|
|
|
|
integer(bit_kind), intent(in) :: key_j(Nint,2)
|
2015-06-29 10:35:29 +02:00
|
|
|
double precision, intent(out) :: s2
|
2015-04-20 16:45:06 +02:00
|
|
|
integer :: exc(0:2,2,2)
|
|
|
|
integer :: degree
|
|
|
|
double precision :: phase_spsm
|
|
|
|
integer :: nup, i
|
|
|
|
|
2015-06-29 10:35:29 +02:00
|
|
|
s2 = 0.d0
|
2015-04-20 16:45:06 +02:00
|
|
|
!$FORCEINLINE
|
|
|
|
call get_excitation_degree(key_i,key_j,degree,Nint)
|
|
|
|
select case (degree)
|
|
|
|
case(2)
|
2015-06-29 10:35:29 +02:00
|
|
|
call get_double_excitation(key_j,key_i,exc,phase_spsm,Nint)
|
2015-04-20 16:45:06 +02:00
|
|
|
if (exc(0,1,1) == 1) then ! Mono alpha + mono-beta
|
|
|
|
if ( (exc(1,1,1) == exc(1,2,2)).and.(exc(1,1,2) == exc(1,2,1)) ) then
|
2015-06-29 10:35:29 +02:00
|
|
|
s2 = -phase_spsm
|
2015-04-20 16:45:06 +02:00
|
|
|
endif
|
|
|
|
endif
|
|
|
|
case(0)
|
|
|
|
nup = 0
|
|
|
|
do i=1,Nint
|
|
|
|
nup += popcnt(iand(xor(key_i(i,1),key_i(i,2)),key_i(i,1)))
|
|
|
|
enddo
|
2015-06-29 10:35:29 +02:00
|
|
|
s2 = dble(nup)
|
2015-04-20 16:45:06 +02:00
|
|
|
end select
|
|
|
|
end
|
|
|
|
|
|
|
|
BEGIN_PROVIDER [ double precision, S_z ]
|
|
|
|
&BEGIN_PROVIDER [ double precision, S_z2_Sz ]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! z component of the Spin
|
|
|
|
END_DOC
|
|
|
|
|
|
|
|
S_z = 0.5d0*dble(elec_alpha_num-elec_beta_num)
|
|
|
|
S_z2_Sz = S_z*(S_z-1.d0)
|
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|
|
|
|
BEGIN_PROVIDER [ double precision, expected_s2]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Expected value of S2 : S*(S+1)
|
|
|
|
END_DOC
|
|
|
|
logical :: has_expected_s2
|
|
|
|
|
|
|
|
call ezfio_has_determinants_expected_s2(has_expected_s2)
|
|
|
|
if (has_expected_s2) then
|
|
|
|
call ezfio_get_determinants_expected_s2(expected_s2)
|
|
|
|
else
|
|
|
|
double precision :: S
|
|
|
|
S = (elec_alpha_num-elec_beta_num)*0.5d0
|
|
|
|
expected_s2 = S * (S+1.d0)
|
|
|
|
! expected_s2 = elec_alpha_num - elec_beta_num + 0.5d0 * ((elec_alpha_num - elec_beta_num)**2*0.5d0 - (elec_alpha_num-elec_beta_num))
|
|
|
|
endif
|
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|
|
|
|
BEGIN_PROVIDER [ double precision, s2_values, (N_states) ]
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! array of the averaged values of the S^2 operator on the various states
|
|
|
|
END_DOC
|
|
|
|
integer :: i
|
|
|
|
double precision :: s2
|
|
|
|
do i = 1, N_states
|
2015-06-29 10:35:29 +02:00
|
|
|
call get_s2_u0(psi_det,psi_coef(1,i),n_det,size(psi_coef,1),s2)
|
2015-04-20 16:45:06 +02:00
|
|
|
s2_values(i) = s2
|
|
|
|
enddo
|
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|
|
|
|
|
|
|
|
subroutine get_s2_u0(psi_keys_tmp,psi_coefs_tmp,n,nmax,s2)
|
|
|
|
implicit none
|
|
|
|
use bitmasks
|
|
|
|
integer(bit_kind), intent(in) :: psi_keys_tmp(N_int,2,nmax)
|
|
|
|
integer, intent(in) :: n,nmax
|
|
|
|
double precision, intent(in) :: psi_coefs_tmp(nmax)
|
|
|
|
double precision, intent(out) :: s2
|
|
|
|
integer :: i,j,l
|
|
|
|
double precision :: s2_tmp
|
2015-05-06 17:01:45 +02:00
|
|
|
s2 = 0.d0
|
2015-06-29 10:35:29 +02:00
|
|
|
!print*,'IN get_s2_u0'
|
|
|
|
!print*,'n,nmax = ',n,nmax
|
|
|
|
double precision :: accu
|
|
|
|
accu = 0.d0
|
|
|
|
do i = 1,n
|
|
|
|
accu += psi_coefs_tmp(i) * psi_coefs_tmp(i)
|
|
|
|
! print*,'psi_coef = ',psi_coefs_tmp(i)
|
|
|
|
enddo
|
|
|
|
!print*,'accu = ',accu
|
|
|
|
!print*,''
|
2015-04-20 16:45:06 +02:00
|
|
|
!$OMP PARALLEL DO DEFAULT(NONE) &
|
2015-06-29 10:35:29 +02:00
|
|
|
!$OMP PRIVATE(i,j,s2_tmp) SHARED(n,psi_coefs_tmp,psi_keys_tmp,N_int) REDUCTION(+:s2) SCHEDULE(dynamic)
|
2015-05-06 15:05:08 +02:00
|
|
|
do i=1,n
|
|
|
|
do j=i+1,n
|
|
|
|
call get_s2(psi_keys_tmp(1,1,i),psi_keys_tmp(1,1,j),s2_tmp,N_int)
|
|
|
|
s2 += (psi_coefs_tmp(i)+psi_coefs_tmp(i))*psi_coefs_tmp(j)*s2_tmp
|
2015-06-29 10:35:29 +02:00
|
|
|
! s2 = s2 + 2.d0 * psi_coefs_tmp(i)*psi_coefs_tmp(j)*s2_tmp
|
2015-05-06 15:05:08 +02:00
|
|
|
enddo
|
2015-04-20 16:45:06 +02:00
|
|
|
enddo
|
|
|
|
!$OMP END PARALLEL DO
|
2015-06-29 10:35:29 +02:00
|
|
|
do i=1,n
|
|
|
|
call get_s2(psi_keys_tmp(1,1,i),psi_keys_tmp(1,1,i),s2_tmp,N_int)
|
|
|
|
s2 += psi_coefs_tmp(i)*psi_coefs_tmp(i)*s2_tmp
|
|
|
|
enddo
|
2015-05-06 17:01:45 +02:00
|
|
|
s2 += S_z2_Sz
|
2015-06-29 10:35:29 +02:00
|
|
|
!print*,'s2 = ',s2
|
|
|
|
!print*,''
|
2015-04-20 16:45:06 +02:00
|
|
|
end
|
|
|
|
|