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

CISD_SC2_selection with n_det_max_cisd_sc2, pt2_max and do_pt2_end

This commit is contained in:
Manu 2014-08-16 16:06:27 +02:00
parent 256efa75b8
commit 6f805d17d5
9 changed files with 151 additions and 82 deletions

View File

@ -1 +1 @@
AOs BiInts Bitmask CISD SC2 CISD_selected Dets Electrons Ezfio_files Hartree_Fock MonoInts MOs Nuclei Output Perturbation Selectors_full SingleRefMethod Utils AOs BiInts Bitmask CISD CISD_selected Dets Electrons Ezfio_files Hartree_Fock MonoInts MOs Nuclei Output Perturbation Selectors_full SingleRefMethod Utils

View File

@ -0,0 +1,4 @@
cisd_sc2_selected
n_det_max_cisd_sc2 integer
pt2_max double precision
do_pt2_end logical

View File

@ -14,8 +14,24 @@ program cisd_sc2_selected
perturbation = "epstein_nesbet_sc2_projected" perturbation = "epstein_nesbet_sc2_projected"
E_old(1) = HF_energy E_old(1) = HF_energy
davidson_threshold = 1.d-6 davidson_threshold = 1.d-6
if (N_det > n_det_max_cisd_sc2) then
call diagonalize_CI_SC2
call save_wavefunction
psi_det = psi_det_sorted
psi_coef = psi_coef_sorted
N_det = n_det_max_cisd_sc2
soft_touch N_det psi_det psi_coef
call diagonalize_CI
call save_wavefunction
print *, 'N_det = ', N_det
print *, 'N_states = ', N_states
print *, 'PT2 = ', pt2
print *, 'E = ', CI_SC2_energy
print *, 'E+PT2 = ', CI_SC2_energy+pt2
print *, '-----'
endif
do while (maxval(abs(pt2(1:N_st))) > 1.d-4) do while (N_det < n_det_max_cisd_sc2.and.maxval(abs(pt2(1:N_st))) > pt2_max)
print*,'----' print*,'----'
print*,'' print*,''
call H_apply_cisd_selection(perturbation,pt2, norm_pert, H_pert_diag, N_st) call H_apply_cisd_selection(perturbation,pt2, norm_pert, H_pert_diag, N_st)
@ -37,37 +53,44 @@ program cisd_sc2_selected
exit exit
endif endif
enddo enddo
pt2 = 0.d0 N_det = min(n_det_max_cisd_sc2,N_det)
call H_apply_PT2(pt2, norm_pert, H_pert_diag, N_st) touch N_det psi_det psi_coef
davidson_threshold = 1.d-10 davidson_threshold = 1.d-10
touch davidson_threshold davidson_criterion touch davidson_threshold davidson_criterion
do i = 1, N_st call diagonalize_CI_SC2
max = 0.d0 pt2 = 0.d0
if(do_pt2_end)then
print*,'' threshold_selectors = 1.d0
print*,'-------------' call H_apply_PT2(pt2, norm_pert, H_pert_diag, N_st)
print*,'for state ',i do i = 1, N_st
print*,'' max = 0.d0
do k = 1, N_det
if(dabs(psi_coef(k,i)).gt.max)then
max = dabs(psi_coef(k,i))
imax = k
endif
enddo
double precision :: max
integer :: imax
print *, 'PT2(SC2) = ', pt2(i)
print *, 'E(SC2) = ', CI_SC2_energy(i)
print *, 'E_before(SC2)+PT2(SC2) = ', (CI_SC2_energy(i)+pt2(i))
if(i==1)then
print *, 'E(SC2)+PT2(projctd)SC2 = ', (CI_SC2_energy(i)+H_pert_diag(i))
endif
print*,'greater coeficient of the state : ',dabs(psi_coef(imax,i))
call get_excitation_degree(ref_bitmask,psi_det(1,1,imax),degree,N_int)
print*,'degree of excitation of such determinant : ',degree
enddo print*,''
print*,'-------------'
print*,'for state ',i
print*,''
print*,'N_det = ',N_det
do k = 1, N_det
if(dabs(psi_coef(k,i)).gt.max)then
max = dabs(psi_coef(k,i))
imax = k
endif
enddo
double precision :: max
integer :: imax
print *, 'PT2(SC2) = ', pt2(i)
print *, 'E(SC2) = ', CI_SC2_energy(i)
print *, 'E_before(SC2)+PT2(SC2) = ', (CI_SC2_energy(i)+pt2(i))
if(i==1)then
print *, 'E(SC2)+PT2(projctd)SC2 = ', (CI_SC2_energy(i)+H_pert_diag(i))
endif
print*,'greater coeficient of the state : ',dabs(psi_coef(imax,i))
call get_excitation_degree(ref_bitmask,psi_det(1,1,imax),degree,N_int)
print*,'degree of excitation of such determinant : ',degree
enddo
print*,'coucou' print*,'coucou'
endif
deallocate(pt2,norm_pert,H_pert_diag) deallocate(pt2,norm_pert,H_pert_diag)
end end

