Octopus
electrons.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2009 X. Andrade
3!! Copyright (C) 2020 M. Oliveira
4!!
5!! This program is free software; you can redistribute it and/or modify
6!! it under the terms of the GNU General Public License as published by
7!! the Free Software Foundation; either version 2, or (at your option)
8!! any later version.
9!!
10!! This program is distributed in the hope that it will be useful,
11!! but WITHOUT ANY WARRANTY; without even the implied warranty of
12!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!! GNU General Public License for more details.
14!!
15!! You should have received a copy of the GNU General Public License
16!! along with this program; if not, write to the Free Software
17!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18!! 02110-1301, USA.
19!!
20
21#include "global.h"
22
23
24module electrons_oct_m
25 use accel_oct_m
33 use debug_oct_m
35 use dipole_oct_m
38 use elf_oct_m
42 use forces_oct_m
44 use global_oct_m
45 use grid_oct_m
55 use ions_oct_m
56 use kick_oct_m
61 use lasers_oct_m
62 use lda_u_oct_m
63 use loct_oct_m
64 use mesh_oct_m
67 use mpi_oct_m
74 use output_oct_m
76 use parser_oct_m
77 use pes_oct_m
92 use rdmft_oct_m
94 use scf_oct_m
95 use space_oct_m
99 use stress_oct_m
100 use sort_oct_m
101 use system_oct_m
102 use td_oct_m
105 use v_ks_oct_m
106 use xc_oct_m
107 use xc_f03_lib_m
108 use xc_oep_oct_m
112
113 implicit none
114
115 private
116 public :: &
118
125 type, extends(system_t) :: electrons_t
126 ! Components are public by default
127 type(electron_space_t) :: space
128 class(ions_t), pointer :: ions => null()
129 type(photons_t), pointer :: photons => null()
130 type(grid_t) :: gr
131 type(states_elec_t) :: st
132 type(v_ks_t) :: ks
133 type(output_t) :: outp
134 type(multicomm_t) :: mc
135 type(hamiltonian_elec_t) :: hm
136 type(td_t) :: td
137 type(current_t) :: current_calculator
138 type(dipole_t) :: dipole
139 type(scf_t) :: scf
140 type(rdm_t) :: rdm
141
142 type(kpoints_t) :: kpoints
143
144 logical :: generate_epot
145
146 type(states_elec_t) :: st_copy
147
148 ! At the moment this is not treated as an external potential
149 class(lasers_t), pointer :: lasers => null()
150 class(gauge_field_t), pointer :: gfield => null()
151
152 ! List with all the external partners
153 ! This will become a list of interactions in the future
154 type(partner_list_t) :: ext_partners
155
156 !TODO: have a list of self interactions
157 type(xc_interaction_t), pointer :: xc_interaction => null()
158
159 logical :: ions_propagated = .false.
160 contains
161 procedure :: init_interaction => electrons_init_interaction
162 procedure :: init_parallelization => electrons_init_parallelization
163 procedure :: new_algorithm => electrons_new_algorithm
164 procedure :: initialize => electrons_initialize
165 procedure :: do_algorithmic_operation => electrons_do_algorithmic_operation
166 procedure :: is_tolerance_reached => electrons_is_tolerance_reached
167 procedure :: update_quantity => electrons_update_quantity
168 procedure :: init_interaction_as_partner => electrons_init_interaction_as_partner
169 procedure :: copy_quantities_to_interaction => electrons_copy_quantities_to_interaction
170 procedure :: output_start => electrons_output_start
171 procedure :: output_write => electrons_output_write
172 procedure :: output_finish => electrons_output_finish
173 procedure :: process_is_slave => electrons_process_is_slave
174 procedure :: restart_write_data => electrons_restart_write_data
175 procedure :: restart_read_data => electrons_restart_read_data
176 procedure :: update_kinetic_energy => electrons_update_kinetic_energy
177 procedure :: algorithm_start => electrons_algorithm_start
178 procedure :: ground_state_run => electrons_ground_state_run_system
179 final :: electrons_finalize
180 end type electrons_t
181
182 interface electrons_t
183 procedure electrons_constructor
184 end interface electrons_t
185
186contains
187
192 function electrons_constructor(namespace, calc_mode_id) result(sys)
193 class(electrons_t), pointer :: sys
194 type(namespace_t), intent(in) :: namespace
195 integer, optional, intent(in) :: calc_mode_id
196
197 integer :: iatom
198 type(lattice_vectors_t) :: latt_inp
199 logical :: has_photons
200
201 push_sub_with_profile(electrons_constructor)
202
203 allocate(sys)
204
205 sys%namespace = namespace
206
207 call messages_obsolete_variable(sys%namespace, 'SystemName')
208
209 sys%space = electron_space_t(sys%namespace)
210 call sys%space%write_info(sys%namespace)
211 if (sys%space%has_mixed_periodicity()) then
212 call messages_experimental('Support for mixed periodicity systems')
213 end if
214
215 sys%ions => ions_t(sys%namespace, latt_inp=latt_inp)
216
217 call grid_init_stage_1(sys%gr, sys%namespace, sys%space, sys%ions%symm, latt_inp, sys%ions%natoms, sys%ions%pos)
218
219 if (sys%space%is_periodic()) then
220 call sys%ions%latt%write_info(sys%namespace)
221 end if
223 ! Sanity check for atomic coordinates
224 do iatom = 1, sys%ions%natoms
225 ! Using the same tolerance of fold_into_cell to avoid problems with mixed periodicities
226 ! as the default tolerance of contains_point is too tight
227 if (.not. sys%gr%box%contains_point(sys%ions%pos(:, iatom), tol=1.0e-6_real64)) then
228 if (sys%space%periodic_dim /= sys%space%dim) then
229 write(message(1), '(a,i5,a)') "Atom ", iatom, " is outside the box."
230 call messages_warning(1, namespace=sys%namespace)
231 end if
232 end if
233 end do
235 ! we need k-points for periodic systems
236 call kpoints_init(sys%kpoints, sys%namespace, sys%gr%symm, sys%space%dim, sys%space%periodic_dim, sys%ions%latt)
238 call states_elec_init(sys%st, sys%namespace, sys%space, sys%ions%val_charge(), sys%kpoints, calc_mode_id)
239 call sys%st%write_info(sys%namespace)
240
241 ! if independent particles in N dimensions are being used, need to initialize them
242 ! after masses are set to 1 in grid_init_stage_1 -> derivatives_init
243 call sys%st%modelmbparticles%copy_masses(sys%gr%der%masses)
245 call elf_init(sys%namespace)
246
247 if (present(calc_mode_id)) then
248 sys%generate_epot = calc_mode_id /= option__calculationmode__dummy
249 else
250 sys%generate_epot = .true.
251 end if
253 call sys%dipole%init(sys%space)
255 sys%supported_interactions_as_partner = [current_to_mxll_field]
257 call sys%quantities%add(quantity_t("wavefunctions", updated_on_demand = .false.))
258 call sys%quantities%add(quantity_t("current", updated_on_demand = .true., parents=["wavefunctions"]))
259 call sys%quantities%add(quantity_t("dipole", updated_on_demand = .true., parents=["wavefunctions"]))
260 call current_init(sys%current_calculator, sys%namespace)
262 !%Variable EnablePhotons
263 !%Type logical
264 !%Default no
265 !%Section Hamiltonian
266 !%Description
267 !% This variable can be used to enable photons in several types of calculations.
268 !% It can be used to activate the one-photon OEP formalism.
269 !% In the case of CalculationMode = casida, it enables photon modes as
270 !% described in ACS Photonics 2019, 6, 11, 2757-2778.
271 !% Finally, if set to yes when solving the ferquency-dependent Sternheimer
272 !% equation, the photons are coupled to the electronic subsystem.
273 !%End
274 call messages_obsolete_variable(namespace, 'OEPPtX', 'EnablePhotons')
275 call parse_variable(namespace, 'EnablePhotons', .false., has_photons)
276 if (has_photons) then
277 call messages_experimental("EnablePhotons = yes")
278 sys%photons => photons_t(sys%namespace)
279 else
280 nullify(sys%photons)
281 end if
282
283 pop_sub_with_profile(electrons_constructor)
284 end function electrons_constructor
285
286 ! ---------------------------------------------------------
287 subroutine electrons_init_interaction(this, interaction)
288 class(electrons_t), target, intent(inout) :: this
289 class(interaction_t), intent(inout) :: interaction
290
291 real(real64) :: dmin
292 integer :: rankmin, depth
293 logical :: mxll_interaction_present
294 logical :: calc_dipole
295
296 push_sub(electrons_init_interactions)
297
298 mxll_interaction_present = .false.
299 select type (interaction)
301 call interaction%init(this%gr, 3)
302 mxll_interaction_present = .true.
303 interaction%type = mxll_field_trans
305 call interaction%init(this%gr, 3)
306 mxll_interaction_present = .true.
308 call interaction%init(this%gr, 3)
309 mxll_interaction_present = .true.
310 class default
311 message(1) = "Trying to initialize an unsupported interaction by the electrons."
312 call messages_fatal(1, namespace=this%namespace)
313 end select
314
315 if (mxll_interaction_present) then
316 calc_dipole = any(this%hm%mxll%coupling_mode == &
318
319 if (calc_dipole) then
320 assert(this%space%periodic_dim == 0)
321 this%hm%mxll%calc_field_dip = .true.
322 this%hm%mxll%center_of_mass(1:3) = this%ions%center_of_mass()
323 this%hm%mxll%center_of_mass_ip = mesh_nearest_point(this%gr, this%hm%mxll%center_of_mass, dmin, rankmin)
324 this%hm%mxll%center_of_mass_rankmin = rankmin
325 end if
326 end if
327
328 ! set interpolation depth for field-transfer interactions
329 select type (interaction)
330 class is (field_transfer_t)
331 ! interpolation depth depends on the propagator
332 select type (algo => this%algo)
333 type is (propagator_exp_mid_t)
334 depth = 3
335 type is (propagator_aetrs_t)
336 depth = 3
337 type is (propagator_bomd_t)
338 depth = 1
339 class default
340 message(1) = "The chosen algorithm does not yet support interaction interpolation"
341 call messages_fatal(1, namespace=this%namespace)
342 end select
343
344 call interaction%init_interpolation(depth, interaction%label)
345 end select
346
347 pop_sub(electrons_init_interactions)
348 end subroutine electrons_init_interaction
349
350 ! ---------------------------------------------------------
351 subroutine electrons_init_parallelization(this, grp)
352 class(electrons_t), intent(inout) :: this
353 type(mpi_grp_t), intent(in) :: grp
354
355 integer(int64) :: index_range(4)
356 real(real64) :: mesh_global, mesh_local, wfns
357 integer :: idir
358 real(real64) :: spiral_q(3), spiral_q_red(3)
359 type(block_t) :: blk
360
362
363 call mpi_grp_copy(this%grp, grp)
364
365 ! store the ranges for these two indices (serves as initial guess
366 ! for parallelization strategy)
367 index_range(1) = this%gr%np_global ! Number of points in mesh
368 index_range(2) = this%st%nst ! Number of states
369 index_range(3) = this%st%nik ! Number of k-points
370 index_range(4) = 100000 ! Some large number
371
372 ! create index and domain communicators
373 call multicomm_init(this%mc, this%namespace, this%grp, calc_mode_par, &
374 mpi_world%size, index_range, (/ 5000, 1, 1, 1 /))
375
376 call this%ions%partition(this%mc)
377 call kpoints_distribute(this%st, this%mc)
378 call states_elec_distribute_nodes(this%st, this%namespace, this%mc)
379
380
381 if (parse_is_defined(this%namespace, 'TDMomentumTransfer') .or. &
382 parse_is_defined(this%namespace, 'TDReducedMomentumTransfer')) then
383 if (parse_block(this%namespace, 'TDMomentumTransfer', blk) == 0) then
384 do idir = 1, this%space%dim
385 call parse_block_float(blk, 0, idir - 1, spiral_q(idir))
386 end do
387 else if(parse_block(this%namespace, 'TDReducedMomentumTransfer', blk) == 0) then
388 do idir = 1, this%space%dim
389 call parse_block_float(blk, 0, idir - 1, spiral_q_red(idir))
390 end do
391 call kpoints_to_absolute(this%kpoints%latt, spiral_q_red(1:this%space%dim), spiral_q(1:this%space%dim))
392 end if
393 call parse_block_end(blk)
394 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc, spiral_q)
395 else
396 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
397 end if
398
399 if (this%st%symmetrize_density) then
400 call mesh_check_symmetries(this%gr, this%gr%symm, this%ions%space%periodic_dim)
401 end if
402
403 call output_init(this%outp, this%namespace, this%space, this%st, this%gr, this%st%nst, this%ks)
404 call states_elec_densities_init(this%st, this%gr)
405 call states_elec_exec_init(this%st, this%namespace, this%mc)
406
407 if (associated(this%photons)) then
408 this%ks%has_photons = .true.
409 end if
410
411 call v_ks_init(this%ks, this%namespace, this%gr, this%st, this%ions, this%mc, this%space, this%kpoints)
412 if (this%ks%theory_level == kohn_sham_dft .or. this%ks%theory_level == generalized_kohn_sham_dft) then
413 this%xc_interaction => xc_interaction_t(this)
414 end if
415
416 ! For the moment the photons are initialized here
417
418 if(this%ks%has_photons) then
419 this%ks%pt => this%photons%modes
420 ! Temporary creation that should go in the system initialization later
421 call photon_mode_set_n_electrons(this%photons%modes, this%st%qtot)
422 write(message(1), '(a,i5,a)') 'Happy to have ', this%photons%modes%nmodes, ' photon modes with us.'
423 call messages_info(1)
424 call mf_init(this%ks%pt_mx, this%gr, this%st, this%ions, this%ks%pt)
425 ! OEP for photons
426 if(bitand(this%ks%xc_family, xc_family_oep) /= 0 .and. this%ks%xc%functional(func_x,1)%id == xc_oep_x) then
427 this%ks%oep_photon%level = this%ks%oep%level
428 call xc_oep_photon_init(this%ks%oep_photon, this%namespace, this%ks%xc_family, this%gr, this%st, this%mc, this%space)
429 else
430 this%ks%oep_photon%level = oep_level_none
431 end if
432
433 end if
434
435
436 ! Temporary place for the initialization of the lasers
437 this%lasers => lasers_t(this%namespace)
438 call lasers_parse_external_fields(this%lasers)
439 call lasers_generate_potentials(this%lasers, this%gr, this%space, this%ions%latt)
440 if(this%lasers%no_lasers > 0) then
441 call this%ext_partners%add(this%lasers)
442 call lasers_check_symmetries(this%lasers, this%kpoints)
443 else
444 deallocate(this%lasers)
445 end if
447 ! Temporary place for the initialization of the gauge field
448 this%gfield => gauge_field_t(this%namespace, this%ions%latt%rcell_volume)
449 if(gauge_field_is_used(this%gfield)) then
450 call this%ext_partners%add(this%gfield)
451 call gauge_field_check_symmetries(this%gfield, this%kpoints)
452 else
453 deallocate(this%gfield)
454 end if
455
456 call hamiltonian_elec_init(this%hm, this%namespace, this%space, this%gr, this%ions, this%ext_partners, &
457 this%st, this%ks%theory_level, this%ks%xc, this%mc, this%kpoints, &
458 need_exchange = output_need_exchange(this%outp) .or. this%ks%oep%level /= oep_level_none, &
459 xc_photons = this%ks%xc_photons )
460
461 if (this%hm%pcm%run_pcm .and. this%mc%par_strategy /= p_strategy_serial .and. this%mc%par_strategy /= p_strategy_states) then
462 call messages_experimental('Parallel in domain calculations with PCM')
463 end if
464
465 ! Print memory requirements
466 call messages_print_with_emphasis(msg='Approximate memory requirements', namespace=this%namespace)
467
468 mesh_global = mesh_global_memory(this%gr)
469 mesh_local = mesh_local_memory(this%gr)
470
471 call messages_write('Mesh')
472 call messages_new_line()
473 call messages_write(' global :')
474 call messages_write(mesh_global, units = unit_megabytes, fmt = '(f10.1)')
475 call messages_new_line()
476 call messages_write(' local :')
477 call messages_write(mesh_local, units = unit_megabytes, fmt = '(f10.1)')
478 call messages_new_line()
479 call messages_write(' total :')
480 call messages_write(mesh_global + mesh_local, units = unit_megabytes, fmt = '(f10.1)')
481 call messages_new_line()
482 call messages_info(namespace=this%namespace)
483
484 wfns = states_elec_wfns_memory(this%st, this%gr)
485 call messages_write('States')
486 call messages_new_line()
487 call messages_write(' real :')
488 call messages_write(wfns, units = unit_megabytes, fmt = '(f10.1)')
489 call messages_write(' (par_kpoints + par_states + par_domains)')
490 call messages_new_line()
491 call messages_write(' complex :')
492 call messages_write(2.0_8*wfns, units = unit_megabytes, fmt = '(f10.1)')
493 call messages_write(' (par_kpoints + par_states + par_domains)')
494 call messages_new_line()
495 call messages_info(namespace=this%namespace)
496
497 call messages_print_with_emphasis(namespace=this%namespace)
498
499 if (this%generate_epot) then
500 message(1) = "Info: Generating external potential"
501 call messages_info(1, namespace=this%namespace)
502 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
503 this%ext_partners, this%st)
504 message(1) = " done."
505 call messages_info(1, namespace=this%namespace)
506 end if
507
508 if (this%ks%theory_level /= independent_particles) then
509 call poisson_async_init(this%hm%psolver, this%mc)
510 ! slave nodes do not call the calculation routine
511 if (multicomm_is_slave(this%mc)) then
512 !for the moment we only have one type of slave
513 call poisson_slave_work(this%hm%psolver, this%namespace)
514 end if
515 end if
516
517 allocate(this%supported_interactions(0))
518 select case (this%hm%mxll%coupling_mode)
520 this%supported_interactions = [this%supported_interactions, mxll_e_field_to_matter]
522 this%supported_interactions = [this%supported_interactions, mxll_vec_pot_to_matter]
523 if (this%hm%mxll%add_zeeman) then
524 this%supported_interactions = [this%supported_interactions, mxll_b_field_to_matter]
525 end if
527 if (this%hm%mxll%add_electric_dip .or. this%hm%mxll%add_electric_quad) then
528 this%supported_interactions = [this%supported_interactions, mxll_e_field_to_matter]
529 end if
530 if (this%hm%mxll%add_magnetic_dip) then
531 this%supported_interactions = [this%supported_interactions, mxll_b_field_to_matter]
532 end if
534 ! Do not initialize any interaction with Maxwell
535 case default
536 message(1) = "Unknown maxwell-matter coupling"
537 call messages_fatal(1, namespace=this%namespace)
538 end select
539
541 end subroutine electrons_init_parallelization
542
543 ! ---------------------------------------------------------
544 subroutine electrons_new_algorithm(this, factory)
545 class(electrons_t), intent(inout) :: this
546 class(algorithm_factory_t), intent(in) :: factory
547
549
550 call system_new_algorithm(this, factory)
551
552 select type (algo => this%algo)
553 class is (propagator_t)
554
555 call td_init(this%td, this%namespace, this%space, this%gr, this%ions, this%st, this%ks, &
556 this%hm, this%ext_partners, this%outp)
557
558 ! this corresponds to the first part of td_init_run
559 call td_allocate_wavefunctions(this%td, this%namespace, this%mc, this%gr, this%ions, this%st, &
560 this%hm, this%space)
561 call td_init_gaugefield(this%td, this%namespace, this%gr, this%st, this%ks, this%hm, &
562 this%ext_partners, this%space)
563
564 class is (minimizer_algorithm_t)
565
566 call electrons_gs_allocate_wavefunctions(this%namespace, this%gr, this%st, this%hm, this%scf, this%ks, &
567 this%ions)
568
569 class default
570 assert(.false.)
571 end select
572
574 end subroutine electrons_new_algorithm
575
576 ! ---------------------------------------------------------
577 subroutine electrons_initialize(this)
578 class(electrons_t), intent(inout) :: this
579
580 push_sub(electrons_initialize)
581
582 if (associated(this%ions)) call this%ions%initialize()
583
584 select type (algo => this%algo)
585 class is (propagator_t)
586 call td_set_from_scratch(this%td, .true.)
587 call td_load_restart_from_gs(this%td, this%namespace, this%space, this%mc, this%gr, &
588 this%ext_partners, this%st, this%ks, this%hm)
589
590 class is (minimizer_algorithm_t)
591
592 call electrons_gs_initialize(this%namespace, this%scf, this%rdm, this%gr, this%mc, this%st, &
593 this%hm, this%ions, this%ks, this%space, this%ext_partners, fromscratch=.true.)
594 class default
595 assert(.false.)
596 end select
597
598 pop_sub(electrons_initialize)
599 end subroutine electrons_initialize
600
601 ! ---------------------------------------------------------
602 subroutine electrons_algorithm_start(this)
603 class(electrons_t), intent(inout) :: this
604
606
607 call system_algorithm_start(this)
608
609 select type (algo => this%algo)
610 class is (propagator_t)
611
612 ! additional initialization needed for electrons
613 call td_init_with_wavefunctions(this%td, this%namespace, this%space, this%mc, this%gr, this%ions, &
614 this%ext_partners, this%st, this%ks, this%hm, this%outp, td_get_from_scratch(this%td))
615
616 end select
617
619 end subroutine electrons_algorithm_start
620
621 ! ---------------------------------------------------------
622 logical function electrons_do_algorithmic_operation(this, operation, updated_quantities) result(done)
623 class(electrons_t), intent(inout) :: this
624 class(algorithmic_operation_t), intent(in) :: operation
625 character(len=:), allocatable, intent(out) :: updated_quantities(:)
626
627 logical :: update_energy_
628 type(gauge_field_t), pointer :: gfield
629 real(real64) :: time
630 integer :: iter
631
633 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
634
635 update_energy_ = .true.
636
637 ! kick at t > 0 still missing!
638
639 done = .true.
640 select type (algo => this%algo)
641 class is (propagator_t)
642 time = algo%iteration%value()
643 select case (operation%id)
644 case (aetrs_first_half)
645 ! propagate half of the time step with H(time)
646 call get_fields_from_interaction(this, time)
647 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
648 this%hm, this%ext_partners, time)
649 call propagation_ops_elec_exp_apply(this%td%tr%te, this%namespace, this%st, this%gr, this%hm, m_half*algo%dt)
650
651 case (aetrs_extrapolate)
652 ! Do the extrapolation of the Hamiltonian
653 ! First the extrapolation of the potential
654 call this%hm%ks_pot%interpolation_new(this%td%tr%vks_old, time+algo%dt, algo%dt)
655
656 !Get the potentials from the interpolator
657 call propagation_ops_elec_interpolate_get(this%hm, this%td%tr%vks_old)
658
659 ! Second the ions and gauge field which later on will be treated as extrapolation
660 ! of interactions, but this is not yet possible
661
662 ! move the ions to time t
663 call propagation_ops_elec_move_ions(this%td%tr%propagation_ops_elec, this%gr, this%hm, &
664 this%st, this%namespace, this%space, this%td%ions_dyn, this%ions, this%ext_partners, &
665 this%mc, time+algo%dt, algo%dt)
666
667 !Propagate gauge field
668 gfield => list_get_gauge_field(this%ext_partners)
669 if(associated(gfield)) then
670 call propagation_ops_elec_propagate_gauge_field(this%td%tr%propagation_ops_elec, gfield, &
671 algo%dt, time+algo%dt)
672 end if
673
674 !Update Hamiltonian and current
675 call get_fields_from_interaction(this, time+algo%dt)
676 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
677 this%hm, this%ext_partners, time+algo%dt)
678
679 case (aetrs_second_half)
680 !Do the time propagation for the second half of the time step
681 call propagation_ops_elec_fuse_density_exp_apply(this%td%tr%te, this%namespace, this%st, &
682 this%gr, this%hm, m_half*algo%dt)
683
684 updated_quantities = ["wavefunctions"]
685
686 case (expmid_extrapolate)
687 ! the half step of this propagator screws with the gauge field kick
688 gfield => list_get_gauge_field(this%ext_partners)
689 if(associated(gfield)) then
690 assert(gauge_field_is_propagated(gfield) .eqv. .false.)
691 end if
692
693 ! Do the extrapolation of the Hamiltonian
694 ! First the extrapolation of the potential
695 call this%hm%ks_pot%interpolation_new(this%td%tr%vks_old, time+algo%dt, algo%dt)
696
697 ! get the potentials from the interpolator
698 call this%hm%ks_pot%interpolate_potentials(this%td%tr%vks_old, 3, time+algo%dt, algo%dt, time + algo%dt/m_two)
699
700 ! Second the ions which later on will be treated as extrapolation of interactions,
701 ! but this is not yet possible
702 ! move the ions to the half step
703 call propagation_ops_elec_move_ions(this%td%tr%propagation_ops_elec, this%gr, this%hm, this%st, &
704 this%namespace, this%space, this%td%ions_dyn, this%ions, this%ext_partners, &
705 this%mc, time + m_half*algo%dt, m_half*algo%dt, save_pos=.true.)
706
707 call get_fields_from_interaction(this, time + m_half*algo%dt)
708 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
709 this%hm, this%ext_partners, time + m_half*algo%dt)
710
711 case (expmid_propagate)
712 ! Do the actual propagation step
713 call propagation_ops_elec_fuse_density_exp_apply(this%td%tr%te, this%namespace, this%st, &
714 this%gr, this%hm, algo%dt)
715
716 ! restore to previous time
717 call propagation_ops_elec_restore_ions(this%td%tr%propagation_ops_elec, this%td%ions_dyn, this%ions)
718
719 updated_quantities = ["wavefunctions"]
720
721 case (bomd_start)
722 call scf_init(this%scf, this%namespace, this%gr, this%ions, this%st, this%mc, this%hm, this%space)
723 ! the ions are propagated inside the propagation step already, so no need to do it at the end
724 this%ions_propagated = .true.
725
726 case (verlet_update_pos)
727 ! move the ions to time t
728 call propagation_ops_elec_propagate_ions_and_cell(this%gr, this%hm, this%st, this%namespace, this%space, &
729 this%td%ions_dyn, this%ions, this%mc, time+algo%dt, algo%dt)
730
731 case (bomd_elec_scf)
732 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
733 this%ext_partners, this%st, time = time+algo%dt)
734 ! now calculate the eigenfunctions
735 call scf_run(this%scf, this%namespace, this%space, this%mc, this%gr, this%ions, &
736 this%ext_partners, this%st, this%ks, this%hm, verbosity = verb_compact)
737 ! TODO: Check if this call is realy needed. - NTD
738 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
739 this%ext_partners, this%st, time = time+algo%dt)
740
741 ! update Hamiltonian and eigenvalues (fermi is *not* called)
742 call v_ks_calc(this%ks, this%namespace, this%space, this%hm, this%st, this%ions, this%ext_partners, &
743 calc_eigenval = .true., time = time+algo%dt, calc_energy = .true.)
744
745 ! Get the energies.
746 call energy_calc_total(this%namespace, this%space, this%hm, this%gr, this%st, this%ext_partners, iunit = -1)
747
748 updated_quantities = ["wavefunctions"]
749
750 case (verlet_compute_acc)
751 ! Do nothing, forces are computing in scf_run
752 if (this%td%ions_dyn%cell_relax()) then
753 assert(this%scf%calc_stress)
754 call this%td%ions_dyn%update_stress(this%ions%space, this%st%stress_tensors%total, &
755 this%ions%latt%rlattice, this%ions%latt%rcell_volume)
756 end if
757
758 case (verlet_compute_vel)
759 call ion_dynamics_propagate_vel(this%td%ions_dyn, this%ions)
760 ! TODO: Check if this call is realy needed. - NTD
761 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
762 this%ext_partners, this%st, time = time+algo%dt)
763 call this%ions%update_kinetic_energy()
764
765 case (expmid_start)
766 this%ions_propagated = .false.
767
768 case (aetrs_start)
769 ! the ions are propagated inside the propagation step already, so no need to do it at the end
770 this%ions_propagated = .true.
771
772 case (iteration_done)
774 done = .false.
775
776 case (bomd_finish)
777 call scf_end(this%scf)
778
780 case default
781 done = .false.
782 end select
783
784 class is(minimizer_algorithm_t)
785
786 ! Clocks starts at 0....
787 iter = nint(this%iteration%value()) + 1
788
789 select case(operation%id)
790 case (gs_scf_start)
791 call scf_start(this%scf, this%namespace, this%gr, this%ions, this%st, this%ks, this%hm, this%outp)
792
793 case (gs_scf_iteration)
794
795 call scf_iter(this%scf, this%namespace, this%space, this%mc, this%gr, this%ions, &
796 this%ext_partners, this%st, this%ks, this%hm, iter, outp = this%outp, &
797 restart_dump=this%scf%restart_dump)
798
799 algo%converged = scf_iter_finish(this%scf, this%namespace, this%space, this%gr, this%ions,&
800 this%st, this%ks, this%hm, iter, outp = this%outp)
801
802 updated_quantities = ["wavefunctions"]
803
804 case (gs_scf_finish)
805
806 ! Here iteration is iter-1, as the clock ticked before SCF_FINISH
807 call scf_finish(this%scf, this%namespace, this%space, this%gr, this%ions, &
808 this%ext_partners, this%st, this%ks, this%hm, iter-1, outp = this%outp)
809
810 case default
811 done = .false.
812 end select
813 class default
814 done = .false.
815 end select
816
817 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
820
821 ! ---------------------------------------------------------
822 logical function electrons_is_tolerance_reached(this, tol) result(converged)
823 class(electrons_t), intent(in) :: this
824 real(real64), intent(in) :: tol
825
827
828 converged = .false.
829
832
833 ! ---------------------------------------------------------
834 subroutine electrons_update_quantity(this, label)
835 class(electrons_t), intent(inout) :: this
836 character(len=*), intent(in) :: label
837
838 class(quantity_t), pointer :: quantity
839
841 call profiling_in(trim(this%namespace%get())//":"//"UPDATE_QUANTITY")
842
843 quantity => this%quantities%get(label)
844 if(associated(quantity)) then
845 assert(quantity%updated_on_demand)
846 endif
847
848 select case (label)
849 case ("current")
850 call states_elec_allocate_current(this%st, this%space, this%gr)
851 call current_calculate(this%current_calculator, this%namespace, this%gr, &
852 this%hm, this%space, this%st)
853 case ("dipole")
854 call this%dipole%calculate(this%gr, this%ions, this%st)
855 case default
856 message(1) = "Incompatible quantity: "//trim(label)//"."
857 call messages_fatal(1, namespace=this%namespace)
858 end select
859
860 call profiling_out(trim(this%namespace%get())//":"//"UPDATE_QUANTITY")
862 end subroutine electrons_update_quantity
863
864 ! ---------------------------------------------------------
865 subroutine electrons_init_interaction_as_partner(partner, interaction)
866 class(electrons_t), intent(in) :: partner
867 class(interaction_surrogate_t), intent(inout) :: interaction
868
870
871 select type (interaction)
873 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
874 class default
875 message(1) = "Unsupported interaction."
876 call messages_fatal(1, namespace=partner%namespace)
877 end select
878
881
882 ! ---------------------------------------------------------
883 subroutine electrons_copy_quantities_to_interaction(partner, interaction)
884 class(electrons_t), intent(inout) :: partner
885 class(interaction_surrogate_t), intent(inout) :: interaction
886
888 call profiling_in(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
889
890 select type (interaction)
892 assert(allocated(partner%st%current))
893 interaction%partner_field(:,:) = partner%st%current(1:partner%gr%np,:,1)
894 call interaction%do_mapping()
895 class default
896 message(1) = "Unsupported interaction."
897 call messages_fatal(1, namespace=partner%namespace)
898 end select
899
900 call profiling_out(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
903
904 ! ---------------------------------------------------------
905 subroutine electrons_output_start(this)
906 class(electrons_t), intent(inout) :: this
907
908 push_sub(electrons_output_start)
909
911 end subroutine electrons_output_start
912
913 ! ---------------------------------------------------------
914 subroutine electrons_output_write(this)
915 class(electrons_t), intent(inout) :: this
916
917 integer :: iter
918
919 push_sub(electrons_output_write)
920 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
921
922 select type (algo => this%algo)
923 class is (propagator_t)
924 iter = this%iteration%counter()
925
926 call td_write_iter(this%td%write_handler, this%namespace, this%space, this%outp, this%gr, &
927 this%st, this%hm, this%ions, this%ext_partners, this%hm%kick, this%ks, algo%dt, iter, this%mc, this%td%recalculate_gs)
928
929 if (this%outp%anything_now(iter)) then ! output
930 call td_write_output(this%namespace, this%space, this%gr, this%st, this%hm, this%ks, &
931 this%outp, this%ions, this%ext_partners, iter, algo%dt)
932 end if
933 end select
934
935 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
937 end subroutine electrons_output_write
938
939 ! ---------------------------------------------------------
940 subroutine electrons_output_finish(this)
941 class(electrons_t), intent(inout) :: this
942
944
946 end subroutine electrons_output_finish
947
948 ! ---------------------------------------------------------
949 logical function electrons_process_is_slave(this) result(is_slave)
950 class(electrons_t), intent(in) :: this
951
953
954 is_slave = multicomm_is_slave(this%mc)
955
957 end function electrons_process_is_slave
958
959 ! ---------------------------------------------------------
961 class(electrons_t), intent(inout) :: this
962 class(propagator_t), intent(in) :: prop
963
964 logical :: stopping
965 logical :: generate
966 logical :: update_energy_
967 integer :: nt
968 real(real64) :: time
969 type(gauge_field_t), pointer :: gfield
970
972 call profiling_in(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
973
974 stopping = .false.
975
976 nt = this%td%iter
977 ! this is the time at the end of the timestep, as required in all routines here
978 time = prop%dt*nt
979 update_energy_ = .true.
980
981 !Apply mask absorbing boundaries
982 if (this%hm%abs_boundaries%abtype == mask_absorbing) call zvmask(this%gr, this%hm, this%st)
983
984 !Photoelectron stuff
985 if (this%td%pesv%calc_spm .or. this%td%pesv%calc_mask .or. this%td%pesv%calc_flux) then
986 call pes_calc(this%td%pesv, this%namespace, this%space, this%gr, this%st, &
987 prop%dt, nt, this%gr%der, this%hm%kpoints, this%ext_partners, stopping)
988 end if
989
990 ! For BOMD, we do not want the lines below to be executed
991 select type(prop)
992 type is(propagator_bomd_t)
993 call profiling_out(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
995 return
996 end select
997
998 ! The propagation of the ions and the gauge field is currently done here.
999 ! TODO: this code is to be moved to their own systems at some point
1000 generate = .false.
1001 if (this%td%ions_dyn%is_active()) then
1002 if (.not. this%ions_propagated) then
1003 call propagation_ops_elec_propagate_ions_and_cell(this%gr, this%hm, this%st, this%namespace, this%space, &
1004 this%td%ions_dyn, this%ions, this%mc, abs(nt*prop%dt), this%td%ions_dyn%ionic_scale*prop%dt)
1005 generate = .true.
1006 end if
1007 end if
1008
1009 gfield => list_get_gauge_field(this%ext_partners)
1010 if(associated(gfield)) then
1011 if (gauge_field_is_propagated(gfield) .and. .not. this%ions_propagated) then
1013 end if
1014 end if
1015
1016 if (generate .or. this%ions%has_time_dependent_species()) then
1017 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
1018 this%ext_partners, this%st, time = abs(nt*prop%dt))
1019 end if
1020
1021 call v_ks_calc(this%ks, this%namespace, this%space, this%hm, this%st, this%ions, this%ext_partners, &
1022 calc_eigenval = update_energy_, time = abs(nt*prop%dt), calc_energy = update_energy_)
1023
1024 if (update_energy_) then
1025 call energy_calc_total(this%namespace, this%space, this%hm, this%gr, this%st, this%ext_partners, iunit = -1)
1026 end if
1027
1028 ! Recalculate forces, update velocities...
1029 if (this%td%ions_dyn%ions_move() .or. this%outp%what(option__output__forces) &
1030 .or. this%td%write_handler%out(out_separate_forces)%write) then
1031 call forces_calculate(this%gr, this%namespace, this%ions, this%hm, this%ext_partners, &
1032 this%st, this%ks, t = abs(nt*prop%dt), dt = prop%dt)
1033 end if
1034
1035 if (this%td%ions_dyn%cell_relax() .or. this%outp%what(option__output__stress)) then
1036 call stress_calculate(this%namespace, this%gr, this%hm, this%st, this%ions, this%ks, this%ext_partners)
1037 end if
1038
1039 if(this%td%ions_dyn%is_active()) then
1040 call ion_dynamics_propagate_vel(this%td%ions_dyn, this%ions, atoms_moved = generate)
1041 call this%ions%update_kinetic_energy()
1042 end if
1043
1044 if(associated(gfield)) then
1045 if(gauge_field_is_propagated(gfield)) then
1046 call gauge_field_get_force(gfield, this%gr, this%st%d%spin_channels, this%st%current, this%ks%xc%lrc)
1048 end if
1049 end if
1050
1051 !We update the occupation matrices
1052 call lda_u_update_occ_matrices(this%hm%lda_u, this%namespace, this%gr, this%st, this%hm%phase, this%hm%energy)
1053
1054 ! this is needed to be compatible with the code in td_*
1055 this%td%iter = this%td%iter + 1
1056
1057 call profiling_out(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
1060
1061 ! ---------------------------------------------------------
1062 subroutine electrons_restart_write_data(this)
1063 class(electrons_t), intent(inout) :: this
1064
1065 integer :: ierr
1066
1068 call profiling_in(trim(this%namespace%get())//":"//"RESTART_WRITE")
1069
1070 select type (algo => this%algo)
1071 class is (propagator_t)
1072 call td_write_data(this%td%write_handler)
1073 call td_dump(this%td, this%namespace, this%space, this%gr, this%st, this%hm, &
1074 this%ks, this%ext_partners, this%iteration%counter(), ierr)
1075 if (ierr /= 0) then
1076 message(1) = "Unable to write time-dependent restart information."
1077 call messages_warning(1, namespace=this%namespace)
1078 end if
1079
1080 ! TODO: this is here because of legacy reasons and should be moved to the output framework
1081 call pes_output(this%td%pesv, this%namespace, this%space, this%gr, this%st, this%iteration%counter(), &
1082 this%outp, algo%dt, this%ions)
1083 end select
1084
1085 call profiling_out(trim(this%namespace%get())//":"//"RESTART_WRITE")
1087 end subroutine electrons_restart_write_data
1088
1089 ! ---------------------------------------------------------
1090 ! this function returns true if restart data could be read
1091 logical function electrons_restart_read_data(this)
1092 class(electrons_t), intent(inout) :: this
1093
1094 logical :: from_scratch
1095
1097 call profiling_in(trim(this%namespace%get())//":"//"RESTART_READ")
1098
1099 select type (algo => this%algo)
1100 class is (propagator_t)
1101 from_scratch = .false.
1102 call td_load_restart_from_td(this%td, this%namespace, this%space, this%mc, this%gr, &
1103 this%ext_partners, this%st, this%ks, this%hm, from_scratch)
1104 call td_set_from_scratch(this%td, from_scratch)
1105
1106 class is (minimizer_algorithm_t)
1107 from_scratch = .false.
1108 call electrons_gs_load_from_restart(this%namespace, this%scf, this%gr, this%mc, this%st, this%hm, &
1109 this%ks, this%space, this%ions, this%ext_partners,from_scratch)
1110
1111 ! electrons_gs_initialize still knows about fromScratch.
1112 assert(.false.)
1113 end select
1114
1115 if (from_scratch) then
1116 ! restart data could not be loaded
1118 else
1119 ! restart data could be loaded
1121 end if
1122
1123 call profiling_out(trim(this%namespace%get())//":"//"RESTART_READ")
1125 end function electrons_restart_read_data
1126
1127 !----------------------------------------------------------
1128 subroutine electrons_update_kinetic_energy(this)
1129 class(electrons_t), intent(inout) :: this
1130
1132
1133 if (states_are_real(this%st)) then
1134 this%kinetic_energy = denergy_calc_electronic(this%namespace, this%hm, this%gr%der, this%st, terms = term_kinetic)
1135 else
1136 this%kinetic_energy = zenergy_calc_electronic(this%namespace, this%hm, this%gr%der, this%st, terms = term_kinetic)
1137 end if
1138
1140
1141 end subroutine electrons_update_kinetic_energy
1142
1143 ! ---------------------------------------------------------
1144 subroutine get_fields_from_interaction(this, time)
1145 class(electrons_t), intent(inout) :: this
1146 real(real64), intent(in) :: time
1147
1148 type(interaction_iterator_t) :: iter
1149 real(real64), allocatable :: field_tmp(:, :)
1150
1152
1153 if (this%hm%mxll%coupling_mode == no_maxwell_coupling) then
1155 return
1156 end if
1158 safe_allocate(field_tmp(1:this%gr%np, 1:this%gr%box%dim))
1159 this%hm%mxll%e_field = m_zero
1160 this%hm%mxll%b_field = m_zero
1161 this%hm%mxll%vec_pot = m_zero
1162
1163 ! interpolate field from interaction
1164 call iter%start(this%interactions)
1165 do while (iter%has_next())
1166 select type (interaction => iter%get_next())
1167 class is (mxll_e_field_to_matter_t)
1168 call interaction%interpolate(time, field_tmp)
1169 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%e_field)
1170 class is (mxll_vec_pot_to_matter_t)
1171 call interaction%interpolate(time, field_tmp)
1172 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%vec_pot)
1173 class is (mxll_b_field_to_matter_t)
1174 call interaction%interpolate(time, field_tmp)
1175 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%b_field)
1176 end select
1177 end do
1178
1179 safe_deallocate_a(field_tmp)
1181
1182 end subroutine get_fields_from_interaction
1183
1184
1186 subroutine electrons_ground_state_run_system(sys, from_scratch)
1187 class(electrons_t), intent(inout) :: sys
1188 logical, intent(inout) :: from_scratch
1189
1191
1192 call electrons_ground_state_run(sys%namespace, sys%mc, sys%gr, sys%ions, &
1193 sys%ext_partners, sys%st, sys%ks, sys%hm, sys%outp, sys%space, from_scratch)
1194
1196
1198
1199
1200 subroutine electrons_finalize(sys)
1201 type(electrons_t), intent(inout) :: sys
1202
1203 type(partner_iterator_t) :: iter
1204 class(interaction_partner_t), pointer :: partner
1205
1206 push_sub(electrons_finalize)
1207
1208 if (associated(sys%algo)) then
1209 select type (algo => sys%algo)
1210 class is (propagator_t)
1211 call td_end_run(sys%td, sys%st, sys%hm)
1212 call td_end(sys%td)
1213 class is(minimizer_algorithm_t)
1214 call electrons_gs_cleanup(sys%ks, sys%scf, sys%rdm, sys%st, sys%hm)
1215 end select
1216 end if
1217
1218 if (sys%ks%theory_level /= independent_particles) then
1219 call poisson_async_end(sys%hm%psolver, sys%mc)
1220 end if
1221
1222 call iter%start(sys%ext_partners)
1223 do while (iter%has_next())
1224 partner => iter%get_next()
1225 safe_deallocate_p(partner)
1226 end do
1227 call sys%ext_partners%empty()
1228
1229 safe_deallocate_p(sys%xc_interaction)
1230
1231 call hamiltonian_elec_end(sys%hm)
1232
1233 nullify(sys%gfield)
1234 nullify(sys%lasers)
1235
1236 call multicomm_end(sys%mc)
1237
1238 call xc_oep_photon_end(sys%ks%oep_photon)
1239 if (sys%ks%has_photons) then
1240 call mf_end(sys%ks%pt_mx)
1241 end if
1242
1243 call sys%dipole%end()
1244
1245 call v_ks_end(sys%ks)
1246
1247 call states_elec_end(sys%st)
1248
1249 deallocate(sys%ions)
1250 safe_deallocate_p(sys%photons)
1251
1252 call kpoints_end(sys%kpoints)
1253
1254 call grid_end(sys%gr)
1255
1256 call system_end(sys)
1257
1258 pop_sub(electrons_finalize)
1259 end subroutine electrons_finalize
1260
1261
1262end module electrons_oct_m
1263
1264!! Local Variables:
1265!! mode: f90
1266!! coding: utf-8
1267!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:173
integer, parameter, public mask_absorbing
This module defines the abstract interfact for algorithm factories.
This module implements the basic elements defining algorithms.
Definition: algorithm.F90:143
character(len=algo_label_len), parameter, public iteration_done
Definition: algorithm.F90:174
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_serial
single domain, all states, k-points on a single processor
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public current_calculate(this, namespace, gr, hm, space, st)
Compute total electronic current density.
Definition: current.F90:373
subroutine, public current_init(this, namespace)
Definition: current.F90:181
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
This modules implements the dipole moment of the matter system.
Definition: dipole.F90:110
A set of subroutines for performing the parts of a ground state calculation with an electrons system....
subroutine, public electrons_ground_state_run(namespace, mc, gr, ions, ext_partners, st, ks, hm, outp, space, fromScratch)
Run a ground state calculation for a system of electrons.
subroutine, public electrons_gs_allocate_wavefunctions(namespace, gr, st, hm, scf, ks, ions)
subroutine, public electrons_gs_initialize(namespace, scf, rdm, gr, mc, st, hm, ions, ks, space, ext_partners, fromScratch)
subroutine, public electrons_gs_load_from_restart(namespace, scf, gr, mc, st, hm, ks, space, ions, ext_partners, fromScratch)
subroutine, public electrons_gs_cleanup(ks, scf, rdm, st, hm)
subroutine electrons_initialize(this)
Definition: electrons.F90:673
logical function electrons_restart_read_data(this)
Definition: electrons.F90:1187
subroutine electrons_init_parallelization(this, grp)
Definition: electrons.F90:447
logical function electrons_process_is_slave(this)
Definition: electrons.F90:1045
subroutine electrons_init_interaction(this, interaction)
Definition: electrons.F90:383
subroutine electrons_finalize(sys)
Definition: electrons.F90:1296
subroutine electrons_exec_end_of_timestep_tasks(this, prop)
Definition: electrons.F90:1056
logical function electrons_do_algorithmic_operation(this, operation, updated_quantities)
Definition: electrons.F90:718
subroutine electrons_new_algorithm(this, factory)
Definition: electrons.F90:640
subroutine electrons_output_write(this)
Definition: electrons.F90:1010
subroutine get_fields_from_interaction(this, time)
Definition: electrons.F90:1240
subroutine electrons_algorithm_start(this)
Definition: electrons.F90:698
subroutine electrons_output_finish(this)
Definition: electrons.F90:1036
subroutine electrons_output_start(this)
Definition: electrons.F90:1001
subroutine electrons_init_interaction_as_partner(partner, interaction)
Definition: electrons.F90:961
subroutine electrons_ground_state_run_system(sys, from_scratch)
Run a ground state calculation for a system of electrons.
Definition: electrons.F90:1282
class(electrons_t) function, pointer electrons_constructor(namespace, calc_mode_id)
@ brief Instantiate an instance of an electrons system
Definition: electrons.F90:288
subroutine electrons_update_kinetic_energy(this)
Definition: electrons.F90:1224
logical function electrons_is_tolerance_reached(this, tol)
Definition: electrons.F90:918
subroutine electrons_copy_quantities_to_interaction(partner, interaction)
Definition: electrons.F90:979
subroutine electrons_restart_write_data(this)
Definition: electrons.F90:1158
subroutine electrons_update_quantity(this, label)
Definition: electrons.F90:930
subroutine, public elf_init(namespace)
Definition: elf.F90:146
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,...
real(real64) function, public zenergy_calc_electronic(namespace, hm, der, st, terms)
real(real64) function, public denergy_calc_electronic(namespace, hm, der, st, terms)
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
This module implements the field transfer.
subroutine, public forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
Definition: forces.F90:340
subroutine, public gauge_field_get_force(this, gr, spin_channels, current, lrc)
subroutine, public gauge_field_do_algorithmic_operation(this, operation, dt, time)
subroutine, public gauge_field_check_symmetries(this, kpoints)
logical pure function, public gauge_field_is_propagated(this)
logical pure function, public gauge_field_is_used(this)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:237
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:237
integer, parameter, public kohn_sham_dft
Definition: global.F90:237
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_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:197
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:473
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:500
integer, parameter, public term_kinetic
subroutine, public zvmask(mesh, hm, st)
subroutine, public hamiltonian_elec_end(hm)
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public hamiltonian_elec_init(hm, namespace, space, gr, ions, ext_partners, st, theory_level, xc, mc, kpoints, need_exchange, xc_photons)
integer, parameter, public mxll_vec_pot_to_matter
integer, parameter, public mxll_b_field_to_matter
integer, parameter, public mxll_e_field_to_matter
integer, parameter, public current_to_mxll_field
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
This module defines classes and functions for interaction partners.
subroutine, public ion_dynamics_propagate_vel(this, ions, atoms_moved)
subroutine, public kpoints_end(this)
Definition: kpoints.F90:1038
subroutine, public kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
Definition: kpoints.F90:334
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1058
A module to handle KS potential, without the external potential.
subroutine, public lasers_check_symmetries(this, kpoints)
Definition: lasers.F90:557
subroutine, public lasers_parse_external_fields(this)
Definition: lasers.F90:246
subroutine, public lasers_generate_potentials(this, mesh, space, latt)
Definition: lasers.F90:445
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:798
System information (time, memory, sysname)
Definition: loct.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public mesh_check_symmetries(mesh, symm, periodic_dim)
Definition: mesh.F90:835
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
Definition: mesh.F90:386
real(real64) pure function, public mesh_global_memory(mesh)
Definition: mesh.F90:791
real(real64) pure function, public mesh_local_memory(mesh)
Definition: mesh.F90:802
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_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module implements the basic minimizer framework.
character(len=algo_label_len), parameter, public gs_scf_start
character(len=algo_label_len), parameter, public gs_scf_finish
character(len=algo_label_len), parameter, public gs_scf_iteration
general module for modelmb particles
subroutine mpi_grp_copy(mpi_grp_out, mpi_grp_in)
MPI_THREAD_FUNNELED allows for calls to MPI from an OMP region if the thread is the team master.
Definition: mpi.F90:383
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
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:697
logical pure function, public multicomm_is_slave(this)
Definition: multicomm.F90:837
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:266
integer, parameter, public mxll_field_trans
integer, parameter, public length_gauge_dipole
integer, parameter, public no_maxwell_coupling
integer, parameter, public velocity_gauge_dipole
integer, parameter, public multipolar_expansion
integer, parameter, public full_minimal_coupling
this module contains the low-level part of the output system
Definition: output_low.F90:117
this module contains the output system
Definition: output.F90:117
logical function, public output_need_exchange(outp)
Definition: output.F90:880
subroutine, public output_init(outp, namespace, space, st, gr, nst, ks)
Definition: output.F90:210
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public pes_calc(pes, namespace, space, mesh, st, dt, iter, der, kpoints, ext_partners, stopping)
Definition: pes.F90:271
subroutine, public pes_output(pes, namespace, space, gr, st, iter, outp, dt, ions)
Definition: pes.F90:297
subroutine, public mf_init(this, gr, st, ions, pt_mode)
subroutine, public mf_end(this)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public poisson_async_init(this, mc)
Definition: poisson.F90:1126
subroutine, public poisson_slave_work(this, namespace)
Definition: poisson.F90:1154
subroutine, public poisson_async_end(this, mc)
Definition: poisson.F90:1138
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
subroutine, public propagation_ops_elec_restore_ions(wo, ions_dyn, ions)
subroutine, public propagation_ops_elec_propagate_gauge_field(wo, gfield, dt, time, save_gf)
subroutine, public propagation_ops_elec_propagate_ions_and_cell(gr, hm, st, namespace, space, ions_dyn, ions, mc, time, dt_ions)
subroutine, public propagation_ops_elec_update_hamiltonian(namespace, space, st, mesh, hm, ext_partners, time)
subroutine, public propagation_ops_elec_interpolate_get(hm, vks_old)
subroutine, public propagation_ops_elec_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, op)
subroutine, public propagation_ops_elec_move_ions(wo, gr, hm, st, namespace, space, ions_dyn, ions, ext_partners, mc, time, dt, save_pos)
subroutine, public propagation_ops_elec_exp_apply(te, namespace, st, mesh, hm, dt, op)
character(len=algo_label_len), parameter, public aetrs_start
character(len=algo_label_len), parameter, public aetrs_finish
character(len=algo_label_len), parameter, public aetrs_extrapolate
character(len=algo_label_len), parameter, public aetrs_first_half
character(len=algo_label_len), parameter, public aetrs_second_half
character(len=algo_label_len), parameter, public bomd_start
character(len=algo_label_len), parameter, public bomd_elec_scf
character(len=algo_label_len), parameter, public bomd_finish
character(len=algo_label_len), parameter, public expmid_extrapolate
character(len=algo_label_len), parameter, public expmid_finish
character(len=algo_label_len), parameter, public expmid_start
character(len=algo_label_len), parameter, public expmid_propagate
This module implements the basic propagator framework.
Definition: propagator.F90:119
character(len=30), parameter, public verlet_compute_acc
type(algorithmic_operation_t), parameter, public op_verlet_compute_acc
character(len=30), parameter, public verlet_update_pos
type(algorithmic_operation_t), parameter, public op_verlet_compute_vel
character(len=30), parameter, public verlet_compute_vel
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
Implementation details for regridding.
Definition: regridding.F90:172
subroutine, public scf_finish(scf, namespace, space, gr, ions, ext_partners, st, ks, hm, iter, outp)
Definition: scf.F90:1290
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_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
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
real(real64) function, public states_elec_wfns_memory(st, mesh)
return the memory usage of a states_elec_t object
subroutine, public states_elec_distribute_nodes(st, namespace, mc)
@Brief. Distribute states over the processes for states parallelization
subroutine, public states_elec_densities_init(st, gr)
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public kpoints_distribute(this, mc)
distribute k-points over the nodes in the corresponding communicator
subroutine, public states_elec_exec_init(st, namespace, mc)
Further initializations.
subroutine, public states_elec_init(st, namespace, space, valence_charge, kpoints, calc_mode_id)
Initialize a new states_elec_t object.
subroutine, public states_elec_allocate_current(st, space, mesh)
This module implements the calculation of the stress tensor.
Definition: stress.F90:120
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
Definition: stress.F90:188
This module implements the abstract system type.
Definition: system.F90:120
subroutine, public system_algorithm_start(this)
Definition: system.F90:1022
subroutine, public system_end(this)
Definition: system.F90:1151
subroutine, public system_new_algorithm(this, factory)
Definition: system.F90:947
Definition: td.F90:116
subroutine, public td_end(td)
Definition: td.F90:641
subroutine, public td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
Definition: td.F90:1226
subroutine, public td_end_run(td, st, hm)
Definition: td.F90:660
subroutine, public td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp)
Definition: td.F90:241
subroutine, public td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
Definition: td.F90:576
logical function, public td_get_from_scratch(td)
Definition: td.F90:1568
subroutine, public td_set_from_scratch(td, from_scratch)
Definition: td.F90:1579
subroutine, public td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
Definition: td.F90:1357
subroutine, public td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
Definition: td.F90:610
subroutine, public td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, from_scratch)
Definition: td.F90:908
subroutine, public td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
Definition: td.F90:1190
subroutine, public td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, dt)
Definition: td_write.F90:1250
subroutine, public td_write_data(writ)
Definition: td_write.F90:1216
subroutine, public td_write_iter(writ, namespace, space, outp, gr, st, hm, ions, ext_partners, kick, ks, dt, iter, mc, recalculate_gs)
Definition: td_write.F90:1043
integer, parameter, public out_separate_forces
Definition: td_write.F90:203
This module defines the unit system, used for input and output.
type(unit_t), public unit_megabytes
For large amounts of data (natural code units are bytes)
subroutine, public v_ks_end(ks)
Definition: v_ks.F90:631
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:752
subroutine, public v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
Definition: v_ks.F90:254
integer, parameter, public xc_oep_x
Exact exchange.
integer, parameter, public func_x
Definition: xc.F90:116
integer, parameter, public oep_level_none
the OEP levels
Definition: xc_oep.F90:174
subroutine, public xc_oep_photon_init(oep, namespace, family, gr, st, mc, space)
subroutine, public xc_oep_photon_end(oep)
Abstract class for the algorithm factories.
Descriptor of one algorithmic operation.
Definition: algorithm.F90:165
Class to transfer a current to a Maxwell field.
Extension of space that contains the knowledge of the spin dimension.
Class describing the electron system.
Definition: electrons.F90:220
class defining the field_transfer interaction
These class extend the list and list iterator to make an interaction list.
abstract interaction class
abstract class for general interaction partners
surrogate interaction class to avoid circular dependencies between modules.
Abstract class implementing minimizers.
This is defined even when running serial.
Definition: mpi.F90:144
class to transfer a Maxwell B field to a matter system
class to transfer a Maxwell field to a medium
class to transfer a Maxwell vector potential to a medium
Implements a propagator for Approximate ETRS.
Implements a propagator for Born-Oppenheimer molecular dynamics.
Implements the explicit exponential midpoint propagator (without predictor-corrector)
Abstract class implementing propagators.
Definition: propagator.F90:144
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:173
Abstract class for systems.
Definition: system.F90:174
int true(void)