9
1
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-12-22 03:23:29 +01:00
qp2/src/iterations/iterations.irp.f
Anthony Scemama 49e9488f62
Develop (#10)
* fixed laplacian of aos

* corrected the laplacians of aos

* added dft_one_e

* added new feature for new dft functionals

* changed the configure to add new functionals

* changed the configure

* added dft_one_e/README.rst

* added README.rst in new_functionals

* added source/programmers_guide/new_ks.rst

* Thesis Yann

* Added gmp installation in configure

* improved qp_e_conv_fci

* Doc

* Typos

* Added variance_max

* Fixed completion in qp_create

* modif TODO

* fixed DFT potential for n_states gt 1

* improved pot pbe

* trying to improve sr PBE

* fixed potential pbe

* fixed the vxc smashed for pbe sr and normal

* Comments in selection

* bug fixed by peter

* Fixed bug with zero beta electrons

* Update README.rst

* Update e_xc_new_func.irp.f

* Update links.rst

* Update quickstart.rst

* Update quickstart.rst

* updated cipsi

* Fixed energies of non-expected s2 (#9)

* Moved diag_algorithm in Davdison
2019-02-22 19:19:58 +01:00

43 lines
1.3 KiB
Fortran

BEGIN_PROVIDER [ double precision, extrapolated_energy, (N_iter,N_states) ]
implicit none
BEGIN_DOC
! Extrapolated energy, using E_var = f(PT2) where PT2=0
END_DOC
integer :: i
do i=1,min(N_states,N_det)
call extrapolate_data(N_iter, &
energy_iterations(i,1:N_iter), &
pt2_iterations(i,1:N_iter), &
extrapolated_energy(1:N_iter,i))
enddo
END_PROVIDER
subroutine save_iterations(e_, pt2_,n_)
implicit none
BEGIN_DOC
! Update the energy in the EZFIO file.
END_DOC
integer, intent(in) :: n_
double precision, intent(in) :: e_(N_states), pt2_(N_states)
integer :: i
if (N_iter == 101) then
do i=2,N_iter-1
energy_iterations(1:N_states,N_iter-1) = energy_iterations(1:N_states,N_iter)
pt2_iterations(1:N_states,N_iter-1) = pt2_iterations(1:N_states,N_iter)
enddo
N_iter = N_iter-1
TOUCH N_iter
endif
energy_iterations(1:N_states,N_iter) = e_(1:N_states)
pt2_iterations(1:N_states,N_iter) = pt2_(1:N_states)
n_det_iterations(N_iter) = n_
call ezfio_set_iterations_N_iter(N_iter)
call ezfio_set_iterations_energy_iterations(energy_iterations)
call ezfio_set_iterations_pt2_iterations(pt2_iterations)
call ezfio_set_iterations_n_det_iterations(n_det_iterations)
end