Octopus
td.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module td_oct_m
27 use debug_oct_m
33 use epot_oct_m
35 use forces_oct_m
37 use global_oct_m
38 use grid_oct_m
41 use io_oct_m
43 use ions_oct_m
44 use kick_oct_m
45 use, intrinsic :: iso_fortran_env
46 use lasers_oct_m
47 use lda_u_oct_m
50 use loct_oct_m
52 use mesh_oct_m
54 use mpi_oct_m
57 use output_oct_m
59 use parser_oct_m
60 use pes_oct_m
70 use scf_oct_m
72 use space_oct_m
76 use stress_oct_m
78 use types_oct_m
79 use unit_oct_m
81 use v_ks_oct_m
84 use xc_oct_m
85
86 implicit none
87
88 private
89 public :: &
90 td_t, &
91 td_run, &
92 td_init, &
94 td_end, &
95 td_end_run, &
98 td_dump, &
106
108 integer, parameter, public :: &
109 EHRENFEST = 1, &
110 bo = 2
111
112 type td_t
113 private
114 type(propagator_base_t), public :: tr
115 type(scf_t), public :: scf
116 type(ion_dynamics_t), public :: ions_dyn
117 real(real64), public :: dt
118 integer, public :: max_iter
119 integer, public :: iter
120 logical, public :: recalculate_gs
121 integer, public :: recalculate_gs_interval
122
123 type(pes_t), public :: pesv
124
125 integer, public :: dynamics
126 integer, public :: energy_update_iter
127 real(real64) :: scissor
128
129 logical :: freeze_occ
130 logical :: freeze_u
131 integer :: freeze_orbitals
132
133 logical :: from_scratch = .false.
134
135 type(td_write_t), public :: write_handler
136 type(restart_t) :: restart_load
137 type(restart_t) :: restart_dump
138
139 type(dmp_t) :: dmp
140 end type td_t
141
142
143contains
144
145 subroutine td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp)
146 type(td_t), intent(inout) :: td
147 type(namespace_t), intent(in) :: namespace
148 class(space_t), intent(in) :: space
149 type(grid_t), intent(in) :: gr
150 type(ions_t), intent(inout) :: ions
151 type(states_elec_t), intent(in) :: st
152 type(v_ks_t), intent(in) :: ks
153 type(hamiltonian_elec_t), intent(in) :: hm
154 type(partner_list_t), intent(in) :: ext_partners
155 type(output_t), intent(in) :: outp
156
157 integer :: default
158 real(real64) :: propagation_time
159 type(lasers_t), pointer :: lasers
160 logical :: symmetrize
161
162 push_sub(td_init)
163
164 if (hm%pcm%run_pcm) call messages_experimental("PCM for CalculationMode = td", namespace=namespace)
165
166 symmetrize = hm%kpoints%use_symmetries .or. st%symmetrize_density
167 call ion_dynamics_init(td%ions_dyn, namespace, ions, symmetrize, gr%symmetrizer)
168
169 if (td%ions_dyn%ions_move()) then
170 if (hm%kpoints%use_symmetries) then
171 message(1) = "KPoints symmetries cannot be used with moving ions."
172 message(2) = "Please set KPointsSymmetries = no."
173 call messages_fatal(2, namespace=namespace)
174 end if
175 if (st%symmetrize_density) then
176 message(1) = "Symmetrization of the density cannot be used with moving ions."
177 message(2) = "Please set SymmetrizeDensity = no."
178 call messages_fatal(2, namespace=namespace)
179 end if
180 end if
181
182 td%iter = 0
183
184 !%Variable TDTimeStep
185 !%Type float
186 !%Section Time-Dependent::Propagation
187 !%Description
188 !% The time-step for the time propagation. For most propagators you
189 !% want to use the largest value that is possible without the
190 !% evolution becoming unstable.
191 !%
192 !% While prior versions of Octopus used to have a default time step, however,
193 !% this now needs to be systematically defined.
194 !%End
195
196 call parse_variable(namespace, 'TDTimeStep', -m_one, td%dt, unit = units_inp%time)
197
198 if (td%dt <= m_zero) then
199 write(message(1),'(a)') 'A positive value for TDTimeStep must be defined in the input file.'
200 call messages_fatal(1, namespace=namespace)
201 end if
202
203 call messages_print_var_value('TDTimeStep', td%dt, unit = units_out%time, namespace=namespace)
204
205
206 if (parse_is_defined(namespace, 'TDMaxSteps') .and. parse_is_defined(namespace, 'TDPropagationTime')) then
207 call messages_write('You cannot set TDMaxSteps and TDPropagationTime at the same time')
208 call messages_fatal(namespace=namespace)
209 end if
211 !%Variable TDPropagationTime
212 !%Type float
213 !%Section Time-Dependent::Propagation
214 !%Description
215 !% The length of the time propagation. You cannot set this variable
216 !% at the same time as <tt>TDMaxSteps</tt>. By default this variable will
217 !% not be used.
218 !%
219 !% The units for this variable are <math>\hbar</math>/Hartree (or <math>\hbar</math>/eV if you
220 !% selected <tt>ev_angstrom</tt> as input units). The approximate conversions to
221 !% femtoseconds are 1 fs = 41.34 <math>\hbar</math>/Hartree = 1.52 <math>\hbar</math>/eV.
222 !%End
223 call parse_variable(namespace, 'TDPropagationTime', -1.0_real64, propagation_time, unit = units_inp%time)
225 call messages_obsolete_variable(namespace, 'TDMaximumIter', 'TDMaxSteps')
227 !%Variable TDMaxSteps
228 !%Type integer
229 !%Default 1500
230 !%Section Time-Dependent::Propagation
231 !%Description
232 !% Number of time-propagation steps that will be performed. You
233 !% cannot use this variable together with <tt>TDPropagationTime</tt>.
234 !%End
235 default = 1500
236 if (propagation_time > m_zero) default = nint(propagation_time/td%dt)
237 call parse_variable(namespace, 'TDMaxSteps', default, td%max_iter)
238
239 if (propagation_time <= m_zero) propagation_time = td%dt*td%max_iter
241 call messages_print_var_value('TDPropagationTime', propagation_time, unit = units_out%time, namespace=namespace)
242 call messages_print_var_value('TDMaxSteps', td%max_iter, namespace=namespace)
243
244 if (td%max_iter < 1) then
245 write(message(1), '(a,i6,a)') "Input: '", td%max_iter, "' is not a valid value for TDMaxSteps."
246 message(2) = '(TDMaxSteps <= 1)'
247 call messages_fatal(2, namespace=namespace)
248 end if
249
250 td%iter = 0
251
252 td%dt = td%dt
253
254 lasers => list_get_lasers(ext_partners)
255
256 ! now the photoelectron stuff
257 call pes_init(td%pesv, namespace, space, gr, gr%box, st, 1, hm%kpoints, &
258 hm%abs_boundaries, ext_partners, td%max_iter, td%dt)
259
260 !%Variable TDDynamics
261 !%Type integer
262 !%Default ehrenfest
263 !%Section Time-Dependent::Propagation
264 !%Description
265 !% Type of dynamics to follow during a time propagation.
266 !% For BO, you must set <tt>MoveIons = yes</tt>.
267 !%Option ehrenfest 1
268 !% Ehrenfest dynamics.
269 !%Option bo 2
270 !% Born-Oppenheimer (Experimental).
271 !%End
272
273 call parse_variable(namespace, 'TDDynamics', ehrenfest, td%dynamics)
274 if (.not. varinfo_valid_option('TDDynamics', td%dynamics)) call messages_input_error(namespace, 'TDDynamics')
275 call messages_print_var_option('TDDynamics', td%dynamics, namespace=namespace)
276 if (td%dynamics .ne. ehrenfest) then
277 if (.not. td%ions_dyn%is_active()) then
278 message(1) = "TDDynamics=bo can only be used if MoveIons=yes or CellDynamics=yes"
279 call messages_fatal(1, namespace=namespace)
280 end if
281 end if
282
283 !%Variable RecalculateGSDuringEvolution
284 !%Type logical
285 !%Default no
286 !%Section Time-Dependent::Propagation
287 !%Description
288 !% In order to calculate some information about the system along the
289 !% evolution (e.g. projection onto the ground-state KS determinant,
290 !% projection of the TDKS spin-orbitals onto the ground-state KS
291 !% spin-orbitals), the ground-state KS orbitals are needed. If the
292 !% ionic potential changes -- that is, the ions move -- one may want
293 !% to recalculate the ground state. You may do this by setting this
294 !% variable.
295 !%
296 !% The recalculation is not done every time step, but only every
297 !% <tt>RecalculateGSInterval</tt> time steps.
298 !%End
299 call parse_variable(namespace, 'RecalculateGSDuringEvolution', .false., td%recalculate_gs)
300 if (hm%lda_u_level /= dft_u_none .and. td%recalculate_gs) then
301 call messages_not_implemented("DFT+U with RecalculateGSDuringEvolution=yes", namespace=namespace)
302 end if
303
304 !%Variable RecalculateGSInterval
305 !%Type integer
306 !%Default 50
307 !%Section Time-Dependent::Propagation
308 !%Description
309 !% This variable controls how often the GS is recalculated for <tt>RecalculateGSDuringEvolution = yes</tt>.
310 !%End
311 call parse_variable(namespace, 'RecalculateGSInterval', 50, td%recalculate_gs_interval)
312
313 !%Variable TDScissor
314 !%Type float
315 !%Default 0.0
316 !%Section Time-Dependent
317 !%Description
318 !% (experimental) If set, a scissor operator will be applied in the
319 !% Hamiltonian, shifting the excitation energies by the amount
320 !% specified. By default, it is not applied.
321 !%End
322 call parse_variable(namespace, 'TDScissor', m_zero, td%scissor)
323 call messages_print_var_value('TDScissor', td%scissor, namespace=namespace)
324
325 call propagator_elec_init(gr, namespace, st, td%tr, hm%ks_pot, td%ions_dyn%is_active() .and.&
326 list_has_gauge_field(ext_partners), family_is_mgga_with_exc(ks%xc), td%ions_dyn%cell_relax())
327
328 if (associated(lasers) .and. mpi_world%is_root()) then
329 call messages_print_with_emphasis(msg="Time-dependent external fields", namespace=namespace)
330 call laser_write_info(lasers%lasers, dt=td%dt, max_iter=td%max_iter, namespace=namespace)
331 call messages_print_with_emphasis(namespace=namespace)
332 end if
333
334 !%Variable TDEnergyUpdateIter
335 !%Type integer
336 !%Section Time-Dependent::Propagation
337 !%Description
338 !% This variable controls after how many iterations Octopus
339 !% updates the total energy during a time-propagation run. For
340 !% iterations where the energy is not updated, the last calculated
341 !% value is reported. If you set this variable to 1, the energy
342 !% will be calculated in each step.
343 !%End
344
345 default = 10
346 call parse_variable(namespace, 'TDEnergyUpdateIter', default, td%energy_update_iter)
347
348 if (gr%der%boundaries%spiralBC .and. hm%ep%reltype == spin_orbit) then
349 message(1) = "Generalized Bloch theorem cannot be used with spin-orbit coupling."
350 call messages_fatal(1, namespace=namespace)
351 end if
352
353 if (gr%der%boundaries%spiralBC) then
354 if (any(abs(hm%kick%easy_axis(1:2)) > m_epsilon)) then
355 message(1) = "Generalized Bloch theorem cannot be used for an easy axis not along the z direction."
356 call messages_fatal(1, namespace=namespace)
357 end if
358 end if
359
360 !%Variable TDFreezeOrbitals
361 !%Type integer
362 !%Default 0
363 !%Section Time-Dependent
364 !%Description
365 !% (Experimental) You have the possibility of "freezing" a number of orbitals during a time-propagation.
366 !% The Hartree and exchange-correlation potential due to these orbitals (which
367 !% will be the lowest-energy ones) will be added during the propagation, but the orbitals
368 !% will not be propagated.
369 !%Option sae -1
370 !% Single-active-electron approximation. This option is only valid for time-dependent
371 !% calculations (<tt>CalculationMode = td</tt>). Also, the nuclei should not move.
372 !% The idea is that all orbitals except the last one are frozen. The orbitals are to
373 !% be read from a previous ground-state calculation. The active orbital is then treated
374 !% as independent (whether it contains one electron or two) -- although it will
375 !% feel the Hartree and exchange-correlation potentials from the ground-state electronic
376 !% configuration.
377 !%
378 !% It is almost equivalent to setting <tt>TDFreezeOrbitals = N-1</tt>, where <tt>N</tt> is the number
379 !% of orbitals, but not completely.
380 !%End
381 call parse_variable(namespace, 'TDFreezeOrbitals', 0, td%freeze_orbitals)
382
383 if (td%freeze_orbitals /= 0) then
384 call messages_experimental('TDFreezeOrbitals', namespace=namespace)
385
386 if (hm%lda_u_level /= dft_u_none) then
387 call messages_not_implemented('TDFreezeOrbitals with DFT+U', namespace=namespace)
388 end if
389 end if
390
391 !%Variable TDDMPropagation
392 !%Type integer
393 !%Default no_propagation
394 !%Section Time-Dependent
395 !%Description
396 !% Decides whether to propagate the density matrix, in the spectral basis,
397 !% together with the KS orbitals.
398 !% The resulting time-depented orbilas are to be considered the natural orbitals
399 !% of the denstiy matrix and with the corresponding eigenvalues the occupation numbers.
400 !% In order to work properly this method requires an number of empty states
401 !% specified with the option "ExtraStates".
402 !% This is an experimental feature.
403 !%Option no_propagation 01
404 !% No density time propagation.
405 !%Option master_equation 02
406 !% Density time propagation with the Lindblad master equation.
407 !%End
408
409 call parse_variable(global_namespace, 'TDDMPropagation', option__tddmpropagation__no_propagation, td%dmp%calculation_mode)
410 if (.not. varinfo_valid_option('TDDMPropagation', td%dmp%calculation_mode)) then
411 call messages_input_error(global_namespace, 'TDDMPropagation')
412 endif
413 call messages_print_var_option('TDDMPropagation', td%dmp%calculation_mode, namespace=namespace)
414
415 if (td%dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
416 call messages_experimental('TDDMPropagation', namespace=namespace)
417 call td%dmp%init(namespace, st)
418 end if
419
420 pop_sub(td_init)
421 nullify(lasers)
422
423 end subroutine td_init
424
425 ! ---------------------------------------------------------
426 subroutine td_init_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, from_scratch)
427 type(td_t), intent(inout) :: td
428 type(namespace_t), intent(in) :: namespace
429 type(multicomm_t), intent(inout) :: mc
430 type(grid_t), intent(inout) :: gr
431 type(ions_t), intent(inout) :: ions
432 type(states_elec_t), intent(inout) :: st
433 type(v_ks_t), intent(inout) :: ks
434 type(hamiltonian_elec_t), intent(inout) :: hm
435 type(partner_list_t), intent(in) :: ext_partners
436 type(output_t), intent(inout) :: outp
437 type(electron_space_t), intent(in) :: space
438 logical, intent(inout) :: from_scratch
439 push_sub(td_init_run)
440
441 ! NOTE: please do not change code in this function, but only in functions
442 ! called from here because the logic of this function is replicated in the
443 ! multisystem framework in different places
444
445 call td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
446 call td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
447
448 td%from_scratch = from_scratch
449
450 if (.not. td%from_scratch) then
451 call td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, td%from_scratch)
452 if (td%from_scratch) then
453 message(1) = "Unable to read time-dependent restart information: Starting from scratch"
454 call messages_warning(1, namespace=namespace)
455 end if
456 end if
457
458 if (td%iter >= td%max_iter) then
459 message(1) = "All requested iterations have already been done. Use FromScratch = yes if you want to redo them."
460 call messages_info(1, namespace=namespace)
462 td%iter = td%iter + 1
463 if ((td%ions_dyn%is_active()) .and. td%recalculate_gs) call td%restart_load%end()
464 pop_sub(td_init_run)
465 return
466 end if
467
468 if (td%from_scratch) then
469 call td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
470 end if
471
472 call td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, td%from_scratch)
473
474 if (td%dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
475 call dm_propagation_init_run(td%dmp%adiabatic_st, namespace, space, gr, st, hm, mc)
476 end if
477
478 pop_sub(td_init_run)
479 end subroutine td_init_run
480
481 ! ---------------------------------------------------------
482 subroutine td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
483 type(td_t), intent(inout) :: td
484 type(namespace_t), intent(in) :: namespace
485 type(multicomm_t), intent(inout) :: mc
486 type(grid_t), intent(inout) :: gr
487 type(ions_t), intent(inout) :: ions
488 type(states_elec_t), intent(inout) :: st
489 type(hamiltonian_elec_t), intent(inout) :: hm
490 class(space_t), intent(in) :: space
491
493
494 ! Allocate wavefunctions during time-propagation
495 if (td%dynamics == ehrenfest) then
496 !Note: this is not really clean to do this
497 if (hm%lda_u_level /= dft_u_none .and. states_are_real(st)) then
498 call lda_u_end(hm%lda_u)
499 !complex wfs are required for Ehrenfest
500 call states_elec_allocate_wfns(st, gr, type_cmplx, packed=.true.)
501 call lda_u_init(hm%lda_u, namespace, space, hm%lda_u_level, gr, ions, st, mc, hm%kpoints)
502 else
503 !complex wfs are required for Ehrenfest
504 call states_elec_allocate_wfns(st, gr, type_cmplx, packed=.true.)
505 end if
506 else
507 call states_elec_allocate_wfns(st, gr, packed=.true.)
508 call scf_init(td%scf, namespace, gr, ions, st, mc, hm, space)
509 end if
510
512 end subroutine td_allocate_wavefunctions
513
514 ! ---------------------------------------------------------
515 subroutine td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
516 type(td_t), intent(inout) :: td
517 type(namespace_t), intent(in) :: namespace
518 type(grid_t), intent(inout) :: gr
519 type(states_elec_t), intent(inout) :: st
520 type(v_ks_t), intent(inout) :: ks
521 type(hamiltonian_elec_t), intent(inout) :: hm
522 type(partner_list_t), intent(in) :: ext_partners
523 class(space_t), intent(in) :: space
524
525 type(gauge_field_t), pointer :: gfield
526
527 push_sub(td_init_gaugefield)
528
529 gfield => list_get_gauge_field(ext_partners)
530 if(associated(gfield)) then
531 if (gauge_field_is_used(gfield)) then
532 !if the gauge field is applied, we need to tell v_ks to calculate the current
533 call v_ks_calculate_current(ks, .true.)
534
535 ! initialize the vector field and update the hamiltonian
536 call gauge_field_init_vec_pot(gfield, st%qtot)
537 call hm%update(gr, namespace, space, ext_partners, time = td%dt*td%iter)
538
539 end if
540 end if
541
542 pop_sub(td_init_gaugefield)
543 end subroutine td_init_gaugefield
544
545 ! ---------------------------------------------------------
546 subroutine td_end(td)
547 type(td_t), intent(inout) :: td
548
549 push_sub(td_end)
550
551 call pes_end(td%pesv)
552 call propagator_elec_end(td%tr) ! clean the evolution method
553 call ion_dynamics_end(td%ions_dyn)
554
555 if (td%dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
556 call states_elec_end(td%dmp%adiabatic_st)
557 endif
558
559 if (td%dynamics == bo) call scf_end(td%scf)
560 pop_sub(td_end)
561
562 end subroutine td_end
563
564 ! ---------------------------------------------------------
565 subroutine td_end_run(td, st, hm)
566 type(td_t), intent(inout) :: td
567 type(states_elec_t), intent(inout) :: st
568 type(hamiltonian_elec_t), intent(inout) :: hm
569
570 push_sub(td_end_run)
571
572 if (st%pack_states .and. hm%apply_packed()) call st%unpack()
573
574 call td%restart_dump%end()
575 call td_write_end(td%write_handler)
576
577 ! free memory
579 if ((td%ions_dyn%is_active()).and. td%recalculate_gs) call td%restart_load%end()
580
581 pop_sub(td_end_run)
582 end subroutine td_end_run
583
584 ! ---------------------------------------------------------
585 subroutine td_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, from_scratch)
586 type(td_t), intent(inout) :: td
587 type(namespace_t), intent(in) :: namespace
588 type(multicomm_t), intent(inout) :: mc
589 type(grid_t), intent(inout) :: gr
590 type(ions_t), intent(inout) :: ions
591 type(states_elec_t), intent(inout) :: st
592 type(v_ks_t), intent(inout) :: ks
593 type(hamiltonian_elec_t), intent(inout) :: hm
594 type(partner_list_t), intent(in) :: ext_partners
595 type(output_t), intent(inout) :: outp
596 type(electron_space_t), intent(in) :: space
597 logical, intent(inout) :: from_scratch
598
599 logical :: stopping
600 integer :: iter, scsteps
601 real(real64) :: etime
602 real(real64) :: wall_time, simulation_time, speed_fs_per_day
603 character(len=20) :: fmt
604
605 push_sub(td_run)
606
607 etime = loct_clock()
608 ! This is the time-propagation loop. It starts at t=0 and finishes
609 ! at td%max_iter*dt. The index i runs from 1 to td%max_iter, and
610 ! step "iter" means propagation from (iter-1)*dt to iter*dt.
611 propagation: do iter = td%iter, td%max_iter
612
613 stopping = clean_stop(mc%master_comm) .or. walltimer_alarm(mc%master_comm)
614
615 call profiling_in("TIME_STEP")
616
617 if (iter > 1) then
618 if (((iter-1)*td%dt <= hm%kick%time) .and. (iter*td%dt > hm%kick%time)) then
619 if (.not. hm%pcm%localf) then
620 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints)
621 else
622 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints, pcm = hm%pcm)
623 end if
624 call td_write_kick(outp, namespace, space, gr, hm%kick, ions, iter)
625 !We activate the sprial BC only after the kick,
626 !to be sure that the first iteration corresponds to the ground state
627 if (gr%der%boundaries%spiralBC) gr%der%boundaries%spiral = .true.
628 end if
629 end if
630
631 ! time iterate the system, one time step.
632 select case (td%dynamics)
633 case (ehrenfest)
634 call propagator_elec_dt(ks, namespace, space, hm, gr, st, td%tr, iter*td%dt, td%dt, iter, td%ions_dyn, &
635 ions, ext_partners, mc, outp, td%write_handler, scsteps = scsteps, &
636 update_energy = (mod(iter, td%energy_update_iter) == 0) .or. (iter == td%max_iter))
637 case (bo)
638 call propagator_elec_dt_bo(td%scf, namespace, space, gr, ks, st, hm, ions, ext_partners, mc, iter, &
639 td%dt, td%ions_dyn, scsteps)
640 end select
642 if (td%dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
643 call dm_propagation_run(td%dmp, namespace, space, gr, ions, st, mc, hm, ks, iter, td%dt, ext_partners, &
644 update_energy = (mod(iter, td%energy_update_iter) == 0) .or. (iter == td%max_iter))
645 end if
646
647 !Apply mask absorbing boundaries
648 if (hm%abs_boundaries%abtype == mask_absorbing) then
649 if (states_are_real(st)) then
650 call dvmask(gr, hm, st)
651 else
652 call zvmask(gr, hm, st)
653 end if
654 end if
655
656 !Photoelectron stuff
657 if (td%pesv%calc_spm .or. td%pesv%calc_mask .or. td%pesv%calc_flux) then
658 call pes_calc(td%pesv, namespace, space, gr, st, td%dt, iter, gr%der, hm%kpoints, ext_partners, stopping)
659 end if
661 call td_write_iter(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
662 hm%kick, ks, td%dt, iter, mc, td%recalculate_gs)
663
664 ! write down data
665 call td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, &
666 iter, scsteps, etime, stopping, from_scratch)
667
668 ! check if debug mode should be enabled or disabled on the fly
669 call io_debug_on_the_fly(namespace)
670
671 call profiling_out("TIME_STEP")
672 if (stopping) exit
673
674 end do propagation
675
676 ! Output propagation speed in fs/day
677 wall_time = loct_clock() - walltimer_get_start_time()
678 simulation_time = td%dt * (iter - td%iter + 1)
679 speed_fs_per_day = units_from_atomic(unit_femtosecond, simulation_time) / (wall_time / 86400.0_real64)
680 if (speed_fs_per_day > 1e4 .or. speed_fs_per_day < 1e-3) then
681 fmt = "ES10.3"
682 else
683 fmt = "F10.3"
684 end if
685 write(message(1), '(a,'//trim(fmt)//',a)') 'Propagation speed: ', speed_fs_per_day, ' fs/day'
686 write(message(2), '(a)') ''
687 call messages_info(2, namespace=namespace)
688
689 pop_sub(td_run)
690 end subroutine td_run
691
692 subroutine td_print_header(namespace)
693 type(namespace_t), intent(in) :: namespace
694
695 push_sub(td_print_header)
696
697 write(message(1), '(a7,1x,a14,a14,a10,a17)') 'Iter ', 'Time ', 'Energy ', 'SC Steps', 'Elapsed Time '
698
699 call messages_info(1, namespace=namespace)
700 call messages_print_with_emphasis(namespace=namespace)
701
702 pop_sub(td_print_header)
703 end subroutine td_print_header
704
705 ! ---------------------------------------------------------
706 subroutine td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, &
707 iter, scsteps, etime, stopping, from_scratch)
708 type(td_t), intent(inout) :: td
709 type(namespace_t), intent(in) :: namespace
710 type(multicomm_t), intent(in) :: mc
711 type(grid_t), intent(inout) :: gr
712 type(ions_t), intent(inout) :: ions
713 type(states_elec_t), intent(inout) :: st
714 type(v_ks_t), intent(inout) :: ks
715 type(hamiltonian_elec_t), intent(inout) :: hm
716 type(partner_list_t), intent(in) :: ext_partners
717 type(output_t), intent(in) :: outp
718 type(electron_space_t), intent(in) :: space
719 integer, intent(in) :: iter
720 integer, intent(in) :: scsteps
721 real(real64), intent(inout) :: etime
722 logical, intent(in) :: stopping
723 logical, intent(inout) :: from_scratch
724
725 integer :: ierr
726
727 push_sub(td_check_point)
728
729 call td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
730
731 if (outp%anything_now(iter)) then ! output
732 call td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, td%dt)
733 end if
734
735 if (restart_walltime_period_alarm(mc%master_comm) .or. iter == td%max_iter .or. stopping) then ! restart
736 !if (iter == td%max_iter) outp%iter = ii - 1
737 call td_write_data(td%write_handler)
738 call td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
739 if (ierr /= 0) then
740 message(1) = "Unable to write time-dependent restart information."
741 call messages_warning(1, namespace=namespace)
742 end if
743
744 call pes_output(td%pesv, namespace, space, gr, st, iter, outp, td%dt, ions)
745 end if
746
747 if (mod(iter, td%recalculate_gs_interval) == 0 .or. iter == td%max_iter .or. stopping) then ! potentially recalculate GS
748 if ((td%ions_dyn%is_active()) .and. td%recalculate_gs) then
749 call messages_print_with_emphasis(msg='Recalculating the ground state.', namespace=namespace)
750 from_scratch = .false.
752 call electrons_ground_state_run(namespace, mc, gr, ions, ext_partners, st, ks, hm, outp, space, from_scratch)
753 call states_elec_allocate_wfns(st, gr, packed=.true.)
754 call td_load(td%restart_load, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
755 if (ierr /= 0) then
756 message(1) = "Unable to load TD states."
757 call messages_fatal(1, namespace=namespace)
758 end if
759 call density_calc(st, gr, st%rho)
760 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
761 calc_eigenval=.true., time = iter*td%dt, calc_energy=.true.)
762 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = iter*td%dt, dt = td%dt)
763 assert(.not. td%ions_dyn%cell_relax())
764 call messages_print_with_emphasis(msg="Time-dependent simulation proceeds", namespace=namespace)
765 call td_print_header(namespace)
766 end if
767 end if
768
769 pop_sub(td_check_point)
770 end subroutine td_check_point
771
772 ! ---------------------------------------------------------
773 subroutine td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
774 type(td_t), intent(inout) :: td
775 type(namespace_t), intent(in) :: namespace
776 type(ions_t), intent(inout) :: ions
777 type(hamiltonian_elec_t), intent(inout) :: hm
778 integer, intent(in) :: iter
779 integer, intent(in) :: scsteps
780 real(real64), intent(inout) :: etime
781
782 push_sub(td_print_message)
783
784 write(message(1), '(i7,1x,2f14.6,i10,f14.3)') iter, units_from_atomic(units_out%time, iter*td%dt), &
785 units_from_atomic(units_out%energy, hm%energy%total + ions%kinetic_energy), scsteps, &
786 loct_clock() - etime
787 call messages_info(1, namespace=namespace)
788 call td_update_elapsed_time(etime)
789
790 pop_sub(td_print_message)
791 end subroutine td_print_message
792
793 ! ---------------------------------------------------------
794 subroutine td_update_elapsed_time(etime)
795 real(real64), intent(inout) :: etime
796
797 push_sub(td_update_elapsed_time)
798
799 etime = loct_clock()
800
802 end subroutine td_update_elapsed_time
803
804 ! ---------------------------------------------------------
805 subroutine td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, from_scratch)
806 type(td_t), intent(inout) :: td
807 type(namespace_t), intent(in) :: namespace
808 type(electron_space_t), intent(in) :: space
809 type(multicomm_t), intent(in) :: mc
810 type(grid_t), intent(inout) :: gr
811 type(ions_t), intent(inout) :: ions
812 type(partner_list_t), intent(in) :: ext_partners
813 type(states_elec_t), target, intent(inout) :: st
814 type(v_ks_t), intent(inout) :: ks
815 type(hamiltonian_elec_t), intent(inout) :: hm
816 type(output_t), intent(inout) :: outp
817 logical, intent(in) :: from_scratch
818
819 integer :: ierr
820 real(real64) :: x
821 real(real64) :: ndinitial(space%dim)
822 logical :: freeze_hxc, freeze_occ, freeze_u
823 type(restart_t) :: restart, restart_frozen
824 type(gauge_field_t), pointer :: gfield
825 type(lasers_t), pointer :: lasers
827
828 !We activate the sprial BC only after the kick,
829 !to be sure that the first iteration corresponds to the ground state
830 if (gr%der%boundaries%spiralBC) then
831 if ((td%iter-1)*td%dt > hm%kick%time) then
832 gr%der%boundaries%spiral = .true.
833 end if
834 hm%vnl%spin => st%spin
835 hm%phase%spin => st%spin
836 !We fill st%spin. In case of restart, we read it in td_load
837 if (from_scratch) call states_elec_fermi(st, namespace, gr)
838 end if
839
840 if (from_scratch) then
841 ! Initialize the occupation matrices and U for DFT+U
842 ! This must be called before parsing TDFreezeOccupations and TDFreezeU
843 ! in order that the code does properly the initialization.
844 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
845 end if
846
847 if (td%freeze_orbitals > 0) then
848 if (from_scratch) then
849 ! In this case, we first freeze the orbitals, then calculate the Hxc potential.
850 call states_elec_freeze_orbitals(st, namespace, space, gr, mc, hm%kpoints, &
851 td%freeze_orbitals, family_is_mgga(ks%xc_family))
852 else
853 call restart%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
854 if (ierr == 0) then
855 call td_load_frozen(namespace, restart, space, gr, st, hm, ierr)
856 end if
857 if (ierr /= 0) then
858 td%iter = 0
859 message(1) = "Unable to read frozen restart information."
860 call messages_fatal(1, namespace=namespace)
861 end if
862 call restart%end()
863 end if
864 write(message(1),'(a,i4,a,i4,a)') 'Info: The lowest', td%freeze_orbitals, &
865 ' orbitals have been frozen.', st%nst, ' will be propagated.'
866 call messages_info(1, namespace=namespace)
868 call density_calc(st, gr, st%rho)
869 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
870 else if (td%freeze_orbitals < 0) then
871 ! This means SAE approximation. We calculate the Hxc first, then freeze all
872 ! orbitals minus one.
873 write(message(1),'(a)') 'Info: The single-active-electron approximation will be used.'
874 call messages_info(1, namespace=namespace)
875 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
876 if (from_scratch) then
877 call states_elec_freeze_orbitals(st, namespace, space, gr, mc, hm%kpoints, st%nst-1, family_is_mgga(ks%xc_family))
878 else
879 call messages_not_implemented("TDFreezeOrbials < 0 with FromScratch=no", namespace=namespace)
880 end if
882 call v_ks_freeze_hxc(ks)
883 call density_calc(st, gr, st%rho)
884 else
885 ! Normal run.
886 call density_calc(st, gr, st%rho)
887 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
888 end if
890 !%Variable TDFreezeHXC
891 !%Type logical
892 !%Default no
893 !%Section Time-Dependent
894 !%Description
895 !% The electrons are evolved as independent particles feeling the Hartree and
896 !% exchange-correlation potentials from the ground-state electronic configuration.
897 !%End
898 call parse_variable(namespace, 'TDFreezeHXC', .false., freeze_hxc)
899 if (freeze_hxc) then
900 write(message(1),'(a)') 'Info: Freezing Hartree and exchange-correlation potentials.'
901 call messages_info(1, namespace=namespace)
902
903 if (.not. from_scratch) then
904
905 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr, exact=.true.)
906 call states_elec_load(restart_frozen, namespace, space, st, gr, hm%kpoints, ierr, label = ": gs")
907 call states_elec_transform(st, namespace, space, restart_frozen, gr, hm%kpoints)
908 call restart_frozen%end()
909
910 call density_calc(st, gr, st%rho)
911 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
912
913 call restart_frozen%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
914 call states_elec_load(restart_frozen, namespace, space, st, gr, hm%kpoints, ierr, iter=td%iter, label = ": td")
915 call restart_frozen%end()
916 call hm%ks_pot%run_zero_iter(td%tr%vks_old)
917
918 end if
919
920 call v_ks_freeze_hxc(ks)
921
922 end if
923
924 x = minval(st%eigenval(st%st_start, :))
925 if (st%parallel_in_states) then
926 call st%mpi_grp%bcast(x, 1, mpi_double_precision, 0)
927 end if
928 call hm%update_span(gr%spacing(1:space%dim), x, namespace)
929 ! initialize Fermi energy
930 call states_elec_fermi(st, namespace, gr, compute_spin = .not. gr%der%boundaries%spiralBC)
931 call energy_calc_total(namespace, space, hm, gr, st, ext_partners)
932
933 !%Variable TDFreezeDFTUOccupations
934 !%Type logical
935 !%Default no
936 !%Section Time-Dependent
937 !%Description
938 !% The occupation matrices than enters in the DFT+U potential
939 !% are not evolved during the time evolution.
940 !%End
941 call parse_variable(namespace, 'TDFreezeDFTUOccupations', .false., freeze_occ)
942 if (freeze_occ) then
943 write(message(1),'(a)') 'Info: Freezing DFT+U occupation matrices that enters in the DFT+U potential.'
944 call messages_info(1, namespace=namespace)
945 call lda_u_freeze_occ(hm%lda_u)
946
947 !In this case we should reload GS wavefunctions
948 if (hm%lda_u_level /= dft_u_none .and. .not. from_scratch) then
949 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr)
950 call lda_u_load(restart_frozen, hm%lda_u, st, hm%energy%dft_u, ierr, occ_only = .true.)
951 call restart_frozen%end()
952 end if
953 end if
954
955 !%Variable TDFreezeU
956 !%Type logical
957 !%Default no
958 !%Section Time-Dependent
959 !%Description
960 !% The effective U of DFT+U is not evolved during the time evolution.
961 !%End
962 call parse_variable(namespace, 'TDFreezeU', .false., freeze_u)
963 if (freeze_u) then
964 write(message(1),'(a)') 'Info: Freezing the effective U of DFT+U.'
965 call messages_info(1, namespace=namespace)
966 call lda_u_freeze_u(hm%lda_u)
967
968 !In this case we should reload GS wavefunctions
969 if (hm%lda_u_level == dft_u_acbn0 .and. .not. from_scratch) then
970 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr)
971 call lda_u_load(restart_frozen, hm%lda_u, st, hm%energy%dft_u, ierr, u_only = .true.)
972 call restart_frozen%end()
973 write(message(1),'(a)') 'Loaded GS effective U of DFT+U'
974 call messages_info(1, namespace=namespace)
975 call lda_u_write_u(hm%lda_u, namespace=namespace)
976 call lda_u_write_v(hm%lda_u, namespace=namespace)
977 end if
978 end if
979
980 ! This needs to be called before the calculation of the forces,
981 ! as we need to test of we output the forces or not
982 call td_write_init(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
983 ks, td%ions_dyn%is_active(), &
984 list_has_gauge_field(ext_partners), hm%kick, td%iter, td%max_iter, td%dt, mc)
985
986 ! Resets the nondipole integration after laser-file has been written.
987 lasers => list_get_lasers(ext_partners)
988 if(associated(lasers)) then
989 if (lasers_with_nondipole_field(lasers)) then
990 ndinitial(1:space%dim)=m_zero
991 call lasers_set_nondipole_parameters(lasers,ndinitial,m_zero)
992 end if
993 end if
994 nullify(lasers)
995
996 call td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
997
998 if (td%scissor > m_epsilon) then
999 call scissor_init(hm%scissor, namespace, space, st, gr, hm%kpoints, hm%phase, td%scissor, mc)
1000 end if
1001
1002 if (td%iter == 0) call td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc)
1003
1004 gfield => list_get_gauge_field(ext_partners)
1005 if(associated(gfield)) then
1006 if (gauge_field_is_propagated(gfield)) then
1007 call gauge_field_get_force(gfield, gr, st%d%spin_channels, st%current, ks%xc%lrc)
1008 if (abs(ks%xc%lrc%alpha) > m_epsilon) then
1009 call messages_experimental('TD-LRC kernel')
1010 end if
1011 endif
1012 end if
1013
1014 !call td_check_trotter(td, sys, h)
1015 td%iter = td%iter + 1
1016
1017 call td%restart_dump%init(namespace, restart_td, restart_type_dump, mc, ierr, mesh=gr)
1018 if (td%ions_dyn%is_active() .and. td%recalculate_gs) then
1019 ! We will also use the TD restart directory as temporary storage during the time propagation
1020 call td%restart_load%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
1021 end if
1022
1023 call messages_print_with_emphasis(msg="Time-Dependent Simulation", namespace=namespace)
1024 call td_print_header(namespace)
1025
1026 if (td%pesv%calc_spm .or. td%pesv%calc_mask .and. from_scratch) then
1027 call pes_init_write(td%pesv,gr,st, namespace)
1028 end if
1029
1030 if (st%pack_states .and. hm%apply_packed()) call st%pack()
1031
1033 end subroutine td_init_with_wavefunctions
1034
1035 ! ---------------------------------------------------------
1036 subroutine td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
1037 type(td_t), intent(inout) :: td
1038 type(namespace_t), intent(in) :: namespace
1039 type(electron_space_t), intent(in) :: space
1040 type(grid_t), intent(inout) :: gr
1041 type(ions_t), intent(inout) :: ions
1042 type(partner_list_t), intent(in) :: ext_partners
1043 type(states_elec_t), target, intent(inout) :: st
1044 type(v_ks_t), intent(inout) :: ks
1045 type(hamiltonian_elec_t), intent(inout) :: hm
1046 type(output_t), intent(inout) :: outp
1047
1048 push_sub(td_init_ions_and_forces)
1049
1050 ! Calculate initial forces and kinetic energy
1051 if (td%ions_dyn%ions_move()) then
1052 if (td%iter > 0) then
1053 call td_read_coordinates(td, namespace, ions)
1054 if (ion_dynamics_drive_ions(td%ions_dyn)) then
1055 call ion_dynamics_propagate(td%ions_dyn, ions, td%iter*td%dt, td%dt, namespace)
1056 end if
1057 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = td%iter*td%dt)
1058 ! recompute potential because the ions have moved
1059 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
1060 end if
1061
1062 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = td%iter*td%dt, dt = td%dt)
1063
1064 call ions%update_kinetic_energy()
1065 else
1066 if (outp%what(option__output__forces) .or. td%write_handler%out(out_separate_forces)%write) then
1067 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = td%iter*td%dt, dt = td%dt)
1068 end if
1069 end if
1070
1071 if (outp%what(option__output__stress) .or. td%ions_dyn%cell_relax()) then
1072 call stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
1073 if (td%ions_dyn%cell_relax()) then
1074 call td%ions_dyn%update_stress(ions%space, st%stress_tensors%total, ions%latt%rlattice, ions%latt%rcell_volume)
1075 end if
1076 end if
1077
1079 end subroutine td_init_ions_and_forces
1080
1081 ! ---------------------------------------------------------
1082 subroutine td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
1083 type(td_t), intent(inout) :: td
1084 type(namespace_t), intent(in) :: namespace
1085 class(space_t), intent(in) :: space
1086 type(multicomm_t), intent(in) :: mc
1087 type(grid_t), intent(inout) :: gr
1088 type(partner_list_t), intent(in) :: ext_partners
1089 type(states_elec_t), target, intent(inout) :: st
1090 type(v_ks_t), intent(inout) :: ks
1091 type(hamiltonian_elec_t), intent(inout) :: hm
1092 logical, intent(inout) :: from_scratch
1093
1094 integer :: ierr
1095 type(restart_t) :: restart
1096
1097 push_sub(td_load_restart_from_td)
1098
1099 !We redistribute the states before the restarting
1100 if (td%freeze_orbitals > 0) then
1101 call states_elec_freeze_redistribute_states(st, namespace, gr, mc, td%freeze_orbitals)
1102 end if
1103
1104 call restart%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
1105 if (ierr == 0) then
1106 call td_load(restart, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
1107 end if
1108 call restart%end()
1109 if (ierr /= 0) then
1110 from_scratch = .true.
1111 td%iter = 0
1112 end if
1113
1115 end subroutine td_load_restart_from_td
1116
1117 ! ---------------------------------------------------------
1118 subroutine td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
1119 type(td_t), intent(inout) :: td
1120 type(namespace_t), intent(in) :: namespace
1121 class(space_t), intent(in) :: space
1122 type(multicomm_t), intent(in) :: mc
1123 type(grid_t), intent(inout) :: gr
1124 type(partner_list_t), intent(in) :: ext_partners
1125 type(states_elec_t), target, intent(inout) :: st
1126 type(v_ks_t), intent(inout) :: ks
1127 type(hamiltonian_elec_t), intent(inout) :: hm
1128
1129 integer :: ierr
1130 type(restart_t) :: restart
1132 push_sub(td_load_restart_from_gs)
1133
1134 call restart%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr, exact=.true.)
1135
1136 if (.not. st%only_userdef_istates) then
1137 if (ierr == 0) then
1138 call states_elec_load(restart, namespace, space, st, gr, hm%kpoints, ierr, label = ": gs")
1139 end if
1140 if (ierr /= 0) then
1141 message(1) = 'Unable to read ground-state wavefunctions.'
1142 call messages_fatal(1, namespace=namespace)
1143 end if
1144 end if
1145
1146 ! check if we should deploy user-defined wavefunctions.
1147 ! according to the settings in the input file the routine
1148 ! overwrites orbitals that were read from restart/gs
1149 if (parse_is_defined(namespace, 'UserDefinedStates')) then
1150 call states_elec_read_user_def_orbitals(gr, namespace, space, st)
1151 end if
1152
1153 call states_elec_transform(st, namespace, space, restart, gr, hm%kpoints)
1154 call restart%end()
1155
1157 end subroutine td_load_restart_from_gs
1158
1159 ! ---------------------------------------------------------
1160 subroutine td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc)
1161 type(td_t), intent(inout) :: td
1162 type(namespace_t), intent(in) :: namespace
1163 type(electron_space_t), intent(in) :: space
1164 type(grid_t), intent(inout) :: gr
1165 type(ions_t), intent(inout) :: ions
1166 type(states_elec_t), intent(inout) :: st
1167 type(v_ks_t), intent(inout) :: ks
1168 type(hamiltonian_elec_t), intent(inout) :: hm
1169 type(partner_list_t), intent(in) :: ext_partners
1170 type(output_t), intent(in) :: outp
1171 type(multicomm_t), intent(in) :: mc
1172
1173 push_sub(td_run_zero_iter)
1174
1175 call td_write_iter(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
1176 hm%kick, ks, td%dt, 0, mc, td%recalculate_gs)
1178 ! I apply the delta electric field *after* td_write_iter, otherwise the
1179 ! dipole matrix elements in write_proj are wrong
1180 if (abs(hm%kick%time) <= m_epsilon) then
1181 if (.not. hm%pcm%localf) then
1182 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints)
1183 else
1184 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints, pcm = hm%pcm)
1185 end if
1186 call td_write_kick(outp, namespace, space, gr, hm%kick, ions, 0)
1187
1188 !We activate the sprial BC only after the kick
1189 if (gr%der%boundaries%spiralBC) then
1190 gr%der%boundaries%spiral = .true.
1191 end if
1192 end if
1193 call hm%ks_pot%run_zero_iter(td%tr%vks_old)
1194
1195 if (any(outp%output_interval > 0)) then
1196 call td_write_data(td%write_handler)
1197 call td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, 0)
1198 end if
1199
1200 pop_sub(td_run_zero_iter)
1201 end subroutine td_run_zero_iter
1202
1203
1204 ! ---------------------------------------------------------
1206 subroutine td_read_coordinates(td, namespace, ions)
1207 type(td_t), intent(in) :: td
1208 type(namespace_t), intent(in) :: namespace
1209 type(ions_t), intent(inout) :: ions
1210
1211 integer :: iatom, iter, iunit
1212
1214
1215 iunit = io_open('td.general/coordinates', namespace, action='read', status='old', die=.false.)
1216 if (iunit == -1) then
1217 message(1) = "Could not open file '"//trim(io_workpath('td.general/coordinates', namespace))//"'."
1218 message(2) = "Starting simulation from initial geometry."
1219 call messages_warning(2, namespace=namespace)
1220 pop_sub(td_read_coordinates)
1221 return
1222 end if
1223
1224 call io_skip_header(iunit)
1225 do iter = 0, td%iter - 1
1226 read(iunit, *) ! skip previous iterations... sorry, but no portable seek in Fortran
1227 end do
1228 read(iunit, '(32x)', advance='no') ! skip the time index.
1229
1230 do iatom = 1, ions%natoms
1231 read(iunit, '(3es24.16)', advance='no') ions%pos(:, iatom)
1232 ions%pos(:, iatom) = units_to_atomic(units_out%length, ions%pos(:, iatom))
1233 end do
1234 do iatom = 1, ions%natoms
1235 read(iunit, '(3es24.16)', advance='no') ions%vel(:, iatom)
1236 ions%vel(:, iatom) = units_to_atomic(units_out%velocity, ions%vel(:, iatom))
1237 end do
1238 do iatom = 1, ions%natoms
1239 read(iunit, '(3es24.16)', advance='no') ions%tot_force(:, iatom)
1240 ions%tot_force(:, iatom) = units_to_atomic(units_out%force, ions%tot_force(:, iatom))
1241 end do
1242
1243 call io_close(iunit)
1244
1245 pop_sub(td_read_coordinates)
1246 end subroutine td_read_coordinates
1247
1248 ! ---------------------------------------------------------
1249 subroutine td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
1250 type(td_t), intent(in) :: td
1251 type(namespace_t), intent(in) :: namespace
1252 class(space_t), intent(in) :: space
1253 type(grid_t), intent(in) :: gr
1254 type(states_elec_t), intent(in) :: st
1255 type(hamiltonian_elec_t), intent(in) :: hm
1256 type(v_ks_t), intent(in) :: ks
1257 type(partner_list_t), intent(in) :: ext_partners
1258 integer, intent(in) :: iter
1259 integer, intent(out) :: ierr
1260
1261 type(gauge_field_t), pointer :: gfield
1262 integer :: err, err2
1263
1264 push_sub(td_dump)
1265
1266 ierr = 0
1267
1268 if (td%restart_dump%skip()) then
1269 pop_sub(td_dump)
1270 return
1271 end if
1272
1273 message(1) = "Debug: Writing td restart."
1274 call messages_info(1, namespace=namespace, debug_only=.true.)
1275
1276 ! first write resume file
1277 call states_elec_dump(td%restart_dump, space, st, gr, hm%kpoints, err, iter=iter)
1278 if (err /= 0) ierr = ierr + 1
1279
1280 call states_elec_dump_rho(td%restart_dump, st, gr, ierr, iter=iter)
1281 if (err /= 0) ierr = ierr + 1
1282
1283 if (hm%lda_u_level /= dft_u_none) then
1284 call lda_u_dump(td%restart_dump, namespace, hm%lda_u, st, gr, ierr)
1285 if (err /= 0) ierr = ierr + 1
1286 end if
1287
1288 call potential_interpolation_dump(td%tr%vks_old, td%restart_dump, gr, st%d%nspin, err2)
1289 if (err2 /= 0) ierr = ierr + 2
1290
1291 call pes_dump(td%pesv, namespace, td%restart_dump, st, gr, err)
1292 if (err /= 0) ierr = ierr + 4
1293
1294 ! Gauge field restart
1295 gfield => list_get_gauge_field(ext_partners)
1296 if(associated(gfield)) then
1297 call gauge_field_dump(td%restart_dump, gfield, ierr)
1298 end if
1299
1300 if (gr%der%boundaries%spiralBC) then
1301 call states_elec_dump_spin(td%restart_dump, st, err)
1302 if (err /= 0) ierr = ierr + 8
1303 end if
1304
1305 if (ks%has_photons) then
1306 call mf_photons_dump(td%restart_dump, ks%pt_mx, gr, td%dt, ks%pt, err)
1307 if (err /= 0) ierr = ierr + 16
1308 end if
1309
1310 if (ks%xc_photon /= 0) then
1311 ! photon-free mean field
1312 call ks%xc_photons%mf_dump(td%restart_dump, err)
1313 if (err /= 0) ierr = ierr + 32
1314 end if
1315
1316 if (allocated(st%frozen_rho)) then
1317 call states_elec_dump_frozen(td%restart_dump, space, st, gr, ierr)
1318 end if
1319 if (err /= 0) ierr = ierr + 64
1320
1321 if (td%ions_dyn%ions_move() .or. td%ions_dyn%cell_relax()) then
1322 call ion_dynamics_dump(td%ions_dyn, td%restart_dump, err)
1323 end if
1324 if (err /= 0) ierr = ierr + 128
1325
1326 message(1) = "Debug: Writing td restart done."
1327 call messages_info(1, namespace=namespace, debug_only=.true.)
1328
1329 pop_sub(td_dump)
1330 end subroutine td_dump
1331
1332 ! ---------------------------------------------------------
1333 subroutine td_load(restart, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
1334 type(restart_t), intent(in) :: restart
1335 type(namespace_t), intent(in) :: namespace
1336 class(space_t), intent(in) :: space
1337 type(grid_t), intent(in) :: gr
1338 type(states_elec_t), intent(inout) :: st
1339 type(hamiltonian_elec_t), intent(inout) :: hm
1340 type(partner_list_t), intent(in) :: ext_partners
1341 type(td_t), intent(inout) :: td
1342 type(v_ks_t), intent(inout) :: ks
1343 integer, intent(out) :: ierr
1345 integer :: err, err2
1346 type(gauge_field_t), pointer :: gfield
1347 push_sub(td_load)
1348
1349 ierr = 0
1350
1351 if (restart%skip()) then
1352 ierr = -1
1353 pop_sub(td_load)
1354 return
1355 end if
1356
1357 message(1) = "Debug: Reading td restart."
1358 call messages_info(1, namespace=namespace, debug_only=.true.)
1359
1360 ! Read states
1361 call states_elec_load(restart, namespace, space, st, gr, hm%kpoints, err, iter=td%iter, label = ": td")
1362 if (err /= 0) then
1363 ierr = ierr + 1
1364 end if
1365
1366 ! read potential from previous interactions
1367 call potential_interpolation_load(td%tr%vks_old, namespace, restart, gr, st%d%nspin, err2)
1368 if (err2 /= 0) ierr = ierr + 2
1369
1370 if (hm%lda_u_level /= dft_u_none) then
1371 call lda_u_load(restart, hm%lda_u, st, hm%energy%dft_u, ierr)
1372 if (err /= 0) ierr = ierr + 1
1373 end if
1374
1375
1376 ! read PES restart
1377 if (td%pesv%calc_spm .or. td%pesv%calc_mask .or. td%pesv%calc_flux) then
1378 call pes_load(td%pesv, namespace, restart, st, err)
1379 if (err /= 0) ierr = ierr + 4
1380 end if
1381
1382 ! Gauge field restart
1383 gfield => list_get_gauge_field(ext_partners)
1384 if (associated(gfield)) then
1385 call gauge_field_load(restart, gfield, err)
1386 if (err /= 0) then
1387 ierr = ierr + 8
1388 else
1389 call hm%update(gr, namespace, space, ext_partners, time = td%dt*td%iter)
1390 end if
1391 end if
1392
1393 ! add photon restart
1394 if (ks%has_photons) then
1395 call mf_photons_load(restart, ks%pt_mx, gr, err)
1396 end if
1397 if (err /= 0) ierr = ierr + 16
1398
1399 if (ks%xc_photon /= 0) then
1400 call ks%xc_photons%mf_load(restart, space, err)
1401 end if
1402 if (err /= 0) ierr = ierr + 32
1403
1404 if (gr%der%boundaries%spiralBC) then
1405 call states_elec_load_spin(restart, st, err)
1406 !To ensure back compatibility, if the file is not present, we use the
1407 !current states to get the spins
1408 if (err /= 0) call states_elec_fermi(st, namespace, gr)
1409 end if
1410
1411 if (td%ions_dyn%is_active()) then
1412 call ion_dynamics_load(td%ions_dyn, restart, err)
1413 end if
1414 if (err /= 0) ierr = ierr + 64
1415
1416 message(1) = "Debug: Reading td restart done."
1417 call messages_info(1, namespace=namespace, debug_only=.true.)
1418
1419 pop_sub(td_load)
1420 end subroutine td_load
1421 ! ---------------------------------------------------------
1422 subroutine td_load_frozen(namespace, restart, space, mesh, st, hm, ierr)
1423 type(namespace_t), intent(in) :: namespace
1424 type(restart_t), intent(in) :: restart
1425 class(space_t), intent(in) :: space
1426 class(mesh_t), intent(in) :: mesh
1427 type(states_elec_t), intent(inout) :: st
1428 type(hamiltonian_elec_t), intent(inout) :: hm
1429 integer, intent(out) :: ierr
1430
1431 push_sub(td_load_frozen)
1432
1433 ierr = 0
1434
1435 if (restart%skip()) then
1436 ierr = -1
1437 pop_sub(td_load_frozen)
1438 return
1439 end if
1440
1441 message(1) = "Debug: Reading td frozen restart."
1442 call messages_info(1, namespace=namespace, debug_only=.true.)
1443
1444 safe_allocate(st%frozen_rho(1:mesh%np, 1:st%d%nspin))
1445 if (family_is_mgga(hm%xc%family)) then
1446 safe_allocate(st%frozen_tau(1:mesh%np, 1:st%d%nspin))
1447 safe_allocate(st%frozen_gdens(1:mesh%np, 1:space%dim, 1:st%d%nspin))
1448 safe_allocate(st%frozen_ldens(1:mesh%np, 1:st%d%nspin))
1449 end if
1450
1451 call states_elec_load_frozen(restart, space, st, mesh, ierr)
1452
1453 message(1) = "Debug: Reading td frozen restart done."
1454 call messages_info(1, namespace=namespace, debug_only=.true.)
1455
1456 pop_sub(td_load_frozen)
1457 end subroutine td_load_frozen
1458
1459 ! ---------------------------------------------------------
1460 logical function td_get_from_scratch(td)
1461 type(td_t), intent(in) :: td
1462
1463 push_sub(td_get_from_scratch)
1464
1465 td_get_from_scratch = td%from_scratch
1466
1467 pop_sub(td_get_from_scratch)
1468 end function td_get_from_scratch
1469
1470 ! ---------------------------------------------------------
1471 subroutine td_set_from_scratch(td, from_scratch)
1472 type(td_t), intent(inout) :: td
1473 logical, intent(in) :: from_scratch
1474
1475 push_sub(td_set_from_scratch)
1476
1477 td%from_scratch = from_scratch
1478
1479 pop_sub(td_set_from_scratch)
1480 end subroutine td_set_from_scratch
1481end module td_oct_m
1482
1483!! Local Variables:
1484!! mode: f90
1485!! coding: utf-8
1486!! End:
Define which routines can be seen from the outside.
Definition: loct.F90:149
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
integer, parameter, public mask_absorbing
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
This module handles the calculation mode.
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public states_elec_freeze_adjust_qtot(st)
Definition: density.F90:814
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
Definition: density.F90:786
subroutine, public states_elec_freeze_orbitals(st, namespace, space, gr, mc, kpoints, n, family_is_mgga)
Calculate partial density for frozen orbitals.
Definition: density.F90:655
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
subroutine, public dm_propagation_init_run(adiabatic_st, namespace, space, gr, st, hm, mc)
Initialise the adiabatic states prior to running TB propagation.
subroutine, public dm_propagation_run(dmp, namespace, space, gr, ions, st, mc, hm, ks, iter, dt, ext_partners, update_energy)
Apply dissipation to a TD run via the Linblad formalism.
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 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,...
integer, parameter, public spin_orbit
Definition: epot.F90:168
logical function, public list_has_gauge_field(partners)
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
subroutine, public forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
Definition: forces.F90:340
subroutine, public gauge_field_load(restart, gfield, ierr)
subroutine, public gauge_field_get_force(this, gr, spin_channels, current, lrc)
subroutine, public gauge_field_dump(restart, gfield, ierr)
logical pure function, public gauge_field_is_propagated(this)
logical pure function, public gauge_field_is_used(this)
subroutine, public gauge_field_init_vec_pot(this, qtot)
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public zvmask(mesh, hm, st)
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public dvmask(mesh, hm, st)
This module defines classes and functions for interaction partners.
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:466
subroutine, public io_skip_header(iunit)
Definition: io.F90:645
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
Definition: io.F90:317
subroutine, public io_debug_on_the_fly(namespace)
check if debug mode should be enabled or disabled on the fly
Definition: io.F90:534
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:401
subroutine, public ion_dynamics_dump(this, restart, ierr)
subroutine, public ion_dynamics_propagate(this, ions, time, dt, namespace)
Interface for the ion/cell dynamics.
subroutine, public ion_dynamics_load(this, restart, ierr)
subroutine, public ion_dynamics_init(this, namespace, ions, symmetrize, symm)
subroutine, public ion_dynamics_end(this)
logical pure function, public ion_dynamics_drive_ions(this)
Is the ion dynamics activated or not.
subroutine, public kick_apply(space, mesh, st, ions_dyn, ions, kick, psolver, kpoints, pcm)
Applies the delta-function electric field where k = kick%delta_strength.
Definition: kick.F90:1137
subroutine, public lasers_set_nondipole_parameters(this, ndfield, nd_integration_time)
Set parameters for nondipole SFA calculation.
Definition: lasers.F90:754
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:741
subroutine, public laser_write_info(lasers, namespace, dt, max_iter, iunit)
Definition: lasers.F90:932
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
Definition: lda_u_io.F90:648
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:532
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:729
subroutine, public lda_u_write_v(this, iunit, namespace)
Definition: lda_u_io.F90:581
integer, parameter, public dft_u_none
Definition: lda_u.F90:203
subroutine, public lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints)
Definition: lda_u.F90:284
subroutine, public lda_u_freeze_occ(this)
Definition: lda_u.F90:887
subroutine, public lda_u_freeze_u(this)
Definition: lda_u.F90:894
subroutine, public lda_u_end(this)
Definition: lda_u.F90:655
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:798
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:203
This module implements fully polymorphic linked lists, and some specializations thereof.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
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
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
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
type(namespace_t), public global_namespace
Definition: namespace.F90:134
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 parse_is_defined(namespace, name)
Definition: parser.F90:455
subroutine, public pes_calc(pes, namespace, space, mesh, st, dt, iter, der, kpoints, ext_partners, stopping)
Definition: pes.F90:271
subroutine, public pes_init(pes, namespace, space, mesh, box, st, save_iter, kpoints, abs_boundaries, ext_partners, max_iter, dt)
Definition: pes.F90:177
subroutine, public pes_output(pes, namespace, space, gr, st, iter, outp, dt, ions)
Definition: pes.F90:297
subroutine, public pes_init_write(pes, mesh, st, namespace)
Definition: pes.F90:401
subroutine, public pes_end(pes)
Definition: pes.F90:257
subroutine, public pes_load(pes, namespace, restart, st, ierr)
Definition: pes.F90:361
subroutine, public pes_dump(pes, namespace, restart, st, mesh, ierr)
Definition: pes.F90:321
subroutine, public mf_photons_load(restart, this, gr, ierr)
subroutine, public mf_photons_dump(restart, this, gr, dt, pt_mode, ierr)
subroutine, public potential_interpolation_load(potential_interpolation, namespace, restart, mesh, nspin, err2)
subroutine, public potential_interpolation_dump(potential_interpolation, restart, mesh, nspin, err2)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:625
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public propagator_elec_dt_bo(scf, namespace, space, gr, ks, st, hm, ions, ext_partners, mc, iter, dt, ions_dyn, scsteps)
subroutine, public propagator_elec_init(gr, namespace, st, tr, ks_pot, have_fields, family_is_mgga_with_exc, relax_cell)
subroutine, public propagator_elec_dt(ks, namespace, space, hm, gr, st, tr, time, dt, nt, ions_dyn, ions, ext_partners, mc, outp, write_handler, scsteps, update_energy, qcchi)
Propagates st from time - dt to t. If dt<0, it propagates backwards from t+|dt| to t.
subroutine, public propagator_elec_end(tr)
This module implements the basic propagator framework.
Definition: propagator.F90:119
logical function, public clean_stop(comm)
returns true if a file named stop exists
Definition: restart.F90:335
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:183
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:256
subroutine, public scf_end(scf)
Definition: scf.F90:526
subroutine, public scissor_init(this, namespace, space, st, mesh, kpoints, phase, gap, mc)
Definition: scissor.F90:164
pure logical function, public states_are_real(st)
subroutine, public states_elec_fermi(st, namespace, mesh, compute_spin)
calculate the Fermi level for the states in this object
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_read_user_def_orbitals(mesh, namespace, space, st)
the routine reads formulas for user-defined wavefunctions from the input file and fills the respectiv...
subroutine, public states_elec_load_frozen(restart, space, st, mesh, ierr)
subroutine, public states_elec_transform(st, namespace, space, restart, mesh, kpoints, prefix)
subroutine, public states_elec_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, verbose)
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public states_elec_load_spin(restart, st, ierr)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
subroutine, public states_elec_dump_frozen(restart, space, st, mesh, ierr)
subroutine, public states_elec_dump_rho(restart, st, mesh, ierr, iter)
subroutine, public states_elec_dump_spin(restart, st, ierr)
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
Definition: td.F90:116
subroutine, public td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, iter, scsteps, etime, stopping, from_scratch)
Definition: td.F90:803
subroutine, public td_end(td)
Definition: td.F90:642
subroutine, public td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
Definition: td.F90:1214
subroutine td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
Definition: td.F90:869
subroutine, public td_end_run(td, st, hm)
Definition: td.F90:661
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:578
logical function, public td_get_from_scratch(td)
Definition: td.F90:1556
subroutine td_load_frozen(namespace, restart, space, mesh, st, hm, ierr)
Definition: td.F90:1518
subroutine td_print_header(namespace)
Definition: td.F90:788
integer, parameter, public bo
Definition: td.F90:203
subroutine td_load(restart, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
Definition: td.F90:1429
subroutine, public td_set_from_scratch(td, from_scratch)
Definition: td.F90:1567
subroutine, public td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
Definition: td.F90:1345
subroutine td_read_coordinates(td, namespace, ions)
reads the pos and vel from coordinates file
Definition: td.F90:1302
subroutine, public td_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, from_scratch)
Definition: td.F90:681
subroutine td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc)
Definition: td.F90:1256
subroutine, public td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
Definition: td.F90:611
subroutine td_update_elapsed_time(etime)
Definition: td.F90:890
subroutine, public td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, from_scratch)
Definition: td.F90:901
subroutine, public td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
Definition: td.F90:1178
subroutine td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
Definition: td.F90:1132
subroutine, public td_init_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, from_scratch)
Definition: td.F90:522
subroutine, public td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, dt)
Definition: td_write.F90:1239
subroutine, public td_write_data(writ)
Definition: td_write.F90:1205
subroutine, public td_write_kick(outp, namespace, space, mesh, kick, ions, iter)
Definition: td_write.F90:337
subroutine, public td_write_end(writ)
Definition: td_write.F90:989
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:1036
subroutine, public td_write_init(writ, namespace, space, outp, gr, st, hm, ions, ext_partners, ks, ions_move, with_gauge_field, kick, iter, max_iter, dt, mc)
Initialize files to write when prograting in time.
Definition: td_write.F90:375
integer, parameter, public out_separate_forces
Definition: td_write.F90:203
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_t), public unit_femtosecond
Time in femtoseconds.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine, public v_ks_freeze_hxc(ks)
Definition: v_ks.F90:1483
subroutine, public v_ks_calculate_current(this, calc_cur)
Definition: v_ks.F90:1494
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:749
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
Definition: walltimer.F90:123
logical function, public walltimer_alarm(comm, print)
indicate whether time is up
Definition: walltimer.F90:333
real(real64) function, public walltimer_get_start_time()
Return the walltimer start time.
Definition: walltimer.F90:407
logical function, public restart_walltime_period_alarm(comm)
Definition: walltimer.F90:375
Definition: xc.F90:116
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:580
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:593
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
output handler class
Definition: output_low.F90:166
The states_elec_t class contains all electronic wave functions.
int true(void)