9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-04 19:03:28 +01:00

openmp for multi state 2rdm work

This commit is contained in:
Emmanuel Giner 2020-03-20 22:20:12 +01:00
parent d7b2714521
commit d96108f772
2 changed files with 46 additions and 33 deletions

View File

@ -14,8 +14,12 @@
! condition for alpha/beta spin
ispin = 1
all_states_act_two_rdm_alpha_alpha_mo = 0.D0
double precision :: wall_1,wall_2
call wall_time(wall_1)
print*,'providing all_states_act_two_rdm_alpha_alpha_mo ...'
call orb_range_all_states_two_rdm(all_states_act_two_rdm_alpha_alpha_mo,n_act_orb,n_act_orb,list_act,list_act_reverse,ispin,psi_coef,size(psi_coef,2),size(psi_coef,1))
call wall_time(wall_2)
print*,'time to provide all_states_act_two_rdm_alpha_alpha_mo',wall_2 - wall_1
END_PROVIDER
@ -32,7 +36,12 @@
! condition for alpha/beta spin
ispin = 2
all_states_act_two_rdm_beta_beta_mo = 0.d0
double precision :: wall_1,wall_2
call wall_time(wall_1)
print*,'providing all_states_act_two_rdm_beta_beta_mo ...'
call orb_range_all_states_two_rdm(all_states_act_two_rdm_beta_beta_mo,n_act_orb,n_act_orb,list_act,list_act_reverse,ispin,psi_coef,size(psi_coef,2),size(psi_coef,1))
call wall_time(wall_2)
print*,'time to provide all_states_act_two_rdm_beta_beta_mo',wall_2 - wall_1
END_PROVIDER
@ -54,10 +63,8 @@
integer :: ispin
double precision :: wall_1, wall_2
! condition for alpha/beta spin
print*,''
print*,''
print*,''
print*,'providint all_states_act_two_rdm_alpha_beta_mo '
call wall_time(wall_1)
print*,'providing all_states_act_two_rdm_alpha_beta_mo ...'
ispin = 3
print*,'ispin = ',ispin
all_states_act_two_rdm_alpha_beta_mo = 0.d0
@ -65,7 +72,7 @@
call orb_range_all_states_two_rdm(all_states_act_two_rdm_alpha_beta_mo,n_act_orb,n_act_orb,list_act,list_act_reverse,ispin,psi_coef,size(psi_coef,2),size(psi_coef,1))
call wall_time(wall_2)
print*,'Wall time to provide all_states_act_two_rdm_alpha_beta_mo',wall_2 - wall_1
print*,'time to provide all_states_act_two_rdm_alpha_beta_mo',wall_2 - wall_1
END_PROVIDER
@ -89,6 +96,11 @@
ispin = 4
all_states_act_two_rdm_spin_trace_mo = 0.d0
double precision :: wall_1,wall_2
call wall_time(wall_1)
print*,'providing all_states_act_two_rdm_spin_trace_mo ...'
call orb_range_all_states_two_rdm(all_states_act_two_rdm_spin_trace_mo,n_act_orb,n_act_orb,list_act,list_act_reverse,ispin,psi_coef,size(psi_coef,2),size(psi_coef,1))
call wall_time(wall_2)
print*,'time to provide all_states_act_two_rdm_spin_trace_mo',wall_2 - wall_1
END_PROVIDER

View File

