Octopus
rdmft.F90
Go to the documentation of this file.
1!! Copyright (C) 2012-2019 I. Theophilou, N. Helbig
2!! Copyright (C) 2019 F. Buchholz, M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module rdmft_oct_m
23 use debug_oct_m
30 use energy_oct_m
32 use global_oct_m
33 use grid_oct_m
37 use io_oct_m
39 use ions_oct_m
40 use, intrinsic :: iso_fortran_env
42 use loct_oct_m
43 use math_oct_m
44 use mesh_oct_m
48 use mpi_oct_m
52 use output_oct_m
55 use parser_oct_m
60 use space_oct_m
65 use unit_oct_m
67 use v_ks_oct_m
68 use xc_oep_oct_m
69
70 implicit none
71
72 private
73 public :: &
74 rdm_t, &
75 rdmft_init, &
76 rdmft_end, &
78
79 type rdm_t
80 private
81 type(eigensolver_t) :: eigens
82 integer :: max_iter
83 integer :: iter
84 integer :: nst
85 integer :: n_twoint !number of unique two electron integrals
86 logical :: do_basis
87 logical :: hf
88 real(real64) :: mu
89 real(real64) :: occsum
90 real(real64) :: qtot
91 real(real64) :: scale_f
92 real(real64) :: toler
93 real(real64) :: conv_ener
94 real(real64) :: maxFO
95 real(real64) :: tolerFO
96
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(:,:,:)
106
107 integer, allocatable :: i_index(:, :)
108 integer, allocatable :: j_index(:, :)
109 integer, allocatable :: k_index(:, :)
110 integer, allocatable :: l_index(:, :)
111 end type rdm_t
112
113 type(rdm_t), pointer :: rdm_ptr
114
115contains
116
117 ! ---------------------------------------------------------
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
126
127 push_sub(rdmft_init)
128
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"
132 call messages_fatal(2, namespace=namespace)
133 end if
134
135 if (states_are_complex(st)) then
136 call messages_not_implemented("Complex states for RDMFT", namespace=namespace)
137 end if
138
139 ! The documentation for the variable is found in scf_init.
140 call parse_variable(namespace, 'MaximumIter', 200, rdm%max_iter)
141
142 !%Variable RDMTolerance
143 !%Type float
144 !%Default 1e-7 Ha
145 !%Section SCF::RDMFT
146 !%Description
147 !% Convergence criterion for stopping the occupation numbers minimization. Minimization is
148 !% stopped when all derivatives of the energy wrt. each occupation number
149 !% are smaller than this criterion. The bisection for finding the correct mu that is needed
150 !% for the occupation number minimization also stops according to this criterion.
151 !%End
152 call parse_variable(namespace, 'RDMTolerance', 1.0e-7_real64, rdm%toler)
153
154 !%Variable RDMToleranceFO
155 !%Type float
156 !%Default 1e-4 Ha
157 !%Section SCF::RDMFT
158 !%Description
159 !% Convergence criterion for stopping the diagonalization of the Fock matrix in the Piris method.
160 !% Orbital minimization is stopped when all off-diagonal ellements of the Fock matrix
161 !% are smaller than this criterion.
162 !%End
163 call parse_variable(namespace, 'RDMToleranceFO', 1.0e-4_real64, rdm%tolerFO)
164
165 !%Variable RDMConvEner
166 !%Type float
167 !%Default 1e-6 Ha
168 !%Section SCF::RDMFT
169 !%Description
170 !% Convergence criterion for stopping the overall minimization of the energy with
171 !% respect to occupation numbers and the orbitals. The minimization of the
172 !% energy stops when the total energy difference between two subsequent
173 !% minimizations of the energy with respect to the occupation numbers and the
174 !% orbitals is smaller than this criterion. It is also used to exit the orbital minimization.
175 !%End
176 call parse_variable(namespace, 'RDMConvEner', 1.0e-7_real64, rdm%conv_ener)
178 !%Variable RDMBasis
179 !%Type logical
180 !%Default yes
181 !%Section SCF::RDMFT
182 !%Description
183 !% If true, all the energy terms and corresponding derivatives involved in RDMFT will
184 !% not be calculated on the grid but on the basis of the initial orbitals
185 !%End
186 call parse_variable(namespace, 'RDMBasis',.true., rdm%do_basis)
188 if (rdm%do_basis .and. fromscratch) then
189 call messages_write("RDMFT calculations with RDMBasis = yes cannot be started FromScratch", new_line=.true.)
190 call messages_write("Run a calculation for independent particles first")
191 call messages_fatal(namespace=namespace)
192 end if
194 !%Variable RDMHartreeFock
195 !%Type logical
196 !%Default no
197 !%Section SCF::RDMFT
198 !%Description
199 !% If true, the code simulates a HF calculation, by omitting the occ.num. optimization
200 !% can be used for test reasons
201 !%End
202 call parse_variable(namespace, 'RDMHartreeFock',.false., rdm%hf)
204 rdm%nst = st%nst
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))
216 rdm%eone_int = m_zero
217 rdm%twoint = m_zero
218 rdm%vecnat(:, :) = diagonal_matrix(rdm%nst, m_one)
219 rdm%i_index = m_zero
220 rdm%j_index = m_zero
221 rdm%k_index = m_zero
222 rdm%l_index = m_zero
223 rdm%Coul = m_zero
224 rdm%Exch = m_zero
225 else
226 ! initialize eigensolver.
227 call eigensolver_init(rdm%eigens, namespace, gr, st, mc, space)
228 if (rdm%eigens%additional_terms) call messages_not_implemented("CG Additional Terms with RDMFT", namespace=namespace)
229 end if
230
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))
235
236 rdm%eone = m_zero
237 rdm%hartree = m_zero
238 rdm%exchange = m_zero
239 rdm%evalues = m_zero
240 rdm%mu = m_two*st%eigenval(max(int(st%qtot*m_half), 1), 1)
241 rdm%qtot = st%qtot
242 rdm%occsum = m_zero
243 rdm%scale_f = 1e-2_real64
244 rdm%maxFO = m_zero
245 rdm%iter = 0
246
247 pop_sub(rdmft_init)
248 end subroutine rdmft_init
249
250 ! ----------------------------------------
251
252 subroutine rdmft_end(rdm)
253 type(rdm_t), intent(inout) :: rdm
254
255 push_sub(rdmft_end)
256
257 safe_deallocate_a(rdm%evalues)
258 safe_deallocate_a(rdm%eone)
259 safe_deallocate_a(rdm%hartree)
260 safe_deallocate_a(rdm%exchange)
261
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)
272 else
273 call eigensolver_end(rdm%eigens)
274 end if
275
276 pop_sub(rdmft_end)
277 end subroutine rdmft_end
278
279 ! ----------------------------------------
280
281 ! scf for the occupation numbers and the natural orbitals
282 subroutine scf_rdmft(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, outp, restart_dump)
283 type(rdm_t), intent(inout) :: rdm
284 type(namespace_t), intent(in) :: namespace
285 type(electron_space_t), intent(in) :: space
286 type(grid_t), intent(in) :: gr
287 type(ions_t), intent(in) :: ions
288 type(partner_list_t), intent(in) :: ext_partners
289 type(states_elec_t), intent(inout) :: st
290 type(v_ks_t), intent(inout) :: ks
291 type(hamiltonian_elec_t), intent(inout) :: hm
292 type(output_t), intent(in) :: outp
293 type(restart_t), intent(in) :: restart_dump
294
295 type(states_elec_t) :: states_save
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(:, :)
300 logical :: conv
301 character(len=MAX_PATH_LEN) :: dirname
302
303 push_sub(scf_rdmft)
304
305 if (hm%d%ispin /= 1) then
306 call messages_not_implemented("RDMFT exchange function not yet implemented for spin_polarized or spinors", &
307 namespace=namespace)
308 end if
309
310 ! problem is about k-points for exchange
311 if (space%is_periodic()) then
312 call messages_not_implemented("Periodic system calculations for RDMFT", namespace=namespace)
313 end if
314
315 ! exchange routine needs all states on each processor currently
316 if(st%parallel_in_states) then
317 call messages_not_implemented("RDMFT parallel in states", namespace=namespace)
318 end if
319
320 call messages_print_with_emphasis(msg='RDMFT Calculation', namespace=namespace)
321 call messages_print_var_value('RDMBasis', rdm%do_basis, namespace=namespace)
322
323 !set initial values
324 energy_old = 1.0e20_real64
325 xpos = m_zero
326 xneg = m_zero
327 energy = m_zero
328 if (.not. rdm%do_basis) then
329 maxcount = 1 !still needs to be checked
330 else
331 maxcount = 50
332 !precalculate matrix elements in basis
333 write(message(1),'(a)') 'Calculating Coulomb and exchange matrix elements in basis'
334 write(message(2),'(a)') '--this may take a while--'
335 call messages_info(2, namespace=namespace)
336
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)
339 call rdm_integrals(rdm, namespace, hm, st, gr)
340 call sum_integrals(rdm)
341 endif
342
343 ! Start the actual minimization, first step is minimization of occupation numbers
344 ! Orbital minimization is according to Piris and Ugalde, Vol. 30, No. 13, J. Comput. Chem. (scf_orb) or
345 ! using conjugated gradient (scf_orb_cg)
346 conv = .false.
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
351 call messages_info(2, namespace=namespace)
352 ! occupation number optimization unless we are doing Hartree-Fock
353 if (rdm%hf) then
354 call scf_occ_no(rdm, namespace, gr, hm, space, st, energy_occ)
355 else
356 call scf_occ(rdm, namespace, gr, hm, space, st, energy_occ)
357 end if
358 ! orbital optimization
359 write(message(1), '(a)') 'Optimization of natural orbitals'
360 call messages_info(1, namespace=namespace)
361 do icount = 1, maxcount
362 if (rdm%do_basis) then
363 call scf_orb(rdm, namespace, gr, st, hm, space, energy)
364 else
365 call scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
366 end if
367 energy_dif = energy - energy_old
368 energy_old = energy
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
372 xneg = xneg + 1
373 else
374 xpos = xpos + 1
375 end if
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
380 end if
381 endif !rdm%do_basis
382 end do !icount
383 xneg = m_zero
384 xpos = m_zero
385
386 rel_ener = abs(energy_occ-energy)/abs(energy)
387
388 write(message(1),'(a,11x,es20.10)') 'Total energy:', units_from_atomic(units_out%energy,energy + hm%ep%eii)
389 write(message(2),'(a,1x,es20.10)') 'Rel. energy difference:', rel_ener
390 call messages_info(2, namespace=namespace)
391
392 if (.not. rdm%hf .and. rdm%do_basis) then
393 write(message(1),'(a,18x,es20.10)') 'Max F0:', rdm%maxFO
394 call messages_info(1, namespace=namespace)
395 end if
396
397
398 if (rdm%do_basis) then
399 conv = (rel_ener < rdm%conv_ener) .and. rdm%maxFO < rdm%tolerFO
400 else
401 conv = rel_ener < rdm%conv_ener
402 endif
403
404 if (rdm%toler > 1e-4_real64) rdm%toler = rdm%toler*1e-1_real64 !Is this still okay or does it restrict the possible convergence? FB: Does this makes sense at all?
405
406 ! save restart information
407 if ((conv .or. (modulo(iter, outp%restart_write_interval) == 0) .or. iter == rdm%max_iter)) then
408 if (rdm%do_basis) then
409 call states_elec_copy(states_save, st)
410 safe_allocate(dpsi(1:gr%np, 1:st%d%dim))
411 safe_allocate(dpsi2(1:gr%np, 1:st%d%dim))
412 do iorb = 1, st%nst
413 dpsi = m_zero
414 do ist = 1, st%nst
415 call states_elec_get_state(st, gr, ist, 1, dpsi2)
416 do ip = 1, gr%np
417 dpsi(ip,1) = dpsi(ip,1) + rdm%vecnat(ist, iorb)*dpsi2(ip,1)
418 end do
419 end do
420 call states_elec_set_state(states_save, gr, iorb, 1, dpsi)
421 end do
422 call density_calc(states_save, gr, states_save%rho)
423 ! if other quantities besides the densities and the states are needed they also have to be recalculated here!
424 call states_elec_dump(restart_dump, space, states_save, gr, hm%kpoints, ierr, iter=iter)
425
426 if (conv .or. iter == rdm%max_iter) then
427 call states_elec_end(st)
428 call states_elec_copy(st, states_save)
429 end if
430
431 call states_elec_end(states_save)
432
433 safe_deallocate_a(dpsi)
434 safe_deallocate_a(dpsi2)
435 else
436 call states_elec_dump(restart_dump, space, st, gr, hm%kpoints, ierr, iter=iter)
437
438 ! calculate maxFO for cg-solver
439 if (.not. rdm%hf) then
440 call calc_maxfo (namespace, hm, st, gr, rdm)
441 write(message(1),'(a,18x,es20.10)') 'Max F0:', rdm%maxFO
442 call messages_info(1, namespace=namespace)
443 end if
444 endif
445
446 if (ierr /= 0) then
447 message(1) = 'Unable to write states wavefunctions.'
448 call messages_warning(1, namespace=namespace)
449 end if
450
451 endif
452
453 ! write output for iterations if requested
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)
460 call scf_write_static(dirname, "info")
461 exit
462 end if
463 end do
464 end if
465
466 if (conv) exit
467 end do
468
469 if(conv) then
470 write(message(1),'(a,i3,a)') 'The calculation converged after ',rdm%iter,' iterations'
471 write(message(2),'(a,9x,es20.10)') 'The total energy is ', units_from_atomic(units_out%energy,energy + hm%ep%eii)
472 call messages_info(2, namespace=namespace)
473 else
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
477 call messages_info(3, namespace=namespace)
478 end if
479
480 call scf_write_static(static_dir, "info")
481 call output_all(outp, namespace, space, static_dir, gr, ions, -1, st, hm, ks)
482 call output_modelmb(outp, namespace, space, static_dir, gr, ions, -1, st)
483
484 pop_sub(scf_rdmft)
485
486 contains
487 ! ---------------------------------------------------------
488 subroutine scf_write_static(dir, fname)
489 character(len=*), intent(in) :: dir, fname
490
491 integer :: iunit, ist
492 real(real64), allocatable :: photon_number_state (:), ekin_state (:), epot_state (:)
493
495
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))
499
500 if(mpi_grp_is_root(mpi_world)) then
501 call io_mkdir(dir, namespace)
502 iunit = io_open(trim(dir) // "/" // trim(fname), namespace, action='write')
503
504 call grid_write_info(gr, iunit=iunit)
505
506 call v_ks_write_info(ks, iunit=iunit)
507
508 if (rdm%do_basis) then
509 write(iunit, '(a)')'Orbital optimization with [basis set]'
510 else
511 write(iunit, '(a)')'Orbital optimization with [conjugated gradients]'
512 end if
513 write(iunit, '(1x)')
514
515 if (rdm%hf) then
516 write(iunit, '(a)')'Hartree Fock calculation'
517 write(iunit, '(1x)')
518 end if
519
520 if (hm%psolver%is_dressed) then
521 write(iunit, '(a)')'Dressed state calculation'
522 call photon_mode_write_info(hm%psolver%photons, iunit=iunit)
523 write(iunit, '(1x)')
524 end if
525
526 ! scf information
527 if(conv) then
528 write(iunit, '(a, i4, a)')'SCF converged in ', iter, ' iterations'
529 else
530 write(iunit, '(a)') 'SCF *not* converged!'
531 end if
532 write(iunit, '(1x)')
533
534 write(iunit, '(3a,es20.10)') 'Total Energy [', trim(units_abbrev(units_out%energy)), ']:', &
535 units_from_atomic(units_out%energy, energy + hm%ep%eii)
536 write(iunit,'(a,1x,f16.12)') 'Sum of occupation numbers:', rdm%occsum
537 else
538 iunit = 0
539 end if
540
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)
543 if(mpi_grp_is_root(mpi_world)) then
544 write(iunit,'(a,1x,f14.12)') 'Total mode occupation:', hm%psolver%photons%number(1)
545 end if
546 end if
547
548 if(mpi_grp_is_root(mpi_world)) then
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
553 write(iunit,'(1x)')
554 end if
555 ! otherwise, these values are uninitialized, and unknown.
556 end if
557
558 if (mpi_grp_is_root(mpi_world)) then
559 ! Write header
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'
564 write(iunit,*)
565
566 ! Write values
567 do ist = 1, st%nst
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)
572 end if
573 write(iunit,*)
574 end do
575 end if
576
577 if (mpi_grp_is_root(mpi_world)) then
578 call io_close(iunit)
579 end if
580
581 safe_deallocate_a(photon_number_state)
582 safe_deallocate_a(ekin_state)
583 safe_deallocate_a(epot_state)
584
586 end subroutine scf_write_static
587 end subroutine scf_rdmft
588
589 ! ---------------------------------------------------------
590 subroutine calc_maxfo (namespace, hm, st, gr, rdm)
591 type(namespace_t), intent(in) :: namespace
592 type(rdm_t), intent(inout) :: rdm
593 type(grid_t), intent(in) :: gr
594 type(hamiltonian_elec_t), intent(inout) :: hm
595 type(states_elec_t), intent(inout) :: st
596
597 real(real64), allocatable :: lambda(:, :), FO(:, :)
598 integer :: ist, jst
599
600 push_sub(calc_maxfo)
601
602 safe_allocate(lambda(1:st%nst,1:st%nst))
603 safe_allocate(fo(1:st%nst, 1:st%nst))
604
605 ! calculate FO operator to check Hermiticity of lagrange multiplier matrix (lambda)
606 lambda = m_zero
607 fo = m_zero
608 call construct_lambda(namespace, hm, st, gr, lambda, rdm)
609
610 !Set up FO matrix to check maxFO
611 do ist = 1, st%nst
612 do jst = 1, ist - 1
613 fo(jst, ist) = - (lambda(jst, ist) - lambda(ist ,jst))
614 end do
615 end do
616 rdm%maxFO = maxval(abs(fo))
617
618 safe_deallocate_a(lambda)
619 safe_deallocate_a(fo)
620
621 pop_sub(calc_maxfo)
622 end subroutine calc_maxfo
623
624 ! ---------------------------------------------------------
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
628 type(states_elec_t), intent(in) :: st
629 type(photon_mode_t), intent(inout) :: photons
630 real(real64), intent(out) :: photon_number_state(:)
631 real(real64), intent(out) :: ekin_state(:)
632 real(real64), intent(out) :: epot_state(:)
633
634 integer :: ist, dim_photon
635 real(real64) :: q2_exp, laplace_exp
636 real(real64), allocatable :: psi(:, :), psi_q2(:), dpsidq(:), d2psidq2(:)
637
638 push_sub(calc_photon_number)
639
640 ! The photon dimension is always the last
641 dim_photon = space%dim
642
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))
647
648 photons%number(1) = m_zero
649
650 do ist = 1, st%nst
651 call states_elec_get_state(st, gr, ist, 1, psi)
652
653 ! <phi(ist)|d^2/dq^2|phi(ist)> ~= <phi(ist)| d/dq (d/dq|phi(ist)>)
654 call dderivatives_partial(gr%der, psi(:, 1), dpsidq(:), dim_photon, ghost_update = .true., set_bc = .true.)
655 call dderivatives_partial(gr%der, dpsidq(1:gr%np_part), d2psidq2(:), dim_photon, ghost_update = .true., set_bc = .true.)
656 laplace_exp = dmf_dotp(gr, psi(:, 1), d2psidq2(:))
657 ekin_state(ist) = -m_half*laplace_exp
658
659 ! <phi(ist)|q^2|psi(ist)>= |q|psi(ist)>|^2
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
663
664 !! N_phot(ist)=( <phi_i|H_ph|phi_i>/omega - 0.5 ) / N_elec
665 !! with <phi_i|H_ph|phi_i>=-0.5* <phi(ist)|d^2/dq^2|phi(ist)> + 0.5*omega <phi(ist)|q^2|psi(ist)>
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
668
669 !! N_phot_total= sum_ist occ_ist*N_phot(ist)
670 photons%number(1) = photons%number(1) + (photon_number_state(ist) + m_half)*st%occ(ist, 1) ! 0.5 must be added again to do the normalization due to the total charge correctly
671 end do
672
673 photons%number(1) = photons%number(1) - st%qtot/m_two
674
675 safe_deallocate_a(psi)
676 safe_deallocate_a(psi_q2)
677 safe_deallocate_a(dpsidq)
678 safe_deallocate_a(d2psidq2)
679
680 pop_sub(calc_photon_number)
681 end subroutine calc_photon_number
682
683 ! ---------------------------------------------------------
684
685 ! reset occ.num. to 2/0
686 subroutine set_occ_pinning(st)
687 type(states_elec_t), intent(inout) :: st
688
689 real(real64), allocatable :: occin(:, :)
690
691 push_sub(set_occ_pinning)
692
693 safe_allocate(occin(1:st%nst, 1:st%nik))
694
695 occin(1:st%nst, 1:st%nik) = st%occ(1:st%nst, 1:st%nik)
696 where(occin(:, :) < m_one) occin(:, :) = m_zero
697 where(occin(:, :) > m_one) occin(:, :) = st%smear%el_per_state
698
699 st%occ(:, :) = occin(:, :)
700
701 safe_deallocate_a(occin)
702
703 pop_sub(set_occ_pinning)
704 end subroutine set_occ_pinning
705
706
707 ! ---------------------------------------------------------
708 ! dummy routine for occupation numbers which only calculates the necessary variables for further use
709 ! used in Hartree-Fock mode
710 subroutine scf_occ_no(rdm, namespace, gr, hm, space, st, energy)
711 type(rdm_t), intent(inout) :: rdm
712 type(namespace_t), intent(in) :: namespace
713 type(grid_t), intent(in) :: gr
714 type(hamiltonian_elec_t), intent(in) :: hm
715 class(space_t), intent(in) :: space
716 type(states_elec_t), intent(inout) :: st
717 real(real64), intent(out) :: energy
719 integer :: ist
720
721 push_sub(scf_occ_no)
722
723 write(message(1),'(a)') 'SKIP Optimization of occupation numbers'
724 call messages_info(1, namespace=namespace)
725
726 call set_occ_pinning(st)
727
728 energy = m_zero
729
730 call rdm_derivatives(rdm, namespace, hm, st, gr, space)
731
732 call total_energy_rdm(rdm, st%occ(:,1), energy)
733
734 rdm%occsum = sum(st%occ(1:st%nst, 1:st%nik))
735
736 write(message(1),'(a4,5x,a12)')'#st','Occupation'
737 call messages_info(1, namespace=namespace)
738
739 do ist = 1, st%nst
740 write(message(1),'(i4,3x,f11.9)') ist, st%occ(ist, 1)
741 call messages_info(1, namespace=namespace)
742 end do
743
744 write(message(1),'(a,1x,f13.9)') 'Sum of occupation numbers', rdm%occsum
745 write(message(2),'(a,es20.10)') 'Total energy occ', units_from_atomic(units_out%energy,energy + hm%ep%eii)
746 call messages_info(2, namespace=namespace)
747
748 pop_sub(scf_occ_no)
749 end subroutine scf_occ_no
750
751 ! scf for the occupation numbers
752 subroutine scf_occ(rdm, namespace, gr, hm, space, st, energy)
753 type(rdm_t), target, intent(inout) :: rdm
754 type(namespace_t), intent(in) :: namespace
755 type(grid_t), intent(in) :: gr
756 type(hamiltonian_elec_t), intent(in) :: hm
757 class(space_t), intent(in) :: space
758 type(states_elec_t), intent(inout) :: st
759 real(real64), intent(out) :: energy
760
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
768
769 push_sub(scf_occ)
770 call profiling_in("SCF_OCC")
771
772 write(message(1),'(a)') 'Optimization of occupation numbers'
773 call messages_info(1, namespace=namespace)
774
775 safe_allocate(occin(1:st%nst, 1:st%nik))
776 safe_allocate(theta(1:st%nst))
777
778 occin = m_zero
779 theta = m_zero
780 energy = m_zero
781
782 ! Defines a threshold on occ nums to avoid numerical instabilities.
783 ! Needs to be changed consistently with the same variable in objective_rdmft
784 thresh_occ = 1e-14_real64
785
786 !Initialize the occin. Smallocc is used for numerical stability
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
790
791 !Renormalize the occupation numbers
792 rdm%occsum = st%qtot
793
794 st%occ(:, :) = occin(:, :)
795
796 call rdm_derivatives(rdm, namespace, hm, st, gr, space)
797
798 !finding the chemical potential mu such that the occupation numbers sum up to the number of electrons
799 !bisection to find the root of rdm%occsum-st%qtot=M_ZERO
800 mu1 = rdm%mu !initial guess for mu
801 mu2 = -1.0e-6_real64
802 dinterv = m_half
804 ! Set pointer to rdm, so that it is available in the functions called by the minimizer
805 rdm_ptr => rdm
806
807 !use n_j=sin^2(2pi*theta_j) to treat pinned states, minimize for both intial mu
808 theta(:) = asin(sqrt(occin(:, 1)/st%smear%el_per_state))*(m_half/m_pi)
809 call minimize_multidim(minmethod_bfgs, st%nst, theta, 0.05_real64, 0.01_real64, &
810 1e-12_real64, 1e-12_real64, 200, objective_rdmft, write_iter_info_rdmft, objective, ierr)
811 sumgi1 = rdm%occsum - st%qtot
812 rdm%mu = mu2
813 theta(:) = asin(sqrt(occin(:, 1)/st%smear%el_per_state))*(m_half/m_pi)
814 call minimize_multidim(minmethod_bfgs, st%nst, theta, 0.05_real64, 0.01_real64, &
815 1e-12_real64, 1e-12_real64, 200, objective_rdmft, write_iter_info_rdmft, objective, ierr)
816 sumgi2 = rdm%occsum - st%qtot
817
818 ! Adjust the interval between the initial mu to include the root of rdm%occsum-st%qtot=M_ZERO
819 do icycle = 1, max_cycle
820 if (sumgi1*sumgi2 <= m_zero) exit
821 if (sumgi2 > m_zero) then
822 mu2 = mu1
823 sumgi2 = sumgi1
824 mu1 = mu1 - dinterv
825 rdm%mu = mu1
826 theta(:) = asin(sqrt(occin(:, 1)/st%smear%el_per_state))*(m_half/m_pi)
827 call minimize_multidim(minmethod_bfgs, st%nst, theta, 0.05_real64, 0.01_real64, &
828 1e-12_real64, 1e-12_real64, 200, objective_rdmft, write_iter_info_rdmft, objective, ierr)
829 sumgi1 = rdm%occsum - st%qtot
830 else
831 mu1 = mu2
832 sumgi1 = sumgi2
833 mu2 = mu2 + dinterv
834 rdm%mu = mu2
835 theta(:) = asin(sqrt(occin(:, 1)/st%smear%el_per_state))*(m_half/m_pi)
836 call minimize_multidim(minmethod_bfgs, st%nst, theta, 0.05_real64, 0.01_real64, &
837 1e-12_real64, 1e-12_real64, 200, objective_rdmft, write_iter_info_rdmft, objective, ierr)
838 sumgi2 = rdm%occsum - st%qtot
839 end if
840 end do
841
842 do icycle = 1, 50
843 mum = (mu1 + mu2)*m_half
844 rdm%mu = mum
845 theta(:) = asin(sqrt(occin(:, 1)/st%smear%el_per_state))*(m_half/m_pi)
846 call minimize_multidim(minmethod_bfgs, st%nst, theta, 0.05_real64, 0.0001_real64, &
847 1e-12_real64, 1e-12_real64, 200, objective_rdmft, write_iter_info_rdmft, objective, ierr)
848 sumgim = rdm%occsum - st%qtot
849
850 if (sumgi1*sumgim < m_zero) then
851 mu2 = mum
852 else
853 mu1 = mum
854 sumgi1 = sumgim
855 end if
856
857 ! check occ.num. threshold again after minimization
858 do ist = 1, st%nst
859 st%occ(ist,1) = m_two*sin(theta(ist)*m_pi*m_two)**2
860 if (st%occ(ist,1) <= thresh_occ ) st%occ(ist,1) = thresh_occ
861 end do
862
863 if (abs(sumgim) < rdm%toler .or. abs((mu1-mu2)*m_half) < rdm%toler) exit
864 end do
865
866 nullify(rdm_ptr)
867
868 if (icycle >= 50) then
869 write(message(1),'(a,1x,f11.4)') 'Bisection ended without finding mu, sum of occupation numbers:', rdm%occsum
870 call messages_fatal(1, namespace=namespace)
871 end if
872
873 do ist = 1, st%nst
874 st%occ(ist, 1) = st%smear%el_per_state*sin(theta(ist)*m_pi*m_two)**2
875 end do
876
877 objective = objective + rdm%mu*(rdm%occsum - rdm%qtot)
878 energy = objective
879
880 write(message(1),'(a4,5x,a12)')'#st','Occupation'
881 call messages_info(1, namespace=namespace)
882
883 do ist = 1, st%nst
884 write(message(1),'(i4,3x,f14.12)') ist, st%occ(ist, 1)
885 call messages_info(1, namespace=namespace)
886 end do
887
888 write(message(1),'(a,3x,f11.9)') 'Sum of occupation numbers: ', rdm%occsum
889 write(message(2),'(a,11x,es20.10)') 'Total energy: ', units_from_atomic(units_out%energy, energy + hm%ep%eii)
890 call messages_info(2, namespace=namespace)
891
892 safe_deallocate_a(occin)
893 safe_deallocate_a(theta)
894
895 call profiling_out("SCF_OCC")
896 pop_sub(scf_occ)
897 end subroutine scf_occ
898
899 ! ---------------------------------------------------------
900 subroutine objective_rdmft(size, theta, objective, getgrad, df)
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)
906
907 integer :: ist
908 real(real64) :: thresh_occ, thresh_theta
909 real(real64), allocatable :: dE_dn(:),occ(:)
910
911 push_sub(objective_rdmft)
912
913 assert(size == rdm_ptr%nst)
914
915 safe_allocate(de_dn(1:size))
916 safe_allocate(occ(1:size))
917
918 occ = m_zero
919
920 ! Defines a threshold on occ nums to avoid numerical instabilities.
921 ! Needs to be changed consistently with the same variable in scf_occ
922 thresh_occ = 1e-14_real64
923 thresh_theta = asin(sqrt(thresh_occ/m_two))*(m_half/m_pi)
924
925 do ist = 1, size
926 occ(ist) = m_two*sin(theta(ist)*m_pi*m_two)**2
927 if (occ(ist) <= thresh_occ ) occ(ist) = thresh_occ
928 end do
929
930 rdm_ptr%occsum = sum(occ(1:size))
931
932 !calculate the total energy without nuclei interaction and the energy
933 !derivatives with respect to the occupation numbers
934
935 call total_energy_rdm(rdm_ptr, occ, objective, de_dn)
936 do ist = 1, size
937 if (occ(ist) <= thresh_occ ) then
938 df(ist) = m_four*m_pi*sin(m_four*thresh_theta*m_pi)*(de_dn(ist) - rdm_ptr%mu)
939 else
940 df(ist) = m_four*m_pi*sin(m_four*theta(ist)*m_pi)*(de_dn(ist) - rdm_ptr%mu)
941 end if
942 end do
943 objective = objective - rdm_ptr%mu*(rdm_ptr%occsum - rdm_ptr%qtot)
944
945 safe_deallocate_a(de_dn)
946 safe_deallocate_a(occ)
947
948 pop_sub(objective_rdmft)
949 end subroutine objective_rdmft
950
951 ! ---------------------------------------------------------
952 subroutine write_iter_info_rdmft(iter, size, energy, maxdr, maxdf, theta)
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)
957
958 push_sub(write_iter_info_rdmft)
959
960 ! Nothing to do.
961
962 pop_sub(write_iter_info_rdmft)
963 end subroutine write_iter_info_rdmft
964
965 ! scf for the natural orbitals
966 subroutine scf_orb(rdm, namespace, gr, st, hm, space, energy)
967 type(rdm_t), intent(inout) :: rdm
968 type(namespace_t), intent(in) :: namespace
969 type(grid_t), intent(in) :: gr
970 type(states_elec_t), intent(inout) :: st
971 type(hamiltonian_elec_t), intent(in) :: hm
972 class(space_t), intent(in) :: space
973 real(real64), intent(out) :: energy
974
975 integer :: ist, jst
976 real(real64), allocatable :: lambda(:, :), fo(:, :)
977
978 push_sub(scf_orb)
979 call profiling_in("SCF_ORB_BASIS")
980
981 !matrix of Lagrange Multipliers from Equation (8), Piris and Ugalde, Vol. 30, No. 13, J. Comput. Chem.
982 safe_allocate(lambda(1:st%nst,1:st%nst))
983 safe_allocate(fo(1:st%nst, 1:st%nst)) !Generalized Fockian Equation (11)
984
985 lambda = m_zero
986 fo = m_zero
987
988 call construct_lambda(namespace, hm, st, gr, lambda, rdm)
989
990 !Set up fo matrix
991 if (rdm%iter==1) then
992 do ist = 1, st%nst
993 do jst = 1, ist
994 fo(ist, jst) = m_half*(lambda(ist, jst) + lambda(jst, ist))
995 fo(jst, ist) = fo(ist, jst)
996 end do
997 end do
998 else
999 do ist = 1, st%nst
1000 do jst = 1, ist - 1
1001 fo(jst, ist) = - ( lambda(jst, ist) - lambda(ist ,jst))
1002 end do
1003 end do
1004 rdm%maxfo = maxval(abs(fo))
1005 do ist = 1, st%nst
1006 fo(ist, ist) = rdm%evalues(ist)
1007 do jst = 1, ist-1
1008 if(abs(fo(jst, ist)) > rdm%scale_f) then
1009 fo(jst, ist) = rdm%scale_f*fo(jst,ist)/abs(fo(jst, ist))
1010 end if
1011 fo(ist, jst) = fo(jst, ist)
1012 end do
1013 end do
1014 end if
1015
1016 call lalg_eigensolve(st%nst, fo, rdm%evalues)
1017 call assign_eigfunctions(rdm, st, fo)
1018 call sum_integrals(rdm) ! to calculate rdm%Coul and rdm%Exch with the new rdm%vecnat
1019 call rdm_derivatives(rdm, namespace, hm, st, gr, space)
1020 call total_energy_rdm(rdm, st%occ(:,1), energy)
1021
1022 safe_deallocate_a(lambda)
1023 safe_deallocate_a(fo)
1024
1025 call profiling_out("SCF_ORB_BASIS")
1026 pop_sub(scf_orb)
1027 end subroutine scf_orb
1028
1029
1030 !-----------------------------------------------------------------
1031 ! Minimize the total energy wrt. an orbital by conjugate gradient
1032 !-----------------------------------------------------------------
1033 subroutine scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
1034 type(rdm_t), intent(inout) :: rdm
1035 type(namespace_t), intent(in) :: namespace
1036 type(electron_space_t), intent(in) :: space
1037 type(grid_t), intent(in) :: gr
1038 type(ions_t), intent(in) :: ions
1039 type(partner_list_t), intent(in) :: ext_partners
1040 type(states_elec_t), intent(inout) :: st
1041 type(v_ks_t), intent(inout) :: ks
1042 type(hamiltonian_elec_t), intent(inout) :: hm
1043 real(real64), intent(out) :: energy
1044
1045 integer :: ik, ist, maxiter
1046
1047
1048 push_sub(scf_orb_cg)
1049 call profiling_in("CG")
1050
1051 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
1052 call hm%update(gr, namespace, space, ext_partners)
1053
1054 rdm%eigens%converged = 0
1055 if(mpi_grp_is_root(mpi_world) .and. .not. debug%info) then
1056 call loct_progress_bar(-1, st%lnst*st%d%kpt%nlocal)
1057 end if
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)
1064
1065 if (.not. rdm%eigens%folded_spectrum) then
1066 ! recheck convergence after subspace diagonalization, since states may have reordered (copied from eigensolver_run)
1067 rdm%eigens%converged(ik) = 0
1068 do ist = 1, st%nst
1069 if(rdm%eigens%diff(ist, ik) < rdm%eigens%tolerance) then
1070 rdm%eigens%converged(ik) = ist
1071 else
1072 exit
1073 end if
1074 end do
1075 end if
1076 end do
1077
1078 if(mpi_grp_is_root(mpi_world) .and. .not. debug%info) then
1079 write(stdout, '(1x)')
1080 end if
1081
1082 ! calculate total energy with new states
1083 call density_calc (st, gr, st%rho)
1084 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
1085 call hm%update(gr, namespace, space, ext_partners)
1086 call rdm_derivatives(rdm, namespace, hm, st, gr, space)
1087
1088 call total_energy_rdm(rdm, st%occ(:,1), energy)
1089
1090 call profiling_out("CG")
1091 pop_sub(scf_orb_cg)
1092 end subroutine scf_orb_cg
1093
1094
1095 ! ----------------------------------------
1096 ! constructs the Lagrange multiplyers needed for the orbital minimization
1097 subroutine construct_lambda(namespace, hm, st, gr, lambda, rdm)
1098 type(namespace_t), intent(in) :: namespace
1099 type(hamiltonian_elec_t), intent(in) :: hm
1100 type(states_elec_t), intent(inout) :: st
1101 type(grid_t), intent(in) :: gr
1102 real(real64), intent(out) :: lambda(:, :)
1103 type(rdm_t), intent(inout) :: rdm
1104
1105 real(real64), allocatable :: hpsi(:, :), hpsi1(:, :), dpsi(:, :), dpsi1(:, :)
1106 real(real64), allocatable :: fock(:,:,:), fvec(:)
1107 integer :: ist, iorb, jorb, jst
1108
1109 push_sub(construct_lambda)
1110
1111 lambda = m_zero
1112
1113 !calculate the Lagrange multiplyer lambda matrix on the grid, Eq. (9), Piris and Ugalde, Vol. 30, No. 13, J. Comput. Chem.
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))
1119
1120 do iorb = 1, st%nst
1121 call states_elec_get_state(st, gr, iorb, 1, dpsi)
1122 call dhamiltonian_elec_apply_single(hm, namespace, gr, dpsi, hpsi, iorb, 1)
1123
1124 do jorb = iorb, st%nst
1125 ! calculate <phi_j|H|phi_i> =lam_ji
1126 call states_elec_get_state(st, gr, jorb, 1, dpsi1)
1127 lambda(jorb, iorb) = dmf_dotp(gr, dpsi1(:,1), hpsi(:,1))
1128
1129 ! calculate <phi_i|H|phi_j>=lam_ij
1130 if (iorb /= jorb ) then
1131 call dhamiltonian_elec_apply_single(hm, namespace, gr, dpsi1, hpsi1, jorb, 1)
1132 lambda(iorb, jorb) = dmf_dotp(gr, dpsi(:,1), hpsi1(:,1))
1133 end if
1134 end do
1135 end do
1136
1137
1138 else ! calculate the same lambda matrix on the basis
1139 !call sum_integrals(rdm)
1140 safe_allocate(fvec(1:st%nst))
1141 safe_allocate(fock(1:st%nst, 1:st%nst, 1:st%nst))
1142 fock = m_zero
1143
1144 do iorb = 1, st%nst
1145 do ist = 1, st%nst
1146 do jst = 1, ist
1147 fock(ist, jst, iorb) = st%occ(iorb, 1)*rdm%eone_int(ist,jst)
1148 do jorb = 1, st%nst
1149 !The coefficient of the Exchange term below is only for the Mueller functional
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)
1152 end do
1153 fock(jst, ist, iorb) = fock(ist, jst, iorb)
1154 end do
1155 end do
1156 end do
1157
1158 do jorb = 1, st%nst
1159 do ist = 1, st%nst
1160 fvec(ist) = m_zero
1161 do jst = 1, st%nst
1162 fvec(ist) = fvec(ist) + fock(ist, jst, jorb)*rdm%vecnat(jst, jorb)
1163 end do
1164 end do
1165 do iorb= 1, st%nst
1166 lambda(iorb, jorb) = m_zero
1167 do ist = 1, st%nst
1168 lambda(iorb, jorb) = lambda(iorb, jorb) + rdm%vecnat(ist, iorb)*fvec(ist)
1169 end do
1170 end do
1171 end do
1172 end if
1173
1174
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)
1180 else
1181 safe_deallocate_a(fvec)
1182 safe_deallocate_a(fock)
1183 end if
1184
1185 pop_sub(construct_lambda)
1186 end subroutine construct_lambda
1187
1188 ! ----------------------------------------
1189
1190 ! finds the new states after the minimization of the orbitals (Piris method)
1191 subroutine assign_eigfunctions(rdm, st, lambda)
1192 type(rdm_t), intent(inout) :: rdm
1193 type(states_elec_t), intent(inout) :: st
1194 real(real64), intent(in) :: lambda(:, :)
1195
1196 integer :: iorb, jorb, ist
1197 real(real64), allocatable :: vecnat_new(:, :)
1198
1199 push_sub(assign_eigenfunctions)
1200
1201 safe_allocate(vecnat_new(1:st%nst, 1:st%nst))
1202 do iorb = 1, st%nst
1203 do ist = 1, st%nst
1204 vecnat_new(ist, iorb) = m_zero
1205 do jorb = 1, st%nst
1206 vecnat_new(ist , iorb) = vecnat_new(ist, iorb) + rdm%vecnat(ist, jorb)*lambda(jorb, iorb)
1207 end do
1208 end do
1209 end do
1210
1211 rdm%vecnat(:, :) = vecnat_new(:, :)
1212
1213 safe_deallocate_a(vecnat_new)
1214
1215 pop_sub(assign_eigenfunctions)
1216 end subroutine assign_eigfunctions
1217
1218 ! --------------------------------------------
1219
1220 ! calculates the total energy when only the occupation numbers are updated
1221 subroutine total_energy_rdm(rdm, occ, energy, dE_dn)
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(:)
1226
1227 integer :: ist, jst
1228 real(real64), allocatable :: V_h(:), V_x(:)
1229
1230 push_sub(total_energy_rdm)
1231
1232 safe_allocate(v_h(1:rdm%nst))
1233 safe_allocate(v_x(1:rdm%nst))
1234
1235 energy = m_zero
1236 v_h = m_zero
1237 v_x = m_zero
1238
1239 !Calculate hartree and exchange contribution
1240 !This is only for the Mueller functional and has to be changed
1241 do ist = 1, rdm%nst
1242 do jst = 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)
1245 end do
1246 v_x(ist) = v_x(ist)*m_half/max(sqrt(occ(ist)), 1.0e-16_real64)
1247 end do
1248
1249
1250 !Calculate the energy derivative with respect to the occupation numbers
1251 if (present(de_dn)) then
1252 de_dn(:) = rdm%eone(:) + v_h(:) + v_x(:)
1253 end if
1254
1255 !Total energy calculation without nuclei interaction
1256 do ist = 1, rdm%nst
1257 energy = energy + occ(ist)*rdm%eone(ist) &
1258 + m_half*occ(ist)*v_h(ist) &
1259 + occ(ist)*v_x(ist)
1260 end do
1261
1262 safe_deallocate_a(v_h)
1263 safe_deallocate_a(v_x)
1264
1265 pop_sub(total_energy_rdm)
1266 end subroutine total_energy_rdm
1267
1268 ! ----------------------------------------
1269 ! calculates the derivatives of the energy terms with respect to the occupation numbers
1270 subroutine rdm_derivatives(rdm, namespace, hm, st, gr, space)
1271 type(rdm_t), intent(inout) :: rdm
1272 type(namespace_t), intent(in) :: namespace
1273 type(hamiltonian_elec_t), intent(in) :: hm
1274 type(states_elec_t), intent(inout) :: st
1275 type(grid_t), intent(in) :: gr
1276 class(space_t), intent(in) :: space
1277
1278
1279 real(real64), allocatable :: hpsi(:, :), rho1(:), rho(:), dpsi(:, :), dpsi2(:, :)
1280 real(real64), allocatable :: v_ij(:,:,:,:,:)
1281 real(real64) :: dd
1282 type(states_elec_t) :: xst
1283
1284 integer :: ist, jst, nspin_, iorb, jorb
1285
1286 push_sub(rdm_derivatives)
1287
1288
1289 nspin_ = min(st%d%nspin, 2)
1290
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))
1298
1299 v_ij = m_zero
1300 rdm%eone = m_zero
1301 rdm%hartree = m_zero
1302 rdm%exchange = m_zero
1303
1304 !derivative of one-electron energy with respect to the natural orbitals occupation number
1305 do ist = 1, st%nst
1306 call states_elec_get_state(st, gr, ist, 1, dpsi)
1307
1308 ! calculate one-body energy
1309 call dhamiltonian_elec_apply_single(hm, namespace, gr, dpsi, hpsi, ist, 1, &
1311 rdm%eone(ist) = dmf_dotp(gr, dpsi(:, 1), hpsi(:, 1))
1312 end do
1313
1314 !integrals used for the hartree and exchange parts of the total energy and their derivatives
1315 ! maybe better to let that be done from the lower level routines like hamiltonian apply?
1316 !
1317 ! only used to calculate total energy
1318 call xst%nullify()
1319 call dexchange_operator_compute_potentials(hm%exxop, namespace, space, gr, st, xst, hm%kpoints, f_out = v_ij)
1320 call states_elec_end(xst)
1321
1322 do ist = 1, st%nst
1323 call states_elec_get_state(st, gr, ist, 1, dpsi)
1324
1325 rho1(1:gr%np) = dpsi(1:gr%np, 1)**2
1326
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)
1330 call states_elec_get_state(st, gr, jst, 1, dpsi2)
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)
1334 end do
1335 end do
1336
1337
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)
1344
1345 else !if energy derivatives are expanded in a basis set
1346
1347 do iorb = 1, st%nst
1348 rdm%eone(iorb) = m_zero
1349 do ist = 1, st%nst
1350 do jst = 1, st%nst
1351 dd = rdm%vecnat(ist, iorb)*rdm%vecnat(jst, iorb)
1352 rdm%eone(iorb) = rdm%eone(iorb) + dd*rdm%eone_int(ist, jst)
1353 end do
1354 end do
1355 end do
1356
1357 do iorb = 1, st%nst
1358 do jorb =1 , iorb
1359 rdm%hartree(iorb ,jorb) = m_zero
1360 rdm%exchange(iorb,jorb) = m_zero
1361 do ist =1, st%nst
1362 do jst =1, st%nst
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
1366 end do
1367 end do
1368 rdm%hartree(jorb, iorb) = rdm%hartree(iorb, jorb)
1369 rdm%exchange(jorb, iorb) = rdm%exchange(iorb, jorb)
1370 end do
1371 end do
1372 end if
1373
1374 pop_sub(rdm_derivatives)
1375 end subroutine rdm_derivatives
1376
1377 ! --------------------------------------------
1378 !calculates the one electron integrals in the basis of the initial orbitals
1379 subroutine rdm_integrals(rdm, namespace, hm, st, mesh)
1380 type(rdm_t), intent(inout) :: rdm
1381 type(namespace_t), intent(in) :: namespace
1382 type(hamiltonian_elec_t), intent(in) :: hm
1383 type(states_elec_t), intent(in) :: st
1384 class(mesh_t), intent(in) :: mesh
1385
1386 real(real64), allocatable :: hpsi(:, :)
1387 real(real64), allocatable :: dpsi(:, :), dpsi2(:, :)
1388 integer :: ist, jst
1389
1390 push_sub(rdm_integrals)
1391
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))
1395
1396 !calculate integrals of the one-electron energy term with respect to the initial orbital basis
1397 do ist = 1, st%nst
1398 call states_elec_get_state(st, mesh, ist, 1, dpsi)
1399 do jst = ist, st%nst
1400 call states_elec_get_state(st, mesh, jst, 1, dpsi2)
1401
1402 ! calculate one-body integrals
1403 call dhamiltonian_elec_apply_single(hm, namespace, mesh, dpsi, hpsi, ist, 1, &
1405 rdm%eone_int(jst, ist) = dmf_dotp(mesh, dpsi2(:, 1), hpsi(:, 1))
1406 rdm%eone_int(ist, jst) = rdm%eone_int(jst, ist)
1407 end do
1408 end do
1409
1410 safe_deallocate_a(hpsi)
1411 safe_deallocate_a(dpsi)
1412 safe_deallocate_a(dpsi2)
1413
1414 pop_sub(rdm_integrals)
1415 end subroutine rdm_integrals
1416
1417 ! --------------------------------------------
1418 ! constructs the Hartree and Exchange part of the RDMFT Fock matrix
1419 subroutine sum_integrals(rdm)
1420 type(rdm_t), intent(inout) :: rdm
1421
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(:,:,:)
1426
1427 push_sub(sum_integrals)
1428
1429 safe_allocate(dm(1:rdm%nst, 1:rdm%nst, 1:rdm%nst))
1430
1431 rdm%Coul = m_zero
1432 rdm%Exch = m_zero
1433 dm = m_zero
1434
1435 do iorb = 1, rdm%nst
1436 do ist = 1, rdm%nst
1437 do jst = 1, ist
1438 dm(ist, jst, iorb) = rdm%vecnat(ist, iorb)*rdm%vecnat(jst, iorb)
1439 dm(jst, ist, iorb) = dm(ist, jst, iorb)
1440 end do
1441 end do
1442 end do
1443
1444 do icount = 1, rdm%n_twoint
1445
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)
1450
1451 two_int = rdm%twoint(icount)
1452
1453 ! create weights of unique integrals
1454 if(ist == jst) then
1455 wij = m_one
1456 else
1457 wij = m_two
1458 endif
1459 if(kst == lst) then
1460 wkl = m_one
1461 else
1462 wkl = m_two
1463 endif
1464
1465 if(ist == kst .and. jst /= lst) then
1466 wik = m_two
1467 else
1468 wik = m_one
1469 endif
1470 if(ist == lst .and. jst /= kst) then
1471 wil = m_two
1472 else
1473 wil = m_one
1474 endif
1475 if(jst == kst .and. ist /= lst) then
1476 wjk = m_two
1477 else
1478 wjk = m_one
1479 endif
1480 if(jst == lst .and. ist /= kst) then
1481 wjl = m_two
1482 else
1483 wjl = m_one
1484 endif
1485
1486 inv_pairs = (ist /= kst .or. jst /= lst)
1487
1488 do iorb = 1, rdm%nst
1489
1490 !the Hartree terms
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
1493
1494 !the exchange terms
1495 !weights are only included if they can differ from one
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
1499 end if
1500 if (ist /= jst) then
1501 if(jst >= kst) then
1502 rdm%Exch(jst, kst, iorb) = rdm%Exch(jst, kst, iorb) + two_int*dm(ist, lst, iorb)*wjk
1503 else
1504 if (inv_pairs) rdm%Exch(kst, jst, iorb) = rdm%Exch(kst, jst, iorb) + two_int*dm(ist, lst, iorb)
1505 end if
1506 end if
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
1510 else
1511 if (inv_pairs) rdm%Exch(lst, jst, iorb) = rdm%Exch(lst, jst, iorb) + two_int*dm(ist, kst, iorb)
1512 end if
1513 end if
1514
1515 end do !iorb
1516 end do !icount
1517
1518 do iorb =1, rdm%nst
1519 do ist = 1, rdm%nst
1520 do jst = 1, ist-1
1521 rdm%Coul(jst, ist, iorb) = rdm%Coul(ist, jst, iorb)
1522 rdm%Exch(jst, ist, iorb) = rdm%Exch(ist, jst, iorb)
1523 end do
1524 end do
1525 end do
1526
1527 safe_deallocate_a(dm)
1528
1529 pop_sub(sum_integrals)
1530 end subroutine sum_integrals
1531
1532end module rdmft_oct_m
1533
1534
1535!! Local Variables:
1536!! mode: f90
1537!! coding: utf-8
1538!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:181
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)
Definition: minimizer.F90:406
type(debug_t), save, public debug
Definition: debug.F90:156
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
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.
Definition: eigen_cg.F90:197
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
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
character(len= *), parameter, public static_dir
Definition: global.F90:252
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public grid_write_info(gr, iunit, namespace)
Definition: grid.F90:522
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.
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
integer minmethod_bfgs
Definition: minimizer.F90:134
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
this module contains the low-level part of the output system
Definition: output_low.F90:115
subroutine, public output_modelmb(outp, namespace, space, dir, gr, ions, iter, st)
this module contains the output system
Definition: output.F90:115
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:494
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.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine scf_occ(rdm, namespace, gr, hm, space, st, energy)
Definition: rdmft.F90:846
subroutine calc_maxfo(namespace, hm, st, gr, rdm)
Definition: rdmft.F90:684
subroutine objective_rdmft(size, theta, objective, getgrad, df)
Definition: rdmft.F90:994
subroutine scf_orb_cg(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, energy)
Definition: rdmft.F90:1127
subroutine, public rdmft_end(rdm)
Definition: rdmft.F90:346
subroutine, public rdmft_init(rdm, namespace, gr, st, mc, space, fromScratch)
Definition: rdmft.F90:212
subroutine write_iter_info_rdmft(iter, size, energy, maxdr, maxdf, theta)
Definition: rdmft.F90:1046
subroutine set_occ_pinning(st)
Definition: rdmft.F90:780
subroutine calc_photon_number(space, gr, st, photons, photon_number_state, ekin_state, epot_state)
Definition: rdmft.F90:719
subroutine assign_eigfunctions(rdm, st, lambda)
Definition: rdmft.F90:1285
subroutine, public scf_rdmft(rdm, namespace, space, gr, ions, ext_partners, st, ks, hm, outp, restart_dump)
Definition: rdmft.F90:376
subroutine construct_lambda(namespace, hm, st, gr, lambda, rdm)
Definition: rdmft.F90:1191
subroutine sum_integrals(rdm)
Definition: rdmft.F90:1513
subroutine rdm_integrals(rdm, namespace, hm, st, mesh)
Definition: rdmft.F90:1473
subroutine scf_orb(rdm, namespace, gr, st, hm, space, energy)
Definition: rdmft.F90:1060
subroutine total_energy_rdm(rdm, occ, energy, dE_dn)
Definition: rdmft.F90:1315
subroutine scf_occ_no(rdm, namespace, gr, hm, space, st, energy)
Definition: rdmft.F90:804
subroutine rdm_derivatives(rdm, namespace, hm, st, gr, space)
Definition: rdmft.F90:1364
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.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
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)
Definition: v_ks.F90:656
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:744
subroutine scf_write_static(dir, fname)
Definition: rdmft.F90:582
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
output handler class
Definition: output_low.F90:164
The states_elec_t class contains all electronic wave functions.
int true(void)