Octopus
scf.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2014 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module scf_oct_m
23 use berry_oct_m
26 use debug_oct_m
34 use epot_oct_m
35 use forces_oct_m
36 use global_oct_m
37 use grid_oct_m
40 use io_oct_m
41 use ions_oct_m
42 use, intrinsic :: iso_fortran_env
46 use lcao_oct_m
47 use lda_u_oct_m
51 use loct_oct_m
53 use math_oct_m
54 use mesh_oct_m
57 use mix_oct_m
59 use mpi_oct_m
62 use output_oct_m
65 use parser_oct_m
69 use smear_oct_m
70 use space_oct_m
75 use stress_oct_m
77 use types_oct_m
78 use unit_oct_m
80 use utils_oct_m
81 use v_ks_oct_m
83 use vdw_ts_oct_m
87 use xc_oct_m
88 use xc_f03_lib_m
91 use xc_oep_oct_m
93
94 implicit none
95
96 private
97 public :: &
98 scf_t, &
99 scf_init, &
101 scf_load, &
102 scf_start, &
103 scf_run, &
104 scf_iter, &
106 scf_finish, &
107 scf_end, &
111
112 integer, public, parameter :: &
113 VERB_NO = 0, &
114 verb_compact = 1, &
115 verb_full = 3
118 type scf_t
119 private
120 integer, public :: max_iter
121
122 real(real64), public :: lmm_r
123
124 ! several convergence criteria
125 logical :: conv_eigen_error
126 logical :: check_conv
127
128 integer :: mix_field
129 logical :: calc_force
130 logical, public :: calc_stress
131 logical :: calc_dipole
132 logical :: calc_partial_charges
133 logical :: calc_orb_moments = .false.
134
135 type(mix_t) :: smix
136 type(mixfield_t), pointer :: mixfield
137 type(eigensolver_t) :: eigens
138 integer :: mixdim1
139 logical :: forced_finish = .false.
140 type(lda_u_mixer_t) :: lda_u_mix
141 type(vtau_mixer_t) :: vtau_mix
142 type(berry_t) :: berry
143 integer :: matvec
144
145 type(restart_t), public :: restart_load, restart_dump
146
147 type(criterion_list_t), public :: criterion_list
148 real(real64) :: energy_in, energy_diff, abs_dens_diff, evsum_in, evsum_out, evsum_diff
149
150 ! Variables needed to store information accross scf_start, scf_run, and scf_finish
151 logical :: converged_current, converged_last
152 integer :: verbosity_
153 type(lcao_t) :: lcao
154 real(real64), allocatable :: rhoout(:,:), rhoin(:,:)
155 real(real64), allocatable :: vhxc_old(:,:)
156 class(wfs_elec_t), allocatable :: psioutb(:, :)
157 logical :: output_forces, calc_current, output_during_scf
158 logical :: finish = .false.
159 end type scf_t
160
161contains
162
163 ! ---------------------------------------------------------
164 subroutine scf_init(scf, namespace, gr, ions, st, mc, hm, space)
165 type(scf_t), intent(inout) :: scf
166 type(grid_t), intent(in) :: gr
167 type(namespace_t), intent(in) :: namespace
168 type(ions_t), intent(in) :: ions
169 type(states_elec_t), intent(in) :: st
170 type(multicomm_t), intent(in) :: mc
171 type(hamiltonian_elec_t), intent(inout) :: hm
172 class(space_t), intent(in) :: space
173
174 real(real64) :: rmin
175 integer :: mixdefault
176 type(type_t) :: mix_type
177 class(convergence_criterion_t), pointer :: crit
178 type(criterion_iterator_t) :: iter
179 logical :: deactivate_oracle
180
181 push_sub(scf_init)
182
183 !%Variable MaximumIter
184 !%Type integer
185 !%Default 200
186 !%Section SCF::Convergence
187 !%Description
188 !% Maximum number of SCF iterations. The code will stop even if convergence
189 !% has not been achieved. -1 means unlimited.
190 !% 0 means just do LCAO (or read from restart), compute the eigenvalues and energy,
191 !% and stop, without updating the wavefunctions or density.
192 !%
193 !% If convergence criteria are set, the SCF loop will only stop once the criteria
194 !% are fulfilled for two consecutive iterations.
195 !%
196 !% Note that this variable is also used in the section Calculation Modes::Unoccupied States,
197 !% where it denotes the maximum number of calls of the eigensolver. In this context, the
198 !% default value is 50.
199 !%End
200 call parse_variable(namespace, 'MaximumIter', 200, scf%max_iter)
201
202 if (allocated(hm%vberry)) then
203 call berry_init(scf%berry, namespace)
204 end if
205
206 !Create the list of convergence criteria
207 call criteria_factory_init(scf%criterion_list, namespace, scf%check_conv)
208 !Setting the pointers
209 call iter%start(scf%criterion_list)
210 do while (iter%has_next())
211 crit => iter%get_next()
212 select type (crit)
214 call crit%set_pointers(scf%energy_diff, scf%energy_in)
216 call crit%set_pointers(scf%abs_dens_diff, st%qtot)
218 call crit%set_pointers(scf%evsum_diff, scf%evsum_out)
219 class default
220 assert(.false.)
221 end select
222 end do
225 if(.not. scf%check_conv .and. scf%max_iter < 0) then
226 call messages_write("All convergence criteria are disabled. Octopus is cowardly refusing")
228 call messages_write("to enter an infinite loop.")
229 call messages_new_line()
231 call messages_write("Please set one of the following variables to a positive value:")
234 call messages_write(" | MaximumIter | ConvEnergy | ConvAbsDens | ConvRelDens |")
236 call messages_write(" | ConvAbsEv | ConvRelEv |")
238 call messages_fatal(namespace=namespace)
239 end if
241 !%Variable ConvEigenError
242 !%Type logical
243 !%Default false
244 !%Section SCF::Convergence
245 !%Description
246 !% If true, the calculation will not be considered converged unless all states have
247 !% individual errors less than <tt>EigensolverTolerance</tt>.
248 !% If <tt>ExtraStatesToConverge</tt> is set, the calculation will stop
249 !% when all occupied states plus <tt>ExtraStatesToConverge</tt> extra states are converged.
250 !%
251 !% If this criterion is used, the SCF loop will only stop once it is
252 !% fulfilled for two consecutive iterations.
253 !%End
254 call parse_variable(namespace, 'ConvEigenError', .false., scf%conv_eigen_error)
255
256 if(scf%max_iter < 0) scf%max_iter = huge(scf%max_iter)
257
258 call messages_obsolete_variable(namespace, 'What2Mix', 'MixField')
260 ! now the eigensolver stuff
261 deactivate_oracle = hm%theory_level == independent_particles
262 call eigensolver_init(scf%eigens, namespace, gr, st, hm, mc, space, deactivate_oracle)
263
264 if(scf%eigens%es_type /= rs_evo) then
265 !%Variable MixField
266 !%Type integer
267 !%Section SCF::Mixing
268 !%Description
269 !% Selects what should be mixed during the SCF cycle. Note that
270 !% currently the exact-exchange part of hybrid functionals is not
271 !% mixed at all, which would require wavefunction-mixing, not yet
272 !% implemented. This may lead to instabilities in the SCF cycle,
273 !% so starting from a converged LDA/GGA calculation is recommended
274 !% for hybrid functionals. The default depends on the <tt>TheoryLevel</tt>
275 !% and the exchange-correlation potential used.
276 !% This is not used in case of imaginary-time evolution.
277 !%Option none 0
278 !% No mixing is done. This is the default for independent
279 !% particles.
280 !%Option potential 1
281 !% The Kohn-Sham potential is mixed. This is the default for other cases.
282 !%Option density 2
283 !% Mix the density.
284 !%Option states 3
285 !% (Experimental) Mix the states. In this case, the mixing is always linear.
286 !%End
287
288 mixdefault = option__mixfield__potential
289 if(hm%theory_level == independent_particles) mixdefault = option__mixfield__none
290
291 call parse_variable(namespace, 'MixField', mixdefault, scf%mix_field)
292 if(.not.varinfo_valid_option('MixField', scf%mix_field)) call messages_input_error(namespace, 'MixField')
293 call messages_print_var_option('MixField', scf%mix_field, "what to mix during SCF cycles", namespace=namespace)
294
295 if (scf%mix_field == option__mixfield__potential .and. hm%theory_level == independent_particles) then
296 call messages_write('Input: Cannot mix the potential for non-interacting particles.')
297 call messages_fatal(namespace=namespace)
298 end if
299
300 if (scf%mix_field == option__mixfield__potential .and. hm%pcm%run_pcm) then
301 call messages_write('Input: You have selected to mix the potential.', new_line = .true.)
302 call messages_write(' This might produce convergence problems for solvated systems.', new_line = .true.)
303 call messages_write(' Mix the Density instead.')
304 call messages_warning(namespace=namespace)
305 end if
306
307 if(scf%mix_field == option__mixfield__density &
308 .and. bitand(hm%xc%family, xc_family_oep + xc_family_mgga + xc_family_hyb_mgga + xc_family_nc_mgga) /= 0) then
309
310 call messages_write('Input: You have selected to mix the density with OEP or MGGA XC functionals.', new_line = .true.)
311 call messages_write(' This might produce convergence problems. Mix the potential instead.')
312 call messages_warning(namespace=namespace)
313 end if
314
315 if(scf%mix_field == option__mixfield__states) then
316 call messages_experimental('MixField = states', namespace=namespace)
317 end if
318
319 ! Handle mixing now...
320 select case(scf%mix_field)
321 case (option__mixfield__potential, option__mixfield__density)
322 scf%mixdim1 = gr%np
323 case(option__mixfield__states)
324 ! we do not really need the mixer, except for the value of the mixing coefficient
325 scf%mixdim1 = 1
326 end select
327
328 mix_type = type_float
329
330 if (scf%mix_field /= option__mixfield__none) then
331 call mix_init(scf%smix, namespace, space, gr%der, scf%mixdim1, st%d%nspin, func_type_ = mix_type)
332 end if
333
334 ! If we use DFT+U, we also have do mix it
335 if (scf%mix_field /= option__mixfield__states .and. scf%mix_field /= option__mixfield__none ) then
336 call lda_u_mixer_init(hm%lda_u, scf%lda_u_mix, st)
337 call lda_u_mixer_init_auxmixer(hm%lda_u, namespace, scf%lda_u_mix, scf%smix, st)
338 end if
339
340 ! If we use tau-dependent MGGA, we need to mix vtau
341 if(scf%mix_field == option__mixfield__potential) then
342 call vtau_mixer_init_auxmixer(namespace, scf%vtau_mix, scf%smix, hm, gr%np, st%d%nspin)
343 end if
344
345 call mix_get_field(scf%smix, scf%mixfield)
346 else
347 scf%mix_field = option__mixfield__none
348 end if
349
350 !%Variable SCFCalculateForces
351 !%Type logical
352 !%Section SCF
353 !%Description
354 !% This variable controls whether the forces on the ions are
355 !% calculated at the end of a self-consistent iteration. The
356 !% default is yes, unless the system only has user-defined
357 !% species.
358 !%End
359 call parse_variable(namespace, 'SCFCalculateForces', .not. ions%only_user_def, scf%calc_force)
360
361 if(scf%calc_force .and. gr%der%boundaries%spiralBC) then
362 message(1) = 'Forces cannot be calculated when using spiral boundary conditions.'
363 write(message(2),'(a)') 'Please use SCFCalculateForces = no.'
364 call messages_fatal(2, namespace=namespace)
365 end if
366 if(scf%calc_force) then
367 if (allocated(hm%ep%b_field) .or. allocated(hm%ep%a_static)) then
368 write(message(1),'(a)') 'The forces are currently not properly calculated if static'
369 write(message(2),'(a)') 'magnetic fields or static vector potentials are present.'
370 write(message(3),'(a)') 'Please use SCFCalculateForces = no.'
371 call messages_fatal(3, namespace=namespace)
372 end if
373 end if
374
375 !%Variable SCFCalculateStress
376 !%Type logical
377 !%Default no
378 !%Section SCF
379 !%Description
380 !% This variable controls whether the stress on the lattice is
381 !% calculated at the end of a self-consistent iteration. The
382 !% default is no.
383 !%End
384 call parse_variable(namespace, 'SCFCalculateStress', .false. , scf%calc_stress)
385
386 !%Variable SCFCalculateDipole
387 !%Type logical
388 !%Section SCF
389 !%Description
390 !% This variable controls whether the dipole is calculated at the
391 !% end of a self-consistent iteration. For finite systems the
392 !% default is yes. For periodic systems the default is no, unless
393 !% an electric field is being applied in a periodic direction.
394 !% The single-point Berry`s phase approximation is used for
395 !% periodic directions. Ref:
396 !% E Yaschenko, L Fu, L Resca, and R Resta, <i>Phys. Rev. B</i> <b>58</b>, 1222-1229 (1998).
397 !%End
398 call parse_variable(namespace, 'SCFCalculateDipole', .not. space%is_periodic(), scf%calc_dipole)
399 if (allocated(hm%vberry)) scf%calc_dipole = .true.
400
401 !%Variable SCFCalculateOrbitalMoments
402 !%Type logical
403 !%Default no
404 !%Section SCF
405 !%Description
406 !% This variable controls whether the local orbital angular moments are
407 !% calculated at the end of a self-consistent iteration. The
408 !% default is no. This is only applicable for spinors with SOC.
409 !%
410 !% This is computed by integrating around atom-centered spheres, and does not include the
411 !% interstitial contribution.
412 !% The same sphere are used as for the spin magnetic moments.
413 !%End
414 call parse_variable(namespace, 'SCFCalculateOrbitalMoments', .false. , scf%calc_orb_moments)
415 if((st%d%ispin /= spinors .or. space%dim /= 3) .and. scf%calc_orb_moments) then
416 message(1) = "Orbital moments are only implemented for spinors and in 3D."
417 call messages_fatal(1, namespace=namespace)
418 end if
419 if (scf%calc_orb_moments .and. .not. (hm%ep%reltype == spin_orbit &
420 .or. hm%ep%reltype == fully_relativistic_zora)) then
421 message(1) = "Orbital moments are only available with SOC."
422 call messages_fatal(1, namespace=namespace)
423 end if
424 if(gr%use_curvilinear .and. scf%calc_orb_moments) then
425 call messages_not_implemented("Orbital angular moments with curvilinear coordinates")
426 end if
427
428 !%Variable SCFCalculatePartialCharges
429 !%Type logical
430 !%Default no
431 !%Section SCF
432 !%Description
433 !% (Experimental) This variable controls whether partial charges
434 !% are calculated at the end of a self-consistent iteration.
435 !%End
436 call parse_variable(namespace, 'SCFCalculatePartialCharges', .false., scf%calc_partial_charges)
437 if (scf%calc_partial_charges) call messages_experimental('SCFCalculatePartialCharges', namespace=namespace)
438
439 rmin = ions%min_distance()
440
441 !%Variable LocalMagneticMomentsSphereRadius
442 !%Type float
443 !%Section Output
444 !%Description
445 !% The local magnetic moments are calculated by integrating the
446 !% magnetization density in spheres centered around each atom.
447 !% This variable controls the radius of the spheres.
448 !% The default is half the minimum distance between two atoms
449 !% in the input coordinates, or 100 a.u. if there is only one atom (for isolated systems).
450 !%End
451 call parse_variable(namespace, 'LocalMagneticMomentsSphereRadius', min(m_half*rmin, lmm_r_single_atom), scf%lmm_r, &
452 unit=units_inp%length)
453 ! this variable is also used in td/td_write.F90
454
455 scf%forced_finish = .false.
456
457 pop_sub(scf_init)
458 end subroutine scf_init
459
460
461 ! ---------------------------------------------------------
462 subroutine scf_end(scf)
463 type(scf_t), intent(inout) :: scf
464
465 class(convergence_criterion_t), pointer :: crit
466 type(criterion_iterator_t) :: iter
467
468 push_sub(scf_end)
469
470 call eigensolver_end(scf%eigens)
471
472 if(scf%mix_field /= option__mixfield__none) call mix_end(scf%smix)
473
474 nullify(scf%mixfield)
475
476 if (scf%mix_field /= option__mixfield__states) then
477 call lda_u_mixer_end(scf%lda_u_mix, scf%smix)
478 call vtau_mixer_end(scf%vtau_mix, scf%smix)
479 end if
480
481 call iter%start(scf%criterion_list)
482 do while (iter%has_next())
483 crit => iter%get_next()
484 safe_deallocate_p(crit)
485 end do
486
487 pop_sub(scf_end)
488 end subroutine scf_end
489
490
491 ! ---------------------------------------------------------
492 subroutine scf_mix_clear(scf)
493 type(scf_t), intent(inout) :: scf
494
495 push_sub(scf_mix_clear)
496
497 call mix_clear(scf%smix)
498
499 if (scf%mix_field /= option__mixfield__states) then
500 call lda_u_mixer_clear(scf%lda_u_mix, scf%smix)
501 call vtau_mixer_clear(scf%vtau_mix, scf%smix)
502 end if
503
504 pop_sub(scf_mix_clear)
505 end subroutine scf_mix_clear
506
507 ! ---------------------------------------------------------
509 subroutine scf_load(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, restart_load)
510 type(scf_t), intent(inout) :: scf
511 type(namespace_t), intent(in) :: namespace
512 type(electron_space_t), intent(in) :: space
513 type(grid_t), intent(inout) :: gr
514 type(ions_t), intent(in) :: ions
515 type(partner_list_t), intent(in) :: ext_partners
516 type(states_elec_t), intent(inout) :: st
517 type(v_ks_t), intent(inout) :: ks
518 type(hamiltonian_elec_t), intent(inout) :: hm
519 type(restart_t), intent(in) :: restart_load
520
521 integer :: ierr, is
522
523 push_sub(scf_load)
524
525 if (restart_load%has_flag(restart_flag_rho)) then
526 ! Load density and used it to recalculated the KS potential.
527 call states_elec_load_rho(restart_load, st, gr, ierr)
528 if (ierr /= 0) then
529 message(1) = 'Unable to read density. Density will be calculated from states.'
530 call messages_warning(1, namespace=namespace)
531 else
532 if (bitand(ks%xc_family, xc_family_oep) == 0) then
533 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
534 else
535 if (.not. restart_load%has_flag(restart_flag_vhxc) .and. ks%oep%level /= oep_level_full) then
536 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
537 end if
538 end if
539 end if
540 end if
541
542 if (restart_load%has_flag(restart_flag_vhxc)) then
543 call hm%ks_pot%load(restart_load, gr, ierr)
544 if (ierr /= 0) then
545 message(1) = 'Unable to read Vhxc. Vhxc will be calculated from states.'
546 call messages_warning(1, namespace=namespace)
547 else
548 call hm%update(gr, namespace, space, ext_partners)
549 if (bitand(ks%xc_family, xc_family_oep) /= 0) then
550 if (ks%oep%level == oep_level_full) then
551 do is = 1, st%d%nspin
552 ks%oep%vxc(1:gr%np, is) = hm%ks_pot%vhxc(1:gr%np, is) - hm%ks_pot%vhartree(1:gr%np)
553 end do
554 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
555 end if
556 end if
557 end if
558 end if
559
560 if (restart_load%has_flag(restart_flag_mix)) then
561 if (scf%mix_field == option__mixfield__density .or. scf%mix_field == option__mixfield__potential) then
562 call mix_load(namespace, restart_load, scf%smix, gr, ierr)
563 end if
564 if (ierr /= 0) then
565 message(1) = "Unable to read mixing information. Mixing will start from scratch."
566 call messages_warning(1, namespace=namespace)
567 end if
568 end if
569
570 if(hm%lda_u_level /= dft_u_none) then
571 call lda_u_load(restart_load, hm%lda_u, st, hm%energy%dft_u, ierr)
572 if (ierr /= 0) then
573 message(1) = "Unable to read DFT+U information. DFT+U data will be calculated from states."
574 call messages_warning(1, namespace=namespace)
575 end if
576
577 ! As v_ks_calc has already been called, we need to update hm%energy%int_dft_u
578 call v_ks_update_dftu_energy(ks, namespace, hm, st, hm%energy%int_dft_u)
579 end if
580
581 !TODO: Create a dedicated routine and call it from the initialize
582
583! if (present(outp) .and. st%system_grp%is_root()) then
584! call io_rm(STATIC_DIR //"info")
585! end if
586! end if
588 pop_sub(scf_load)
589 end subroutine scf_load
590
591 ! ---------------------------------------------------------
593 subroutine scf_start(scf, namespace, gr, ions, st, ks, hm, outp, verbosity)
594 type(scf_t), intent(inout) :: scf
595 type(namespace_t), intent(in) :: namespace
596 type(grid_t), intent(inout) :: gr
597 type(ions_t), intent(inout) :: ions
598 type(states_elec_t), intent(inout) :: st
599 type(v_ks_t), intent(inout) :: ks
600 type(hamiltonian_elec_t), intent(inout) :: hm
601 type(output_t), optional, intent(in) :: outp
602 integer, optional, intent(in) :: verbosity
603
604 integer :: ib, iqn
605
606 push_sub(scf_start)
607
608 if(scf%forced_finish) then
609 message(1) = "Previous clean stop, not doing SCF and quitting."
610 call messages_fatal(1, only_root_writes = .true., namespace=namespace)
611 end if
612
613 if (.not. hm%is_hermitian()) then
614 message(1) = "Trying to run a SCF calculation for a non-hermitian Hamiltonian. This is not supported."
615 call messages_fatal(1, namespace=namespace)
616 end if
617
618 scf%verbosity_ = optional_default(verbosity, verb_full)
619
620 scf%output_during_scf = .false.
621 scf%output_forces = .false.
622 scf%calc_current = .false.
623
624 if (present(outp)) then
625 ! if the user has activated output=stress but not SCFCalculateStress,
626 ! we assume that is implied
627 if (outp%what(option__output__stress)) then
628 scf%calc_stress = .true.
629 end if
630
631 scf%output_during_scf = outp%duringscf
632 scf%calc_current = output_needs_current(outp, states_are_real(st))
633
634 if (outp%duringscf .and. outp%what(option__output__forces)) then
635 scf%output_forces = .true.
636 end if
637 end if
638
639 safe_allocate(scf%rhoout(1:gr%np, 1:st%d%nspin))
640 safe_allocate(scf%rhoin (1:gr%np, 1:st%d%nspin))
641
642 call lalg_copy(gr%np, st%d%nspin, st%rho, scf%rhoin)
643 scf%rhoout = m_zero
644
645 if (scf%calc_force .or. scf%output_forces) then
646 !We store the Hxc potential for the contribution to the forces
647 safe_allocate(scf%vhxc_old(1:gr%np, 1:st%d%nspin))
648 call lalg_copy(gr%np, st%d%nspin, hm%ks_pot%vhxc, scf%vhxc_old)
649 end if
650
651
652 select case(scf%mix_field)
653 case(option__mixfield__potential)
654 call mixfield_set_vin(scf%mixfield, hm%ks_pot%vhxc)
655 call vtau_mixer_set_vin(scf%vtau_mix, hm)
656 case(option__mixfield__density)
657 call mixfield_set_vin(scf%mixfield, scf%rhoin)
658
659 case(option__mixfield__states)
660
661 ! There is a ICE with foss2022a-serial. I am changing to allocate - NTD
662 allocate(wfs_elec_t::scf%psioutb (st%group%block_start:st%group%block_end, st%d%kpt%start:st%d%kpt%end))
663
664 do iqn = st%d%kpt%start, st%d%kpt%end
665 do ib = st%group%block_start, st%group%block_end
666 call st%group%psib(ib, iqn)%copy_to(scf%psioutb(ib, iqn))
667 end do
668 end do
669
670 end select
671
672 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
673 ! If we use DFT+U, we also have do mix it
674 if (scf%mix_field /= option__mixfield__states) call lda_u_mixer_set_vin(hm%lda_u, scf%lda_u_mix)
675
676 call create_convergence_file(static_dir, "convergence")
677
678 if ( scf%verbosity_ /= verb_no ) then
679 if(scf%max_iter > 0) then
680 write(message(1),'(a)') 'Info: Starting SCF iteration.'
681 else
682 write(message(1),'(a)') 'Info: No SCF iterations will be done.'
683 ! we cannot tell whether it is converged.
684 scf%finish = .false.
685 end if
686 call messages_info(1, namespace=namespace)
687 end if
689 scf%converged_current = .false.
690 scf%matvec = 0
691
692 pop_sub(scf_start)
693
694 contains
695
696 ! -----------------------------------------------------
697
698 subroutine create_convergence_file(dir, fname)
699 character(len=*), intent(in) :: dir
700 character(len=*), intent(in) :: fname
701
702 integer :: iunit
703 character(len=12) :: label
704 if(st%system_grp%is_root()) then
705 call io_mkdir(dir, namespace)
706 iunit = io_open(trim(dir) // "/" // trim(fname), namespace, action='write')
707 write(iunit, '(a)', advance = 'no') '#iter energy '
708 label = 'energy_diff'
709 write(iunit, '(1x,a)', advance = 'no') label
710 label = 'abs_dens'
711 write(iunit, '(1x,a)', advance = 'no') label
712 label = 'rel_dens'
713 write(iunit, '(1x,a)', advance = 'no') label
714 label = 'abs_ev'
715 write(iunit, '(1x,a)', advance = 'no') label
716 label = 'rel_ev'
717 write(iunit, '(1x,a)', advance = 'no') label
718 if (bitand(ks%xc_family, xc_family_oep) /= 0 .and. ks%theory_level /= hartree_fock &
719 .and. ks%theory_level /= generalized_kohn_sham_dft) then
720 if (ks%oep%level == oep_level_full) then
721 label = 'OEP norm2ss'
722 write(iunit, '(1x,a)', advance = 'no') label
723 end if
724 end if
725 write(iunit,'(a)') ''
726 call io_close(iunit)
727 end if
728
729 end subroutine create_convergence_file
730
731 end subroutine scf_start
732
733 ! ---------------------------------------------------------
735 subroutine scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, &
736 verbosity, iters_done, restart_dump)
737 type(scf_t), intent(inout) :: scf
738 type(namespace_t), intent(in) :: namespace
739 type(electron_space_t), intent(in) :: space
740 type(multicomm_t), intent(in) :: mc
741 type(grid_t), intent(inout) :: gr
742 type(ions_t), intent(inout) :: ions
743 type(partner_list_t), intent(in) :: ext_partners
744 type(states_elec_t), intent(inout) :: st
745 type(v_ks_t), intent(inout) :: ks
746 type(hamiltonian_elec_t), intent(inout) :: hm
747 type(output_t), optional, intent(in) :: outp
748 integer, optional, intent(in) :: verbosity
749 integer, optional, intent(out) :: iters_done
750 type(restart_t), optional, intent(in) :: restart_dump
751
752 integer :: iter
753 logical :: completed
754
755 push_sub(scf_run)
756
757 call scf_start(scf, namespace, gr, ions, st, ks, hm, outp, verbosity)
758
759 ! SCF cycle
760 do iter = 1, scf%max_iter
761
762 call scf_iter(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, iter, outp, &
763 restart_dump)
764
765 completed = scf_iter_finish(scf, namespace, space, gr, ions, st, ks, hm, iter, outp, iters_done)
766
767 if(scf%forced_finish .or. completed) then
768 exit
769 end if
770 end do
771
772 if (.not.scf%forced_finish) then
773 ! this is only executed if the computation has converged
774 call scf_finish(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, iter, outp)
775 end if
776
777 pop_sub(scf_run)
778 end subroutine scf_run
779
780 ! ---------------------------------------------------------
781 subroutine scf_iter(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, iter, outp, &
782 restart_dump)
783 type(scf_t), intent(inout) :: scf
784 type(namespace_t), intent(in) :: namespace
785 type(electron_space_t), intent(in) :: space
786 type(multicomm_t), intent(in) :: mc
787 type(grid_t), intent(inout) :: gr
788 type(ions_t), intent(inout) :: ions
789 type(partner_list_t), intent(in) :: ext_partners
790 type(states_elec_t), intent(inout) :: st
791 type(v_ks_t), intent(inout) :: ks
792 type(hamiltonian_elec_t), intent(inout) :: hm
793 integer, intent(in) :: iter
794 type(output_t), optional, intent(in) :: outp
795 type(restart_t), optional, intent(in) :: restart_dump
796
797 integer :: iqn, ib, ierr
798 class(convergence_criterion_t), pointer :: crit
799 type(criterion_iterator_t) :: iterator
800 logical :: is_crit_conv
801 real(real64) :: etime, itime
802
803 push_sub(scf_iter)
804
805 call profiling_in("SCF_CYCLE")
806 itime = loct_clock()
807
808 ! this initialization seems redundant but avoids improper optimization at -O3 by PGI 7 on chum,
809 ! which would cause a failure of testsuite/linear_response/04-vib_modes.03-vib_modes_fd.inp
810 scf%eigens%converged = 0
811
812 ! keep the information about the spectrum up to date, needed e.g. for Chebyshev expansion for imaginary time
813 call hm%update_span(gr%spacing(1:space%dim), minval(st%eigenval(:, :)), namespace)
814
815 !We update the quantities at the begining of the scf cycle
816 if (iter == 1) then
817 scf%evsum_in = states_elec_eigenvalues_sum(st)
818 end if
819 call iterator%start(scf%criterion_list)
820 do while (iterator%has_next())
821 crit => iterator%get_next()
822 call scf_update_initial_quantity(scf, hm, crit)
823 end do
824
825 if (scf%calc_force .or. scf%output_forces) then
826 !Used for computing the imperfect convegence contribution to the forces
827 scf%vhxc_old(1:gr%np, 1:st%d%nspin) = hm%ks_pot%vhxc(1:gr%np, 1:st%d%nspin)
828 end if
829
830 !We check if the system is coupled with a partner that requires self-consistency
831 ! if(hamiltonian_has_scf_partner(hm)) then
832 if (allocated(hm%vberry)) then
833 !In this case, v_Hxc is frozen and we do an internal SCF loop over the
834 ! partners that require SCF
835 ks%frozen_hxc = .true.
836 ! call perform_scf_partners()
837 call berry_perform_internal_scf(scf%berry, namespace, space, scf%eigens, gr, st, hm, iter, ks, ions, ext_partners)
838 !and we unfreeze the potential once finished
839 ks%frozen_hxc = .false.
840 else
841 scf%eigens%converged = 0
842 call scf%eigens%run(namespace, gr, st, hm, space, ext_partners, iter)
843 end if
844
845 scf%matvec = scf%matvec + scf%eigens%matvec
846
847 ! occupations
848 call states_elec_fermi(st, namespace, gr)
849 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
850
851 ! compute output density, potential (if needed) and eigenvalues sum
852 call density_calc(st, gr, st%rho)
853
854 call lalg_copy(gr%np, st%d%nspin, st%rho, scf%rhoout)
855
856 select case (scf%mix_field)
857 case (option__mixfield__potential)
858 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_current=scf%output_during_scf)
859 call mixfield_set_vout(scf%mixfield, hm%ks_pot%vhxc)
860 call vtau_mixer_set_vout(scf%vtau_mix, hm)
861 case (option__mixfield__density)
862 call mixfield_set_vout(scf%mixfield, scf%rhoout)
863 case(option__mixfield__states)
864 do iqn = st%d%kpt%start, st%d%kpt%end
865 do ib = st%group%block_start, st%group%block_end
866 call st%group%psib(ib, iqn)%copy_data_to(gr%np, scf%psioutb(ib, iqn))
867 end do
868 end do
869 end select
870
871 if (scf%mix_field /= option__mixfield__states .and. scf%mix_field /= option__mixfield__none) then
872 call lda_u_mixer_set_vout(hm%lda_u, scf%lda_u_mix)
873 endif
874
875 ! recalculate total energy
876 call energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit = -1)
877
878 if (present(outp)) then
879 ! compute forces only if requested
880 if (outp%duringscf .and. outp%what_now(option__output__forces, iter)) then
881 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old=scf%vhxc_old)
882 end if
883 end if
884
885 !We update the quantities at the end of the scf cycle
886 call iterator%start(scf%criterion_list)
887 do while (iterator%has_next())
888 crit => iterator%get_next()
889 call scf_update_diff_quantity(scf, hm, st, gr, scf%rhoout, scf%rhoin, crit)
890 end do
891
892 ! are we finished?
893 scf%converged_last = scf%converged_current
894
895 scf%converged_current = scf%check_conv .and. &
896 (.not. scf%conv_eigen_error .or. all(scf%eigens%converged >= st%nst_conv))
897 !Loop over the different criteria
898 call iterator%start(scf%criterion_list)
899 do while (iterator%has_next())
900 crit => iterator%get_next()
901 call crit%is_converged(is_crit_conv)
902 scf%converged_current = scf%converged_current .and. is_crit_conv
903 end do
904
905 ! only finish if the convergence criteria are fulfilled in two
906 ! consecutive iterations
907 scf%finish = scf%converged_last .and. scf%converged_current
908
909 etime = loct_clock() - itime
910 call scf_write_iter(namespace)
911
912 ! mixing
913 select case (scf%mix_field)
914 case (option__mixfield__density)
915 ! mix input and output densities and compute new potential
916 call mixing(namespace, scf%smix)
917 call mixfield_get_vnew(scf%mixfield, st%rho)
918 ! for spinors, having components 3 or 4 be negative is not unphysical
919 if (minval(st%rho(1:gr%np, 1:st%d%spin_channels)) < -1e-6_real64) then
920 write(message(1),*) 'Negative density after mixing. Minimum value = ', &
921 minval(st%rho(1:gr%np, 1:st%d%spin_channels))
922 call messages_warning(1, namespace=namespace)
923 end if
924 call lda_u_mixer_get_vnew(hm%lda_u, scf%lda_u_mix, st)
925 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_current=scf%output_during_scf)
926
927 case (option__mixfield__potential)
928 ! mix input and output potentials
929 call mixing(namespace, scf%smix)
930 call mixfield_get_vnew(scf%mixfield, hm%ks_pot%vhxc)
931 call lda_u_mixer_get_vnew(hm%lda_u, scf%lda_u_mix, st)
932 call vtau_mixer_get_vnew(scf%vtau_mix, hm)
933 call hamiltonian_elec_update_pot(hm, gr)
934
935 case(option__mixfield__states)
936 do iqn = st%d%kpt%start, st%d%kpt%end
937 do ib = st%group%block_start, st%group%block_end
938 call batch_axpby(gr%np, mix_coefficient(scf%smix), scf%psioutb(ib, iqn), &
939 m_one - mix_coefficient(scf%smix), st%group%psib(ib, iqn))
940 end do
941 end do
942 call density_calc(st, gr, st%rho)
943 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_current=scf%output_during_scf)
944
945 case (option__mixfield__none)
946 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_current=scf%output_during_scf)
947 end select
948
949 ! Are we asked to stop? (Whenever Fortran is ready for signals, this should go away)
950 scf%forced_finish = clean_stop(mc%master_comm) .or. walltimer_alarm(mc%master_comm)
951
952 if (scf%finish .and. st%modelmbparticles%nparticle > 0) then
953 call modelmb_sym_all_states(space, gr, st)
954 end if
955
956 if (present(outp) .and. present(restart_dump)) then
957 ! save restart information
958
959 if ( (scf%finish .or. restart_walltime_period_alarm(mc%master_comm) &
960 .or. iter == scf%max_iter .or. scf%forced_finish) ) then
961
962 call states_elec_dump(scf%restart_dump, space, st, gr, hm%kpoints, ierr, iter=iter)
963 if (ierr /= 0) then
964 message(1) = 'Unable to write states wavefunctions.'
965 call messages_warning(1, namespace=namespace)
966 end if
967
968 call states_elec_dump_rho(scf%restart_dump, st, gr, ierr, iter=iter)
969 if (ierr /= 0) then
970 message(1) = 'Unable to write density.'
971 call messages_warning(1, namespace=namespace)
972 end if
973
974 if(hm%lda_u_level /= dft_u_none) then
975 call lda_u_dump(scf%restart_dump, namespace, hm%lda_u, st, gr, ierr)
976 if (ierr /= 0) then
977 message(1) = 'Unable to write DFT+U information.'
978 call messages_warning(1, namespace=namespace)
979 end if
980 end if
981
982 select case (scf%mix_field)
983 case (option__mixfield__density)
984 call mix_dump(namespace, scf%restart_dump, scf%smix, gr, ierr)
985 if (ierr /= 0) then
986 message(1) = 'Unable to write mixing information.'
987 call messages_warning(1, namespace=namespace)
988 end if
989 case (option__mixfield__potential)
990 call hm%ks_pot%dump(scf%restart_dump, gr, ierr)
991 if (ierr /= 0) then
992 message(1) = 'Unable to write Vhxc.'
993 call messages_warning(1, namespace=namespace)
994 end if
995
996 call mix_dump(namespace, scf%restart_dump, scf%smix, gr, ierr)
997 if (ierr /= 0) then
998 message(1) = 'Unable to write mixing information.'
999 call messages_warning(1, namespace=namespace)
1000 end if
1001 end select
1002
1003 end if
1004 end if
1005
1006 call write_convergence_file(static_dir, "convergence")
1007
1008 call profiling_out("SCF_CYCLE")
1009
1010 pop_sub(scf_iter)
1011 contains
1012
1013 ! ---------------------------------------------------------
1014 subroutine scf_write_iter(namespace)
1015 type(namespace_t), intent(in) :: namespace
1016
1017 character(len=50) :: str
1018 real(real64) :: dipole(1:space%dim)
1019
1020 push_sub(scf_run.scf_write_iter)
1021
1022 if ( scf%verbosity_ == verb_full ) then
1023
1024 write(str, '(a,i5)') 'SCF CYCLE ITER #' ,iter
1025 call messages_print_with_emphasis(msg=trim(str), namespace=namespace)
1026 write(message(1),'(a,es15.8,2(a,es9.2))') ' etot = ', units_from_atomic(units_out%energy, hm%energy%total), &
1027 ' abs_ev = ', units_from_atomic(units_out%energy, scf%evsum_diff), &
1028 ' rel_ev = ', scf%evsum_diff/(abs(scf%evsum_out)+1e-20)
1029 write(message(2),'(a,es15.2,2(a,es9.2))') &
1030 ' ediff = ', scf%energy_diff, ' abs_dens = ', scf%abs_dens_diff, &
1031 ' rel_dens = ', scf%abs_dens_diff/st%qtot
1032 call messages_info(2, namespace=namespace)
1033
1034 write(message(1),'(a,i0)') 'Matrix vector products: ', scf%eigens%matvec
1035 write(message(2),'(a,i0)') 'Converged eigenvectors: ', sum(scf%eigens%converged(1:st%nik))
1036 call messages_info(2, namespace=namespace)
1037 call states_elec_write_eigenvalues(st%nst, st, space, hm%kpoints, scf%eigens%diff, compact = .true., namespace=namespace)
1038
1039 if (allocated(hm%vberry)) then
1040 call calc_dipole(dipole, space, gr, st, ions)
1041 call write_dipole(st, hm, space, dipole, namespace=namespace)
1042 end if
1043
1044 if(st%d%ispin > unpolarized) then
1045 call compute_and_write_magnetic_moments(gr, st, hm%phase, hm%ep, ions, scf%lmm_r, namespace=namespace)
1046 end if
1047
1048 if(hm%lda_u_level == dft_u_acbn0) then
1049 call lda_u_write_u(hm%lda_u, namespace=namespace)
1050 call lda_u_write_v(hm%lda_u, namespace=namespace)
1051 end if
1052
1053 write(message(1),'(a)') ''
1054 write(message(2),'(a,i5,a,f14.2)') 'Elapsed time for SCF step ', iter,':', etime
1055 call messages_info(2, namespace=namespace)
1056
1057 call scf_print_mem_use(namespace)
1058
1059 call messages_print_with_emphasis(namespace=namespace)
1060
1061 end if
1062
1063 if ( scf%verbosity_ == verb_compact ) then
1064 write(message(1),'(a,i4,a,es15.8, a,es9.2, a, f7.1, a)') &
1065 'iter ', iter, &
1066 ' : etot ', units_from_atomic(units_out%energy, hm%energy%total), &
1067 ' : abs_dens', scf%abs_dens_diff, &
1068 ' : etime ', etime, 's'
1069 call messages_info(1, namespace=namespace)
1070 end if
1071
1072 pop_sub(scf_run.scf_write_iter)
1073 end subroutine scf_write_iter
1074
1075
1076 ! -----------------------------------------------------
1077 subroutine write_convergence_file(dir, fname)
1078 character(len=*), intent(in) :: dir
1079 character(len=*), intent(in) :: fname
1080
1081 integer :: iunit
1082
1083 if(st%system_grp%is_root()) then ! this the absolute master writes
1084 call io_mkdir(dir, namespace)
1085 iunit = io_open(trim(dir) // "/" // trim(fname), namespace, action='write', position='append')
1086 write(iunit, '(i5,es18.8)', advance = 'no') iter, units_from_atomic(units_out%energy, hm%energy%total)
1087 call iterator%start(scf%criterion_list)
1088 do while (iterator%has_next())
1089 crit => iterator%get_next()
1090 select type (crit)
1091 type is (energy_criterion_t)
1092 write(iunit, '(es13.5)', advance = 'no') units_from_atomic(units_out%energy, crit%val_abs)
1093 type is (density_criterion_t)
1094 write(iunit, '(2es13.5)', advance = 'no') crit%val_abs, crit%val_rel
1095 type is (eigenval_criterion_t)
1096 write(iunit, '(es13.5)', advance = 'no') units_from_atomic(units_out%energy, crit%val_abs)
1097 write(iunit, '(es13.5)', advance = 'no') crit%val_rel
1098 class default
1099 assert(.false.)
1100 end select
1101 end do
1102 if (bitand(ks%xc_family, xc_family_oep) /= 0 .and. ks%theory_level /= hartree_fock &
1103 .and. ks%theory_level /= generalized_kohn_sham_dft) then
1104 if (ks%oep%level == oep_level_full) then
1105 write(iunit, '(es13.5)', advance = 'no') ks%oep%norm2ss
1106 end if
1107 end if
1108 write(iunit,'(a)') ''
1109 call io_close(iunit)
1110 end if
1111 end subroutine write_convergence_file
1112
1113 end subroutine scf_iter
1114
1115 logical function scf_iter_finish(scf, namespace, space, gr, ions, st, ks, hm, iter, outp, &
1116 iters_done) result(completed)
1117 type(scf_t), intent(inout) :: scf
1118 type(namespace_t), intent(in) :: namespace
1119 type(electron_space_t), intent(in) :: space
1120 type(grid_t), intent(inout) :: gr
1121 type(ions_t), intent(inout) :: ions
1122 type(states_elec_t), intent(inout) :: st
1123 type(v_ks_t), intent(inout) :: ks
1124 type(hamiltonian_elec_t), intent(inout) :: hm
1125 integer, intent(in) :: iter
1126 type(output_t), optional, intent(in) :: outp
1127 integer, optional, intent(out) :: iters_done
1128
1129 character(len=MAX_PATH_LEN) :: dirname
1130 integer(int64) :: what_i
1131
1132 push_sub(scf_iter_finish)
1133
1134 completed = .false.
1135
1136 if(scf%finish) then
1137 if(present(iters_done)) iters_done = iter
1138 if(scf%verbosity_ >= verb_compact) then
1139 write(message(1), '(a, i4, a)') 'Info: SCF converged in ', iter, ' iterations'
1140 write(message(2), '(a)') ''
1141 call messages_info(2, namespace=namespace)
1142 end if
1143 completed = .true.
1144 pop_sub(scf_iter_finish)
1145 return
1146 end if
1147 if (present(outp)) then
1148 if (any(outp%what) .and. outp%duringscf) then
1149 do what_i = lbound(outp%what, 1), ubound(outp%what, 1)
1150 if (outp%what_now(what_i, iter)) then
1151 write(dirname,'(a,a,i4.4)') trim(outp%iter_dir),"scf.", iter
1152 call output_all(outp, namespace, space, dirname, gr, ions, iter, st, hm, ks)
1153 call output_modelmb(outp, namespace, space, dirname, gr, ions, iter, st)
1154 exit
1155 end if
1156 end do
1157 end if
1158 end if
1159
1160 ! save information for the next iteration
1161 call lalg_copy(gr%np, st%d%nspin, st%rho, scf%rhoin)
1162
1163 ! restart mixing
1164 if (scf%mix_field /= option__mixfield__none) then
1165 if (scf%smix%ns_restart > 0) then
1166 if (mix_scheme(scf%smix) /= option__mixingscheme__broyden_adaptive .and. &
1167 mod(iter, scf%smix%ns_restart) == 0) then
1168 message(1) = "Info: restarting mixing."
1169 call messages_info(1, namespace=namespace)
1170 call scf_mix_clear(scf)
1171 end if
1172 end if
1173 end if
1174
1175 select case(scf%mix_field)
1176 case(option__mixfield__potential)
1177 call mixfield_set_vin(scf%mixfield, hm%ks_pot%vhxc(1:gr%np, 1:st%d%nspin))
1178 call vtau_mixer_set_vin(scf%vtau_mix, hm)
1179 case (option__mixfield__density)
1180 call mixfield_set_vin(scf%mixfield, scf%rhoin)
1181 end select
1182
1183 !If we use LDA+U, we also have do mix it
1184 if (scf%mix_field /= option__mixfield__states) then
1185 call lda_u_mixer_set_vin(hm%lda_u, scf%lda_u_mix)
1186 end if
1187
1188 ! check if debug mode should be enabled or disabled on the fly
1189 call io_debug_on_the_fly(namespace)
1190
1191 pop_sub(scf_iter_finish)
1192 end function scf_iter_finish
1193
1194 ! ---------------------------------------------------------
1195 subroutine scf_finish(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, iter, outp)
1196 type(scf_t), intent(inout) :: scf
1197 type(namespace_t), intent(in) :: namespace
1198 type(electron_space_t), intent(in) :: space
1199 type(grid_t), intent(inout) :: gr
1200 type(ions_t), intent(inout) :: ions
1201 type(partner_list_t), intent(in) :: ext_partners
1202 type(states_elec_t), intent(inout) :: st
1203 type(v_ks_t), intent(inout) :: ks
1204 type(hamiltonian_elec_t), intent(inout) :: hm
1205 integer, intent(in) :: iter
1206 type(output_t), optional, intent(in) :: outp
1207
1208 integer :: iqn, ib
1209 class(convergence_criterion_t), pointer :: crit
1210 type(criterion_iterator_t) :: iterator
1211
1212
1213 push_sub(scf_finish)
1214
1215 ! Compute the KS potential corresponding to the final density
1216 ! This is critical for getting consistent TD calculations
1217 if ((scf%max_iter > 0 .and. scf%mix_field == option__mixfield__potential) .or. scf%calc_current) then
1218 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
1219 calc_current=scf%calc_current)
1220 end if
1221
1222 select case(scf%mix_field)
1223 case(option__mixfield__states)
1224
1225 do iqn = st%d%kpt%start, st%d%kpt%end
1226 do ib = st%group%block_start, st%group%block_end
1227 call scf%psioutb(ib, iqn)%end()
1228 end do
1229 end do
1230
1231 ! There is a ICE with foss2022a-serial. I am changing to deallocate - NTD
1232 deallocate(scf%psioutb)
1233 end select
1234
1235 safe_deallocate_a(scf%rhoout)
1236 safe_deallocate_a(scf%rhoin)
1237
1238 if (scf%max_iter > 0 .and. any(scf%eigens%converged < st%nst)) then
1239 write(message(1),'(a)') 'Some of the states are not fully converged!'
1240 if (all(scf%eigens%converged >= st%nst_conv)) then
1241 write(message(2),'(a)') 'But all requested states to converge are converged.'
1242 call messages_info(2, namespace=namespace)
1243 else
1244 if(scf%eigens%es_type == rs_chebyshev) then
1245 write(message(2),'(a)') 'With the Chebyshev filtering eigensolver, it usually helps to'
1246 write(message(3),'(a)') 'increase ExtraStates and set ExtraStatesToConverge to the number'
1247 write(message(4),'(a)') 'of states to be converged.'
1248 call messages_warning(4, namespace=namespace)
1249 else
1250 call messages_warning(1, namespace=namespace)
1251 end if
1252 end if
1253 end if
1254
1255 if (.not.scf%finish) then
1256 write(message(1), '(a,i4,a)') 'SCF *not* converged after ', iter - 1, ' iterations.'
1257 if(scf%eigens%es_type == rs_chebyshev) then
1258 write(message(2),'(a)') 'With the Chebyshev filtering eigensolver, it usually helps to'
1259 write(message(3),'(a)') 'increase ExtraStates to improve convergence.'
1260 call messages_warning(3, namespace=namespace)
1261 else
1262 call messages_warning(1, namespace=namespace)
1263 end if
1264 end if
1265
1266 write(message(1), '(a,i10)') 'Info: Number of matrix-vector products: ', scf%matvec
1267 call messages_info(1)
1268
1269 if (scf%calc_force) then
1270 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old=scf%vhxc_old)
1271 end if
1272
1273 if (scf%calc_stress) call stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
1274
1275 ! Update the eigenvalues, to match the KS potential that just got recomputed
1276 if (scf%mix_field == option__mixfield__potential) then
1277 call energy_calc_eigenvalues(namespace, hm, gr%der, st)
1278 call states_elec_fermi(st, namespace, gr)
1279 end if
1280
1281 if(present(outp)) then
1282 ! output final information
1283 call scf_write_static(static_dir, "info")
1284 call output_all(outp, namespace, space, static_dir, gr, ions, -1, st, hm, ks)
1285 call output_modelmb(outp, namespace, space, static_dir, gr, ions, -1, st)
1286 end if
1287
1288 if (space%is_periodic() .and. st%nik > st%d%nspin) then
1289 if (bitand(hm%kpoints%method, kpoints_path) /= 0) then
1290 call states_elec_write_bandstructure(static_dir, namespace, st%nst, st, &
1291 ions, gr, hm%kpoints, hm%phase, vec_pot = hm%hm_base%uniform_vector_potential, &
1292 vec_pot_var = hm%hm_base%vector_potential)
1293 end if
1294 end if
1295
1296 if (ks%vdw%vdw_correction == option__vdwcorrection__vdw_ts) then
1297 call vdw_ts_write_c6ab(ks%vdw%vdw_ts, ions, static_dir, 'c6ab_eff', namespace)
1298 end if
1299
1300 safe_deallocate_a(scf%vhxc_old)
1301
1302 pop_sub(scf_finish)
1303
1304 contains
1305
1306 ! ---------------------------------------------------------
1307 subroutine scf_write_static(dir, fname)
1308 character(len=*), intent(in) :: dir, fname
1309
1310 integer :: iunit
1311 real(real64) :: dipole(1:space%dim)
1312 real(real64) :: ex_virial
1313
1314 push_sub(scf_run.scf_write_static)
1315
1316 if(st%system_grp%is_root()) then ! this the absolute master writes
1317 call io_mkdir(dir, namespace)
1318 iunit = io_open(trim(dir) // "/" // trim(fname), namespace, action='write')
1319
1320 call grid_write_info(gr, iunit=iunit)
1321
1322 call symmetries_write_info(gr%symm, space, iunit=iunit)
1323
1324 if (space%is_periodic()) then
1325 call hm%kpoints%write_info(iunit=iunit)
1326 write(iunit,'(1x)')
1327 end if
1328
1329 call v_ks_write_info(ks, iunit=iunit)
1330
1331 ! scf information
1332 if(scf%finish) then
1333 write(iunit, '(a, i4, a)')'SCF converged in ', iter, ' iterations'
1334 else
1335 write(iunit, '(a)') 'SCF *not* converged!'
1336 end if
1337 write(iunit, '(1x)')
1338
1339 if(any(scf%eigens%converged < st%nst)) then
1340 write(iunit,'(a)') 'Some of the states are not fully converged!'
1341 if (all(scf%eigens%converged >= st%nst_conv)) then
1342 write(iunit,'(a)') 'But all requested states to converge are converged.'
1343 end if
1344 end if
1345
1346 call states_elec_write_eigenvalues(st%nst, st, space, hm%kpoints, iunit=iunit)
1347 write(iunit, '(1x)')
1348
1349 if (space%is_periodic()) then
1350 call states_elec_write_gaps(iunit, st, space)
1351 write(iunit, '(1x)')
1352 end if
1353
1354 write(iunit, '(3a)') 'Energy [', trim(units_abbrev(units_out%energy)), ']:'
1355 else
1356 iunit = -1
1357 end if
1358
1359 call energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full = .true.)
1360
1361 if(st%system_grp%is_root()) write(iunit, '(1x)')
1362 if(st%d%ispin > unpolarized) then
1363 call compute_and_write_magnetic_moments(gr, st, hm%phase, hm%ep, ions, scf%lmm_r, iunit=iunit, &
1364 calc_orb_moments=scf%calc_orb_moments)
1365 if (st%system_grp%is_root()) write(iunit, '(1x)')
1366 end if
1367
1368 if(st%d%ispin == spinors .and. space%dim == 3 .and. &
1369 (ks%theory_level == kohn_sham_dft .or. ks%theory_level == generalized_kohn_sham_dft) ) then
1370 call write_total_xc_torque(iunit, gr, hm%ks_pot%vxc, st)
1371 if(st%system_grp%is_root()) write(iunit, '(1x)')
1372 end if
1373
1374 if(hm%lda_u_level == dft_u_acbn0) then
1375 call lda_u_write_u(hm%lda_u, iunit=iunit)
1376 call lda_u_write_v(hm%lda_u, iunit=iunit)
1377 if(st%system_grp%is_root()) write(iunit, '(1x)')
1378 end if
1379
1380 if(scf%calc_dipole) then
1381 call calc_dipole(dipole, space, gr, st, ions)
1382 call write_dipole(st, hm, space, dipole, iunit=iunit)
1383 end if
1384
1385 ! This only works when we do not have a correlation part
1386 if(ks%theory_level == kohn_sham_dft .and. &
1387 hm%xc%functional(func_c,1)%family == xc_family_none .and. st%d%ispin /= spinors &
1388 .and. .not. space%is_periodic()) then
1389 call energy_calc_virial_ex(gr%der, hm%ks_pot%vxc, st, ex_virial)
1390
1391 if (st%system_grp%is_root()) then
1392 write(iunit, '(3a)') 'Virial relation for exchange [', trim(units_abbrev(units_out%energy)), ']:'
1393 write(iunit,'(a,es14.6)') "Energy from the orbitals ", units_from_atomic(units_out%energy, hm%energy%exchange)
1394 write(iunit,'(a,es14.6)') "Energy from the potential (virial) ", units_from_atomic(units_out%energy, ex_virial)
1395 write(iunit, '(1x)')
1396 end if
1397 end if
1398
1399 if(st%system_grp%is_root()) then
1400 if(scf%max_iter > 0) then
1401 write(iunit, '(a)') 'Convergence:'
1402 call iterator%start(scf%criterion_list)
1403 do while (iterator%has_next())
1404 crit => iterator%get_next()
1405 call crit%write_info(iunit)
1406 end do
1407 write(iunit,'(1x)')
1408 end if
1409 ! otherwise, these values are uninitialized, and unknown.
1410
1411 if (bitand(ks%xc_family, xc_family_oep) /= 0 .and. ks%theory_level /= hartree_fock &
1412 .and. ks%theory_level /= generalized_kohn_sham_dft) then
1413 if ((ks%oep_photon%level == oep_level_full) .or. (ks%oep_photon%level == oep_level_kli)) then
1414 write(iunit, '(a)') 'Photon observables:'
1415 write(iunit, '(6x, a, es15.8,a,es15.8,a)') 'Photon number = ', ks%oep_photon%pt%number(1)
1416 write(iunit, '(6x, a, es15.8,a,es15.8,a)') 'Photon ex. = ', ks%oep_photon%pt%ex
1417 write(iunit,'(1x)')
1418 end if
1419 end if
1420
1421 if (scf%calc_force) call forces_write_info(iunit, ions, dir, namespace)
1422
1423 if (scf%calc_stress) then
1424 call output_stress(iunit, space%periodic_dim, st%stress_tensors, all_terms=.false.)
1425 call output_pressure(iunit, space%periodic_dim, st%stress_tensors%total)
1426 end if
1427
1428 end if
1429
1430 if(scf%calc_partial_charges) then
1431 call partial_charges_compute_and_print_charges(gr, st, ions, iunit)
1432 end if
1433
1434 if(st%system_grp%is_root()) then
1435 call io_close(iunit)
1436 end if
1437
1438 pop_sub(scf_run.scf_write_static)
1439 end subroutine scf_write_static
1440
1441 end subroutine scf_finish
1442
1443 ! ---------------------------------------------------------
1444 subroutine scf_state_info(namespace, st)
1445 type(namespace_t), intent(in) :: namespace
1446 class(states_abst_t), intent(in) :: st
1447
1448 push_sub(scf_state_info)
1449
1450 if (states_are_real(st)) then
1451 call messages_write('Info: SCF using real wavefunctions.')
1452 else
1453 call messages_write('Info: SCF using complex wavefunctions.')
1454 end if
1455 call messages_info(namespace=namespace)
1456
1457 pop_sub(scf_state_info)
1458
1459 end subroutine scf_state_info
1460
1461 ! ---------------------------------------------------------
1462 subroutine scf_print_mem_use(namespace)
1463 type(namespace_t), intent(in) :: namespace
1464 real(real64) :: mem
1465 real(real64) :: mem_tmp
1466
1467 push_sub(scf_print_mem_use)
1468
1469 if(conf%report_memory) then
1470 mem = loct_get_memory_usage()/(1024.0_real64**2)
1471 call mpi_world%allreduce(mem, mem_tmp, 1, mpi_double_precision, mpi_sum)
1472 mem = mem_tmp
1473 write(message(1),'(a,f14.2)') 'Memory usage [Mbytes] :', mem
1474 call messages_info(1, namespace=namespace)
1475 end if
1476
1477 pop_sub(scf_print_mem_use)
1478 end subroutine scf_print_mem_use
1479
1480 ! --------------------------------------------------------
1482 subroutine scf_update_initial_quantity(scf, hm, criterion)
1483 type(scf_t), intent(inout) :: scf
1484 type(hamiltonian_elec_t), intent(in) :: hm
1485 class(convergence_criterion_t), intent(in) :: criterion
1486
1488
1489 select type (criterion)
1490 type is (energy_criterion_t)
1491 scf%energy_in = hm%energy%total
1492 type is (density_criterion_t)
1493 !Do nothing here
1494 type is (eigenval_criterion_t)
1495 !Setting of the value is done in the scf_update_diff_quantity routine
1496 class default
1497 assert(.false.)
1498 end select
1499
1501 end subroutine scf_update_initial_quantity
1502
1503 ! --------------------------------------------------------
1505 subroutine scf_update_diff_quantity(scf, hm, st, gr, rhoout, rhoin, criterion)
1506 type(scf_t), intent(inout) :: scf
1507 type(hamiltonian_elec_t), intent(in) :: hm
1508 type(states_elec_t), intent(in) :: st
1509 type(grid_t), intent(in) :: gr
1510 real(real64), intent(in) :: rhoout(:,:), rhoin(:,:)
1511 class(convergence_criterion_t), intent(in) :: criterion
1512
1513 integer :: is
1514 real(real64), allocatable :: tmp(:)
1515
1516 push_sub(scf_update_diff_quantity)
1517
1518 select type (criterion)
1519 type is (energy_criterion_t)
1520 scf%energy_diff = abs(hm%energy%total - scf%energy_in)
1521
1522 type is (density_criterion_t)
1523 scf%abs_dens_diff = m_zero
1524 safe_allocate(tmp(1:gr%np))
1525 do is = 1, st%d%nspin
1526 tmp(:) = abs(rhoin(1:gr%np, is) - rhoout(1:gr%np, is))
1527 scf%abs_dens_diff = scf%abs_dens_diff + dmf_integrate(gr, tmp)
1528 end do
1529 safe_deallocate_a(tmp)
1530
1531 type is (eigenval_criterion_t)
1532 scf%evsum_out = states_elec_eigenvalues_sum(st)
1533 scf%evsum_diff = abs(scf%evsum_out - scf%evsum_in)
1534 scf%evsum_in = scf%evsum_out
1535
1536 class default
1537 assert(.false.)
1538 end select
1541 end subroutine scf_update_diff_quantity
1542
1543 ! ---------------------------------------------------------
1544 subroutine write_dipole(st, hm, space, dipole, iunit, namespace)
1545 type(states_elec_t), intent(in) :: st
1546 type(hamiltonian_elec_t), intent(in) :: hm
1547 type(electron_space_t), intent(in) :: space
1548 real(real64), intent(in) :: dipole(:)
1549 integer, optional, intent(in) :: iunit
1550 type(namespace_t), optional, intent(in) :: namespace
1551
1552 push_sub(write_dipole)
1553
1554 if(st%system_grp%is_root()) then
1555 call output_dipole(dipole, space%dim, iunit=iunit, namespace=namespace)
1556
1557 if (space%is_periodic()) then
1558 message(1) = "Defined only up to quantum of polarization (e * lattice vector)."
1559 message(2) = "Single-point Berry's phase method only accurate for large supercells."
1560 call messages_info(2, iunit=iunit, namespace=namespace)
1561
1562 if (hm%kpoints%full%npoints > 1) then
1563 message(1) = &
1564 "WARNING: Single-point Berry's phase method for dipole should not be used when there is more than one k-point."
1565 message(2) = "Instead, finite differences on k-points (not yet implemented) are needed."
1566 call messages_info(2, iunit=iunit, namespace=namespace)
1567 end if
1568
1569 if(.not. smear_is_semiconducting(st%smear)) then
1570 message(1) = "Single-point Berry's phase dipole calculation not correct without integer occupations."
1571 call messages_info(1, iunit=iunit, namespace=namespace)
1572 end if
1573 end if
1574
1575 call messages_info(iunit=iunit, namespace=namespace)
1576 end if
1578 pop_sub(write_dipole)
1579 end subroutine write_dipole
1580
1582 subroutine scf_set_lower_bound_is_known(scf, known_lower_bound)
1583 type(scf_t), intent(inout) :: scf
1584 logical, intent(in) :: known_lower_bound
1585
1586 call scf%eigens%set_lower_bound_is_known(known_lower_bound)
1587 end subroutine scf_set_lower_bound_is_known
1588
1589end module scf_oct_m
1590
1591
1592!! Local Variables:
1593!! mode: f90
1594!! coding: utf-8
1595!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:188
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
subroutine, public berry_perform_internal_scf(this, namespace, space, eigensolver, gr, st, hm, iter, ks, ions, ext_partners)
Definition: berry.F90:186
subroutine, public berry_init(this, namespace)
Definition: berry.F90:161
subroutine, public calc_dipole(dipole, space, mesh, st, ions)
Definition: berry.F90:252
subroutine, public criteria_factory_init(list, namespace, check_conv)
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:619
integer, parameter, public rs_evo
subroutine, public eigensolver_init(eigens, namespace, gr, st, hm, mc, space, deactivate_oracle)
integer, parameter, public rs_chebyshev
subroutine, public eigensolver_end(eigens)
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
subroutine, public energy_calc_virial_ex(der, vxc, st, ex)
subroutine, public energy_calc_eigenvalues(namespace, hm, der, st)
integer, parameter, public spin_orbit
Definition: epot.F90:169
integer, parameter, public fully_relativistic_zora
Definition: epot.F90:169
subroutine, public forces_write_info(iunit, ions, dir, namespace)
Definition: forces.F90:594
subroutine, public forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
Definition: forces.F90:340
real(real64), parameter, public m_zero
Definition: global.F90:200
integer, parameter, public hartree_fock
Definition: global.F90:250
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:250
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:250
real(real64), parameter, public lmm_r_single_atom
Default local magnetic moments sphere radius for an isolated system.
Definition: global.F90:234
integer, parameter, public kohn_sham_dft
Definition: global.F90:250
type(conf_t), public conf
Global instance of Octopus configuration.
Definition: global.F90:190
character(len= *), parameter, public static_dir
Definition: global.F90:279
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public grid_write_info(gr, iunit, namespace)
Definition: grid.F90:529
subroutine, public hamiltonian_elec_update_pot(this, mesh, accumulate)
Update the KS potential of the electronic Hamiltonian.
This module defines classes and functions for interaction partners.
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
subroutine, public io_debug_on_the_fly(namespace)
check if debug mode should be enabled or disabled on the fly
Definition: io.F90:535
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:361
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
integer, parameter, public kpoints_path
Definition: kpoints.F90:223
A module to handle KS potential, without the external potential.
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
Definition: lda_u_io.F90:642
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:528
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:723
subroutine, public lda_u_write_v(this, iunit, namespace)
Definition: lda_u_io.F90:576
subroutine, public lda_u_mixer_set_vin(this, mixer)
subroutine, public lda_u_mixer_init(this, mixer, st)
subroutine, public lda_u_mixer_clear(mixer, smix)
subroutine, public lda_u_mixer_init_auxmixer(this, namespace, mixer, smix, st)
subroutine, public lda_u_mixer_get_vnew(this, mixer, st)
subroutine, public lda_u_mixer_set_vout(this, mixer)
subroutine, public lda_u_mixer_end(mixer, smix)
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:892
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:205
System information (time, memory, sysname)
Definition: loct.F90:117
subroutine, public compute_and_write_magnetic_moments(gr, st, phase, ep, ions, lmm_r, calc_orb_moments, iunit, namespace)
Computes and prints the global and local magnetic moments.
Definition: magnetic.F90:214
subroutine, public write_total_xc_torque(iunit, mesh, vxc, st)
Definition: magnetic.F90:567
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
subroutine, public messages_new_line()
Definition: messages.F90:1089
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
integer pure function, public mix_scheme(this)
Definition: mix.F90:827
real(real64) pure function, public mix_coefficient(this)
Definition: mix.F90:821
subroutine, public mixing(namespace, smix)
Main entry-point to SCF mixer.
Definition: mix.F90:847
subroutine, public mix_get_field(this, mixfield)
Definition: mix.F90:839
subroutine, public mix_dump(namespace, restart, smix, mesh, ierr)
Definition: mix.F90:592
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
Definition: mix.F90:270
subroutine, public mix_load(namespace, restart, smix, mesh, ierr)
Definition: mix.F90:691
subroutine, public mix_end(smix)
Definition: mix.F90:569
subroutine, public mix_clear(smix)
Definition: mix.F90:553
subroutine, public modelmb_sym_all_states(space, mesh, st)
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
this module contains the low-level part of the output system
Definition: output_low.F90:117
subroutine, public output_modelmb(outp, namespace, space, dir, gr, ions, iter, st)
this module contains the output system
Definition: output.F90:117
logical function, public output_needs_current(outp, states_are_real)
Definition: output.F90:984
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:498
subroutine, public partial_charges_compute_and_print_charges(mesh, st, ions, iunit)
Computes and write partial charges to a file.
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:335
integer, parameter, public restart_flag_mix
Definition: restart.F90:188
integer, parameter, public restart_flag_rho
Definition: restart.F90:188
integer, parameter, public restart_flag_vhxc
Definition: restart.F90:188
subroutine, public scf_finish(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, iter, outp)
Definition: scf.F90:1291
subroutine, public scf_set_lower_bound_is_known(scf, known_lower_bound)
Set the flag lower_bound_is_known.
Definition: scf.F90:1678
subroutine, public scf_load(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, restart_load)
Loading of restarting data of the SCF cycle.
Definition: scf.F90:605
subroutine write_dipole(st, hm, space, dipole, iunit, namespace)
Definition: scf.F90:1640
subroutine scf_update_initial_quantity(scf, hm, criterion)
Update the quantity at the begining of a SCF cycle.
Definition: scf.F90:1578
subroutine scf_update_diff_quantity(scf, hm, st, gr, rhoout, rhoin, criterion)
Update the quantity at the begining of a SCF cycle.
Definition: scf.F90:1601
subroutine, public scf_state_info(namespace, st)
Definition: scf.F90:1540
subroutine, public scf_print_mem_use(namespace)
Definition: scf.F90:1558
subroutine, public scf_mix_clear(scf)
Definition: scf.F90:588
subroutine, public scf_start(scf, namespace, gr, ions, st, ks, hm, outp, verbosity)
Preparation of the SCF cycle.
Definition: scf.F90:689
integer, parameter, public verb_full
Definition: scf.F90:207
integer, parameter, public verb_compact
Definition: scf.F90:207
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:260
subroutine, public scf_end(scf)
Definition: scf.F90:558
subroutine, public scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, verbosity, iters_done, restart_dump)
Legacy version of the SCF code.
Definition: scf.F90:832
subroutine, public scf_iter(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, iter, outp, restart_dump)
Definition: scf.F90:878
logical function, public scf_iter_finish(scf, namespace, space, gr, ions, st, ks, hm, iter, outp, iters_done)
Definition: scf.F90:1212
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:1042
pure logical function, public states_are_real(st)
This module defines routines to write information about states.
subroutine, public states_elec_write_eigenvalues(nst, st, space, kpoints, error, st_start, compact, iunit, namespace)
write the eigenvalues for some states to a file.
subroutine, public states_elec_write_gaps(iunit, st, space)
calculate gaps and write to a file.
subroutine, public states_elec_write_bandstructure(dir, namespace, nst, st, ions, mesh, kpoints, phase, vec_pot, vec_pot_var)
calculate and write the bandstructure
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
real(real64) function, public states_elec_eigenvalues_sum(st, alt_eig)
function to calculate the eigenvalues sum using occupations as weights
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, verbose)
subroutine, public states_elec_load_rho(restart, st, mesh, ierr)
subroutine, public states_elec_dump_rho(restart, st, mesh, ierr, iter)
This module implements the calculation of the stress tensor.
Definition: stress.F90:120
subroutine, public output_pressure(iunit, space_dim, total_stress_tensor)
Definition: stress.F90:1145
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
Definition: stress.F90:188
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
Definition: stress.F90:1080
subroutine, public symmetries_write_info(this, space, iunit, namespace)
Definition: symmetries.F90:615
type(type_t), public type_float
Definition: types.F90:135
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
subroutine, public output_dipole(dipole, ndim, iunit, namespace)
Definition: utils.F90:281
subroutine, public v_ks_write_info(ks, iunit, namespace)
Definition: v_ks.F90:665
subroutine, public v_ks_update_dftu_energy(ks, namespace, hm, st, int_dft_u)
Update the value of <\psi | V_U | \psi>, where V_U is the DFT+U potential.
Definition: v_ks.F90:1511
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:754
Tkatchenko-Scheffler pairwise method for van der Waals (vdW, dispersion) interactions.
Definition: vdw_ts.F90:121
subroutine, public vdw_ts_write_c6ab(this, ions, dir, fname, namespace)
Definition: vdw_ts.F90:549
subroutine, public vtau_mixer_end(mixer, smix)
Definition: vtau_mixer.F90:189
subroutine, public vtau_mixer_init_auxmixer(namespace, mixer, smix, hm, np, nspin)
Definition: vtau_mixer.F90:152
subroutine, public vtau_mixer_set_vout(mixer, hm)
Definition: vtau_mixer.F90:202
subroutine, public vtau_mixer_get_vnew(mixer, hm)
Definition: vtau_mixer.F90:228
subroutine, public vtau_mixer_clear(mixer, smix)
Definition: vtau_mixer.F90:176
subroutine, public vtau_mixer_set_vin(mixer, hm)
Definition: vtau_mixer.F90:215
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
Definition: walltimer.F90:123
logical function, public walltimer_alarm(comm, print)
indicate whether time is up
Definition: walltimer.F90:333
logical function, public restart_walltime_period_alarm(comm)
Definition: walltimer.F90:375
integer, parameter, public xc_family_nc_mgga
integer, parameter, public func_c
Definition: xc.F90:116
integer, parameter, public oep_level_full
Definition: xc_oep.F90:175
integer, parameter, public oep_level_kli
Definition: xc_oep.F90:175
subroutine scf_write_static(dir, fname)
Definition: rdmft.F90:587
subroutine create_convergence_file(dir, fname)
Definition: scf.F90:794
subroutine scf_write_iter(namespace)
Definition: scf.F90:1110
subroutine write_convergence_file(dir, fname)
Definition: scf.F90:1173
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:171
Stores all communicators and groups.
Definition: multicomm.F90:208
output handler class
Definition: output_low.F90:166
some variables used for the SCF cycle
Definition: scf.F90:213
abstract class for states
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)