@ -16,10 +16,10 @@ subroutine routine_active_only
double precision :: vijkl,rdmaa,get_two_e_integral,rdmab,rdmbb,rdmtot
double precision :: accu_aa(N_states),accu_bb(N_states),accu_ab(N_states),accu_tot(N_states)
double precision :: accu_ab_omp,rdmab_omp
double precision :: accu_bb_omp,rdmbb_omp
double precision :: accu_aa_omp,rdmaa_omp
double precision :: accu_tot_omp,rdmtot_omp
double precision :: accu_ab_omp(N_states),rdmab_omp
double precision :: accu_bb_omp(N_states),rdmbb_omp
double precision :: accu_aa_omp(N_states),rdmaa_omp
double precision :: accu_tot_omp(N_states),rdmtot_omp
accu_ab_omp = 0.d0
accu_bb_omp = 0.d0
accu_aa_omp = 0.d0
@ -48,40 +48,41 @@ subroutine routine_active_only
! rdmtot_omp = state_av_act_two_rdm_spin_trace_mo(l,k,j,i)
rdmab_omp = all_states_openmp_act_two_rdm_alpha_beta_mo(l,k,j,i,istate)
rdmaa_omp = all_states_openmp_act_two_rdm_alpha_alpha_mo(l,k,j,i,istate)
rdmbb_omp = all_states_openmp_act_two_rdm_beta_beta_mo(l,k,j,i,istate)
! rdmab_omp = all_states_openmp_act_two_rdm_alpha_beta_mo(l,k,j,i,istate)
! rdmaa_omp = all_states_openmp_act_two_rdm_alpha_alpha_mo(l,k,j,i,istate)
! rdmbb_omp = all_states_openmp_act_two_rdm_beta_beta_mo(l,k,j,i,istate)
rdmtot_omp = all_states_openmp_act_two_rdm_spin_trace_mo(l,k,j,i,istate)
rdmaa = all_states_act_two_rdm_alpha_alpha_mo(l,k,j,i,istate)
rdmbb = all_states_act_two_rdm_beta_beta_mo(l,k,j,i,istate)
rdmab = all_states_act_two_rdm_alpha_beta_mo(l,k,j,i,istate)
! rdmtot = all_states_act_two_rdm_spin_trace_mo(l,k,j,i,istate)
! rdmab = all_states_act_two_rdm_alpha_beta_mo(l,k,j,i,istate)
! rdmaa = all_states_act_two_rdm_alpha_alpha_mo(l,k,j,i,istate)
! rdmbb = all_states_act_two_rdm_beta_beta_mo(l,k,j,i,istate)
rdmtot = all_states_act_two_rdm_spin_trace_mo(l,k,j,i,istate)
accu_ab_omp += vijkl * rdmab_omp
accu_aa_omp += vijkl * rdmaa_omp
accu_bb_omp += vijkl * rdmbb_omp
! accu_tot_omp += vijkl * rdmtot_omp
! accu_ab_omp(istate) += vijkl * rdmab_omp
! accu_aa_omp(istate) += vijkl * rdmaa_omp
! accu_bb_omp(istate) += vijkl * rdmbb_omp
accu_tot_omp(istate) += vijkl * rdmtot_omp
accu_ab(istate) += vijkl * rdmab
accu_aa(istate) += vijkl * rdmaa
accu_bb(istate) += vijkl * rdmbb
! accu_tot(istate) += vijkl * rdmtot
! accu_ab(istate) += vijkl * rdmab
! accu_aa(istate) += vijkl * rdmaa
! accu_bb(istate) += vijkl * rdmbb
accu_tot(istate) += vijkl * rdmtot
enddo
enddo
enddo
enddo
print*,''
print*,'Active space only energy '
print*,'accu_aa(istate) = ',accu_aa(istate)
print*,'accu_aa_omp = ',accu_aa_omp
print*,'accu_bb(istate) = ',accu_bb(istate)
print*,'accu_bb_omp = ',accu_bb_omp
print*,'accu_ab(istate) = ',accu_ab(istate)
print*,'accu_ab_omp = ',accu_ab_omp
! print*,'accu_aa(istate) = ',accu_aa(istate)
! print*,'accu_aa_omp = ',accu_aa_omp(istate)
! print*,'accu_bb(istate) = ',accu_bb(istate)
! print*,'accu_bb_omp = ',accu_bb_omp(istate)
! print*,'accu_ab(istate) = ',accu_ab(istate)
! print*,'accu_ab_omp = ',accu_ab_omp(istate)
print*,''
! print*,'sum (istate) = ',accu_aa(istate) + accu_bb(istate) + accu_ab(istate)
! print*,'accu_tot(istate) = ',accu_tot(istate)
! print*,'accu_tot_omp = ',accu_tot_omp
print*,'accu_tot(istate) = ',accu_tot(istate)
print*,'accu_tot_omp = ',accu_tot_omp(istate)
print*,'psi_energy_two_e(istate) = ',psi_energy_two_e(istate)
enddo