View File

@ -0,0 +1,51 @@
BEGIN_PROVIDER [ integer, n_det_max_cisd_sc2 ]
implicit none
BEGIN_DOC
! Get n_det_max_cisd_sc2 from EZFIO file
END_DOC
logical :: has_n_det_max_cisd_sc2
PROVIDE ezfio_filename
call ezfio_has_cisd_sc2_selected_n_det_max_cisd_sc2(has_n_det_max_cisd_sc2)
if (has_n_det_max_cisd_sc2) then
call ezfio_get_cisd_sc2_selected_n_det_max_cisd_sc2(n_det_max_cisd_sc2)
else
n_det_max_cisd_sc2 = 1000
call ezfio_set_cisd_sc2_selected_n_det_max_cisd_sc2(n_det_max_cisd_sc2)
endif
print*,'n_det_max_cisd_sc2 = ',n_det_max_cisd_sc2
END_PROVIDER
BEGIN_PROVIDER [ double precision , pt2_max ]
implicit none
BEGIN_DOC
! Get pt2_max from EZFIO file
END_DOC
logical :: has_pt2_max
PROVIDE ezfio_filename
call ezfio_has_cisd_sc2_selected_pt2_max(has_pt2_max)
if (has_pt2_max) then
call ezfio_get_cisd_sc2_selected_pt2_max(pt2_max)
else
pt2_max = 1.d-3
call ezfio_set_cisd_sc2_selected_pt2_max(pt2_max)
endif
print*,'pt2_max = ',pt2_max
END_PROVIDER
BEGIN_PROVIDER [ logical, do_pt2_end ]
implicit none
BEGIN_DOC
! Get do_pt2_end from EZFIO file
END_DOC
logical :: has_do_pt2_end
PROVIDE ezfio_filename
call ezfio_has_cisd_sc2_selected_do_pt2_end(has_do_pt2_end)
if (has_do_pt2_end) then
call ezfio_get_cisd_sc2_selected_do_pt2_end(do_pt2_end)
else
do_pt2_end = .True.
call ezfio_set_cisd_sc2_selected_do_pt2_end(do_pt2_end)
endif
print*,'do_pt2_end = ',do_pt2_end
END_PROVIDER

View File

