40 use,
intrinsic :: iso_fortran_env
81 type(eigensolver_t) :: eigens
89 real(real64) :: occsum
91 real(real64) :: scale_f
93 real(real64) :: conv_ener
95 real(real64) :: tolerFO
97 real(real64),
allocatable :: eone(:)
98 real(real64),
allocatable :: eone_int(:, :)
99 real(real64),
allocatable :: twoint(:)
100 real(real64),
allocatable :: hartree(:, :)
101 real(real64),
allocatable :: exchange(:, :)
102 real(real64),
allocatable :: evalues(:)
103 real(real64),
allocatable :: vecnat(:, :)
104 real(real64),
allocatable :: Coul(:,:,:)
105 real(real64),
allocatable :: Exch(:,:,:)
107 integer,
allocatable :: i_index(:, :)
108 integer,
allocatable :: j_index(:, :)
109 integer,
allocatable :: k_index(:, :)
110 integer,
allocatable :: l_index(:, :)
113 type(rdm_t),
pointer :: rdm_ptr
118 subroutine rdmft_init(rdm, namespace, gr, st, mc, space, fromScratch)
119 type(rdm_t),
intent(out) :: rdm
120 type(namespace_t),
intent(in) :: namespace
121 type(grid_t),
intent(inout) :: gr
122 type(states_elec_t),
intent(in) :: st
123 type(multicomm_t),
intent(in) :: mc
124 class(space_t),
intent(in) :: space
125 logical,
intent(in) :: fromScratch
129 if(st%nst < st%qtot + 1)
then
130 message(1) =
"Too few states to run RDMFT calculation"
131 message(2) =
"Number of states should be at least the number of electrons plus one"
152 call parse_variable(namespace,
'RDMTolerance', 1.0e-7_real64, rdm%toler)
163 call parse_variable(namespace,
'RDMToleranceFO', 1.0e-4_real64, rdm%tolerFO)
188 if (rdm%do_basis .and. fromscratch)
then
189 call messages_write(
"RDMFT calculations with RDMBasis = yes cannot be started FromScratch", new_line=.
true.)
205 if (rdm%do_basis)
then
206 rdm%n_twoint = rdm%nst*(rdm%nst + 1)*(rdm%nst**2 + rdm%nst + 2)/8
207 safe_allocate(rdm%eone_int(1:rdm%nst, 1:rdm%nst))
208 safe_allocate(rdm%twoint(1:rdm%n_twoint))
209 safe_allocate(rdm%i_index(1:2,1:rdm%n_twoint))
210 safe_allocate(rdm%j_index(1:2,1:rdm%n_twoint))
211 safe_allocate(rdm%k_index(1:2,1:rdm%n_twoint))
212 safe_allocate(rdm%l_index(1:2,1:rdm%n_twoint))
213 safe_allocate(rdm%vecnat(1:rdm%nst, 1:rdm%nst))
214 safe_allocate(rdm%Coul(1:rdm%nst, 1:rdm%nst, 1:rdm%nst))
215 safe_allocate(rdm%Exch(1:rdm%nst, 1:rdm%nst, 1:rdm%nst))
228 if (rdm%eigens%additional_terms)
call messages_not_implemented(
"CG Additional Terms with RDMFT", namespace=namespace)
231 safe_allocate(rdm%eone(1:rdm%nst))
232 safe_allocate(rdm%hartree(1:rdm%nst, 1:rdm%nst))
233 safe_allocate(rdm%exchange(1:rdm%nst, 1:rdm%nst))
234 safe_allocate(rdm%evalues(1:rdm%nst))
240 rdm%mu =
m_two*st%eigenval(max(int(st%qtot*
m_half), 1), 1)
243 rdm%scale_f = 1e-2_real64
253 type(
rdm_t),
intent(inout) :: rdm
257 safe_deallocate_a(rdm%evalues)
258 safe_deallocate_a(rdm%eone)
259 safe_deallocate_a(rdm%hartree)
260 safe_deallocate_a(rdm%exchange)
262 if (rdm%do_basis)
then
263 safe_deallocate_a(rdm%eone_int)
264 safe_deallocate_a(rdm%twoint)
265 safe_deallocate_a(rdm%i_index)
266 safe_deallocate_a(rdm%j_index)
267 safe_deallocate_a(rdm%k_index)
268 safe_deallocate_a(rdm%l_index)
269 safe_deallocate_a(rdm%vecnat)
270 safe_deallocate_a(rdm%Coul)
271 safe_deallocate_a(rdm%Exch)
282 subroutine scf_rdmft(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, outp, restart_dump)
283 type(
rdm_t),
intent(inout) :: rdm
286 type(
grid_t),
intent(in) :: gr
287 type(
ions_t),
intent(in) :: ions
290 type(
v_ks_t),
intent(inout) :: ks
293 type(
restart_t),
intent(in) :: restart_dump
296 integer :: iter, icount, ip, ist, ierr, maxcount, iorb
297 integer(int64) :: what_i
298 real(real64) :: energy, energy_dif, energy_old, energy_occ, xpos, xneg, rel_ener
299 real(real64),
allocatable :: dpsi(:, :), dpsi2(:, :)
301 character(len=MAX_PATH_LEN) :: dirname
305 if (hm%d%ispin /= 1)
then
311 if (space%is_periodic())
then
316 if(st%parallel_in_states)
then
324 energy_old = 1.0e20_real64
328 if (.not. rdm%do_basis)
then
333 write(
message(1),
'(a)')
'Calculating Coulomb and exchange matrix elements in basis'
334 write(
message(2),
'(a)')
'--this may take a while--'
337 call two_body_me(gr, st, space, namespace, hm%kpoints, hm%exxop%psolver, 1, st%nst, rdm%i_index, rdm%j_index, rdm%k_index, &
338 rdm%l_index, rdm%twoint)
347 do iter = 1, rdm%max_iter
348 rdm%iter = rdm%iter + 1
349 write(
message(1),
'(a)')
'**********************************************************************'
350 write(
message(2),
'(a, i4)')
'Iteration:', iter
354 call scf_occ_no(rdm, namespace, gr, hm, space, st, energy_occ)
356 call scf_occ(rdm, namespace, gr, hm, space, st, energy_occ)
359 write(
message(1),
'(a)')
'Optimization of natural orbitals'
361 do icount = 1, maxcount
362 if (rdm%do_basis)
then
363 call scf_orb(rdm, namespace, gr, st, hm, space, energy)
365 call scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
367 energy_dif = energy - energy_old
369 if (rdm%do_basis)
then
370 if (abs(energy_dif)/abs(energy) < rdm%conv_ener .and. rdm%maxFO < rdm%tolerFO)
exit
371 if (energy_dif <
m_zero)
then
376 if (xneg > 1.5e0_real64*xpos)
then
377 rdm%scale_f = 1.01_real64*rdm%scale_f
378 elseif (xneg < 1.1e0_real64*xpos)
then
379 rdm%scale_f = 0.95_real64* rdm%scale_f
386 rel_ener = abs(energy_occ-energy)/abs(energy)
389 write(
message(2),
'(a,1x,es20.10)')
'Rel. energy difference:', rel_ener
392 if (.not. rdm%hf .and. rdm%do_basis)
then
393 write(
message(1),
'(a,18x,es20.10)')
'Max F0:', rdm%maxFO
398 if (rdm%do_basis)
then
399 conv = (rel_ener < rdm%conv_ener) .and. rdm%maxFO < rdm%tolerFO
401 conv = rel_ener < rdm%conv_ener
404 if (rdm%toler > 1e-4_real64) rdm%toler = rdm%toler*1e-1_real64
407 if ((conv .or. (modulo(iter, outp%restart_write_interval) == 0) .or. iter == rdm%max_iter))
then
408 if (rdm%do_basis)
then
410 safe_allocate(dpsi(1:gr%np, 1:st%d%dim))
411 safe_allocate(dpsi2(1:gr%np, 1:st%d%dim))
417 dpsi(ip,1) = dpsi(ip,1) + rdm%vecnat(ist, iorb)*dpsi2(ip,1)
424 call states_elec_dump(restart_dump, space, states_save, gr, hm%kpoints, ierr, iter=iter)
426 if (conv .or. iter == rdm%max_iter)
then
433 safe_deallocate_a(dpsi)
434 safe_deallocate_a(dpsi2)
436 call states_elec_dump(restart_dump, space, st, gr, hm%kpoints, ierr, iter=iter)
439 if (.not. rdm%hf)
then
441 write(
message(1),
'(a,18x,es20.10)')
'Max F0:', rdm%maxFO
447 message(1) =
'Unable to write states wavefunctions.'
454 if (any(outp%what) .and. outp%duringscf)
then
455 do what_i = lbound(outp%what, 1), ubound(outp%what, 1)
456 if (outp%what_now(what_i, iter))
then
457 write(dirname,
'(a,a,i4.4)') trim(outp%iter_dir),
"scf.", iter
458 call output_all(outp, namespace, space, dirname, gr, ions, iter, st, hm, ks)
459 call output_modelmb(outp, namespace, space, dirname, gr, ions, iter, st)
470 write(
message(1),
'(a,i3,a)')
'The calculation converged after ',rdm%iter,
' iterations'
474 write(
message(1),
'(a,i3,a)')
'The calculation did not converge after ', iter-1,
' iterations '
475 write(
message(2),
'(a,es15.5)')
'Relative energy difference between the last two iterations ', rel_ener
476 write(
message(3),
'(a,es15.5)')
'The maximal non-diagonal element of the Hermitian matrix F is ', rdm%maxFO
489 character(len=*),
intent(in) :: dir, fname
491 integer :: iunit, ist
492 real(real64),
allocatable :: photon_number_state (:), ekin_state (:), epot_state (:)
496 safe_allocate(photon_number_state(1:st%nst))
497 safe_allocate(ekin_state(1:st%nst))
498 safe_allocate(epot_state(1:st%nst))
502 iunit =
io_open(trim(dir) //
"/" // trim(fname), namespace, action=
'write')
508 if (rdm%do_basis)
then
509 write(iunit,
'(a)')
'Orbital optimization with [basis set]'
511 write(iunit,
'(a)')
'Orbital optimization with [conjugated gradients]'
516 write(iunit,
'(a)')
'Hartree Fock calculation'
520 if (hm%psolver%is_dressed)
then
521 write(iunit,
'(a)')
'Dressed state calculation'
528 write(iunit,
'(a, i4, a)')
'SCF converged in ', iter,
' iterations'
530 write(iunit,
'(a)')
'SCF *not* converged!'
536 write(iunit,
'(a,1x,f16.12)')
'Sum of occupation numbers:', rdm%occsum
541 if (hm%psolver%is_dressed)
then
542 call calc_photon_number(space, gr, st, hm%psolver%photons, photon_number_state, ekin_state, epot_state)
544 write(iunit,
'(a,1x,f14.12)')
'Total mode occupation:', hm%psolver%photons%number(1)
549 if (rdm%max_iter > 0)
then
550 write(iunit,
'(a)')
'Convergence:'
551 write(iunit,
'(6x, a, es15.8)')
'maxFO = ', rdm%maxFO
552 write(iunit,
'(6x, a, es15.8)')
'rel_ener = ', rel_ener
560 write(iunit,
'(a)')
'Natural occupation numbers:'
561 write(iunit,
'(a4,5x,a12)', advance=
'no')
'#st',
'Occupation'
562 if (.not. rdm%do_basis)
write(iunit,
'(5x,a12)', advance=
'no')
'conv'
563 if (hm%psolver%is_dressed)
write(iunit,
'(3(5x,a12))', advance=
'no')
'Mode Occ.',
'-1/2d^2/dq^2',
'1/2w^2q^2'
568 write(iunit,
'(i4,3x,f14.12)', advance=
'no') ist, st%occ(ist, 1)
569 if (.not. rdm%do_basis)
write(iunit,
'(3x,f14.12)', advance=
'no') rdm%eigens%diff(ist, 1)
570 if (hm%psolver%is_dressed)
then
571 write(iunit,
'(3(3x,f14.12))', advance=
'no') photon_number_state(ist), ekin_state(ist), epot_state(ist)
581 safe_deallocate_a(photon_number_state)
582 safe_deallocate_a(ekin_state)
583 safe_deallocate_a(epot_state)
590 subroutine calc_maxfo (namespace, hm, st, gr, rdm)
592 type(
rdm_t),
intent(inout) :: rdm
593 type(
grid_t),
intent(in) :: gr
597 real(real64),
allocatable :: lambda(:, :), FO(:, :)
602 safe_allocate(lambda(1:st%nst,1:st%nst))
603 safe_allocate(fo(1:st%nst, 1:st%nst))
613 fo(jst, ist) = - (lambda(jst, ist) - lambda(ist ,jst))
616 rdm%maxFO = maxval(abs(fo))
618 safe_deallocate_a(lambda)
619 safe_deallocate_a(fo)
625 subroutine calc_photon_number(space, gr, st, photons, photon_number_state, ekin_state, epot_state)
626 class(
space_t),
intent(in) :: space
627 type(
grid_t),
intent(in) :: gr
630 real(real64),
intent(out) :: photon_number_state(:)
631 real(real64),
intent(out) :: ekin_state(:)
632 real(real64),
intent(out) :: epot_state(:)
634 integer :: ist, dim_photon
635 real(real64) :: q2_exp, laplace_exp
636 real(real64),
allocatable :: psi(:, :), psi_q2(:), dpsidq(:), d2psidq2(:)
641 dim_photon = space%dim
643 safe_allocate(psi(1:gr%np_part, 1))
644 safe_allocate(psi_q2(1:gr%np))
645 safe_allocate(dpsidq(1:gr%np_part))
646 safe_allocate(d2psidq2(1:gr%np))
648 photons%number(1) =
m_zero
656 laplace_exp =
dmf_dotp(gr, psi(:, 1), d2psidq2(:))
657 ekin_state(ist) = -
m_half*laplace_exp
660 psi_q2(1:gr%np) = psi(1:gr%np, 1) * gr%x(1:gr%np, dim_photon)**2
661 q2_exp =
dmf_dotp(gr, psi(:, 1), psi_q2(:))
662 epot_state(ist) =
m_half * photons%omega(1)**2 * q2_exp
666 photon_number_state(ist) = -
m_half*laplace_exp / photons%omega(1) +
m_half * photons%omega(1) * q2_exp
667 photon_number_state(ist) = photon_number_state(ist) -
m_half
670 photons%number(1) = photons%number(1) + (photon_number_state(ist) +
m_half)*st%occ(ist, 1)
673 photons%number(1) = photons%number(1) - st%qtot/
m_two
675 safe_deallocate_a(psi)
676 safe_deallocate_a(psi_q2)
677 safe_deallocate_a(dpsidq)
678 safe_deallocate_a(d2psidq2)
689 real(real64),
allocatable :: occin(:, :)
693 safe_allocate(occin(1:st%nst, 1:st%nik))
695 occin(1:st%nst, 1:st%nik) = st%occ(1:st%nst, 1:st%nik)
697 where(occin(:, :) >
m_one) occin(:, :) = st%smear%el_per_state
699 st%occ(:, :) = occin(:, :)
701 safe_deallocate_a(occin)
710 subroutine scf_occ_no(rdm, namespace, gr, hm, space, st, energy)
711 type(
rdm_t),
intent(inout) :: rdm
713 type(
grid_t),
intent(in) :: gr
715 class(
space_t),
intent(in) :: space
717 real(real64),
intent(out) :: energy
723 write(
message(1),
'(a)')
'SKIP Optimization of occupation numbers'
734 rdm%occsum = sum(st%occ(1:st%nst, 1:st%nik))
736 write(
message(1),
'(a4,5x,a12)')
'#st',
'Occupation'
740 write(
message(1),
'(i4,3x,f11.9)') ist, st%occ(ist, 1)
744 write(
message(1),
'(a,1x,f13.9)')
'Sum of occupation numbers', rdm%occsum
752 subroutine scf_occ(rdm, namespace, gr, hm, space, st, energy)
753 type(
rdm_t),
target,
intent(inout) :: rdm
755 type(
grid_t),
intent(in) :: gr
757 class(
space_t),
intent(in) :: space
759 real(real64),
intent(out) :: energy
761 integer :: ist, icycle, ierr
762 real(real64) :: sumgi1, sumgi2, sumgim, mu1, mu2, mum, dinterv, thresh_occ
763 real(real64),
allocatable :: occin(:, :)
764 real(real64),
parameter :: smallocc = 0.00001_real64
765 real(real64),
allocatable :: theta(:)
766 real(real64) :: objective
767 integer,
parameter :: max_cycle = 200
772 write(
message(1),
'(a)')
'Optimization of occupation numbers'
775 safe_allocate(occin(1:st%nst, 1:st%nik))
776 safe_allocate(theta(1:st%nst))
784 thresh_occ = 1e-14_real64
787 occin(1:st%nst, 1:st%nik) = st%occ(1:st%nst, 1:st%nik)
788 where(occin(:, :) < smallocc) occin(:, :) = smallocc
789 where(occin(:, :) > st%smear%el_per_state - smallocc) occin(:, :) = st%smear%el_per_state - smallocc
794 st%occ(:, :) = occin(:, :)
811 sumgi1 = rdm%occsum - st%qtot
816 sumgi2 = rdm%occsum - st%qtot
819 do icycle = 1, max_cycle
820 if (sumgi1*sumgi2 <=
m_zero)
exit
829 sumgi1 = rdm%occsum - st%qtot
838 sumgi2 = rdm%occsum - st%qtot
848 sumgim = rdm%occsum - st%qtot
850 if (sumgi1*sumgim <
m_zero)
then
860 if (st%occ(ist,1) <= thresh_occ ) st%occ(ist,1) = thresh_occ
863 if (abs(sumgim) < rdm%toler .or. abs((mu1-mu2)*
m_half) < rdm%toler)
exit
868 if (icycle >= 50)
then
869 write(
message(1),
'(a,1x,f11.4)')
'Bisection ended without finding mu, sum of occupation numbers:', rdm%occsum
874 st%occ(ist, 1) = st%smear%el_per_state*
sin(theta(ist)*
m_pi*
m_two)**2
877 objective = objective + rdm%mu*(rdm%occsum - rdm%qtot)
880 write(
message(1),
'(a4,5x,a12)')
'#st',
'Occupation'
884 write(
message(1),
'(i4,3x,f14.12)') ist, st%occ(ist, 1)
888 write(
message(1),
'(a,3x,f11.9)')
'Sum of occupation numbers: ', rdm%occsum
892 safe_deallocate_a(occin)
893 safe_deallocate_a(theta)
901 integer,
intent(in) :: size
902 real(real64),
intent(in) :: theta(size)
903 real(real64),
intent(inout) :: objective
904 integer,
intent(in) :: getgrad
905 real(real64),
intent(inout) :: df(size)
908 real(real64) :: thresh_occ, thresh_theta
909 real(real64),
allocatable :: dE_dn(:),occ(:)
913 assert(
size == rdm_ptr%nst)
915 safe_allocate(de_dn(1:size))
916 safe_allocate(occ(1:size))
922 thresh_occ = 1e-14_real64
927 if (occ(ist) <= thresh_occ ) occ(ist) = thresh_occ
930 rdm_ptr%occsum = sum(occ(1:size))
937 if (occ(ist) <= thresh_occ )
then
943 objective = objective - rdm_ptr%mu*(rdm_ptr%occsum - rdm_ptr%qtot)
945 safe_deallocate_a(de_dn)
946 safe_deallocate_a(occ)
953 integer,
intent(in) :: iter
954 integer,
intent(in) :: size
955 real(real64),
intent(in) :: energy, maxdr, maxdf
956 real(real64),
intent(in) :: theta(size)
966 subroutine scf_orb(rdm, namespace, gr, st, hm, space, energy)
967 type(
rdm_t),
intent(inout) :: rdm
969 type(
grid_t),
intent(in) :: gr
972 class(
space_t),
intent(in) :: space
973 real(real64),
intent(out) :: energy
976 real(real64),
allocatable :: lambda(:, :), fo(:, :)
982 safe_allocate(lambda(1:st%nst,1:st%nst))
983 safe_allocate(fo(1:st%nst, 1:st%nst))
991 if (rdm%iter==1)
then
994 fo(ist, jst) =
m_half*(lambda(ist, jst) + lambda(jst, ist))
995 fo(jst, ist) = fo(ist, jst)
1001 fo(jst, ist) = - ( lambda(jst, ist) - lambda(ist ,jst))
1004 rdm%maxfo = maxval(abs(fo))
1006 fo(ist, ist) = rdm%evalues(ist)
1008 if(abs(fo(jst, ist)) > rdm%scale_f)
then
1009 fo(jst, ist) = rdm%scale_f*fo(jst,ist)/abs(fo(jst, ist))
1011 fo(ist, jst) = fo(jst, ist)
1022 safe_deallocate_a(lambda)
1023 safe_deallocate_a(fo)
1033 subroutine scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
1034 type(
rdm_t),
intent(inout) :: rdm
1037 type(
grid_t),
intent(in) :: gr
1038 type(
ions_t),
intent(in) :: ions
1041 type(
v_ks_t),
intent(inout) :: ks
1043 real(real64),
intent(out) :: energy
1045 integer :: ik, ist, maxiter
1051 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
1052 call hm%update(gr, namespace, space, ext_partners)
1054 rdm%eigens%converged = 0
1058 do ik = st%d%kpt%start, st%d%kpt%end
1059 rdm%eigens%matvec = 0
1060 maxiter = rdm%eigens%es_maxiter
1061 call deigensolver_cg(namespace, gr, st, hm, hm%xc, rdm%eigens%pre, rdm%eigens%tolerance, maxiter, &
1062 rdm%eigens%converged(ik), ik, rdm%eigens%diff(:, ik), rdm%eigens%energy_change_threshold, &
1063 rdm%eigens%orthogonalize_to_all, rdm%eigens%conjugate_direction, rdm%eigens%additional_terms)
1065 if (.not. rdm%eigens%folded_spectrum)
then
1067 rdm%eigens%converged(ik) = 0
1069 if(rdm%eigens%diff(ist, ik) < rdm%eigens%tolerance)
then
1070 rdm%eigens%converged(ik) = ist
1079 write(stdout,
'(1x)')
1084 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
1085 call hm%update(gr, namespace, space, ext_partners)
1101 type(
grid_t),
intent(in) :: gr
1102 real(real64),
intent(out) :: lambda(:, :)
1103 type(
rdm_t),
intent(inout) :: rdm
1105 real(real64),
allocatable :: hpsi(:, :), hpsi1(:, :), dpsi(:, :), dpsi1(:, :)
1106 real(real64),
allocatable :: fock(:,:,:), fvec(:)
1107 integer :: ist, iorb, jorb, jst
1114 if (.not. rdm%do_basis)
then
1115 safe_allocate(hpsi(1:gr%np,1:st%d%dim))
1116 safe_allocate(hpsi1(1:gr%np,1:st%d%dim))
1117 safe_allocate(dpsi(1:gr%np_part ,1:st%d%dim))
1118 safe_allocate(dpsi1(1:gr%np_part ,1:st%d%dim))
1124 do jorb = iorb, st%nst
1127 lambda(jorb, iorb) =
dmf_dotp(gr, dpsi1(:,1), hpsi(:,1))
1130 if (iorb /= jorb )
then
1132 lambda(iorb, jorb) =
dmf_dotp(gr, dpsi(:,1), hpsi1(:,1))
1140 safe_allocate(fvec(1:st%nst))
1141 safe_allocate(fock(1:st%nst, 1:st%nst, 1:st%nst))
1147 fock(ist, jst, iorb) = st%occ(iorb, 1)*rdm%eone_int(ist,jst)
1150 fock(ist ,jst, iorb) = fock(ist, jst, iorb) + st%occ(iorb, 1)*st%occ(jorb, 1)*rdm%Coul(ist, jst, jorb) &
1151 -
sqrt(st%occ(iorb, 1))*
sqrt(st%occ(jorb, 1))*rdm%Exch(ist, jst, jorb)
1153 fock(jst, ist, iorb) = fock(ist, jst, iorb)
1162 fvec(ist) = fvec(ist) + fock(ist, jst, jorb)*rdm%vecnat(jst, jorb)
1166 lambda(iorb, jorb) =
m_zero
1168 lambda(iorb, jorb) = lambda(iorb, jorb) + rdm%vecnat(ist, iorb)*fvec(ist)
1175 if (.not. rdm%do_basis)
then
1176 safe_deallocate_a(hpsi)
1177 safe_deallocate_a(hpsi1)
1178 safe_deallocate_a(dpsi)
1179 safe_deallocate_a(dpsi1)
1181 safe_deallocate_a(fvec)
1182 safe_deallocate_a(fock)
1192 type(
rdm_t),
intent(inout) :: rdm
1194 real(real64),
intent(in) :: lambda(:, :)
1196 integer :: iorb, jorb, ist
1197 real(real64),
allocatable :: vecnat_new(:, :)
1199 push_sub(assign_eigenfunctions)
1201 safe_allocate(vecnat_new(1:st%nst, 1:st%nst))
1204 vecnat_new(ist, iorb) =
m_zero
1206 vecnat_new(ist , iorb) = vecnat_new(ist, iorb) + rdm%vecnat(ist, jorb)*lambda(jorb, iorb)
1211 rdm%vecnat(:, :) = vecnat_new(:, :)
1213 safe_deallocate_a(vecnat_new)
1215 pop_sub(assign_eigenfunctions)
1222 type(
rdm_t),
intent(in) :: rdm
1223 real(real64),
intent(in) :: occ(:)
1224 real(real64),
intent(out) :: energy
1225 real(real64),
optional,
intent(out) :: dE_dn(:)
1228 real(real64),
allocatable :: V_h(:), V_x(:)
1232 safe_allocate(v_h(1:rdm%nst))
1233 safe_allocate(v_x(1:rdm%nst))
1243 v_h(ist) = v_h(ist) + occ(jst)*rdm%hartree(ist, jst)
1244 v_x(ist) = v_x(ist) -
sqrt(occ(jst))*rdm%exchange(ist, jst)
1246 v_x(ist) = v_x(ist)*
m_half/max(
sqrt(occ(ist)), 1.0e-16_real64)
1251 if (
present(de_dn))
then
1252 de_dn(:) = rdm%eone(:) + v_h(:) + v_x(:)
1257 energy = energy + occ(ist)*rdm%eone(ist) &
1258 +
m_half*occ(ist)*v_h(ist) &
1262 safe_deallocate_a(v_h)
1263 safe_deallocate_a(v_x)
1271 type(
rdm_t),
intent(inout) :: rdm
1275 type(
grid_t),
intent(in) :: gr
1276 class(
space_t),
intent(in) :: space
1279 real(real64),
allocatable :: hpsi(:, :), rho1(:), rho(:), dpsi(:, :), dpsi2(:, :)
1280 real(real64),
allocatable :: v_ij(:,:,:,:,:)
1284 integer :: ist, jst, nspin_, iorb, jorb
1289 nspin_ = min(st%d%nspin, 2)
1291 if (.not. rdm%do_basis)
then
1292 safe_allocate(hpsi(1:gr%np, 1:st%d%dim))
1293 safe_allocate(rho1(1:gr%np))
1294 safe_allocate(rho(1:gr%np))
1295 safe_allocate(dpsi(1:gr%np_part, 1:st%d%dim))
1296 safe_allocate(dpsi2(1:gr%np, 1:st%d%dim))
1297 safe_allocate(v_ij(1:gr%np, 1:st%nst, 1:st%nst, 1:st%nik, 1:st%nik))
1311 rdm%eone(ist) =
dmf_dotp(gr, dpsi(:, 1), hpsi(:, 1))
1325 rho1(1:gr%np) = dpsi(1:gr%np, 1)**2
1327 do jst = ist, st%nst
1328 rdm%hartree(ist, jst) =
dmf_dotp(gr, rho1, v_ij(:,jst, jst, 1, 1))
1329 rdm%hartree(jst, ist) = rdm%hartree(ist, jst)
1331 rho(1:gr%np) = dpsi2(1:gr%np, 1)*dpsi(1:gr%np, 1)
1332 rdm%exchange(ist, jst) =
dmf_dotp(gr, rho, v_ij(:, ist, jst, 1, 1))
1333 rdm%exchange(jst, ist) = rdm%exchange(ist, jst)
1338 safe_deallocate_a(hpsi)
1339 safe_deallocate_a(rho)
1340 safe_deallocate_a(rho1)
1341 safe_deallocate_a(dpsi)
1342 safe_deallocate_a(dpsi2)
1343 safe_deallocate_a(v_ij)
1351 dd = rdm%vecnat(ist, iorb)*rdm%vecnat(jst, iorb)
1352 rdm%eone(iorb) = rdm%eone(iorb) + dd*rdm%eone_int(ist, jst)
1359 rdm%hartree(iorb ,jorb) =
m_zero
1360 rdm%exchange(iorb,jorb) =
m_zero
1363 dd = rdm%vecnat(ist, iorb)*rdm%vecnat(jst, iorb)
1364 rdm%hartree(iorb ,jorb) = rdm%hartree(iorb ,jorb)+rdm%Coul(ist,jst, jorb)*dd
1365 rdm%exchange(iorb ,jorb) = rdm%exchange(iorb ,jorb)+rdm%Exch(ist,jst, jorb)*dd
1368 rdm%hartree(jorb, iorb) = rdm%hartree(iorb, jorb)
1369 rdm%exchange(jorb, iorb) = rdm%exchange(iorb, jorb)
1380 type(
rdm_t),
intent(inout) :: rdm
1384 class(
mesh_t),
intent(in) :: mesh
1386 real(real64),
allocatable :: hpsi(:, :)
1387 real(real64),
allocatable :: dpsi(:, :), dpsi2(:, :)
1392 safe_allocate(dpsi(1:mesh%np_part, 1:st%d%dim))
1393 safe_allocate(dpsi2(1:mesh%np, 1:st%d%dim))
1394 safe_allocate(hpsi(1:mesh%np, 1:st%d%dim))
1399 do jst = ist, st%nst
1405 rdm%eone_int(jst, ist) =
dmf_dotp(mesh, dpsi2(:, 1), hpsi(:, 1))
1406 rdm%eone_int(ist, jst) = rdm%eone_int(jst, ist)
1410 safe_deallocate_a(hpsi)
1411 safe_deallocate_a(dpsi)
1412 safe_deallocate_a(dpsi2)
1420 type(
rdm_t),
intent(inout) :: rdm
1422 integer :: ist, jst, kst, lst, iorb, icount
1423 logical :: inv_pairs
1424 real(real64) :: two_int, wij, wik, wil, wjk, wjl, wkl
1425 real(real64),
allocatable :: dm(:,:,:)
1429 safe_allocate(dm(1:rdm%nst, 1:rdm%nst, 1:rdm%nst))
1435 do iorb = 1, rdm%nst
1438 dm(ist, jst, iorb) = rdm%vecnat(ist, iorb)*rdm%vecnat(jst, iorb)
1439 dm(jst, ist, iorb) = dm(ist, jst, iorb)
1444 do icount = 1, rdm%n_twoint
1446 ist = rdm%i_index(1,icount)
1447 jst = rdm%j_index(1,icount)
1448 kst = rdm%k_index(1,icount)
1449 lst = rdm%l_index(1,icount)
1451 two_int = rdm%twoint(icount)
1465 if(ist == kst .and. jst /= lst)
then
1470 if(ist == lst .and. jst /= kst)
then
1475 if(jst == kst .and. ist /= lst)
then
1480 if(jst == lst .and. ist /= kst)
then
1486 inv_pairs = (ist /= kst .or. jst /= lst)
1488 do iorb = 1, rdm%nst
1491 rdm%Coul(ist, jst, iorb) = rdm%Coul(ist, jst, iorb) + dm(kst, lst, iorb)*wkl*two_int
1492 if (inv_pairs) rdm%Coul(kst, lst, iorb) = rdm%Coul(kst, lst, iorb) + dm(ist, jst, iorb)*wij*two_int
1496 rdm%Exch(ist, kst, iorb) = rdm%Exch(ist, kst, iorb) + two_int*dm(jst, lst, iorb)*wik
1497 if (kst /= lst)
then
1498 rdm%Exch(ist, lst, iorb) = rdm%Exch(ist, lst, iorb) + two_int*dm(jst, kst, iorb)*wil
1500 if (ist /= jst)
then
1502 rdm%Exch(jst, kst, iorb) = rdm%Exch(jst, kst, iorb) + two_int*dm(ist, lst, iorb)*wjk
1504 if (inv_pairs) rdm%Exch(kst, jst, iorb) = rdm%Exch(kst, jst, iorb) + two_int*dm(ist, lst, iorb)
1507 if (ist /=jst .and. kst /= lst)
then
1508 if (jst >= lst)
then
1509 rdm%Exch(jst, lst, iorb) = rdm%Exch(jst, lst, iorb) + two_int*dm(ist, kst, iorb)*wjl
1511 if (inv_pairs) rdm%Exch(lst, jst, iorb) = rdm%Exch(lst, jst, iorb) + two_int*dm(ist, kst, iorb)
1521 rdm%Coul(jst, ist, iorb) = rdm%Coul(ist, jst, iorb)
1522 rdm%Exch(jst, ist, iorb) = rdm%Exch(ist, jst, iorb)
1527 safe_deallocate_a(dm)
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
double sin(double __x) __attribute__((__nothrow__
double asin(double __x) __attribute__((__nothrow__
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
type(debug_t), save, public debug
This module implements a calculator for the density and defines related functions.
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_partial(der, ff, op_ff, dir, ghost_update, set_bc)
apply the partial derivative along dir to a mesh function
subroutine, public deigensolver_cg(namespace, mesh, st, hm, xc, pre, tol, niter, converged, ik, diff, energy_change_threshold, orthogonalize_to_all, conjugate_direction, additional_terms, shift)
conjugate-gradients method.
subroutine, public eigensolver_init(eigens, namespace, gr, st, mc, space, deactivate_oracle)
subroutine, public eigensolver_end(eigens)
subroutine, public dexchange_operator_compute_potentials(this, namespace, space, gr, st, xst, kpoints, ex, F_out)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
character(len= *), parameter, public static_dir
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public grid_write_info(gr, iunit, namespace)
integer, parameter, public term_local_external
integer, parameter, public term_non_local_potential
integer, parameter, public term_kinetic
subroutine, public dhamiltonian_elec_apply_single(hm, namespace, mesh, psi, hpsi, ist, ik, terms, set_bc, set_phase)
This module defines classes and functions for interaction partners.
subroutine, public io_close(iunit, grp)
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module is intended to contain "only mathematical" functions and procedures.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
subroutine, public messages_not_implemented(feature, namespace)
character(len=512), private msg
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module contains some common usage patterns of MPI routines.
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
this module contains the low-level part of the output system
subroutine, public output_modelmb(outp, namespace, space, dir, gr, ions, iter, st)
this module contains the output system
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
subroutine, public photon_mode_write_info(this, iunit, namespace)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine scf_occ(rdm, namespace, gr, hm, space, st, energy)
subroutine calc_maxfo(namespace, hm, st, gr, rdm)
subroutine objective_rdmft(size, theta, objective, getgrad, df)
subroutine scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
subroutine, public rdmft_end(rdm)
subroutine, public rdmft_init(rdm, namespace, gr, st, mc, space, fromScratch)
subroutine write_iter_info_rdmft(iter, size, energy, maxdr, maxdf, theta)
subroutine set_occ_pinning(st)
subroutine calc_photon_number(space, gr, st, photons, photon_number_state, ekin_state, epot_state)
subroutine assign_eigfunctions(rdm, st, lambda)
subroutine, public scf_rdmft(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, outp, restart_dump)
subroutine construct_lambda(namespace, hm, st, gr, lambda, rdm)
subroutine sum_integrals(rdm)
subroutine rdm_integrals(rdm, namespace, hm, st, mesh)
subroutine scf_orb(rdm, namespace, gr, st, hm, space, energy)
subroutine total_energy_rdm(rdm, occ, energy, dE_dn)
subroutine scf_occ_no(rdm, namespace, gr, hm, space, st, energy)
subroutine rdm_derivatives(rdm, namespace, hm, st, gr, space)
pure logical function, public states_are_complex(st)
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public v_ks_write_info(ks, iunit, namespace)
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
subroutine scf_write_static(dir, fname)
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.