10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-07-03 09:55:59 +02:00

Introduced ib_jb pairs in H_apply

This commit is contained in:
Anthony Scemama 2014-06-06 14:28:47 +02:00
parent 0f8f18497e
commit 7e0b254c48
3 changed files with 78 additions and 56 deletions

View File

@ -41,7 +41,7 @@ class H_apply(object):
s["omp_parallel"] = """!$OMP PARALLEL DEFAULT(SHARED) & s["omp_parallel"] = """!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(i,j,k,l,keys_out,hole,particle, & !$OMP PRIVATE(i,j,k,l,keys_out,hole,particle, &
!$OMP occ_particle,occ_hole,j_a,k_a,other_spin, & !$OMP occ_particle,occ_hole,j_a,k_a,other_spin, &
!$OMP hole_save,ispin,jj,l_a, & !$OMP hole_save,ispin,jj,l_a,ib_jb_pairs, &
!$OMP accu,i_a,hole_tmp,particle_tmp,occ_particle_tmp, & !$OMP accu,i_a,hole_tmp,particle_tmp,occ_particle_tmp, &
!$OMP occ_hole_tmp,key_idx,i_b,j_b,key,N_elec_in_key_part_1,& !$OMP occ_hole_tmp,key_idx,i_b,j_b,key,N_elec_in_key_part_1,&
!$OMP N_elec_in_key_hole_1,N_elec_in_key_part_2, & !$OMP N_elec_in_key_hole_1,N_elec_in_key_part_2, &

View File

@ -26,6 +26,7 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2, i_gene
double precision :: mo_bielec_integral double precision :: mo_bielec_integral
integer, allocatable :: ia_ja_pairs(:,:,:) integer, allocatable :: ia_ja_pairs(:,:,:)
integer, allocatable :: ib_jb_pairs(:,:)
double precision :: diag_H_mat_elem double precision :: diag_H_mat_elem
integer :: iproc integer :: iproc
integer(omp_lock_kind), save :: lck, ifirst=0 integer(omp_lock_kind), save :: lck, ifirst=0
@ -60,7 +61,8 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2, i_gene
call bitstring_to_list(particle(1,2),occ_particle(1,2),N_elec_in_key_part_1(2),N_int) call bitstring_to_list(particle(1,2),occ_particle(1,2),N_elec_in_key_part_1(2),N_int)
call bitstring_to_list(hole(1,1),occ_hole(1,1),N_elec_in_key_hole_1(1),N_int) call bitstring_to_list(hole(1,1),occ_hole(1,1),N_elec_in_key_hole_1(1),N_int)
call bitstring_to_list(hole(1,2),occ_hole(1,2),N_elec_in_key_hole_1(2),N_int) call bitstring_to_list(hole(1,2),occ_hole(1,2),N_elec_in_key_hole_1(2),N_int)
allocate (ia_ja_pairs(2,0:(elec_alpha_num)*mo_tot_num,2)) allocate (ia_ja_pairs(2,0:(elec_alpha_num)*mo_tot_num,2), &
ib_jb_pairs(2,0:(elec_alpha_num)*mo_tot_num))
do ispin=1,2 do ispin=1,2
i=0 i=0
@ -129,18 +131,29 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2, i_gene
if (ispin == 1) then if (ispin == 1) then
integer :: jjj integer :: jjj
i=0
do kk = 1,N_elec_in_key_hole_2(other_spin) do kk = 1,N_elec_in_key_hole_2(other_spin)
hole = hole_save
i_b = occ_hole_tmp(kk,other_spin) i_b = occ_hole_tmp(kk,other_spin)
ASSERT (i_b > 0) ASSERT (i_b > 0)
ASSERT (i_b <= mo_tot_num) ASSERT (i_b <= mo_tot_num)
k = ishft(i_b-1,-bit_kind_shift)+1
j = i_b-ishft(k-1,bit_kind_shift)-1
hole(k,other_spin) = ibclr(hole(k,other_spin),j)
do jjj=1,N_elec_in_key_part_2(other_spin) ! particule do jjj=1,N_elec_in_key_part_2(other_spin) ! particule
j_b = occ_particle_tmp(jjj,other_spin) j_b = occ_particle_tmp(jjj,other_spin)
ASSERT (j_b > 0) ASSERT (j_b > 0)
ASSERT (j_b <= mo_tot_num) ASSERT (j_b <= mo_tot_num)
i+= 1
ib_jb_pairs(1,i) = i_b
ib_jb_pairs(2,i) = j_b
enddo
enddo
ib_jb_pairs(1,0) = i
do kk = 1,ib_jb_pairs(1,0)
hole = hole_save
i_b = ib_jb_pairs(1,kk)
j_b = ib_jb_pairs(2,kk)
k = ishft(i_b-1,-bit_kind_shift)+1
j = i_b-ishft(k-1,bit_kind_shift)-1
hole(k,other_spin) = ibclr(hole(k,other_spin),j)
key = hole key = hole
k = ishft(j_b-1,-bit_kind_shift)+1 k = ishft(j_b-1,-bit_kind_shift)+1
l = j_b-ishft(k-1,bit_kind_shift)-1 l = j_b-ishft(k-1,bit_kind_shift)-1
@ -155,34 +168,42 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2, i_gene
$keys_work $keys_work
key_idx = 0 key_idx = 0
endif endif
! endif
enddo
if (abort_here) then if (abort_here) then
exit exit
endif endif
enddo enddo
endif endif
! does all the mono excitations of the same spin ! does all the mono excitations of the same spin
i=0
do kk = 1,N_elec_in_key_hole_2(ispin) do kk = 1,N_elec_in_key_hole_2(ispin)
i_b = occ_hole_tmp(kk,ispin) i_b = occ_hole_tmp(kk,ispin)
if (i_b <= i_a.or.i_b == j_a) cycle
ASSERT (i_b > 0) ASSERT (i_b > 0)
ASSERT (i_b <= mo_tot_num) ASSERT (i_b <= mo_tot_num)
if (i_b <= i_a.or.i_b == j_a) cycle do jjj=1,N_elec_in_key_part_2(ispin) ! particule
hole = hole_save
k = ishft(i_b-1,-bit_kind_shift)+1
j = i_b-ishft(k-1,bit_kind_shift)-1
hole(k,ispin) = ibclr(hole(k,ispin),j)
do jjj=1,N_elec_in_key_part_2(ispin)
j_b = occ_particle_tmp(jjj,ispin) j_b = occ_particle_tmp(jjj,ispin)
ASSERT (j_b > 0) ASSERT (j_b > 0)
ASSERT (j_b <= mo_tot_num) ASSERT (j_b <= mo_tot_num)
if (j_b <= j_a) cycle if (j_b <= j_a) cycle
i+= 1
ib_jb_pairs(1,i) = i_b
ib_jb_pairs(2,i) = j_b
enddo
enddo
ib_jb_pairs(1,0) = i
do kk = 1,ib_jb_pairs(1,0)
hole = hole_save
i_b = ib_jb_pairs(1,kk)
j_b = ib_jb_pairs(2,kk)
k = ishft(i_b-1,-bit_kind_shift)+1
j = i_b-ishft(k-1,bit_kind_shift)-1
hole(k,ispin) = ibclr(hole(k,ispin),j)
key = hole key = hole
k = ishft(j_b-1,-bit_kind_shift)+1 k = ishft(j_b-1,-bit_kind_shift)+1
l = j_b-ishft(k-1,bit_kind_shift)-1 l = j_b-ishft(k-1,bit_kind_shift)-1
key(k,ispin) = ibset(key(k,ispin),l) key(k,ispin) = ibset(key(k,ispin),l)
!! a^((+)_j_b(ispin) a_i_b(ispin) : mono exc :: orb(i_b,ispin) --> orb(j_b,ispin)
key_idx += 1 key_idx += 1
do k=1,N_int do k=1,N_int
keys_out(k,1,key_idx) = key(k,1) keys_out(k,1,key_idx) = key(k,1)
@ -196,14 +217,14 @@ subroutine $subroutine_diexc(key_in, hole_1,particl_1, hole_2, particl_2, i_gene
if (abort_here) then if (abort_here) then
exit exit
endif endif
enddo enddo ! kk
enddo! kk
enddo ! ii enddo ! ii
$omp_enddo $omp_enddo
enddo ! ispin enddo ! ispin
$keys_work $keys_work
$deinit_thread $deinit_thread
deallocate (ia_ja_pairs, & deallocate (ia_ja_pairs, ib_jb_pairs, &
keys_out, hole_save, & keys_out, hole_save, &
key,hole, particle, hole_tmp,& key,hole, particle, hole_tmp,&
particle_tmp, occ_particle, & particle_tmp, occ_particle, &
@ -235,6 +256,7 @@ subroutine $subroutine_monoexc(key_in, hole_1,particl_1,i_generator $parameters
integer :: ii,i,jj,j,k,ispin,l integer :: ii,i,jj,j,k,ispin,l
integer,allocatable :: occ_particle(:,:), occ_hole(:,:) integer,allocatable :: occ_particle(:,:), occ_hole(:,:)
integer,allocatable :: occ_particle_tmp(:,:), occ_hole_tmp(:,:) integer,allocatable :: occ_particle_tmp(:,:), occ_hole_tmp(:,:)
integer,allocatable :: ib_jb_pairs(:,:)
integer :: kk,pp,other_spin,key_idx integer :: kk,pp,other_spin,key_idx
integer :: N_elec_in_key_hole_1(2),N_elec_in_key_part_1(2) integer :: N_elec_in_key_hole_1(2),N_elec_in_key_part_1(2)
integer :: N_elec_in_key_hole_2(2),N_elec_in_key_part_2(2) integer :: N_elec_in_key_hole_2(2),N_elec_in_key_part_2(2)

View File

@ -83,7 +83,7 @@ Documentation
.br .br
Initial guess vectors are not necessarily orthonormal Initial guess vectors are not necessarily orthonormal
`repeat_excitation <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/SC2.irp.f#L215>`_ `repeat_excitation <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/SC2.irp.f#L220>`_
Undocumented Undocumented
`connected_to_ref <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/connected_to_ref.irp.f#L95>`_ `connected_to_ref <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/connected_to_ref.irp.f#L95>`_
@ -152,7 +152,7 @@ Documentation
`davidson_threshold <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/davidson.irp.f#L374>`_ `davidson_threshold <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/davidson.irp.f#L374>`_
Can be : [ energy | residual | both | wall_time | cpu_time | iterations ] Can be : [ energy | residual | both | wall_time | cpu_time | iterations ]
`det_search_key <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L180>`_ `det_search_key <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L195>`_
Return an integer*8 corresponding to a determinant index for searching Return an integer*8 corresponding to a determinant index for searching
`n_det <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L20>`_ `n_det <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L20>`_
@ -164,35 +164,35 @@ Documentation
`n_states <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L3>`_ `n_states <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L3>`_
Number of states to consider Number of states to consider
`psi_average_norm_contrib <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L93>`_ `psi_average_norm_contrib <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L108>`_
Contribution of determinants to the state-averaged density Contribution of determinants to the state-averaged density
`psi_average_norm_contrib_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L114>`_ `psi_average_norm_contrib_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L129>`_
Wave function sorted by determinants contribution to the norm (state-averaged) Wave function sorted by determinants contribution to the norm (state-averaged)
`psi_coef <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L65>`_ `psi_coef <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L86>`_
The wave function. Initialized with Hartree-Fock if the EZFIO file The wave function coefficients. Initialized with Hartree-Fock if the EZFIO file
is empty is empty
`psi_coef_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L113>`_ `psi_coef_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L128>`_
Wave function sorted by determinants contribution to the norm (state-averaged) Wave function sorted by determinants contribution to the norm (state-averaged)
`psi_coef_sorted_bit <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L144>`_ `psi_coef_sorted_bit <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L159>`_
Determinants on which we apply <i|H|psi> for perturbation. Determinants on which we apply <i|H|psi> for perturbation.
o They are sorted by determinants interpreted as integers. Useful o They are sorted by determinants interpreted as integers. Useful
to accelerate the search of a determinant to accelerate the search of a determinant
`psi_det <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L64>`_ `psi_det <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L64>`_
The wave function. Initialized with Hartree-Fock if the EZFIO file The wave function determinants. Initialized with Hartree-Fock if the EZFIO file
is empty is empty
`psi_det_size <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L56>`_ `psi_det_size <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L56>`_
Size of the psi_det/psi_coef arrays Size of the psi_det/psi_coef arrays
`psi_det_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L112>`_ `psi_det_sorted <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L127>`_
Wave function sorted by determinants contribution to the norm (state-averaged) Wave function sorted by determinants contribution to the norm (state-averaged)
`psi_det_sorted_bit <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L143>`_ `psi_det_sorted_bit <http://github.com/LCPQ/quantum_package/tree/master/src/Dets/determinants.irp.f#L158>`_
Determinants on which we apply <i|H|psi> for perturbation. Determinants on which we apply <i|H|psi> for perturbation.
o They are sorted by determinants interpreted as integers. Useful o They are sorted by determinants interpreted as integers. Useful
to accelerate the search of a determinant to accelerate the search of a determinant