@ -115,3 +115,46 @@ double precision function repeat_all_e_corr(key_in)
repeat_all_e_corr = accu repeat_all_e_corr = accu
end end
subroutine pt2_epstein_nesbet_sc2(det_pert,c_pert,e_2_pert,H_pert_diag,Nint,ndet,N_st)
use bitmasks
implicit none
integer, intent(in) :: Nint,ndet,N_st
integer(bit_kind), intent(in) :: det_pert(Nint,2)
double precision , intent(out) :: c_pert(N_st),e_2_pert(N_st),H_pert_diag(N_st)
double precision :: i_H_psi_array(N_st)
BEGIN_DOC
! compute the standard Epstein-Nesbet perturbative first order coefficient and second order energetic contribution
!
! for the various N_st states, but with the CISD_SC2 energies and coefficients
!
! c_pert(i) = <psi(i)|H|det_pert>/( E(i) - <det_pert|H|det_pert> )
!
! e_2_pert(i) = <psi(i)|H|det_pert>^2/( E(i) - <det_pert|H|det_pert> )
!
END_DOC
integer :: i,j
double precision :: diag_H_mat_elem, h
ASSERT (Nint == N_int)
ASSERT (Nint > 0)
call i_H_psi(det_pert,psi_selectors,psi_selectors_coef,Nint,N_det_selectors,psi_selectors_size,N_st,i_H_psi_array)
h = diag_H_mat_elem(det_pert,Nint)
do i =1,N_st
if(CI_SC2_electronic_energy(i)>h.and.CI_SC2_electronic_energy(i).ne.0.d0)then
c_pert(i) = -1.d0
e_2_pert(i) = selection_criterion*selection_criterion_factor*2.d0
else if (dabs(CI_SC2_electronic_energy(i) - h) > 1.d-6) then
c_pert(i) = i_H_psi_array(i) / (CI_SC2_electronic_energy(i) - h)
H_pert_diag(i) = h*c_pert(i)*c_pert(i)
e_2_pert(i) = c_pert(i) * i_H_psi_array(i)
else
c_pert(i) = -1.d0
e_2_pert(i) = -dabs(i_H_psi_array(i))
H_pert_diag(i) = h
endif
enddo
end

View File

@ -1 +0,0 @@
AOs BiInts Bitmask CISD Dets Electrons Ezfio_files Hartree_Fock MonoInts MOs Nuclei Output SingleRefMethod Utils Selectors_full

View File

@ -1,37 +0,0 @@
===============
CISD_SC2 Module
===============
Documentation
=============
.. Do not edit this section. It was auto-generated from the
.. NEEDED_MODULES file.
`cisd <http://github.com/LCPQ/quantum_package/tree/master/src/SC2/cisd_SC2.irp.f#L1>`_
Undocumented
Needed Modules
==============
.. Do not edit this section. It was auto-generated from the
.. NEEDED_MODULES file.
* `AOs <http://github.com/LCPQ/quantum_package/tree/master/src/AOs>`_
* `BiInts <http://github.com/LCPQ/quantum_package/tree/master/src/BiInts>`_
* `Bitmask <http://github.com/LCPQ/quantum_package/tree/master/src/Bitmask>`_
* `CISD <http://github.com/LCPQ/quantum_package/tree/master/src/CISD>`_
* `Dets <http://github.com/LCPQ/quantum_package/tree/master/src/Dets>`_
* `Electrons <http://github.com/LCPQ/quantum_package/tree/master/src/Electrons>`_
* `Ezfio_files <http://github.com/LCPQ/quantum_package/tree/master/src/Ezfio_files>`_
* `Hartree_Fock <http://github.com/LCPQ/quantum_package/tree/master/src/Hartree_Fock>`_
* `MonoInts <http://github.com/LCPQ/quantum_package/tree/master/src/MonoInts>`_
* `MOs <http://github.com/LCPQ/quantum_package/tree/master/src/MOs>`_
* `Nuclei <http://github.com/LCPQ/quantum_package/tree/master/src/Nuclei>`_
* `Output <http://github.com/LCPQ/quantum_package/tree/master/src/Output>`_
* `SingleRefMethod <http://github.com/LCPQ/quantum_package/tree/master/src/SingleRefMethod>`_
* `Utils <http://github.com/LCPQ/quantum_package/tree/master/src/Utils>`_
* `Selectors_full <http://github.com/LCPQ/quantum_package/tree/master/src/Selectors_full>`_

View File

@ -1,14 +0,0 @@
program cisd
implicit none
integer :: i
print *, ' HF = ', HF_energy
print *, 'N_states = ', N_states
call H_apply_cisd
print *, 'N_det = ', N_det
do i = 1,N_states
print *, 'energy = ',CI_SC2_energy(i)
print *, 'E_corr = ',CI_SC2_electronic_energy(i) - ref_bitmask_energy
enddo
end