2015-04-09 23:59:06 +02:00
|
|
|
use omp_lib
|
2015-10-21 15:50:15 +02:00
|
|
|
use bitmasks
|
2015-04-09 23:59:06 +02:00
|
|
|
|
2015-07-13 18:00:38 +02:00
|
|
|
BEGIN_PROVIDER [ integer(omp_lock_kind), psi_ref_lock, (psi_det_size) ]
|
2015-04-09 23:59:06 +02:00
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
2015-07-13 18:00:38 +02:00
|
|
|
! Locks on ref determinants to fill delta_ij
|
2015-04-09 23:59:06 +02:00
|
|
|
END_DOC
|
|
|
|
integer :: i
|
|
|
|
do i=1,psi_det_size
|
2015-07-13 18:00:38 +02:00
|
|
|
call omp_init_lock( psi_ref_lock(i) )
|
2015-04-09 23:59:06 +02:00
|
|
|
enddo
|
|
|
|
|
|
|
|
END_PROVIDER
|
|
|
|
|
2015-10-16 13:06:19 +02:00
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
subroutine mrcc_dress(delta_ij_, delta_ii_, Nstates, Ndet_non_ref, Ndet_ref,i_generator,n_selected,det_buffer,Nint,iproc,key_mask)
|
2015-04-09 21:46:37 +02:00
|
|
|
use bitmasks
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer, intent(in) :: i_generator,n_selected, Nint, iproc
|
2016-03-29 23:18:26 +02:00
|
|
|
integer, intent(in) :: Nstates, Ndet_ref, Ndet_non_ref
|
|
|
|
double precision, intent(inout) :: delta_ij_(Nstates,Ndet_non_ref,Ndet_ref)
|
|
|
|
double precision, intent(inout) :: delta_ii_(Nstates,Ndet_ref)
|
2015-04-09 21:46:37 +02:00
|
|
|
|
|
|
|
integer(bit_kind), intent(in) :: det_buffer(Nint,2,n_selected)
|
2016-03-30 17:37:22 +02:00
|
|
|
integer :: i,j,k,l,m
|
2016-03-31 10:40:39 +02:00
|
|
|
integer :: degree_alpha(psi_det_size)
|
|
|
|
integer :: idx_alpha(0:psi_det_size)
|
2015-11-19 20:57:44 +01:00
|
|
|
logical :: good, fullMatch
|
2015-04-09 21:46:37 +02:00
|
|
|
|
|
|
|
integer(bit_kind) :: tq(Nint,2,n_selected)
|
|
|
|
integer :: N_tq, c_ref ,degree
|
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
double precision :: hIk, hla, hIl, dIk(Nstates), dka(Nstates), dIa(Nstates)
|
2015-04-09 23:59:06 +02:00
|
|
|
double precision, allocatable :: dIa_hla(:,:)
|
|
|
|
double precision :: haj, phase, phase2
|
2016-03-29 23:18:26 +02:00
|
|
|
double precision :: f(Nstates), ci_inv(Nstates)
|
2015-04-09 23:59:06 +02:00
|
|
|
integer :: exc(0:2,2,2)
|
|
|
|
integer :: h1,h2,p1,p2,s1,s2
|
|
|
|
integer(bit_kind) :: tmp_det(Nint,2)
|
|
|
|
integer :: iint, ipos
|
|
|
|
integer :: i_state, k_sd, l_sd, i_I, i_alpha
|
2015-10-16 13:06:19 +02:00
|
|
|
|
2015-11-19 20:57:44 +01:00
|
|
|
integer(bit_kind),allocatable :: miniList(:,:,:)
|
2015-10-29 11:39:26 +01:00
|
|
|
integer(bit_kind),intent(in) :: key_mask(Nint, 2)
|
|
|
|
integer,allocatable :: idx_miniList(:)
|
2015-11-19 20:57:44 +01:00
|
|
|
integer :: N_miniList, ni, leng
|
2016-03-29 23:18:26 +02:00
|
|
|
double precision, allocatable :: hij_cache(:)
|
2015-10-16 13:06:19 +02:00
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
integer(bit_kind), allocatable :: microlist(:,:,:), microlist_zero(:,:,:)
|
|
|
|
integer, allocatable :: idx_microlist(:), N_microlist(:), ptr_microlist(:), idx_microlist_zero(:)
|
|
|
|
integer :: mobiles(2), smallerlist
|
2016-05-20 09:44:22 +02:00
|
|
|
logical, external :: is_generable
|
2016-03-30 17:37:22 +02:00
|
|
|
|
2015-11-05 15:05:19 +01:00
|
|
|
leng = max(N_det_generators, N_det_non_ref)
|
2016-03-30 17:37:22 +02:00
|
|
|
allocate(miniList(Nint, 2, leng), idx_minilist(leng), hij_cache(N_det_non_ref))
|
2015-10-29 11:39:26 +01:00
|
|
|
|
2015-11-19 20:57:44 +01:00
|
|
|
!create_minilist_find_previous(key_mask, fullList, miniList, N_fullList, N_miniList, fullMatch, Nint)
|
|
|
|
call create_minilist_find_previous(key_mask, psi_det_generators, miniList, i_generator-1, N_miniList, fullMatch, Nint)
|
2015-11-05 15:05:19 +01:00
|
|
|
|
2015-11-19 20:57:44 +01:00
|
|
|
if(fullMatch) then
|
|
|
|
return
|
2015-10-29 11:39:26 +01:00
|
|
|
end if
|
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
allocate(ptr_microlist(0:mo_tot_num*2+1), &
|
|
|
|
N_microlist(0:mo_tot_num*2) )
|
|
|
|
allocate( microlist(Nint,2,N_minilist*4), &
|
|
|
|
idx_microlist(N_minilist*4))
|
|
|
|
|
2016-05-20 09:44:22 +02:00
|
|
|
if(key_mask(1,1) /= 0_8) then
|
2016-03-30 17:37:22 +02:00
|
|
|
call create_microlist(miniList, N_minilist, key_mask, microlist, idx_microlist, N_microlist, ptr_microlist, Nint)
|
|
|
|
call find_triples_and_quadruples_micro(i_generator,n_selected,det_buffer,Nint,tq,N_tq,microlist,ptr_microlist,N_microlist,key_mask)
|
|
|
|
else
|
|
|
|
call find_triples_and_quadruples(i_generator,n_selected,det_buffer,Nint,tq,N_tq,miniList,N_minilist)
|
|
|
|
end if
|
|
|
|
|
|
|
|
|
2015-10-23 15:58:43 +02:00
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
deallocate(microlist, idx_microlist)
|
2016-03-29 23:18:26 +02:00
|
|
|
|
|
|
|
allocate (dIa_hla(Nstates,Ndet_non_ref))
|
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
! |I>
|
2016-03-29 23:18:26 +02:00
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
! |alpha>
|
2015-10-16 13:06:19 +02:00
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
if(N_tq > 0) then
|
2016-05-20 09:44:22 +02:00
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
call create_minilist(key_mask, psi_non_ref, miniList, idx_minilist, N_det_non_ref, N_minilist, Nint)
|
|
|
|
if(N_minilist == 0) return
|
|
|
|
|
|
|
|
|
|
|
|
if(key_mask(1,1) /= 0) then !!!!!!!!!!! PAS GENERAL !!!!!!!!!
|
|
|
|
allocate(microlist_zero(Nint,2,N_minilist), idx_microlist_zero(N_minilist))
|
|
|
|
|
|
|
|
allocate( microlist(Nint,2,N_minilist*4), &
|
|
|
|
idx_microlist(N_minilist*4))
|
|
|
|
call create_microlist(miniList, N_minilist, key_mask, microlist, idx_microlist, N_microlist, ptr_microlist, Nint)
|
|
|
|
|
|
|
|
|
|
|
|
do i=0,mo_tot_num*2
|
|
|
|
do k=ptr_microlist(i),ptr_microlist(i+1)-1
|
|
|
|
idx_microlist(k) = idx_minilist(idx_microlist(k))
|
|
|
|
end do
|
|
|
|
end do
|
|
|
|
|
|
|
|
do l=1,N_microlist(0)
|
|
|
|
do k=1,Nint
|
|
|
|
microlist_zero(k,1,l) = microlist(k,1,l)
|
|
|
|
microlist_zero(k,2,l) = microlist(k,2,l)
|
|
|
|
enddo
|
|
|
|
idx_microlist_zero(l) = idx_microlist(l)
|
|
|
|
enddo
|
|
|
|
end if
|
2016-03-29 23:18:26 +02:00
|
|
|
end if
|
2015-10-16 13:06:19 +02:00
|
|
|
|
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
|
2015-04-09 22:32:25 +02:00
|
|
|
do i_alpha=1,N_tq
|
2016-05-20 09:44:22 +02:00
|
|
|
! ok = .false.
|
|
|
|
! do i=N_det_generators, 1, -1
|
|
|
|
! if(is_generable(psi_det_generators(1,1,i), tq(1,1,i_alpha), Nint)) then
|
|
|
|
! ok = .true.
|
|
|
|
! exit
|
|
|
|
! end if
|
|
|
|
! end do
|
|
|
|
! if(.not. ok) then
|
|
|
|
! cycle
|
|
|
|
! end if
|
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
if(key_mask(1,1) /= 0) then
|
|
|
|
call getMobiles(tq(1,1,i_alpha), key_mask, mobiles, Nint)
|
|
|
|
|
|
|
|
if(N_microlist(mobiles(1)) < N_microlist(mobiles(2))) then
|
|
|
|
smallerlist = mobiles(1)
|
|
|
|
else
|
|
|
|
smallerlist = mobiles(2)
|
|
|
|
end if
|
|
|
|
|
2015-10-16 13:06:19 +02:00
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
do l=0,N_microlist(smallerlist)-1
|
|
|
|
microlist_zero(:,:,ptr_microlist(1) + l) = microlist(:,:,ptr_microlist(smallerlist) + l)
|
|
|
|
idx_microlist_zero(ptr_microlist(1) + l) = idx_microlist(ptr_microlist(smallerlist) + l)
|
|
|
|
end do
|
|
|
|
|
2016-03-31 10:40:39 +02:00
|
|
|
call get_excitation_degree_vector(microlist_zero,tq(1,1,i_alpha),degree_alpha,Nint,N_microlist(smallerlist)+N_microlist(0),idx_alpha)
|
|
|
|
do j=1,idx_alpha(0)
|
|
|
|
idx_alpha(j) = idx_microlist_zero(idx_alpha(j))
|
2016-03-30 17:37:22 +02:00
|
|
|
end do
|
|
|
|
else
|
|
|
|
call get_excitation_degree_vector(miniList,tq(1,1,i_alpha),degree_alpha,Nint,N_minilist,idx_alpha)
|
|
|
|
do j=1,idx_alpha(0)
|
|
|
|
idx_alpha(j) = idx_miniList(idx_alpha(j))
|
|
|
|
end do
|
|
|
|
end if
|
|
|
|
|
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
do l_sd=1,idx_alpha(0)
|
|
|
|
k_sd = idx_alpha(l_sd)
|
|
|
|
call i_h_j(tq(1,1,i_alpha),psi_non_ref(1,1,idx_alpha(l_sd)),Nint,hij_cache(k_sd))
|
|
|
|
enddo
|
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
! |I>
|
2015-07-13 18:00:38 +02:00
|
|
|
do i_I=1,N_det_ref
|
2016-03-29 23:18:26 +02:00
|
|
|
! Find triples and quadruple grand parents
|
|
|
|
call get_excitation_degree(tq(1,1,i_alpha),psi_ref(1,1,i_I),degree,Nint)
|
|
|
|
if (degree > 4) then
|
|
|
|
cycle
|
|
|
|
endif
|
|
|
|
|
|
|
|
do i_state=1,Nstates
|
|
|
|
dIa(i_state) = 0.d0
|
|
|
|
enddo
|
|
|
|
|
|
|
|
! <I| <> |alpha>
|
|
|
|
do k_sd=1,idx_alpha(0)
|
|
|
|
|
|
|
|
! Loop if lambda == 0
|
|
|
|
logical :: loop
|
|
|
|
loop = .True.
|
|
|
|
do i_state=1,Nstates
|
|
|
|
if (lambda_mrcc(i_state,idx_alpha(k_sd)) /= 0.d0) then
|
|
|
|
loop = .False.
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
if (loop) then
|
|
|
|
cycle
|
|
|
|
endif
|
|
|
|
|
|
|
|
call get_excitation_degree(psi_ref(1,1,i_I),psi_non_ref(1,1,idx_alpha(k_sd)),degree,Nint)
|
|
|
|
if (degree > 2) then
|
|
|
|
cycle
|
|
|
|
endif
|
|
|
|
|
|
|
|
! <I| /k\ |alpha>
|
|
|
|
! <I|H|k>
|
|
|
|
hIk = hij_mrcc(idx_alpha(k_sd),i_I)
|
|
|
|
! call i_h_j(psi_ref(1,1,i_I),psi_non_ref(1,1,idx_alpha(k_sd)),Nint,hIk)
|
|
|
|
do i_state=1,Nstates
|
|
|
|
dIk(i_state) = hIk * lambda_mrcc(i_state,idx_alpha(k_sd))
|
|
|
|
enddo
|
|
|
|
! |l> = Exc(k -> alpha) |I>
|
|
|
|
call get_excitation(psi_non_ref(1,1,idx_alpha(k_sd)),tq(1,1,i_alpha),exc,degree,phase,Nint)
|
|
|
|
call decode_exc(exc,degree,h1,p1,h2,p2,s1,s2)
|
|
|
|
do k=1,N_int
|
|
|
|
tmp_det(k,1) = psi_ref(k,1,i_I)
|
|
|
|
tmp_det(k,2) = psi_ref(k,2,i_I)
|
|
|
|
enddo
|
|
|
|
|
2016-03-31 10:40:39 +02:00
|
|
|
logical :: ok
|
|
|
|
call apply_excitation(psi_ref(1,1,i_I), exc, tmp_det, ok, Nint)
|
|
|
|
if(.not. ok) cycle
|
2016-03-29 23:18:26 +02:00
|
|
|
|
|
|
|
! <I| \l/ |alpha>
|
|
|
|
do i_state=1,Nstates
|
|
|
|
dka(i_state) = 0.d0
|
|
|
|
enddo
|
|
|
|
do l_sd=k_sd+1,idx_alpha(0)
|
2016-03-31 10:40:39 +02:00
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
call get_excitation_degree(tmp_det,psi_non_ref(1,1,idx_alpha(l_sd)),degree,Nint)
|
|
|
|
if (degree == 0) then
|
|
|
|
|
|
|
|
loop = .True.
|
|
|
|
do i_state=1,Nstates
|
|
|
|
if (lambda_mrcc(i_state,idx_alpha(l_sd)) /= 0.d0) then
|
|
|
|
loop = .False.
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
if (.not.loop) then
|
|
|
|
call get_excitation(psi_ref(1,1,i_I),psi_non_ref(1,1,idx_alpha(l_sd)),exc,degree,phase2,Nint)
|
|
|
|
hIl = hij_mrcc(idx_alpha(l_sd),i_I)
|
2016-03-31 10:40:39 +02:00
|
|
|
! call i_h_j(psi_ref(1,1,i_I),psi_non_ref(1,1,idx_alpha(l_sd)),Nint,hIl)
|
2016-03-29 23:18:26 +02:00
|
|
|
do i_state=1,Nstates
|
|
|
|
dka(i_state) = hIl * lambda_mrcc(i_state,idx_alpha(l_sd)) * phase * phase2
|
|
|
|
enddo
|
|
|
|
endif
|
2016-03-31 10:40:39 +02:00
|
|
|
|
2016-03-29 23:18:26 +02:00
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
do i_state=1,Nstates
|
|
|
|
dIa(i_state) = dIa(i_state) + dIk(i_state) * dka(i_state)
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
|
|
|
|
do i_state=1,Nstates
|
2016-05-02 17:06:18 +02:00
|
|
|
ci_inv(i_state) = psi_ref_coef_inv(i_I,i_state)
|
2016-03-29 23:18:26 +02:00
|
|
|
enddo
|
|
|
|
do l_sd=1,idx_alpha(0)
|
|
|
|
k_sd = idx_alpha(l_sd)
|
|
|
|
hla = hij_cache(k_sd)
|
|
|
|
! call i_h_j(tq(1,1,i_alpha),psi_non_ref(1,1,idx_alpha(l_sd)),Nint,hla)
|
|
|
|
do i_state=1,Nstates
|
|
|
|
dIa_hla(i_state,k_sd) = dIa(i_state) * hla
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
call omp_set_lock( psi_ref_lock(i_I) )
|
2016-05-20 09:44:22 +02:00
|
|
|
|
|
|
|
|
2016-05-02 17:06:18 +02:00
|
|
|
do i_state=1,Nstates
|
|
|
|
if(dabs(psi_ref_coef(i_I,i_state)).ge.5.d-5)then
|
|
|
|
do l_sd=1,idx_alpha(0)
|
|
|
|
k_sd = idx_alpha(l_sd)
|
2016-05-20 09:44:22 +02:00
|
|
|
delta_ij_(i_state,k_sd,i_I) = delta_ij_(i_state,k_sd,i_I) + dIa_hla(i_state,k_sd)
|
|
|
|
delta_ii_(i_state,i_I) = delta_ii_(i_state,i_I) - dIa_hla(i_state,k_sd) * ci_inv(i_state) * psi_non_ref_coef_transp(i_state,k_sd)
|
2016-05-02 17:06:18 +02:00
|
|
|
enddo
|
|
|
|
else
|
2016-05-20 09:44:22 +02:00
|
|
|
!delta_ii_(i_state,i_I) = 0.d0
|
2016-05-02 17:06:18 +02:00
|
|
|
do l_sd=1,idx_alpha(0)
|
|
|
|
k_sd = idx_alpha(l_sd)
|
2016-10-18 23:07:03 +02:00
|
|
|
delta_ij_(i_state,k_sd,i_I) = delta_ij_(i_state,k_sd,i_I) + 0.5d0 * dIa_hla(i_state,k_sd)
|
2016-05-02 17:06:18 +02:00
|
|
|
enddo
|
|
|
|
endif
|
2016-03-29 23:18:26 +02:00
|
|
|
enddo
|
|
|
|
call omp_unset_lock( psi_ref_lock(i_I) )
|
2015-04-09 21:46:37 +02:00
|
|
|
enddo
|
|
|
|
enddo
|
2016-05-20 09:44:22 +02:00
|
|
|
deallocate (dIa_hla,hij_cache)
|
|
|
|
deallocate(miniList, idx_miniList)
|
2015-04-09 21:46:37 +02:00
|
|
|
end
|
|
|
|
|
|
|
|
|
2015-10-23 15:58:43 +02:00
|
|
|
subroutine find_triples_and_quadruples(i_generator,n_selected,det_buffer,Nint,tq,N_tq,miniList,N_miniList)
|
2015-10-26 08:31:27 +01:00
|
|
|
|
2015-04-03 14:26:14 +02:00
|
|
|
use bitmasks
|
2015-04-02 10:13:33 +02:00
|
|
|
implicit none
|
2015-04-03 14:26:14 +02:00
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
integer, intent(in) :: i_generator,n_selected, Nint
|
2015-04-03 14:26:14 +02:00
|
|
|
|
|
|
|
integer(bit_kind), intent(in) :: det_buffer(Nint,2,n_selected)
|
|
|
|
integer :: i,j,k,m
|
|
|
|
logical :: is_in_wavefunction
|
|
|
|
integer :: degree(psi_det_size)
|
|
|
|
integer :: idx(0:psi_det_size)
|
|
|
|
logical :: good
|
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
integer(bit_kind), intent(out) :: tq(Nint,2,n_selected)
|
|
|
|
integer, intent(out) :: N_tq
|
2015-10-21 15:50:15 +02:00
|
|
|
|
2015-10-23 15:58:43 +02:00
|
|
|
|
2015-10-29 11:39:26 +01:00
|
|
|
integer :: nt,ni
|
2015-11-19 21:20:43 +01:00
|
|
|
logical, external :: is_connected_to
|
2015-10-23 15:58:43 +02:00
|
|
|
|
|
|
|
|
2015-10-29 11:39:26 +01:00
|
|
|
integer(bit_kind),intent(in) :: miniList(Nint,2,N_det_generators)
|
2015-10-23 15:58:43 +02:00
|
|
|
integer,intent(in) :: N_miniList
|
2015-11-05 15:05:19 +01:00
|
|
|
|
2015-10-23 15:58:43 +02:00
|
|
|
|
|
|
|
|
2015-04-03 14:26:14 +02:00
|
|
|
N_tq = 0
|
2015-11-05 15:05:19 +01:00
|
|
|
|
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
|
2015-10-29 11:39:26 +01:00
|
|
|
i_loop : do i=1,N_selected
|
2015-12-02 11:59:01 +01:00
|
|
|
if(is_connected_to(det_buffer(1,1,i), miniList, Nint, N_miniList)) then
|
2015-11-19 21:20:43 +01:00
|
|
|
cycle
|
|
|
|
end if
|
2015-10-23 15:58:43 +02:00
|
|
|
|
2015-04-03 14:26:14 +02:00
|
|
|
! Select determinants that are triple or quadruple excitations
|
2015-07-13 18:00:38 +02:00
|
|
|
! from the ref
|
2015-04-03 14:26:14 +02:00
|
|
|
good = .True.
|
2015-10-16 13:06:19 +02:00
|
|
|
call get_excitation_degree_vector(psi_ref,det_buffer(1,1,i),degree,Nint,N_det_ref,idx)
|
|
|
|
!good=(idx(0) == 0) tant que degree > 2 pas retourné par get_excitation_degree_vector
|
2015-04-03 14:26:14 +02:00
|
|
|
do k=1,idx(0)
|
|
|
|
if (degree(k) < 3) then
|
|
|
|
good = .False.
|
2015-04-02 10:13:33 +02:00
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
2015-04-03 14:26:14 +02:00
|
|
|
if (good) then
|
2016-05-20 09:44:22 +02:00
|
|
|
if (.not. is_in_wavefunction(det_buffer(1,1,i),Nint)) then
|
2015-04-03 14:26:14 +02:00
|
|
|
N_tq += 1
|
|
|
|
do k=1,N_int
|
|
|
|
tq(k,1,N_tq) = det_buffer(k,1,i)
|
|
|
|
tq(k,2,N_tq) = det_buffer(k,2,i)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
endif
|
2015-10-29 11:39:26 +01:00
|
|
|
enddo i_loop
|
2015-04-03 14:26:14 +02:00
|
|
|
end
|
2015-04-02 10:13:33 +02:00
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
|
2016-03-30 17:37:22 +02:00
|
|
|
subroutine find_triples_and_quadruples_micro(i_generator,n_selected,det_buffer,Nint,tq,N_tq,microlist,ptr_microlist,N_microlist,key_mask)
|
|
|
|
|
|
|
|
use bitmasks
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
integer, intent(in) :: i_generator,n_selected, Nint
|
|
|
|
|
|
|
|
integer(bit_kind), intent(in) :: det_buffer(Nint,2,n_selected)
|
|
|
|
integer :: i,j,k,m
|
|
|
|
logical :: is_in_wavefunction
|
|
|
|
integer :: degree(psi_det_size)
|
|
|
|
integer :: idx(0:psi_det_size)
|
|
|
|
logical :: good
|
|
|
|
|
|
|
|
integer(bit_kind), intent(out) :: tq(Nint,2,n_selected)
|
|
|
|
integer, intent(out) :: N_tq
|
|
|
|
|
|
|
|
|
|
|
|
integer :: nt,ni
|
|
|
|
logical, external :: is_connected_to
|
|
|
|
|
|
|
|
|
|
|
|
integer(bit_kind),intent(in) :: microlist(Nint,2,*)
|
|
|
|
integer,intent(in) :: ptr_microlist(0:*)
|
|
|
|
integer,intent(in) :: N_microlist(0:*)
|
|
|
|
integer(bit_kind),intent(in) :: key_mask(Nint, 2)
|
|
|
|
|
|
|
|
integer :: mobiles(2), smallerlist
|
|
|
|
|
|
|
|
N_tq = 0
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
i_loop : do i=1,N_selected
|
|
|
|
call getMobiles(det_buffer(1,1,i), key_mask, mobiles, Nint)
|
|
|
|
if(N_microlist(mobiles(1)) < N_microlist(mobiles(2))) then
|
|
|
|
smallerlist = mobiles(1)
|
|
|
|
else
|
|
|
|
smallerlist = mobiles(2)
|
|
|
|
end if
|
|
|
|
|
|
|
|
if(N_microlist(smallerlist) > 0) then
|
|
|
|
if(is_connected_to(det_buffer(1,1,i), microlist(1,1,ptr_microlist(smallerlist)), Nint, N_microlist(smallerlist))) then
|
|
|
|
cycle
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
if(N_microlist(0) > 0) then
|
|
|
|
if(is_connected_to(det_buffer(1,1,i), microlist, Nint, N_microlist(0))) then
|
|
|
|
cycle
|
|
|
|
end if
|
|
|
|
end if
|
|
|
|
|
|
|
|
! Select determinants that are triple or quadruple excitations
|
|
|
|
! from the ref
|
|
|
|
good = .True.
|
|
|
|
call get_excitation_degree_vector(psi_ref,det_buffer(1,1,i),degree,Nint,N_det_ref,idx)
|
|
|
|
!good=(idx(0) == 0) tant que degree > 2 pas retourné par get_excitation_degree_vector
|
|
|
|
do k=1,idx(0)
|
|
|
|
if (degree(k) < 3) then
|
|
|
|
good = .False.
|
|
|
|
exit
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
if (good) then
|
2016-05-20 09:44:22 +02:00
|
|
|
if (.not. is_in_wavefunction(det_buffer(1,1,i),Nint)) then
|
2016-03-30 17:37:22 +02:00
|
|
|
N_tq += 1
|
|
|
|
do k=1,N_int
|
|
|
|
tq(k,1,N_tq) = det_buffer(k,1,i)
|
|
|
|
tq(k,2,N_tq) = det_buffer(k,2,i)
|
|
|
|
enddo
|
|
|
|
endif
|
|
|
|
endif
|
|
|
|
enddo i_loop
|
|
|
|
end
|
|
|
|
|
2015-04-09 21:46:37 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|