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
119
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 final :: electrons_finalize
179 end type electrons_t
180
181 interface electrons_t
182 procedure electrons_constructor
183 end interface electrons_t
184
185contains
186
187 !----------------------------------------------------------
188 function electrons_constructor(namespace, generate_epot) result(sys)
189 class(electrons_t), pointer :: sys
190 type(namespace_t), intent(in) :: namespace
191 logical, optional, intent(in) :: generate_epot
192
193 integer :: iatom
194 type(lattice_vectors_t) :: latt_inp
195 logical :: has_photons
196
197 push_sub_with_profile(electrons_constructor)
198
199 allocate(sys)
200
201 sys%namespace = namespace
202
203 call messages_obsolete_variable(sys%namespace, 'SystemName')
204
205 sys%space = electron_space_t(sys%namespace)
206 call sys%space%write_info(sys%namespace)
207 if (sys%space%has_mixed_periodicity()) then
208 call messages_experimental('Support for mixed periodicity systems')
209 end if
210
211 sys%ions => ions_t(sys%namespace, latt_inp=latt_inp)
212
213 call grid_init_stage_1(sys%gr, sys%namespace, sys%space, sys%ions%symm, latt_inp, sys%ions%natoms, sys%ions%pos)
214
215 if (sys%space%is_periodic()) then
216 call sys%ions%latt%write_info(sys%namespace)
217 end if
219 ! Sanity check for atomic coordinates
220 do iatom = 1, sys%ions%natoms
221 ! Using the same tolerance of fold_into_cell to avoid problems with mixed periodicities
222 ! as the default tolerance of contains_point is too tight
223 if (.not. sys%gr%box%contains_point(sys%ions%pos(:, iatom), tol=1.0e-6_real64)) then
224 if (sys%space%periodic_dim /= sys%space%dim) then
225 write(message(1), '(a,i5,a)') "Atom ", iatom, " is outside the box."
226 call messages_warning(1, namespace=sys%namespace)
227 end if
228 end if
229 end do
231 ! we need k-points for periodic systems
232 call kpoints_init(sys%kpoints, sys%namespace, sys%gr%symm, sys%space%dim, sys%space%periodic_dim, sys%ions%latt)
234 call states_elec_init(sys%st, sys%namespace, sys%space, sys%ions%val_charge(), sys%kpoints)
235 call sys%st%write_info(sys%namespace)
236
237 ! if independent particles in N dimensions are being used, need to initialize them
238 ! after masses are set to 1 in grid_init_stage_1 -> derivatives_init
239 call sys%st%modelmbparticles%copy_masses(sys%gr%der%masses)
240
241 call elf_init(sys%namespace)
243 sys%generate_epot = optional_default(generate_epot, .true.)
244
245 call sys%dipole%init(sys%space)
246
247 sys%supported_interactions_as_partner = [current_to_mxll_field]
248
249 call sys%quantities%add(quantity_t("wavefunctions", updated_on_demand = .false.))
250 call sys%quantities%add(quantity_t("current", updated_on_demand = .true., parents=["wavefunctions"]))
251 call sys%quantities%add(quantity_t("dipole", updated_on_demand = .true., parents=["wavefunctions"]))
252 call current_init(sys%current_calculator, sys%namespace)
253
254 !%Variable EnablePhotons
255 !%Type logical
256 !%Default no
257 !%Section Hamiltonian
258 !%Description
259 !% This variable can be used to enable photons in several types of calculations.
260 !% It can be used to activate the one-photon OEP formalism.
261 !% In the case of CalculationMode = casida, it enables photon modes as
262 !% described in ACS Photonics 2019, 6, 11, 2757-2778.
263 !% Finally, if set to yes when solving the ferquency-dependent Sternheimer
264 !% equation, the photons are coupled to the electronic subsystem.
265 !%End
266 call messages_obsolete_variable(namespace, 'OEPPtX', 'EnablePhotons')
267 call parse_variable(namespace, 'EnablePhotons', .false., has_photons)
268 if (has_photons) then
269 call messages_experimental("EnablePhotons = yes")
270 sys%photons => photons_t(sys%namespace)
271 else
272 nullify(sys%photons)
273 end if
274
275 pop_sub_with_profile(electrons_constructor)
276 end function electrons_constructor
277
278 ! ---------------------------------------------------------
279 subroutine electrons_init_interaction(this, interaction)
280 class(electrons_t), target, intent(inout) :: this
281 class(interaction_t), intent(inout) :: interaction
282
283 real(real64) :: dmin
284 integer :: rankmin, depth
285 logical :: mxll_interaction_present = .false.
286 logical :: calc_dipole
287
288 push_sub(electrons_init_interactions)
289
290 select type (interaction)
292 call interaction%init(this%gr, 3)
293 mxll_interaction_present = .true.
294 interaction%type = mxll_field_trans
296 call interaction%init(this%gr, 3)
297 mxll_interaction_present = .true.
299 call interaction%init(this%gr, 3)
300 mxll_interaction_present = .true.
301 class default
302 message(1) = "Trying to initialize an unsupported interaction by the electrons."
303 call messages_fatal(1, namespace=this%namespace)
304 end select
305
306 if (mxll_interaction_present) then
307 calc_dipole = any(this%hm%mxll%coupling_mode == &
309
310 if (calc_dipole) then
311 assert(this%space%periodic_dim == 0)
312 this%hm%mxll%calc_field_dip = .true.
313 this%hm%mxll%center_of_mass(1:3) = this%ions%center_of_mass()
314 this%hm%mxll%center_of_mass_ip = mesh_nearest_point(this%gr, this%hm%mxll%center_of_mass, dmin, rankmin)
315 this%hm%mxll%center_of_mass_rankmin = rankmin
316 end if
317 end if
318
319 ! set interpolation depth for field-transfer interactions
320 select type (interaction)
321 class is (field_transfer_t)
322 ! interpolation depth depends on the propagator
323 select type (algo => this%algo)
324 type is (propagator_exp_mid_t)
325 depth = 3
326 type is (propagator_aetrs_t)
327 depth = 3
328 type is (propagator_bomd_t)
329 depth = 1
330 class default
331 message(1) = "The chosen algorithm does not yet support interaction interpolation"
332 call messages_fatal(1, namespace=this%namespace)
333 end select
334
335 call interaction%init_interpolation(depth, interaction%label)
336 end select
337
338 pop_sub(electrons_init_interactions)
339 end subroutine electrons_init_interaction
340
341 ! ---------------------------------------------------------
342 subroutine electrons_init_parallelization(this, grp)
343 class(electrons_t), intent(inout) :: this
344 type(mpi_grp_t), intent(in) :: grp
345
346 integer(int64) :: index_range(4)
347 real(real64) :: mesh_global, mesh_local, wfns
348 integer :: idir
349 real(real64) :: spiral_q(3), spiral_q_red(3)
350 type(block_t) :: blk
351
353
354 call mpi_grp_copy(this%grp, grp)
355
356 ! store the ranges for these two indices (serves as initial guess
357 ! for parallelization strategy)
358 index_range(1) = this%gr%np_global ! Number of points in mesh
359 index_range(2) = this%st%nst ! Number of states
360 index_range(3) = this%st%nik ! Number of k-points
361 index_range(4) = 100000 ! Some large number
362
363 ! create index and domain communicators
364 call multicomm_init(this%mc, this%namespace, this%grp, calc_mode_par, &
365 mpi_world%size, index_range, (/ 5000, 1, 1, 1 /))
366
367 call this%ions%partition(this%mc)
368 call kpoints_distribute(this%st, this%mc)
369 call states_elec_distribute_nodes(this%st, this%namespace, this%mc)
370
371
372 if (parse_is_defined(this%namespace, 'TDMomentumTransfer') .or. &
373 parse_is_defined(this%namespace, 'TDReducedMomentumTransfer')) then
374 spiral_q = m_zero
375 if (parse_block(this%namespace, 'TDMomentumTransfer', blk) == 0) then
376 do idir = 1, this%space%dim
377 call parse_block_float(blk, 0, idir - 1, spiral_q(idir))
378 end do
379 else if(parse_block(this%namespace, 'TDReducedMomentumTransfer', blk) == 0) then
380 do idir = 1, this%space%dim
381 call parse_block_float(blk, 0, idir - 1, spiral_q_red(idir))
382 end do
383 call kpoints_to_absolute(this%kpoints%latt, spiral_q_red(1:this%space%dim), spiral_q(1:this%space%dim))
384 end if
385 call parse_block_end(blk)
386 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc, spiral_q)
387 else
388 call grid_init_stage_2(this%gr, this%namespace, this%space, this%mc)
389 end if
390
391 if (this%st%symmetrize_density) then
392 call mesh_check_symmetries(this%gr, this%gr%symm, this%ions%space%periodic_dim)
393 end if
394
395 call output_init(this%outp, this%namespace, this%space, this%st, this%gr, this%st%nst, this%ks)
396 call states_elec_densities_init(this%st, this%gr)
397 call states_elec_exec_init(this%st, this%namespace, this%mc)
398
399 if (associated(this%photons)) then
400 this%ks%has_photons = .true.
401 end if
402
403 call v_ks_init(this%ks, this%namespace, this%gr, this%st, this%ions, this%mc, this%space, this%kpoints)
404 if (this%ks%theory_level == kohn_sham_dft .or. this%ks%theory_level == generalized_kohn_sham_dft) then
405 this%xc_interaction => xc_interaction_t(this)
406 end if
407
408 ! For the moment the photons are initialized here
409
410 if(this%ks%has_photons) then
411 this%ks%pt => this%photons%modes
412 ! Temporary creation that should go in the system initialization later
413 call photon_mode_set_n_electrons(this%photons%modes, this%st%qtot)
414 write(message(1), '(a,i5,a)') 'Happy to have ', this%photons%modes%nmodes, ' photon modes with us.'
415 call messages_info(1)
416 call mf_init(this%ks%pt_mx, this%gr, this%st, this%ions, this%ks%pt)
417 ! OEP for photons
418 if(bitand(this%ks%xc_family, xc_family_oep) /= 0 .and. this%ks%xc%functional(func_x,1)%id == xc_oep_x) then
419 this%ks%oep_photon%level = this%ks%oep%level
420 call xc_oep_photon_init(this%ks%oep_photon, this%namespace, this%ks%xc_family, this%gr, this%st, this%mc, this%space)
421 else
422 this%ks%oep_photon%level = oep_level_none
423 end if
424
425 end if
426
427
428 ! Temporary place for the initialization of the lasers
429 this%lasers => lasers_t(this%namespace)
430 call lasers_parse_external_fields(this%lasers)
431 call lasers_generate_potentials(this%lasers, this%gr, this%space, this%ions%latt)
432 if(this%lasers%no_lasers > 0) then
433 call this%ext_partners%add(this%lasers)
434 call lasers_check_symmetries(this%lasers, this%kpoints)
435 else
436 deallocate(this%lasers)
437 end if
438
439 ! Temporary place for the initialization of the gauge field
440 this%gfield => gauge_field_t(this%namespace, this%ions%latt%rcell_volume)
441 if(gauge_field_is_used(this%gfield)) then
442 call this%ext_partners%add(this%gfield)
443 call gauge_field_check_symmetries(this%gfield, this%kpoints)
444 else
445 deallocate(this%gfield)
446 end if
447
448 call hamiltonian_elec_init(this%hm, this%namespace, this%space, this%gr, this%ions, this%ext_partners, &
449 this%st, this%ks%theory_level, this%ks%xc, this%mc, this%kpoints, &
450 need_exchange = output_need_exchange(this%outp) .or. this%ks%oep%level /= oep_level_none, &
451 xc_photons = this%ks%xc_photons )
452
453 if (this%hm%pcm%run_pcm .and. this%mc%par_strategy /= p_strategy_serial .and. this%mc%par_strategy /= p_strategy_states) then
454 call messages_experimental('Parallel in domain calculations with PCM')
455 end if
456
457 ! Print memory requirements
458 call messages_print_with_emphasis(msg='Approximate memory requirements', namespace=this%namespace)
459
460 mesh_global = mesh_global_memory(this%gr)
461 mesh_local = mesh_local_memory(this%gr)
462
463 call messages_write('Mesh')
464 call messages_new_line()
465 call messages_write(' global :')
466 call messages_write(mesh_global, units = unit_megabytes, fmt = '(f10.1)')
467 call messages_new_line()
468 call messages_write(' local :')
469 call messages_write(mesh_local, units = unit_megabytes, fmt = '(f10.1)')
470 call messages_new_line()
471 call messages_write(' total :')
472 call messages_write(mesh_global + mesh_local, units = unit_megabytes, fmt = '(f10.1)')
473 call messages_new_line()
474 call messages_info(namespace=this%namespace)
475
476 wfns = states_elec_wfns_memory(this%st, this%gr)
477 call messages_write('States')
478 call messages_new_line()
479 call messages_write(' real :')
480 call messages_write(wfns, units = unit_megabytes, fmt = '(f10.1)')
481 call messages_write(' (par_kpoints + par_states + par_domains)')
482 call messages_new_line()
483 call messages_write(' complex :')
484 call messages_write(2.0_8*wfns, units = unit_megabytes, fmt = '(f10.1)')
485 call messages_write(' (par_kpoints + par_states + par_domains)')
486 call messages_new_line()
487 call messages_info(namespace=this%namespace)
488
489 call messages_print_with_emphasis(namespace=this%namespace)
490
491 if (this%generate_epot) then
492 message(1) = "Info: Generating external potential"
493 call messages_info(1, namespace=this%namespace)
494 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
495 this%ext_partners, this%st)
496 message(1) = " done."
497 call messages_info(1, namespace=this%namespace)
498 end if
499
500 if (this%ks%theory_level /= independent_particles) then
501 call poisson_async_init(this%hm%psolver, this%mc)
502 ! slave nodes do not call the calculation routine
503 if (multicomm_is_slave(this%mc)) then
504 !for the moment we only have one type of slave
505 call poisson_slave_work(this%hm%psolver, this%namespace)
506 end if
507 end if
508
509 allocate(this%supported_interactions(0))
510 select case (this%hm%mxll%coupling_mode)
512 this%supported_interactions = [this%supported_interactions, mxll_e_field_to_matter]
514 this%supported_interactions = [this%supported_interactions, mxll_vec_pot_to_matter]
515 if (this%hm%mxll%add_zeeman) then
516 this%supported_interactions = [this%supported_interactions, mxll_b_field_to_matter]
517 end if
519 if (this%hm%mxll%add_electric_dip .or. this%hm%mxll%add_electric_quad) then
520 this%supported_interactions = [this%supported_interactions, mxll_e_field_to_matter]
521 end if
522 if (this%hm%mxll%add_magnetic_dip) then
523 this%supported_interactions = [this%supported_interactions, mxll_b_field_to_matter]
524 end if
526 ! Do not initialize any interaction with Maxwell
527 case default
528 message(1) = "Unknown maxwell-matter coupling"
529 call messages_fatal(1, namespace=this%namespace)
530 end select
531
533 end subroutine electrons_init_parallelization
534
535 ! ---------------------------------------------------------
536 subroutine electrons_new_algorithm(this, factory)
537 class(electrons_t), intent(inout) :: this
538 class(algorithm_factory_t), intent(in) :: factory
539
541
542 call system_new_algorithm(this, factory)
543
544 select type (algo => this%algo)
545 class is (propagator_t)
546
547 call td_init(this%td, this%namespace, this%space, this%gr, this%ions, this%st, this%ks, &
548 this%hm, this%ext_partners, this%outp)
549
550 ! this corresponds to the first part of td_init_run
551 call td_allocate_wavefunctions(this%td, this%namespace, this%mc, this%gr, this%ions, this%st, &
552 this%hm, this%space)
553 call td_init_gaugefield(this%td, this%namespace, this%gr, this%st, this%ks, this%hm, &
554 this%ext_partners, this%space)
555
556 class is (minimizer_algorithm_t)
557
558 call gs_allocate_wavefunctions(this%namespace, this%gr, this%st, this%hm, this%scf, this%ks, &
559 this%ions, this%mc, this%space)
560
561 class default
562 assert(.false.)
563 end select
564
566 end subroutine electrons_new_algorithm
567
568 ! ---------------------------------------------------------
569 subroutine electrons_initialize(this)
570 class(electrons_t), intent(inout) :: this
571
572 push_sub(electrons_initialize)
573
574 select type (algo => this%algo)
575 class is (propagator_t)
576 call td_set_from_scratch(this%td, .true.)
577 call td_load_restart_from_gs(this%td, this%namespace, this%space, this%mc, this%gr, &
578 this%ext_partners, this%st, this%ks, this%hm)
579
580 class is (minimizer_algorithm_t)
581
582 call gs_initialize(this%namespace, this%scf, this%rdm, this%gr, this%mc, this%st, &
583 this%hm, this%ions, this%ks, this%space, this%ext_partners, fromscratch=.true.)
584 class default
585 assert(.false.)
586 end select
587
588 pop_sub(electrons_initialize)
589 end subroutine electrons_initialize
590
591 ! ---------------------------------------------------------
592 subroutine electrons_algorithm_start(this)
593 class(electrons_t), intent(inout) :: this
594
596
597 call system_algorithm_start(this)
598
599 select type (algo => this%algo)
600 class is (propagator_t)
601
602 ! additional initialization needed for electrons
603 call td_init_with_wavefunctions(this%td, this%namespace, this%space, this%mc, this%gr, this%ions, &
604 this%ext_partners, this%st, this%ks, this%hm, this%outp, td_get_from_scratch(this%td))
605
606 end select
607
609 end subroutine electrons_algorithm_start
610
611 ! ---------------------------------------------------------
612 logical function electrons_do_algorithmic_operation(this, operation, updated_quantities) result(done)
613 class(electrons_t), intent(inout) :: this
614 class(algorithmic_operation_t), intent(in) :: operation
615 character(len=:), allocatable, intent(out) :: updated_quantities(:)
616
617 logical :: update_energy_
618 type(gauge_field_t), pointer :: gfield
619 real(real64) :: time
620 integer :: iter
621
623 call profiling_in(trim(this%namespace%get())//":"//trim(operation%id))
624
625 update_energy_ = .true.
626
627 ! kick at t > 0 still missing!
628
629 done = .true.
630 select type (algo => this%algo)
631 class is (propagator_t)
632 time = algo%iteration%value()
633 select case (operation%id)
634 case (aetrs_first_half)
635 ! propagate half of the time step with H(time)
636 call get_fields_from_interaction(this, time)
637 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
638 this%hm, this%ext_partners, time)
639 call propagation_ops_elec_exp_apply(this%td%tr%te, this%namespace, this%st, this%gr, this%hm, m_half*algo%dt)
640
641 case (aetrs_extrapolate)
642 ! Do the extrapolation of the Hamiltonian
643 ! First the extrapolation of the potential
644 call this%hm%ks_pot%interpolation_new(this%td%tr%vks_old, time+algo%dt, algo%dt)
645
646 !Get the potentials from the interpolator
647 call propagation_ops_elec_interpolate_get(this%hm, this%gr, this%td%tr%vks_old)
648
649 ! Second the ions and gauge field which later on will be treated as extrapolation
650 ! of interactions, but this is not yet possible
651
652 ! move the ions to time t
653 call propagation_ops_elec_move_ions(this%td%tr%propagation_ops_elec, this%gr, this%hm, &
654 this%st, this%namespace, this%space, this%td%ions_dyn, this%ions, this%ext_partners, &
655 this%mc, time+algo%dt, algo%dt)
656
657 !Propagate gauge field
658 gfield => list_get_gauge_field(this%ext_partners)
659 if(associated(gfield)) then
660 call propagation_ops_elec_propagate_gauge_field(this%td%tr%propagation_ops_elec, gfield, &
661 algo%dt, time+algo%dt)
662 end if
663
664 !Update Hamiltonian and current
665 call get_fields_from_interaction(this, time+algo%dt)
666 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
667 this%hm, this%ext_partners, time+algo%dt)
668
669 case (aetrs_second_half)
670 !Do the time propagation for the second half of the time step
671 call propagation_ops_elec_fuse_density_exp_apply(this%td%tr%te, this%namespace, this%st, &
672 this%gr, this%hm, m_half*algo%dt)
673
674 updated_quantities = ["wavefunctions"]
675
676 case (expmid_extrapolate)
677 ! the half step of this propagator screws with the gauge field kick
678 gfield => list_get_gauge_field(this%ext_partners)
679 if(associated(gfield)) then
680 assert(gauge_field_is_propagated(gfield) .eqv. .false.)
681 end if
682
683 ! Do the extrapolation of the Hamiltonian
684 ! First the extrapolation of the potential
685 call this%hm%ks_pot%interpolation_new(this%td%tr%vks_old, time+algo%dt, algo%dt)
686
687 ! get the potentials from the interpolator
688 call this%hm%ks_pot%interpolate_potentials(this%td%tr%vks_old, 3, time+algo%dt, algo%dt, time + algo%dt/m_two)
689
690 ! Second the ions which later on will be treated as extrapolation of interactions,
691 ! but this is not yet possible
692 ! move the ions to the half step
693 call propagation_ops_elec_move_ions(this%td%tr%propagation_ops_elec, this%gr, this%hm, this%st, &
694 this%namespace, this%space, this%td%ions_dyn, this%ions, this%ext_partners, &
695 this%mc, time + m_half*algo%dt, m_half*algo%dt, save_pos=.true.)
696
697 call get_fields_from_interaction(this, time + m_half*algo%dt)
698 call propagation_ops_elec_update_hamiltonian(this%namespace, this%space, this%st, this%gr, &
699 this%hm, this%ext_partners, time + m_half*algo%dt)
700
701 case (expmid_propagate)
702 ! Do the actual propagation step
703 call propagation_ops_elec_fuse_density_exp_apply(this%td%tr%te, this%namespace, this%st, &
704 this%gr, this%hm, algo%dt)
706 ! restore to previous time
707 call propagation_ops_elec_restore_ions(this%td%tr%propagation_ops_elec, this%td%ions_dyn, this%ions)
708
709 updated_quantities = ["wavefunctions"]
710
711 case (bomd_start)
712 call scf_init(this%scf, this%namespace, this%gr, this%ions, this%st, this%mc, this%hm, this%space)
713 ! the ions are propagated inside the propagation step already, so no need to do it at the end
714 this%ions_propagated = .true.
715
716 case (verlet_update_pos)
717 ! move the ions to time t
718 call propagation_ops_elec_propagate_ions_and_cell(this%gr, this%hm, this%st, this%namespace, this%space, &
719 this%td%ions_dyn, this%ions, this%mc, time+algo%dt, algo%dt)
720
721 case (bomd_elec_scf)
722 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
723 this%ext_partners, this%st, time = time+algo%dt)
724 ! now calculate the eigenfunctions
725 call scf_run(this%scf, this%namespace, this%space, this%mc, this%gr, this%ions, &
726 this%ext_partners, this%st, this%ks, this%hm, verbosity = verb_compact)
727 ! TODO: Check if this call is realy needed. - NTD
728 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
729 this%ext_partners, this%st, time = time+algo%dt)
730
731 ! update Hamiltonian and eigenvalues (fermi is *not* called)
732 call v_ks_calc(this%ks, this%namespace, this%space, this%hm, this%st, this%ions, this%ext_partners, &
733 calc_eigenval = .true., time = time+algo%dt, calc_energy = .true.)
734
735 ! Get the energies.
736 call energy_calc_total(this%namespace, this%space, this%hm, this%gr, this%st, this%ext_partners, iunit = -1)
737
738 updated_quantities = ["wavefunctions"]
739
740 case (verlet_compute_acc)
741 ! Do nothing, forces are computing in scf_run
742 if (this%td%ions_dyn%cell_relax()) then
743 assert(this%scf%calc_stress)
744 call this%td%ions_dyn%update_stress(this%ions%space, this%st%stress_tensors%total, &
745 this%ions%latt%rlattice, this%ions%latt%rcell_volume)
746 end if
747
748 case (verlet_compute_vel)
749 call ion_dynamics_propagate_vel(this%td%ions_dyn, this%ions)
750 ! TODO: Check if this call is realy needed. - NTD
751 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
752 this%ext_partners, this%st, time = time+algo%dt)
753 call this%ions%update_kinetic_energy()
754
755 case (expmid_start)
756 this%ions_propagated = .false.
757
758 case (aetrs_start)
759 ! the ions are propagated inside the propagation step already, so no need to do it at the end
760 this%ions_propagated = .true.
761
762 case (iteration_done)
764 done = .false.
765
766 case (bomd_finish)
767 call scf_end(this%scf)
768
770 case default
771 done = .false.
772 end select
773
774 class is(minimizer_algorithm_t)
775
776 ! Clocks starts at 0....
777 iter = nint(this%iteration%value()) + 1
778
779 select case(operation%id)
780 case (gs_scf_start)
781 call scf_start(this%scf, this%namespace, this%space, this%gr, this%ions, this%st, this%ks, &
782 this%hm, this%outp)
783
784 case (gs_scf_iteration)
785
786 call scf_iter(this%scf, this%namespace, this%space, this%mc, this%gr, this%ions, &
787 this%ext_partners, this%st, this%ks, this%hm, iter, outp = this%outp, &
788 restart_dump=this%scf%restart_dump)
789
790 algo%converged = scf_iter_finish(this%scf, this%namespace, this%space, this%gr, this%ions,&
791 this%st, this%ks, this%hm, iter, outp = this%outp)
792
793 updated_quantities = ["wavefunctions"]
794
795 case (gs_scf_finish)
796
797 ! Here iteration is iter-1, as the clock ticked before SCF_FINISH
798 call scf_finish(this%scf, this%namespace, this%space, this%mc, this%gr, this%ions, &
799 this%ext_partners, this%st, this%ks, this%hm, iter-1, &
800 outp = this%outp, restart_dump=this%scf%restart_dump)
801
802 case default
803 done = .false.
804 end select
805 class default
806 done = .false.
807 end select
808
809 call profiling_out(trim(this%namespace%get())//":"//trim(operation%id))
812
813 ! ---------------------------------------------------------
814 logical function electrons_is_tolerance_reached(this, tol) result(converged)
815 class(electrons_t), intent(in) :: this
816 real(real64), intent(in) :: tol
817
819
820 converged = .false.
821
824
825 ! ---------------------------------------------------------
826 subroutine electrons_update_quantity(this, label)
827 class(electrons_t), intent(inout) :: this
828 character(len=*), intent(in) :: label
829
830 class(quantity_t), pointer :: quantity
831
833 call profiling_in(trim(this%namespace%get())//":"//"UPDATE_QUANTITY")
834
835 quantity => this%quantities%get(label)
836 if(associated(quantity)) then
837 assert(quantity%updated_on_demand)
838 endif
839
840 select case (label)
841 case ("current")
842 call states_elec_allocate_current(this%st, this%space, this%gr)
843 call current_calculate(this%current_calculator, this%namespace, this%gr, &
844 this%hm, this%space, this%st)
845 case ("dipole")
846 call this%dipole%calculate(this%gr, this%ions, this%st)
847 case default
848 message(1) = "Incompatible quantity: "//trim(label)//"."
849 call messages_fatal(1, namespace=this%namespace)
850 end select
851
852 call profiling_out(trim(this%namespace%get())//":"//"UPDATE_QUANTITY")
854 end subroutine electrons_update_quantity
855
856 ! ---------------------------------------------------------
857 subroutine electrons_init_interaction_as_partner(partner, interaction)
858 class(electrons_t), intent(in) :: partner
859 class(interaction_surrogate_t), intent(inout) :: interaction
860
862
863 select type (interaction)
865 call interaction%init_from_partner(partner%gr, partner%space, partner%namespace)
866 class default
867 message(1) = "Unsupported interaction."
868 call messages_fatal(1, namespace=partner%namespace)
869 end select
870
873
874 ! ---------------------------------------------------------
875 subroutine electrons_copy_quantities_to_interaction(partner, interaction)
876 class(electrons_t), intent(inout) :: partner
877 class(interaction_surrogate_t), intent(inout) :: interaction
878
880 call profiling_in(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
881
882 select type (interaction)
884 assert(allocated(partner%st%current))
885 interaction%partner_field(:,:) = partner%st%current(1:partner%gr%np,:,1)
886 call interaction%do_mapping()
887 class default
888 message(1) = "Unsupported interaction."
889 call messages_fatal(1, namespace=partner%namespace)
890 end select
891
892 call profiling_out(trim(partner%namespace%get())//":"//"COPY_QUANTITY_INTER")
895
896 ! ---------------------------------------------------------
897 subroutine electrons_output_start(this)
898 class(electrons_t), intent(inout) :: this
899
900 push_sub(electrons_output_start)
901
903 end subroutine electrons_output_start
904
905 ! ---------------------------------------------------------
906 subroutine electrons_output_write(this)
907 class(electrons_t), intent(inout) :: this
908
909 integer :: iter
910
911 push_sub(electrons_output_write)
912 call profiling_in(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
913
914 select type (algo => this%algo)
915 class is (propagator_t)
916 iter = this%iteration%counter()
917
918 call td_write_iter(this%td%write_handler, this%namespace, this%space, this%outp, this%gr, &
919 this%st, this%hm, this%ions, this%ext_partners, this%hm%kick, this%ks, algo%dt, iter, this%mc, this%td%recalculate_gs)
920
921 if (this%outp%anything_now(iter)) then ! output
922 call td_write_output(this%namespace, this%space, this%gr, this%st, this%hm, this%ks, &
923 this%outp, this%ions, this%ext_partners, iter, algo%dt)
924 end if
925 end select
926
927 call profiling_out(trim(this%namespace%get())//":"//"OUTPUT_WRITE")
929 end subroutine electrons_output_write
930
931 ! ---------------------------------------------------------
932 subroutine electrons_output_finish(this)
933 class(electrons_t), intent(inout) :: this
934
936
938 end subroutine electrons_output_finish
939
940 ! ---------------------------------------------------------
941 logical function electrons_process_is_slave(this) result(is_slave)
942 class(electrons_t), intent(in) :: this
943
945
946 is_slave = multicomm_is_slave(this%mc)
947
949 end function electrons_process_is_slave
951 ! ---------------------------------------------------------
952 subroutine electrons_exec_end_of_timestep_tasks(this, prop)
953 class(electrons_t), intent(inout) :: this
954 class(propagator_t), intent(in) :: prop
955
956 logical :: stopping
957 logical :: generate
958 logical :: update_energy_
959 integer :: nt
960 real(real64) :: time
961 type(gauge_field_t), pointer :: gfield
962
964 call profiling_in(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
965
966 stopping = .false.
967
968 nt = this%td%iter
969 ! this is the time at the end of the timestep, as required in all routines here
970 time = prop%dt*nt
971 update_energy_ = .true.
972
973 !Apply mask absorbing boundaries
974 if (this%hm%abs_boundaries%abtype == mask_absorbing) call zvmask(this%gr, this%hm, this%st)
975
976 !Photoelectron stuff
977 if (this%td%pesv%calc_spm .or. this%td%pesv%calc_mask .or. this%td%pesv%calc_flux) then
978 call pes_calc(this%td%pesv, this%namespace, this%space, this%gr, this%st, &
979 prop%dt, nt, this%gr%der, this%hm%kpoints, this%ext_partners, stopping)
980 end if
981
982 ! For BOMD, we do not want the lines below to be executed
983 select type(prop)
984 type is(propagator_bomd_t)
985 call profiling_out(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
987 return
988 end select
989
990 ! The propagation of the ions and the gauge field is currently done here.
991 ! TODO: this code is to be moved to their own systems at some point
992 generate = .false.
993 if (this%td%ions_dyn%is_active()) then
994 if (.not. this%ions_propagated) then
995 call propagation_ops_elec_propagate_ions_and_cell(this%gr, this%hm, this%st, this%namespace, this%space, &
996 this%td%ions_dyn, this%ions, this%mc, abs(nt*prop%dt), this%td%ions_dyn%ionic_scale*prop%dt)
997 generate = .true.
998 end if
999 end if
1000
1001 gfield => list_get_gauge_field(this%ext_partners)
1002 if(associated(gfield)) then
1003 if (gauge_field_is_propagated(gfield) .and. .not. this%ions_propagated) then
1005 end if
1006 end if
1007
1008 if (generate .or. this%ions%has_time_dependent_species()) then
1009 call hamiltonian_elec_epot_generate(this%hm, this%namespace, this%space, this%gr, this%ions, &
1010 this%ext_partners, this%st, time = abs(nt*prop%dt))
1011 end if
1012
1013 call v_ks_calc(this%ks, this%namespace, this%space, this%hm, this%st, this%ions, this%ext_partners, &
1014 calc_eigenval = update_energy_, time = abs(nt*prop%dt), calc_energy = update_energy_)
1015
1016 if (update_energy_) then
1017 call energy_calc_total(this%namespace, this%space, this%hm, this%gr, this%st, this%ext_partners, iunit = -1)
1018 end if
1019
1020 ! Recalculate forces, update velocities...
1021 if (this%td%ions_dyn%ions_move() .or. this%outp%what(option__output__forces) &
1022 .or. this%td%write_handler%out(out_separate_forces)%write) then
1023 call forces_calculate(this%gr, this%namespace, this%ions, this%hm, this%ext_partners, &
1024 this%st, this%ks, t = abs(nt*prop%dt), dt = prop%dt)
1025 end if
1026
1027 if (this%td%ions_dyn%cell_relax() .or. this%outp%what(option__output__stress)) then
1028 call stress_calculate(this%namespace, this%gr, this%hm, this%st, this%ions, this%ks, this%ext_partners)
1029 end if
1030
1031 if(this%td%ions_dyn%is_active()) then
1032 call ion_dynamics_propagate_vel(this%td%ions_dyn, this%ions, atoms_moved = generate)
1033 call this%ions%update_kinetic_energy()
1034 end if
1035
1036 if(associated(gfield)) then
1037 if(gauge_field_is_propagated(gfield)) then
1038 if(this%ks%xc%kernel_lrc_alpha>m_epsilon) then
1039 call gauge_field_get_force(gfield, this%gr, this%st%d%spin_channels, this%st%current, this%ks%xc%kernel_lrc_alpha)
1040 else
1041 call gauge_field_get_force(gfield, this%gr, this%st%d%spin_channels, this%st%current)
1042 endif
1044 end if
1045 end if
1046
1047 !We update the occupation matrices
1048 call lda_u_update_occ_matrices(this%hm%lda_u, this%namespace, this%gr, this%st, this%hm%hm_base, this%hm%phase, this%hm%energy)
1049
1050 ! this is needed to be compatible with the code in td_*
1051 this%td%iter = this%td%iter + 1
1052
1053 call profiling_out(trim(this%namespace%get())//":"//"END_OF_TIMESTEP")
1056
1057 ! ---------------------------------------------------------
1058 subroutine electrons_restart_write_data(this)
1059 class(electrons_t), intent(inout) :: this
1060
1061 integer :: ierr
1062
1064 call profiling_in(trim(this%namespace%get())//":"//"RESTART_WRITE")
1065
1066 select type (algo => this%algo)
1067 class is (propagator_t)
1068 call td_write_data(this%td%write_handler)
1069 call td_dump(this%td, this%namespace, this%space, this%gr, this%st, this%hm, &
1070 this%ks, this%ext_partners, this%iteration%counter(), ierr)
1071 if (ierr /= 0) then
1072 message(1) = "Unable to write time-dependent restart information."
1073 call messages_warning(1, namespace=this%namespace)
1074 end if
1075
1076 ! TODO: this is here because of legacy reasons and should be moved to the output framework
1077 call pes_output(this%td%pesv, this%namespace, this%space, this%gr, this%st, this%iteration%counter(), &
1078 this%outp, algo%dt, this%ions)
1079 end select
1080
1081 call profiling_out(trim(this%namespace%get())//":"//"RESTART_WRITE")
1083 end subroutine electrons_restart_write_data
1084
1085 ! ---------------------------------------------------------
1086 ! this function returns true if restart data could be read
1087 logical function electrons_restart_read_data(this)
1088 class(electrons_t), intent(inout) :: this
1089
1090 logical :: from_scratch
1091
1093 call profiling_in(trim(this%namespace%get())//":"//"RESTART_READ")
1094
1095 select type (algo => this%algo)
1096 class is (propagator_t)
1097 from_scratch = .false.
1098 call td_load_restart_from_td(this%td, this%namespace, this%space, this%mc, this%gr, &
1099 this%ext_partners, this%st, this%ks, this%hm, from_scratch)
1100 call td_set_from_scratch(this%td, from_scratch)
1101
1102 class is (minimizer_algorithm_t)
1103 from_scratch = .false.
1104 call gs_load_from_restart(this%namespace, this%scf, this%gr, this%mc, this%st, this%hm, &
1105 this%ks, this%space, this%ions, this%ext_partners,from_scratch)
1106
1107 ! gs_initialize still knows about fromScratch.
1108 assert(.false.)
1109 end select
1110
1111 if (from_scratch) then
1112 ! restart data could not be loaded
1114 else
1115 ! restart data could be loaded
1117 end if
1118
1119 call profiling_out(trim(this%namespace%get())//":"//"RESTART_READ")
1121 end function electrons_restart_read_data
1122
1123 !----------------------------------------------------------
1124 subroutine electrons_update_kinetic_energy(this)
1125 class(electrons_t), intent(inout) :: this
1126
1128
1129 if (states_are_real(this%st)) then
1130 this%kinetic_energy = denergy_calc_electronic(this%namespace, this%hm, this%gr%der, this%st, terms = term_kinetic)
1131 else
1132 this%kinetic_energy = zenergy_calc_electronic(this%namespace, this%hm, this%gr%der, this%st, terms = term_kinetic)
1133 end if
1134
1136
1137 end subroutine electrons_update_kinetic_energy
1138
1139 ! ---------------------------------------------------------
1140 subroutine get_fields_from_interaction(this, time)
1141 class(electrons_t), intent(inout) :: this
1142 real(real64), intent(in) :: time
1143
1144 type(interaction_iterator_t) :: iter
1145 real(real64), allocatable :: field_tmp(:, :)
1146
1148
1149 if (this%hm%mxll%coupling_mode == no_maxwell_coupling) then
1151 return
1152 end if
1153
1154 assert(this%gr%box%dim == 3)
1155
1156 safe_allocate(field_tmp(1:this%gr%np, 1:this%gr%box%dim))
1157 this%hm%mxll%e_field = m_zero
1158 this%hm%mxll%b_field = m_zero
1159 this%hm%mxll%vec_pot = m_zero
1160
1161 ! interpolate field from interaction
1162 call iter%start(this%interactions)
1163 do while (iter%has_next())
1164 select type (interaction => iter%get_next())
1165 class is (mxll_e_field_to_matter_t)
1166 call interaction%interpolate(time, field_tmp)
1167 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%e_field)
1168 class is (mxll_vec_pot_to_matter_t)
1169 call interaction%interpolate(time, field_tmp)
1170 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%vec_pot)
1171 class is (mxll_b_field_to_matter_t)
1172 call interaction%interpolate(time, field_tmp)
1173 call lalg_axpy(this%gr%np, 3, m_one, field_tmp, this%hm%mxll%b_field)
1174 end select
1175 end do
1176
1177 safe_deallocate_a(field_tmp)
1179
1181
1182 !----------------------------------------------------------
1183 subroutine electrons_finalize(sys)
1184 type(electrons_t), intent(inout) :: sys
1185
1186 type(partner_iterator_t) :: iter
1187 class(interaction_partner_t), pointer :: partner
1188
1189 push_sub(electrons_finalize)
1190
1191 if (associated(sys%algo)) then
1192 select type (algo => sys%algo)
1193 class is (propagator_t)
1194 call td_end_run(sys%td, sys%st, sys%hm)
1195 call td_end(sys%td)
1196 class is(minimizer_algorithm_t)
1197 call gs_cleanup(sys%ks, sys%scf, sys%rdm, sys%st, sys%hm)
1198 end select
1199 end if
1200
1201 if (sys%ks%theory_level /= independent_particles) then
1202 call poisson_async_end(sys%hm%psolver, sys%mc)
1203 end if
1204
1205 call iter%start(sys%ext_partners)
1206 do while (iter%has_next())
1207 partner => iter%get_next()
1208 safe_deallocate_p(partner)
1209 end do
1210 call sys%ext_partners%empty()
1211
1212 safe_deallocate_p(sys%xc_interaction)
1213
1214 call hamiltonian_elec_end(sys%hm)
1215
1216 nullify(sys%gfield)
1217 nullify(sys%lasers)
1218
1219 call multicomm_end(sys%mc)
1220
1221 call xc_oep_photon_end(sys%ks%oep_photon)
1222 if (sys%ks%has_photons) then
1223 call mf_end(sys%ks%pt_mx)
1224 end if
1225
1226 call sys%dipole%end()
1227
1228 call v_ks_end(sys%ks)
1229
1230 call states_elec_end(sys%st)
1231
1232 deallocate(sys%ions)
1233 safe_deallocate_p(sys%photons)
1234
1235 call kpoints_end(sys%kpoints)
1236
1237 call grid_end(sys%gr)
1238
1239 call system_end(sys)
1240
1241 pop_sub(electrons_finalize)
1242 end subroutine electrons_finalize
1243
1244end module electrons_oct_m
1245
1246!! Local Variables:
1247!! mode: f90
1248!! coding: utf-8
1249!! End:
constant times a vector plus a vector
Definition: lalg_basic.F90:171
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:141
character(len=algo_label_len), parameter, public iteration_done
Definition: algorithm.F90:172
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:369
subroutine, public current_init(this, namespace)
Definition: current.F90:179
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
This modules implements the dipole moment of the matter system.
Definition: dipole.F90:108
subroutine, public gs_cleanup(ks, scf, rdm, st, hm)
subroutine, public gs_initialize(namespace, scf, rdm, gr, mc, st, hm, ions, ks, space, ext_partners, fromScratch)
subroutine, public gs_allocate_wavefunctions(namespace, gr, st, hm, scf, ks, ions, mc, space)
subroutine, public gs_load_from_restart(namespace, scf, gr, mc, st, hm, ks, space, ions, ext_partners, fromScratch)
subroutine electrons_initialize(this)
Definition: electrons.F90:663
logical function electrons_restart_read_data(this)
Definition: electrons.F90:1181
subroutine electrons_init_parallelization(this, grp)
Definition: electrons.F90:436
logical function electrons_process_is_slave(this)
Definition: electrons.F90:1035
subroutine electrons_init_interaction(this, interaction)
Definition: electrons.F90:373
subroutine electrons_finalize(sys)
Definition: electrons.F90:1277
subroutine electrons_exec_end_of_timestep_tasks(this, prop)
Definition: electrons.F90:1046
logical function electrons_do_algorithmic_operation(this, operation, updated_quantities)
Definition: electrons.F90:706
subroutine electrons_new_algorithm(this, factory)
Definition: electrons.F90:630
subroutine electrons_output_write(this)
Definition: electrons.F90:1000
subroutine get_fields_from_interaction(this, time)
Definition: electrons.F90:1234
subroutine electrons_algorithm_start(this)
Definition: electrons.F90:686
subroutine electrons_output_finish(this)
Definition: electrons.F90:1026
subroutine electrons_output_start(this)
Definition: electrons.F90:991
subroutine electrons_init_interaction_as_partner(partner, interaction)
Definition: electrons.F90:951
class(electrons_t) function, pointer electrons_constructor(namespace, generate_epot)
Definition: electrons.F90:282
subroutine electrons_update_kinetic_energy(this)
Definition: electrons.F90:1218
logical function electrons_is_tolerance_reached(this, tol)
Definition: electrons.F90:908
subroutine electrons_copy_quantities_to_interaction(partner, interaction)
Definition: electrons.F90:969
subroutine electrons_restart_write_data(this)
Definition: electrons.F90:1152
subroutine electrons_update_quantity(this, label)
Definition: electrons.F90:920
subroutine, public elf_init(namespace)
Definition: elf.F90:144
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:339
subroutine, public gauge_field_do_algorithmic_operation(this, operation, dt, time)
subroutine, public gauge_field_check_symmetries(this, kpoints)
subroutine, public gauge_field_get_force(this, gr, spin_channels, current, lrc_alpha)
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:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public grid_init_stage_1(gr, namespace, space, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:194
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:465
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:493
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:1012
subroutine, public kpoints_init(this, namespace, symm, dim, periodic_dim, latt)
Definition: kpoints.F90:322
subroutine, public kpoints_to_absolute(latt, kin, kout)
Definition: kpoints.F90:1031
A module to handle KS potential, without the external potential.
integer, parameter, public independent_particles
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
subroutine, public lasers_check_symmetries(this, kpoints)
Definition: lasers.F90:555
subroutine, public lasers_parse_external_fields(this)
Definition: lasers.F90:244
subroutine, public lasers_generate_potentials(this, mesh, space, latt)
Definition: lasers.F90:444
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, hm_base, phase, energy)
Definition: lda_u.F90:797
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public mesh_check_symmetries(mesh, symm, periodic_dim)
Definition: mesh.F90:836
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:380
real(real64) pure function, public mesh_global_memory(mesh)
Definition: mesh.F90:783
real(real64) pure function, public mesh_local_memory(mesh)
Definition: mesh.F90:794
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:921
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
subroutine, public messages_new_line()
Definition: messages.F90:1135
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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)
Definition: mpi.F90:403
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public multicomm_end(mc)
Definition: multicomm.F90:897
logical pure function, public multicomm_is_slave(this)
Definition: multicomm.F90:1141
subroutine, public multicomm_init(mc, namespace, base_grp, mode_para, n_node, index_range, min_range)
create index and domain communicators
Definition: multicomm.F90:268
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:115
this module contains the output system
Definition: output.F90:115
logical function, public output_need_exchange(outp)
Definition: output.F90:874
subroutine, public output_init(outp, namespace, space, st, gr, nst, ks)
Definition: output.F90:208
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public pes_calc(pes, namespace, space, mesh, st, dt, iter, der, kpoints, ext_partners, stopping)
Definition: pes.F90:269
subroutine, public pes_output(pes, namespace, space, gr, st, iter, outp, dt, ions)
Definition: pes.F90:295
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:1143
subroutine, public poisson_slave_work(this, namespace)
Definition: poisson.F90:1195
subroutine, public poisson_async_end(this, mc)
Definition: poisson.F90:1168
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine, 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_exp_apply(te, namespace, st, mesh, hm, dt)
subroutine, public propagation_ops_elec_interpolate_get(hm, mesh, vks_old)
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_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, vmagnus)
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:117
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:138
Implementation details for regridding.
Definition: regridding.F90:170
logical function, public scf_iter_finish(scf, namespace, space, gr, ions, st, ks, hm, iter, outp, verbosity, iters_done)
Definition: scf.F90:1226
subroutine, public scf_finish(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, iter, outp, verbosity, iters_done, restart_dump)
Definition: scf.F90:1306
integer, parameter, public verb_compact
Definition: scf.F90:203
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:254
subroutine, public scf_end(scf)
Definition: scf.F90:547
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:833
subroutine, public scf_iter(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, iter, outp, verbosity, iters_done, restart_dump)
Definition: scf.F90:877
subroutine, public scf_start(scf, namespace, space, gr, ions, st, ks, hm, outp, verbosity)
Preparation of the SCF cycle.
Definition: scf.F90:679
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
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 states_elec_init(st, namespace, space, valence_charge, kpoints)
Initialize a new 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_allocate_current(st, space, mesh)
This module implements the calculation of the stress tensor.
Definition: stress.F90:118
subroutine, public stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
This computes the total stress on the lattice.
Definition: stress.F90:186
This module implements the abstract system type.
Definition: system.F90:118
subroutine, public system_algorithm_start(this)
Definition: system.F90:1018
subroutine, public system_end(this)
Definition: system.F90:1145
subroutine, public system_new_algorithm(this, factory)
Definition: system.F90:943
Definition: td.F90:114
subroutine, public td_end(td)
Definition: td.F90:609
subroutine, public td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
Definition: td.F90:1157
subroutine, public td_end_run(td, st, hm)
Definition: td.F90:624
subroutine, public td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp)
Definition: td.F90:247
subroutine, public td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
Definition: td.F90:544
logical function, public td_get_from_scratch(td)
Definition: td.F90:1503
subroutine, public td_set_from_scratch(td, from_scratch)
Definition: td.F90:1514
subroutine, public td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
Definition: td.F90:1289
subroutine, public td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
Definition: td.F90:578
subroutine, public td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, from_scratch)
Definition: td.F90:842
subroutine, public td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
Definition: td.F90:1121
subroutine, public td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, dt)
Definition: td_write.F90:1242
subroutine, public td_write_data(writ)
Definition: td_write.F90:1208
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:1038
integer, parameter, public out_separate_forces
Definition: td_write.F90:201
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:624
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:744
subroutine, public v_ks_init(ks, namespace, gr, st, ions, mc, space, kpoints)
Definition: v_ks.F90:250
integer, parameter, public xc_oep_x
Exact exchange.
integer, parameter, public func_x
Definition: xc.F90:114
integer, parameter, public oep_level_none
the OEP levels
Definition: xc_oep.F90:171
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:163
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:218
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:142
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:138
Systems (system_t) can expose quantities that can be used to calculate interactions with other system...
Definition: quantity.F90:171
Abstract class for systems.
Definition: system.F90:172
int true(void)