program printSQ_ij_T_kl implicit none BEGIN_DOC ! perturbative approach to build psi_postsvd END_DOC read_wf = .True. TOUCH read_wf PROVIDE N_int call run() end subroutine run USE OMP_LIB USE bitmasks implicit none integer(bit_kind) :: det1(N_int,2), det2(N_int,2) integer :: degree, i_state double precision :: h12 integer :: i, j, k, l double precision, allocatable :: T(:,:,:,:) double precision :: ti, tf integer :: nb_taches !$OMP PARALLEL nb_taches = OMP_GET_NUM_THREADS() !$OMP END PARALLEL call wall_time(ti) i_state = 1 det1(:,1) = psi_det_alpha_unique(:,1) det2(:,1) = psi_det_alpha_unique(:,1) det1(:,2) = psi_det_beta_unique(:,1) det2(:,2) = psi_det_beta_unique(:,1) call get_excitation_degree_spin(det1(1,1),det2(1,1),degree,N_int) call get_excitation_degree(det1,det2,degree,N_int) call i_H_j(det1, det2, N_int, h12) ! --------------------------------------------------------------------------------------- ! construct the initial CI matrix print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~' print *, ' CI matrix :', n_det_alpha_unique,'x',n_det_beta_unique print *, ' N det :', N_det print *, ' ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~ ~' ! --------------------------------------------------------------------------------------- ! --------------------------------------------------------------------------------------- allocate( T(n_det_alpha_unique,n_det_beta_unique,n_det_alpha_unique,n_det_beta_unique) ) call const_2b(T) open(UNIT=11, FILE="ij_T_kl.dat", ACTION="WRITE") do i = 1, n_det_alpha_unique do j = 1, n_det_beta_unique do k = 1, n_det_alpha_unique do l = 1, n_det_beta_unique write(11, '(E15.7)') T(i,j,k,l) !write(11, '(4(I5,2X), 5X, E15.7)') i, j, k, l, T(i,j,k,l) enddo enddo enddo enddo close(11) deallocate( T ) call wall_time(tf) print *, ' ___________________________________________________________________' print *, ' ' !print *, " Execution avec ", nb_taches, " threads" print *, " Execution avec 1 threads" print *, " total elapsed time (min) = ", (tf-ti)/60.d0 print *, ' ___________________________________________________________________' end !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !---! !---! !---! !---! !---! !---! !---! !---! !---! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !---! !---! !---! !---! !---! !---! !---! !---! !---! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !---! !---! !---! !---! !---! !---! !---! !---! !---! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !/-\! !---! !---! !---! !---! !---! !---! !---! !---! !---! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! ! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! !___! subroutine const_2b(T) USE bitmasks implicit none double precision, external :: get_two_e_integral double precision, intent(out) :: T(n_det_alpha_unique,n_det_beta_unique,n_det_alpha_unique,n_det_beta_unique ) integer :: na, nb integer :: i, k, j, l integer(bit_kind) :: psi_ij(N_int,2), psi_kl(N_int,2) double precision :: phase integer :: degree, h1, h2, p1, p2, s1, s2, e1, e2 integer :: ii, jj integer :: exc(0:2,2,2) integer :: occ(N_int*bit_kind_size,2), n_occ_alpha double precision :: two_body_fact na = n_det_alpha_unique nb = n_det_beta_unique T(:,:,:,:) = 0.d0 ! ----------------------------------------------------------------------------------------------------------------- do i = 1, na psi_ij(1:N_int,1) = psi_det_alpha_unique(1:N_int,i) do j = 1, nb psi_ij(1:N_int,2) = psi_det_beta_unique(1:N_int,j) call bitstring_to_list(psi_ij(1,1), occ(1,1), n_occ_alpha, N_int) call bitstring_to_list(psi_ij(1,2), occ(1,2), n_occ_alpha, N_int) do k = 1, na psi_kl(1:N_int,1) = psi_det_alpha_unique(1:N_int,k) do l = 1, nb psi_kl(1:N_int,2) = psi_det_beta_unique(1:N_int,l) call get_excitation_degree(psi_ij, psi_kl, degree, N_int) two_body_fact = 0.d0 if(degree .eq. 2) then call get_double_excitation(psi_ij, psi_kl, exc, phase, N_int) call decode_exc(exc, degree, h1, p1, h2, p2, s1, s2) select case(s1+s2) case(2,4) two_body_fact += phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact -= phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map) case(3) two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map) end select else if(degree .eq. 1) then call get_single_excitation(psi_ij, psi_kl, exc, phase, N_int) call decode_exc(exc, degree, h1, p1, h2, p2, s1, s2) select case(s1) case(1) do ii = 1, elec_alpha_num p2 = occ(ii,1) h2 = p2 two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact -= 0.5d0 * phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map) two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map) two_body_fact -= 0.5d0 * phase * get_two_e_integral(h2, h1, p1, p2, mo_integrals_map) enddo do ii = 1, elec_beta_num p2 = occ(ii,2) h2 = p2 two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map) enddo case(2) do ii = 1, elec_alpha_num p2 = occ(ii,1) h2 = p2 two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map) enddo do ii = 1, elec_beta_num p2 = occ(ii,2) h2 = p2 two_body_fact += 0.5d0 * phase * get_two_e_integral(h1, h2, p1, p2, mo_integrals_map) two_body_fact -= 0.5d0 * phase * get_two_e_integral(h1, h2, p2, p1, mo_integrals_map) two_body_fact += 0.5d0 * phase * get_two_e_integral(h2, h1, p2, p1, mo_integrals_map) two_body_fact -= 0.5d0 * phase * get_two_e_integral(h2, h1, p1, p2, mo_integrals_map) enddo end select else if(degree .eq. 0) then do ii = 1, elec_alpha_num e1 = occ(ii,1) do jj = 1, elec_alpha_num e2 = occ(jj,1) two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map) two_body_fact -= 0.5d0 * get_two_e_integral(e1, e2, e2, e1, mo_integrals_map) enddo do jj = 1, elec_beta_num e2 = occ(jj,2) two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map) two_body_fact += 0.5d0 * get_two_e_integral(e2, e1, e2, e1, mo_integrals_map) enddo enddo do ii = 1, elec_beta_num e1 = occ(ii,2) do jj = 1, elec_beta_num e2 = occ(jj,2) two_body_fact += 0.5d0 * get_two_e_integral(e1, e2, e1, e2, mo_integrals_map) two_body_fact -= 0.5d0 * get_two_e_integral(e1, e2, e2, e1, mo_integrals_map) enddo enddo end if T(i,j,k,l) = two_body_fact enddo enddo enddo enddo ! ----------------------------------------------------------------------------------------------------------------- return end subroutine const_2b