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