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 end type td_t
140
141
142contains
143
144 subroutine td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, dmp)
145 type(td_t), intent(inout) :: td
146 type(namespace_t), intent(in) :: namespace
147 class(space_t), intent(in) :: space
148 type(grid_t), intent(in) :: gr
149 type(ions_t), intent(inout) :: ions
150 type(states_elec_t), intent(in) :: st
151 type(v_ks_t), intent(in) :: ks
152 type(hamiltonian_elec_t), intent(in) :: hm
153 type(partner_list_t), intent(in) :: ext_partners
154 type(output_t), intent(in) :: outp
155 type(dmp_t), intent(inout) :: dmp
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
240
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 lasers => list_get_lasers(ext_partners)
253
254 ! now the photoelectron stuff
255 call pes_init(td%pesv, namespace, space, gr, gr%box, st, 1, hm%kpoints, &
256 hm%abs_boundaries, ext_partners, td%max_iter, td%dt)
257
258 !%Variable TDDynamics
259 !%Type integer
260 !%Default ehrenfest
261 !%Section Time-Dependent::Propagation
262 !%Description
263 !% Type of dynamics to follow during a time propagation.
264 !% For BO, you must set <tt>MoveIons = yes</tt>.
265 !%Option ehrenfest 1
266 !% Ehrenfest dynamics.
267 !%Option bo 2
268 !% Born-Oppenheimer (Experimental).
269 !%End
270
271 call parse_variable(namespace, 'TDDynamics', ehrenfest, td%dynamics)
272 if (.not. varinfo_valid_option('TDDynamics', td%dynamics)) call messages_input_error(namespace, 'TDDynamics')
273 call messages_print_var_option('TDDynamics', td%dynamics, namespace=namespace)
274 if (td%dynamics .ne. ehrenfest) then
275 if (.not. td%ions_dyn%is_active()) then
276 message(1) = "TDDynamics=bo can only be used if MoveIons=yes or CellDynamics=yes"
277 call messages_fatal(1, namespace=namespace)
278 end if
279 end if
280
281 !%Variable RecalculateGSDuringEvolution
282 !%Type logical
283 !%Default no
284 !%Section Time-Dependent::Propagation
285 !%Description
286 !% In order to calculate some information about the system along the
287 !% evolution (e.g. projection onto the ground-state KS determinant,
288 !% projection of the TDKS spin-orbitals onto the ground-state KS
289 !% spin-orbitals), the ground-state KS orbitals are needed. If the
290 !% ionic potential changes -- that is, the ions move -- one may want
291 !% to recalculate the ground state. You may do this by setting this
292 !% variable.
293 !%
294 !% The recalculation is not done every time step, but only every
295 !% <tt>RecalculateGSInterval</tt> time steps.
296 !%End
297 call parse_variable(namespace, 'RecalculateGSDuringEvolution', .false., td%recalculate_gs)
298 if (hm%lda_u_level /= dft_u_none .and. td%recalculate_gs) then
299 call messages_not_implemented("DFT+U with RecalculateGSDuringEvolution=yes", namespace=namespace)
300 end if
301
302 !%Variable RecalculateGSInterval
303 !%Type integer
304 !%Default 50
305 !%Section Time-Dependent::Propagation
306 !%Description
307 !% This variable controls how often the GS is recalculated for
308 !% <tt>RecalculateGSDuringEvolution = yes</tt>. It is the number of
309 !% time-propagation steps between two ground-state recalculations.
310 !% A smaller value follows the instantaneous Hamiltonian more closely,
311 !% while a larger value reduces the cost of outputs that need the
312 !% instantaneous ground-state basis, such as <tt>TDOutput = n_excited_el</tt>.
313 !%End
314 call parse_variable(namespace, 'RecalculateGSInterval', 50, td%recalculate_gs_interval)
315
316 !%Variable TDScissor
317 !%Type float
318 !%Default 0.0
319 !%Section Time-Dependent
320 !%Description
321 !% (experimental) If set, a scissor operator will be applied in the
322 !% Hamiltonian, shifting the excitation energies by the amount
323 !% specified. By default, it is not applied.
324 !%End
325 call parse_variable(namespace, 'TDScissor', m_zero, td%scissor)
326 call messages_print_var_value('TDScissor', td%scissor, namespace=namespace)
327
328 call propagator_elec_init(gr, namespace, st, td%tr, hm%ks_pot, td%ions_dyn%is_active() .and.&
329 list_has_gauge_field(ext_partners), family_is_mgga_with_exc(ks%xc), td%ions_dyn%cell_relax())
330
331 if (associated(lasers) .and. st%system_grp%is_root()) then
332 call messages_print_with_emphasis(msg="Time-dependent external fields", namespace=namespace)
333 call laser_write_info(lasers%lasers, dt=td%dt, max_iter=td%max_iter, namespace=namespace)
334 call messages_print_with_emphasis(namespace=namespace)
335 end if
336
337 !%Variable TDEnergyUpdateIter
338 !%Type integer
339 !%Section Time-Dependent::Propagation
340 !%Description
341 !% This variable controls after how many iterations Octopus
342 !% updates the total energy during a time-propagation run. For
343 !% iterations where the energy is not updated, the last calculated
344 !% value is reported. If you set this variable to 1, the energy
345 !% will be calculated in each step.
346 !%End
347
348 default = 10
349 call parse_variable(namespace, 'TDEnergyUpdateIter', default, td%energy_update_iter)
350 if (td%energy_update_iter < 1) then
351 call messages_input_error(namespace, 'TDEnergyUpdateIter')
352 end if
353
354 if (gr%der%boundaries%spiralBC .and. hm%ep%reltype == spin_orbit) then
355 message(1) = "Generalized Bloch theorem cannot be used with spin-orbit coupling."
356 call messages_fatal(1, namespace=namespace)
357 end if
358
359 if (gr%der%boundaries%spiralBC) then
360 if (any(abs(hm%kick%easy_axis(1:2)) > m_epsilon)) then
361 message(1) = "Generalized Bloch theorem cannot be used for an easy axis not along the z direction."
362 call messages_fatal(1, namespace=namespace)
363 end if
364 end if
365
366 !%Variable TDFreezeOrbitals
367 !%Type integer
368 !%Default 0
369 !%Section Time-Dependent
370 !%Description
371 !% (Experimental) You have the possibility of "freezing" a number of orbitals during a time-propagation.
372 !% The Hartree and exchange-correlation potential due to these orbitals (which
373 !% will be the lowest-energy ones) will be added during the propagation, but the orbitals
374 !% will not be propagated.
375 !%Option sae -1
376 !% Single-active-electron approximation. This option is only valid for time-dependent
377 !% calculations (<tt>CalculationMode = td</tt>). Also, the nuclei should not move.
378 !% The idea is that all orbitals except the last one are frozen. The orbitals are to
379 !% be read from a previous ground-state calculation. The active orbital is then treated
380 !% as independent (whether it contains one electron or two) -- although it will
381 !% feel the Hartree and exchange-correlation potentials from the ground-state electronic
382 !% configuration.
383 !%
384 !% It is almost equivalent to setting <tt>TDFreezeOrbitals = N-1</tt>, where <tt>N</tt> is the number
385 !% of orbitals, but not completely.
386 !%End
387 call parse_variable(namespace, 'TDFreezeOrbitals', 0, td%freeze_orbitals)
388
389 if (td%freeze_orbitals /= 0) then
390 call messages_experimental('TDFreezeOrbitals', namespace=namespace)
391
392 if (hm%lda_u_level /= dft_u_none) then
393 call messages_not_implemented('TDFreezeOrbitals with DFT+U', namespace=namespace)
394 end if
395 end if
396
397 !%Variable TDDMPropagation
398 !%Type integer
399 !%Default 1
400 !%Section Time-Dependent
401 !%Description
402 !% Decides whether to propagate the density matrix, in the spectral basis,
403 !% together with the KS orbitals.
404 !% The resulting time-depented orbilas are to be considered the natural orbitals
405 !% of the denstiy matrix and with the corresponding eigenvalues the occupation numbers.
406 !% In order to work properly this method requires an number of empty states
407 !% specified with the option "ExtraStates".
408 !% This is an experimental feature.
409 !%Option no_propagation 1
410 !% No density time propagation.
411 !%Option master_equation 2
412 !% Density time propagation with the Lindblad master equation.
413 !%Option collision_integral 3
414 !% Density time propagation with the collosion integral.
415 !%End
416
417 call parse_variable(global_namespace, 'TDDMPropagation', option__tddmpropagation__no_propagation, dmp%calculation_mode)
418 if (.not. varinfo_valid_option('TDDMPropagation', dmp%calculation_mode)) then
419 call messages_input_error(global_namespace, 'TDDMPropagation')
420 endif
421 call messages_print_var_option('TDDMPropagation', dmp%calculation_mode, namespace=namespace)
422
423 if (dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
424 call messages_experimental('TDDMPropagation', namespace=namespace)
425 call dmp%init(namespace, st, space, hm)
426 end if
427
428 pop_sub(td_init)
429 nullify(lasers)
430
431 end subroutine td_init
432
433 ! ---------------------------------------------------------
434 subroutine td_init_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, dmp, from_scratch)
435 type(td_t), intent(inout) :: td
436 type(namespace_t), intent(in) :: namespace
437 type(multicomm_t), intent(inout) :: mc
438 type(grid_t), intent(inout) :: gr
439 type(ions_t), intent(inout) :: ions
440 type(states_elec_t), intent(inout) :: st
441 type(v_ks_t), intent(inout) :: ks
442 type(hamiltonian_elec_t), intent(inout) :: hm
443 type(partner_list_t), intent(in) :: ext_partners
444 type(output_t), intent(inout) :: outp
445 type(electron_space_t), intent(in) :: space
446 type(dmp_t), intent(inout) :: dmp
447 logical, intent(inout) :: from_scratch
448 push_sub(td_init_run)
449
450 ! NOTE: please do not change code in this function, but only in functions
451 ! called from here because the logic of this function is replicated in the
452 ! multisystem framework in different places
453
454 call td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
455 call td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
456
457 call ions%initialize()
458
459 td%from_scratch = from_scratch
460
461 if (.not. td%from_scratch) then
462 call td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, td%from_scratch)
463 if (td%from_scratch) then
464 message(1) = "Unable to read time-dependent restart information: Starting from scratch"
465 call messages_warning(1, namespace=namespace)
466 end if
467 end if
468
469 if (td%iter >= td%max_iter) then
470 message(1) = "All requested iterations have already been done. Use FromScratch = yes if you want to redo them."
471 call messages_info(1, namespace=namespace)
473 td%iter = td%iter + 1
474 if ((td%ions_dyn%is_active()) .and. td%recalculate_gs) call td%restart_load%end()
475 pop_sub(td_init_run)
476 return
477 end if
478
479 if (td%from_scratch) then
480 call td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
481 end if
482
483 call td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, dmp, td%from_scratch)
484
485 pop_sub(td_init_run)
486 end subroutine td_init_run
487
488 ! ---------------------------------------------------------
489 subroutine td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
490 type(td_t), intent(inout) :: td
491 type(namespace_t), intent(in) :: namespace
492 type(multicomm_t), intent(inout) :: mc
493 type(grid_t), intent(inout) :: gr
494 type(ions_t), intent(inout) :: ions
495 type(states_elec_t), intent(inout) :: st
496 type(hamiltonian_elec_t), intent(inout) :: hm
497 class(space_t), intent(in) :: space
498
500
501 ! Allocate wavefunctions during time-propagation
502 if (td%dynamics == ehrenfest) then
503 !Note: this is not really clean to do this
504 if (hm%lda_u_level /= dft_u_none .and. states_are_real(st)) then
505 call lda_u_end(hm%lda_u)
506 !complex wfs are required for Ehrenfest
507 call states_elec_allocate_wfns(st, gr, type_cmplx, packed=.true.)
508 call lda_u_init(hm%lda_u, namespace, space, hm%lda_u_level, gr, ions, st, mc, hm%kpoints)
509 else
510 !complex wfs are required for Ehrenfest
511 call states_elec_allocate_wfns(st, gr, type_cmplx, packed=.true.)
512 end if
513 else
514 call states_elec_allocate_wfns(st, gr, packed=.true.)
515 call scf_init(td%scf, namespace, gr, ions, st, mc, hm, space)
516 call scf_set_lower_bound_is_known(td%scf, .true.)
517 end if
518
520 end subroutine td_allocate_wavefunctions
521
522 ! ---------------------------------------------------------
523 subroutine td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
524 type(td_t), intent(inout) :: td
525 type(namespace_t), intent(in) :: namespace
526 type(grid_t), intent(inout) :: gr
527 type(states_elec_t), intent(inout) :: st
528 type(v_ks_t), intent(inout) :: ks
529 type(hamiltonian_elec_t), intent(inout) :: hm
530 type(partner_list_t), intent(in) :: ext_partners
531 class(space_t), intent(in) :: space
532
533 type(gauge_field_t), pointer :: gfield
534
535 push_sub(td_init_gaugefield)
536
537 gfield => list_get_gauge_field(ext_partners)
538 if(associated(gfield)) then
539 if (gauge_field_is_used(gfield)) then
540 !if the gauge field is applied, we need to tell v_ks to calculate the current
541 call v_ks_calculate_current(ks, .true.)
542
543 ! initialize the vector field and update the hamiltonian
544 call gauge_field_init_vec_pot(gfield, st%qtot)
545 call hm%update(gr, namespace, space, ext_partners, time = td%dt*td%iter)
546
547 end if
548 end if
549
550 pop_sub(td_init_gaugefield)
551 end subroutine td_init_gaugefield
552
553 ! ---------------------------------------------------------
554 subroutine td_end(td)
555 type(td_t), intent(inout) :: td
556
557 push_sub(td_end)
558
559 call pes_end(td%pesv)
560 call propagator_elec_end(td%tr) ! clean the evolution method
561 call ion_dynamics_end(td%ions_dyn)
562
563 if (td%dynamics == bo) call scf_end(td%scf)
564 pop_sub(td_end)
565
566 end subroutine td_end
567
568 ! ---------------------------------------------------------
569 subroutine td_end_run(td, st, hm, dmp)
570 type(td_t), intent(inout) :: td
571 type(states_elec_t), intent(inout) :: st
572 type(hamiltonian_elec_t), intent(inout) :: hm
573 type(dmp_t), intent(inout) :: dmp
574
575 push_sub(td_end_run)
576
577 if (st%pack_states .and. hm%apply_packed()) call st%unpack()
578
579 call td%restart_dump%end()
580 call td_write_end(td%write_handler)
581
582 ! free memory
584 if ((td%ions_dyn%is_active()).and. td%recalculate_gs) call td%restart_load%end()
585
586 if (dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
587 call dm_end_run(st%system_grp, dmp)
588 endif
589
590 pop_sub(td_end_run)
591 end subroutine td_end_run
592
593 ! ---------------------------------------------------------
594 subroutine td_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, dmp, from_scratch)
595 type(td_t), intent(inout) :: td
596 type(namespace_t), intent(in) :: namespace
597 type(multicomm_t), intent(inout) :: mc
598 type(grid_t), intent(inout) :: gr
599 type(ions_t), intent(inout) :: ions
600 type(states_elec_t), intent(inout) :: st
601 type(v_ks_t), intent(inout) :: ks
602 type(hamiltonian_elec_t), intent(inout) :: hm
603 type(partner_list_t), intent(in) :: ext_partners
604 type(output_t), intent(inout) :: outp
605 type(electron_space_t), intent(in) :: space
606 type(dmp_t), intent(inout) :: dmp
607 logical, intent(inout) :: from_scratch
608
609 logical :: stopping
610 integer :: iter, scsteps
611 real(real64) :: etime
612 real(real64) :: wall_time, simulation_time, speed_fs_per_day
613 character(len=20) :: fmt
614
615 push_sub(td_run)
616
617 etime = loct_clock()
618 ! This is the time-propagation loop. It starts at t=0 and finishes
619 ! at td%max_iter*dt. The index i runs from 1 to td%max_iter, and
620 ! step "iter" means propagation from (iter-1)*dt to iter*dt.
621 propagation: do iter = td%iter, td%max_iter
622
623 stopping = clean_stop(mc%master_comm) .or. walltimer_alarm(mc%master_comm)
624
625 call profiling_in("TIME_STEP")
626
627 if (iter > 1) then
628 if (((iter-1)*td%dt <= hm%kick%time) .and. (iter*td%dt > hm%kick%time)) then
629 if (.not. hm%pcm%localf) then
630 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints)
631 else
632 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints, pcm = hm%pcm)
633 end if
634 call td_write_kick(outp, namespace, space, gr, hm%kick, ions, iter)
635 !We activate the sprial BC only after the kick,
636 !to be sure that the first iteration corresponds to the ground state
637 if (gr%der%boundaries%spiralBC) gr%der%boundaries%spiral = .true.
638 end if
639 end if
640
641 ! time iterate the system, one time step.
642 select case (td%dynamics)
643 case (ehrenfest)
644 call propagator_elec_dt(ks, namespace, space, hm, gr, st, td%tr, iter*td%dt, td%dt, iter, td%ions_dyn, &
645 ions, ext_partners, mc, outp, td%write_handler, scsteps = scsteps, &
646 update_energy = (mod(iter, td%energy_update_iter) == 0) .or. (iter == td%max_iter))
647 case (bo)
648 call propagator_elec_dt_bo(td%scf, namespace, space, gr, ks, st, hm, ions, ext_partners, mc, iter, &
649 td%dt, td%ions_dyn, scsteps)
650 end select
651
652 if (dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
653 call dm_propagation_run(dmp, namespace, space, gr, ions, st, mc, hm, ks, iter, td%dt, ext_partners, &
654 update_energy = (mod(iter, td%energy_update_iter) == 0) .or. (iter == td%max_iter))
655 end if
656
657 !Apply mask absorbing boundaries
658 if (hm%abs_boundaries%abtype == mask_absorbing) then
659 if (states_are_real(st)) then
660 call dvmask(gr, hm, st)
661 else
662 call zvmask(gr, hm, st)
663 end if
664 end if
665
666 !Photoelectron stuff
667 if (td%pesv%calc_spm .or. td%pesv%calc_mask .or. td%pesv%calc_flux) then
668 call pes_calc(td%pesv, namespace, space, gr, st, td%dt, iter, gr%der, hm%kpoints, ext_partners, stopping)
669 end if
670
671 call td_write_iter(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
672 hm%kick, ks, td%dt, iter, mc, td%recalculate_gs, dmp%adiabatic_st)
673
674 ! write down data
675 call td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, &
676 iter, scsteps, etime, dmp, stopping, from_scratch)
677
678 ! check if debug mode should be enabled or disabled on the fly
679 call io_debug_on_the_fly(namespace)
680
681 call profiling_out("TIME_STEP")
682 if (stopping) exit
683
684 end do propagation
685
686 ! Output propagation speed in fs/day
687 wall_time = loct_clock() - walltimer_get_start_time()
688 simulation_time = td%dt * (iter - td%iter + 1)
689 speed_fs_per_day = units_from_atomic(unit_femtosecond, simulation_time) / (wall_time / 86400.0_real64)
690 if (speed_fs_per_day > 1e4 .or. speed_fs_per_day < 1e-3) then
691 fmt = "ES10.3"
692 else
693 fmt = "F10.3"
694 end if
695 write(message(1), '(a,'//trim(fmt)//',a)') 'Propagation speed: ', speed_fs_per_day, ' fs/day'
696 write(message(2), '(a)') ''
697 call messages_info(2, namespace=namespace)
698
699 pop_sub(td_run)
700 end subroutine td_run
701
702 subroutine td_print_header(namespace)
703 type(namespace_t), intent(in) :: namespace
704
705 push_sub(td_print_header)
706
707 write(message(1), '(a7,1x,a14,a14,a10,a17)') 'Iter ', 'Time ', 'Energy ', 'SC Steps', 'Elapsed Time '
708
709 call messages_info(1, namespace=namespace)
710 call messages_print_with_emphasis(namespace=namespace)
711
712 pop_sub(td_print_header)
713 end subroutine td_print_header
714
715 ! ---------------------------------------------------------
716 subroutine td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, &
717 iter, scsteps, etime, dmp, stopping, from_scratch)
718 type(td_t), intent(inout) :: td
719 type(namespace_t), intent(in) :: namespace
720 type(multicomm_t), intent(in) :: mc
721 type(grid_t), intent(inout) :: gr
722 type(ions_t), intent(inout) :: ions
723 type(states_elec_t), intent(inout) :: st
724 type(v_ks_t), intent(inout) :: ks
725 type(hamiltonian_elec_t), intent(inout) :: hm
726 type(partner_list_t), intent(in) :: ext_partners
727 type(output_t), intent(in) :: outp
728 type(electron_space_t), intent(in) :: space
729 integer, intent(in) :: iter
730 integer, intent(in) :: scsteps
731 real(real64), intent(inout) :: etime
732 type(dmp_t), intent(inout) :: dmp
733 logical, intent(in) :: stopping
734 logical, intent(inout) :: from_scratch
735
736 integer :: ierr
737
738 push_sub(td_check_point)
739
740 call td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
741
742 if (outp%anything_now(iter)) then ! output
743 call td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, td%dt)
744 end if
745
746 if (restart_walltime_period_alarm(mc%master_comm) .or. iter == td%max_iter .or. stopping) then ! restart
747 !if (iter == td%max_iter) outp%iter = ii - 1
748 call td_write_data(td%write_handler)
749 call td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
750 if (ierr /= 0) then
751 message(1) = "Unable to write time-dependent restart information."
752 call messages_warning(1, namespace=namespace)
753 end if
754 ! dm_propagation adiabatic states
755 if (dmp%calculation_mode /= option__tddmpropagation__no_propagation .and. &
756 dmp%basis == option__tddmpropagationbasis__adiabatic) then
757 call states_elec_dump(dmp%restart_dump, space, dmp%adiabatic_st, gr, hm%kpoints, ierr, iter=iter)
758 if (ierr /= 0) then
759 message(1) = "Unable to write dm-adiabatic restart information."
760 call messages_warning(1, namespace=namespace)
761 end if
762 end if
763
764 call pes_output(td%pesv, namespace, space, gr, st, iter, outp, td%dt, ions)
765 end if
766
767 if (mod(iter, td%recalculate_gs_interval) == 0 .or. iter == td%max_iter .or. stopping) then ! potentially recalculate GS
768 if ((td%ions_dyn%is_active()) .and. td%recalculate_gs) then
769 ! write TD restart data
770 call td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
771 if (ierr /= 0) then
772 message(1) = "Unable to write time-dependent restart information."
773 call messages_warning(1, namespace=namespace)
774 end if
775 ! run GS
776 call messages_print_with_emphasis(msg='Recalculating the ground state.', namespace=namespace)
777 from_scratch = .false.
779 call electrons_ground_state_run(namespace, mc, gr, ions, ext_partners, st, ks, hm, outp, space, from_scratch)
780 call states_elec_allocate_wfns(st, gr, packed=.true.)
781 ! load TD restart data again
782 call td_load(td%restart_load, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
783 if (ierr /= 0) then
784 message(1) = "Unable to load TD states."
785 call messages_fatal(1, namespace=namespace)
786 end if
787 call density_calc(st, gr, st%rho)
788 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
789 calc_eigenval=.true., time = iter*td%dt, calc_energy=.true.)
790 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = iter*td%dt, dt = td%dt)
791 assert(.not. td%ions_dyn%cell_relax())
792 call messages_print_with_emphasis(msg="Time-dependent simulation proceeds", namespace=namespace)
793 call td_print_header(namespace)
794 end if
795 end if
796
798 end subroutine td_check_point
799
800 ! ---------------------------------------------------------
801 subroutine td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
802 type(td_t), intent(inout) :: td
803 type(namespace_t), intent(in) :: namespace
804 type(ions_t), intent(inout) :: ions
805 type(hamiltonian_elec_t), intent(inout) :: hm
806 integer, intent(in) :: iter
807 integer, intent(in) :: scsteps
808 real(real64), intent(inout) :: etime
809
810 push_sub(td_print_message)
812 write(message(1), '(i7,1x,2f14.6,i10,f14.3)') iter, units_from_atomic(units_out%time, iter*td%dt), &
813 units_from_atomic(units_out%energy, hm%energy%total + ions%kinetic_energy), scsteps, &
814 loct_clock() - etime
815 call messages_info(1, namespace=namespace)
816 call td_update_elapsed_time(etime)
817
818 pop_sub(td_print_message)
819 end subroutine td_print_message
820
821 ! ---------------------------------------------------------
822 subroutine td_update_elapsed_time(etime)
823 real(real64), intent(inout) :: etime
824
825 push_sub(td_update_elapsed_time)
826
827 etime = loct_clock()
828
830 end subroutine td_update_elapsed_time
831
832 ! ---------------------------------------------------------
833 subroutine td_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, dmp, from_scratch)
834 type(td_t), intent(inout) :: td
835 type(namespace_t), intent(in) :: namespace
836 type(electron_space_t), intent(in) :: space
837 type(multicomm_t), intent(in) :: mc
838 type(grid_t), intent(inout) :: gr
839 type(ions_t), intent(inout) :: ions
840 type(partner_list_t), intent(in) :: ext_partners
841 type(states_elec_t), target, intent(inout) :: st
842 type(v_ks_t), intent(inout) :: ks
843 type(hamiltonian_elec_t), intent(inout) :: hm
844 type(output_t), intent(inout) :: outp
845 type(dmp_t), intent(inout) :: dmp
846 logical, intent(in) :: from_scratch
847
848 integer :: ierr
849 real(real64) :: x
850 real(real64) :: ndinitial(space%dim)
851 logical :: freeze_hxc, freeze_occ, freeze_u
852 type(restart_t) :: restart, restart_frozen
853 type(gauge_field_t), pointer :: gfield
854 type(lasers_t), pointer :: lasers
856
857 !We activate the sprial BC only after the kick,
858 !to be sure that the first iteration corresponds to the ground state
859 if (gr%der%boundaries%spiralBC) then
860 if ((td%iter-1)*td%dt > hm%kick%time) then
861 gr%der%boundaries%spiral = .true.
862 end if
863 hm%vnl%spin => st%spin
864 hm%phase%spin => st%spin
865 !We fill st%spin. In case of restart, we read it in td_load
866 if (from_scratch) call states_elec_fermi(st, namespace, gr)
867 end if
868
869 if (from_scratch) then
870 ! Initialize the occupation matrices and U for DFT+U
871 ! This must be called before parsing TDFreezeOccupations and TDFreezeU
872 ! in order that the code does properly the initialization.
873 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
874 end if
875
876 if (td%freeze_orbitals > 0) then
877 if (from_scratch) then
878 ! In this case, we first freeze the orbitals, then calculate the Hxc potential.
879 call states_elec_freeze_orbitals(st, namespace, space, gr, mc, hm%kpoints, &
880 td%freeze_orbitals, family_is_mgga(ks%xc_family))
881 else
882 call restart%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
883 if (ierr == 0) then
884 call td_load_frozen(namespace, restart, space, gr, st, hm, ierr)
885 end if
886 if (ierr /= 0) then
887 td%iter = 0
888 message(1) = "Unable to read frozen restart information."
889 call messages_fatal(1, namespace=namespace)
890 end if
891 call restart%end()
892 end if
893 write(message(1),'(a,i4,a,i4,a)') 'Info: The lowest', td%freeze_orbitals, &
894 ' orbitals have been frozen.', st%nst, ' will be propagated.'
895 call messages_info(1, namespace=namespace)
897 call density_calc(st, gr, st%rho)
898 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
899 else if (td%freeze_orbitals < 0) then
900 ! This means SAE approximation. We calculate the Hxc first, then freeze all
901 ! orbitals minus one.
902 write(message(1),'(a)') 'Info: The single-active-electron approximation will be used.'
903 call messages_info(1, namespace=namespace)
904 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
905 if (from_scratch) then
906 call states_elec_freeze_orbitals(st, namespace, space, gr, mc, hm%kpoints, st%nst-1, family_is_mgga(ks%xc_family))
907 else
908 call messages_not_implemented("TDFreezeOrbials < 0 with FromScratch=no", namespace=namespace)
909 end if
911 call v_ks_freeze_hxc(ks)
912 call density_calc(st, gr, st%rho)
913 else
914 ! Normal run.
915 call density_calc(st, gr, st%rho)
916 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
917 end if
918
919 !%Variable TDFreezeHXC
920 !%Type logical
921 !%Default no
922 !%Section Time-Dependent
923 !%Description
924 !% The electrons are evolved as independent particles feeling the Hartree and
925 !% exchange-correlation potentials from the ground-state electronic configuration.
926 !%End
927 call parse_variable(namespace, 'TDFreezeHXC', .false., freeze_hxc)
928 if (freeze_hxc) then
929 write(message(1),'(a)') 'Info: Freezing Hartree and exchange-correlation potentials.'
930 call messages_info(1, namespace=namespace)
931
932 if (.not. from_scratch) then
933
934 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr, exact=.true.)
935 call states_elec_load(restart_frozen, namespace, space, st, gr, hm%kpoints, fixed_occ=.true., ierr=ierr, label = ": gs")
936 call states_elec_transform(st, namespace, space, restart_frozen, gr, hm%kpoints)
937 call restart_frozen%end()
938
939 call density_calc(st, gr, st%rho)
940 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
941
942 call restart_frozen%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
943 call states_elec_load(restart_frozen, namespace, space, st, gr, hm%kpoints, st%restart_fixed_occ, &
944 ierr=ierr, iter=td%iter, label = ": td")
945 call restart_frozen%end()
946 call hm%ks_pot%run_zero_iter(td%tr%vks_old)
947
948 end if
949
950 call v_ks_freeze_hxc(ks)
951
952 end if
953
954 x = minval(st%eigenval(st%st_start, :))
955 if (st%parallel_in_states) then
956 call st%mpi_grp%bcast(x, 1, mpi_double_precision, 0)
957 end if
958 call hm%update_span(gr%spacing(1:space%dim), x, namespace)
959 ! initialize Fermi energy
960 call states_elec_fermi(st, namespace, gr, compute_spin = .not. gr%der%boundaries%spiralBC)
961 call energy_calc_total(namespace, space, hm, gr, st, ext_partners)
962
963 !%Variable TDFreezeDFTUOccupations
964 !%Type logical
965 !%Default no
966 !%Section Time-Dependent
967 !%Description
968 !% The occupation matrices than enters in the DFT+U potential
969 !% are not evolved during the time evolution.
970 !%End
971 call parse_variable(namespace, 'TDFreezeDFTUOccupations', .false., freeze_occ)
972 if (freeze_occ) then
973 write(message(1),'(a)') 'Info: Freezing DFT+U occupation matrices that enters in the DFT+U potential.'
974 call messages_info(1, namespace=namespace)
975 call lda_u_freeze_occ(hm%lda_u)
976
977 !In this case we should reload GS wavefunctions
978 if (hm%lda_u_level /= dft_u_none .and. .not. from_scratch) then
979 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr)
980 call lda_u_load(restart_frozen, hm%lda_u, st, hm%energy%dft_u, ierr, occ_only = .true.)
981 call restart_frozen%end()
982 end if
983 end if
984
985 !%Variable TDFreezeU
986 !%Type logical
987 !%Default no
988 !%Section Time-Dependent
989 !%Description
990 !% The effective U of DFT+U is not evolved during the time evolution.
991 !%End
992 call parse_variable(namespace, 'TDFreezeU', .false., freeze_u)
993 if (freeze_u) then
994 write(message(1),'(a)') 'Info: Freezing the effective U of DFT+U.'
995 call messages_info(1, namespace=namespace)
996 call lda_u_freeze_u(hm%lda_u)
997
998 !In this case we should reload GS wavefunctions
999 if (hm%lda_u_level == dft_u_acbn0 .and. .not. from_scratch) then
1000 call restart_frozen%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr)
1001 call lda_u_load(restart_frozen, hm%lda_u, st, hm%energy%dft_u, ierr, u_only = .true.)
1002 call restart_frozen%end()
1003 write(message(1),'(a)') 'Loaded GS effective U of DFT+U'
1004 call messages_info(1, namespace=namespace)
1005 call lda_u_write_u(hm%lda_u, namespace=namespace)
1006 call lda_u_write_v(hm%lda_u, namespace=namespace)
1007 end if
1008 end if
1009
1010 ! This needs to be called before the calculation of the forces,
1011 ! as we need to test of we output the forces or not
1012 call td_write_init(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
1013 ks, td%ions_dyn%is_active(), &
1014 list_has_gauge_field(ext_partners), hm%kick, td%iter, td%max_iter, td%dt, mc, dmp)
1015
1016 ! Resets the nondipole integration after laser-file has been written.
1017 lasers => list_get_lasers(ext_partners)
1018 if(associated(lasers)) then
1019 if (lasers_with_nondipole_field(lasers)) then
1020 ndinitial(1:space%dim)=m_zero
1021 call lasers_set_nondipole_parameters(lasers,ndinitial,m_zero)
1022 end if
1023 end if
1024 nullify(lasers)
1025
1026 call td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
1027
1028 if (td%scissor > m_epsilon) then
1029 call scissor_init(hm%scissor, namespace, space, st, gr, hm%kpoints, hm%phase, td%scissor, mc)
1030 end if
1031
1032 ! Density matrix-based dissipation. initialize ad_st before td_run_zero_iter
1033 if (dmp%calculation_mode /= option__tddmpropagation__no_propagation) then
1034 call dm_propagation_init_run(dmp, namespace, space, gr, ions, st, hm, mc, td%from_scratch)
1035 end if
1036
1037 if (td%iter == 0) then
1038 call td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc, dmp%adiabatic_st)
1039 end if
1040
1041 gfield => list_get_gauge_field(ext_partners)
1042 if(associated(gfield)) then
1043 if (gauge_field_is_propagated(gfield)) then
1044 call gauge_field_get_force(gfield, gr, st%d%spin_channels, st%current, ks%xc%lrc)
1045 if (abs(ks%xc%lrc%alpha) > m_epsilon) then
1046 call messages_experimental('TD-LRC kernel')
1047 end if
1048 endif
1049 end if
1050
1051 !call td_check_trotter(td, sys, h)
1052 td%iter = td%iter + 1
1053
1054 call td%restart_dump%init(namespace, restart_td, restart_type_dump, mc, ierr, mesh=gr)
1055 if (td%ions_dyn%is_active() .and. td%recalculate_gs) then
1056 ! We will also use the TD restart directory as temporary storage during the time propagation
1057 call td%restart_load%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
1058 end if
1059
1060 call messages_print_with_emphasis(msg="Time-Dependent Simulation", namespace=namespace)
1061 call td_print_header(namespace)
1062
1063 if ((td%pesv%calc_spm .or. td%pesv%calc_mask) .and. from_scratch) then
1064 call pes_init_write(td%pesv,gr,st, namespace)
1065 end if
1066
1067 if (st%pack_states .and. hm%apply_packed()) call st%pack()
1068
1070 end subroutine td_init_with_wavefunctions
1071
1072 ! ---------------------------------------------------------
1073 subroutine td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
1074 type(td_t), intent(inout) :: td
1075 type(namespace_t), intent(in) :: namespace
1076 type(electron_space_t), intent(in) :: space
1077 type(grid_t), intent(inout) :: gr
1078 type(ions_t), intent(inout) :: ions
1079 type(partner_list_t), intent(in) :: ext_partners
1080 type(states_elec_t), target, intent(inout) :: st
1081 type(v_ks_t), intent(inout) :: ks
1082 type(hamiltonian_elec_t), intent(inout) :: hm
1083 type(output_t), intent(inout) :: outp
1084
1085 push_sub(td_init_ions_and_forces)
1086
1087 ! Calculate initial forces and kinetic energy
1088 if (td%ions_dyn%ions_move()) then
1089 if (td%iter > 0) then
1090 call td_read_coordinates(td, namespace, ions)
1091 if (ion_dynamics_drive_ions(td%ions_dyn)) then
1092 call ion_dynamics_propagate(td%ions_dyn, ions, td%iter*td%dt, td%dt, namespace)
1093 end if
1094 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = td%iter*td%dt)
1095 ! recompute potential because the ions have moved
1096 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval=.true., time = td%iter*td%dt)
1097 end if
1098
1099 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = td%iter*td%dt, dt = td%dt)
1100
1101 call ions%update_kinetic_energy()
1102 else
1103 if (outp%what(option__output__forces) .or. td%write_handler%out(out_separate_forces)%write) then
1104 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = td%iter*td%dt, dt = td%dt)
1105 end if
1106 end if
1107
1108 if (outp%what(option__output__stress) .or. td%ions_dyn%cell_relax()) then
1109 call stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
1110 if (td%ions_dyn%cell_relax()) then
1111 call td%ions_dyn%update_stress(ions%space, st%stress_tensors%total, ions%latt%rlattice, ions%latt%rcell_volume)
1112 end if
1113 end if
1114
1116 end subroutine td_init_ions_and_forces
1117
1118 ! ---------------------------------------------------------
1119 subroutine td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
1120 type(td_t), intent(inout) :: td
1121 type(namespace_t), intent(in) :: namespace
1122 class(space_t), intent(in) :: space
1123 type(multicomm_t), intent(in) :: mc
1124 type(grid_t), intent(inout) :: gr
1125 type(partner_list_t), intent(in) :: ext_partners
1126 type(states_elec_t), target, intent(inout) :: st
1127 type(v_ks_t), intent(inout) :: ks
1128 type(hamiltonian_elec_t), intent(inout) :: hm
1129 logical, intent(inout) :: from_scratch
1130
1131 integer :: ierr
1132 type(restart_t) :: restart
1133
1134 push_sub(td_load_restart_from_td)
1135
1136 !We redistribute the states before the restarting
1137 if (td%freeze_orbitals > 0) then
1138 call states_elec_freeze_redistribute_states(st, namespace, gr, mc, td%freeze_orbitals)
1139 end if
1140
1141 call restart%init(namespace, restart_td, restart_type_load, mc, ierr, mesh=gr)
1142 if (ierr == 0) then
1143 call td_load(restart, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
1144 end if
1145 call restart%end()
1146 if (ierr /= 0) then
1147 from_scratch = .true.
1148 td%iter = 0
1149 end if
1150
1152 end subroutine td_load_restart_from_td
1153
1154 ! ---------------------------------------------------------
1155 subroutine td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
1156 type(td_t), intent(inout) :: td
1157 type(namespace_t), intent(in) :: namespace
1158 class(space_t), intent(in) :: space
1159 type(multicomm_t), intent(in) :: mc
1160 type(grid_t), intent(inout) :: gr
1161 type(partner_list_t), intent(in) :: ext_partners
1162 type(states_elec_t), target, intent(inout) :: st
1163 type(v_ks_t), intent(inout) :: ks
1164 type(hamiltonian_elec_t), intent(inout) :: hm
1165
1166 integer :: ierr
1167 type(restart_t) :: restart
1169 push_sub(td_load_restart_from_gs)
1170
1171 call restart%init(namespace, restart_gs, restart_type_load, mc, ierr, mesh=gr, exact=.true.)
1172
1173 if (.not. st%only_userdef_istates) then
1174 if (ierr == 0) then
1175 call states_elec_load(restart, namespace, space, st, gr, hm%kpoints, &
1176 st%restart_fixed_occ, ierr, label = ": gs")
1177 end if
1178 if (ierr /= 0) then
1179 message(1) = 'Unable to read ground-state wavefunctions.'
1180 call messages_fatal(1, namespace=namespace)
1181 end if
1182 end if
1183
1184 ! check if we should deploy user-defined wavefunctions.
1185 ! according to the settings in the input file the routine
1186 ! overwrites orbitals that were read from restart/gs
1187 if (parse_is_defined(namespace, 'UserDefinedStates')) then
1188 call states_elec_read_user_def_orbitals(gr, namespace, space, st)
1189 end if
1190
1191 call states_elec_transform(st, namespace, space, restart, gr, hm%kpoints)
1192 call restart%end()
1193
1195 end subroutine td_load_restart_from_gs
1196
1197 ! ---------------------------------------------------------
1198 subroutine td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc, dmp_st)
1199 type(td_t), intent(inout) :: td
1200 type(namespace_t), intent(in) :: namespace
1201 type(electron_space_t), intent(in) :: space
1202 type(grid_t), intent(inout) :: gr
1203 type(ions_t), intent(inout) :: ions
1204 type(states_elec_t), intent(inout) :: st
1205 type(v_ks_t), intent(inout) :: ks
1206 type(hamiltonian_elec_t), intent(inout) :: hm
1207 type(partner_list_t), intent(in) :: ext_partners
1208 type(output_t), intent(in) :: outp
1209 type(multicomm_t), intent(in) :: mc
1210 type(states_elec_t), intent(in) :: dmp_st
1211
1212 push_sub(td_run_zero_iter)
1213
1214 call td_write_iter(td%write_handler, namespace, space, outp, gr, st, hm, ions, ext_partners, &
1215 hm%kick, ks, td%dt, 0, mc, td%recalculate_gs, dmp_st)
1216
1217 ! I apply the delta electric field *after* td_write_iter, otherwise the
1218 ! dipole matrix elements in write_proj are wrong
1219 if (abs(hm%kick%time) <= m_epsilon) then
1220 if (.not. hm%pcm%localf) then
1221 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints)
1222 else
1223 call kick_apply(space, gr, st, td%ions_dyn, ions, hm%kick, hm%psolver, hm%kpoints, pcm = hm%pcm)
1224 end if
1225 call td_write_kick(outp, namespace, space, gr, hm%kick, ions, 0)
1226
1227 !We activate the sprial BC only after the kick
1228 if (gr%der%boundaries%spiralBC) then
1229 gr%der%boundaries%spiral = .true.
1230 end if
1231 end if
1232 call hm%ks_pot%run_zero_iter(td%tr%vks_old)
1233
1234 if (any(outp%output_interval > 0)) then
1235 call td_write_data(td%write_handler)
1236 call td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, 0)
1237 end if
1238
1239 pop_sub(td_run_zero_iter)
1240 end subroutine td_run_zero_iter
1241
1242
1243 ! ---------------------------------------------------------
1245 subroutine td_read_coordinates(td, namespace, ions)
1246 type(td_t), intent(in) :: td
1247 type(namespace_t), intent(in) :: namespace
1248 type(ions_t), intent(inout) :: ions
1249
1250 integer :: iatom, iter, iunit
1251
1252 push_sub(td_read_coordinates)
1253
1254 iunit = io_open('td.general/coordinates', namespace, action='read', status='old', die=.false.)
1255 if (iunit == -1) then
1256 message(1) = "Could not open file '"//trim(io_workpath('td.general/coordinates', namespace))//"'."
1257 message(2) = "Starting simulation from initial geometry."
1258 call messages_warning(2, namespace=namespace)
1259 pop_sub(td_read_coordinates)
1260 return
1261 end if
1262
1263 call io_skip_header(iunit)
1264 do iter = 0, td%iter - 1
1265 read(iunit, *) ! skip previous iterations... sorry, but no portable seek in Fortran
1266 end do
1267 read(iunit, '(32x)', advance='no') ! skip the time index.
1268
1269 do iatom = 1, ions%natoms
1270 read(iunit, '(3es24.16)', advance='no') ions%pos(:, iatom)
1271 ions%pos(:, iatom) = units_to_atomic(units_out%length, ions%pos(:, iatom))
1272 end do
1273 do iatom = 1, ions%natoms
1274 read(iunit, '(3es24.16)', advance='no') ions%vel(:, iatom)
1275 ions%vel(:, iatom) = units_to_atomic(units_out%velocity, ions%vel(:, iatom))
1276 end do
1277 do iatom = 1, ions%natoms
1278 read(iunit, '(3es24.16)', advance='no') ions%tot_force(:, iatom)
1279 ions%tot_force(:, iatom) = units_to_atomic(units_out%force, ions%tot_force(:, iatom))
1280 end do
1281
1282 call io_close(iunit)
1283
1284 pop_sub(td_read_coordinates)
1285 end subroutine td_read_coordinates
1286
1287 ! ---------------------------------------------------------
1288 subroutine td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
1289 type(td_t), intent(in) :: td
1290 type(namespace_t), intent(in) :: namespace
1291 class(space_t), intent(in) :: space
1292 type(grid_t), intent(in) :: gr
1293 type(states_elec_t), intent(in) :: st
1294 type(hamiltonian_elec_t), intent(in) :: hm
1295 type(v_ks_t), intent(in) :: ks
1296 type(partner_list_t), intent(in) :: ext_partners
1297 integer, intent(in) :: iter
1298 integer, intent(out) :: ierr
1299
1300 type(gauge_field_t), pointer :: gfield
1301 integer :: err, err2
1302
1303 push_sub(td_dump)
1304
1305 ierr = 0
1306
1307 if (td%restart_dump%skip()) then
1308 pop_sub(td_dump)
1309 return
1310 end if
1311
1312 message(1) = "Debug: Writing td restart."
1313 call messages_info(1, namespace=namespace, debug_only=.true.)
1314
1315 ! first write resume file
1316 call states_elec_dump(td%restart_dump, space, st, gr, hm%kpoints, err, iter=iter)
1317 if (err /= 0) ierr = ierr + 1
1318
1319 call states_elec_dump_rho(td%restart_dump, st, gr, ierr, iter=iter)
1320 if (err /= 0) ierr = ierr + 1
1321
1322 if (hm%lda_u_level /= dft_u_none) then
1323 call lda_u_dump(td%restart_dump, namespace, hm%lda_u, st, gr, err)
1324 if (err /= 0) ierr = ierr + 1
1325 end if
1326
1327 call potential_interpolation_dump(td%tr%vks_old, td%restart_dump, gr, st%d%nspin, err2)
1328 if (err2 /= 0) ierr = ierr + 2
1329
1330 call pes_dump(td%pesv, namespace, td%restart_dump, st, gr, err)
1331 if (err /= 0) ierr = ierr + 4
1332
1333 ! Gauge field restart
1334 gfield => list_get_gauge_field(ext_partners)
1335 if(associated(gfield)) then
1336 call gauge_field_dump(td%restart_dump, gfield, err)
1337 if (err /= 0) ierr = ierr + 8
1338 end if
1339
1340 if (gr%der%boundaries%spiralBC) then
1341 call states_elec_dump_spin(td%restart_dump, st, err)
1342 if (err /= 0) ierr = ierr + 8
1343 end if
1344
1345 if (ks%has_photons) then
1346 call mf_photons_dump(td%restart_dump, ks%pt_mx, gr, td%dt, ks%pt, err)
1347 if (err /= 0) ierr = ierr + 16
1348 end if
1349
1350 if (ks%xc_photon /= 0) then
1351 ! photon-free mean field
1352 call ks%xc_photons%mf_dump(td%restart_dump, err)
1353 if (err /= 0) ierr = ierr + 32
1354 end if
1355
1356 if (allocated(st%frozen_rho)) then
1357 call states_elec_dump_frozen(td%restart_dump, space, st, gr, err)
1358 end if
1359 if (err /= 0) ierr = ierr + 64
1360
1361 if (td%ions_dyn%ions_move() .or. td%ions_dyn%cell_relax()) then
1362 call ion_dynamics_dump(td%ions_dyn, td%restart_dump, err)
1363 end if
1364 if (err /= 0) ierr = ierr + 128
1365
1366 message(1) = "Debug: Writing td restart done."
1367 call messages_info(1, namespace=namespace, debug_only=.true.)
1368
1369 pop_sub(td_dump)
1370 end subroutine td_dump
1371
1372 ! ---------------------------------------------------------
1373 subroutine td_load(restart, namespace, space, gr, st, hm, ext_partners, td, ks, ierr)
1374 type(restart_t), intent(in) :: restart
1375 type(namespace_t), intent(in) :: namespace
1376 class(space_t), intent(in) :: space
1377 type(grid_t), intent(in) :: gr
1378 type(states_elec_t), intent(inout) :: st
1379 type(hamiltonian_elec_t), intent(inout) :: hm
1380 type(partner_list_t), intent(in) :: ext_partners
1381 type(td_t), intent(inout) :: td
1382 type(v_ks_t), intent(inout) :: ks
1383 integer, intent(out) :: ierr
1384
1385 integer :: err, err2
1386 type(gauge_field_t), pointer :: gfield
1387 push_sub(td_load)
1388
1389 ierr = 0
1390
1391 if (restart%skip()) then
1392 ierr = -1
1393 pop_sub(td_load)
1394 return
1395 end if
1396
1397 message(1) = "Debug: Reading td restart."
1398 call messages_info(1, namespace=namespace, debug_only=.true.)
1399
1400 ! Read states
1401 call states_elec_load(restart, namespace, space, st, gr, hm%kpoints, st%restart_fixed_occ, &
1402 err, iter=td%iter, label = ": td")
1403 if (err /= 0) then
1404 ierr = ierr + 1
1405 end if
1406
1407 ! read potential from previous interactions
1408 call potential_interpolation_load(td%tr%vks_old, namespace, restart, gr, st%d%nspin, err2)
1409 if (err2 /= 0) ierr = ierr + 2
1410
1411 if (hm%lda_u_level /= dft_u_none) then
1412 call lda_u_load(restart, hm%lda_u, st, hm%energy%dft_u, err)
1413 if (err /= 0) ierr = ierr + 1
1414 end if
1415
1416
1417 ! read PES restart
1418 if (td%pesv%calc_spm .or. td%pesv%calc_mask .or. td%pesv%calc_flux) then
1419 call pes_load(td%pesv, namespace, restart, st, err)
1420 if (err /= 0) ierr = ierr + 4
1421 end if
1422
1423 ! Gauge field restart
1424 gfield => list_get_gauge_field(ext_partners)
1425 if (associated(gfield)) then
1426 call gauge_field_load(restart, gfield, err)
1427 if (err /= 0) then
1428 ierr = ierr + 8
1429 else
1430 call hm%update(gr, namespace, space, ext_partners, time = td%dt*td%iter)
1431 end if
1432 end if
1433
1434 ! add photon restart
1435 if (ks%has_photons) then
1436 call mf_photons_load(restart, ks%pt_mx, gr, err)
1437 if (err /= 0) ierr = ierr + 16
1438 end if
1439
1440 if (ks%xc_photon /= 0) then
1441 call ks%xc_photons%mf_load(restart, space, err)
1442 if (err /= 0) ierr = ierr + 32
1443 end if
1444
1445 if (gr%der%boundaries%spiralBC) then
1446 call states_elec_load_spin(restart, st, err)
1447 !To ensure back compatibility, if the file is not present, we use the
1448 !current states to get the spins
1449 if (err /= 0) call states_elec_fermi(st, namespace, gr)
1450 end if
1451
1452 if (td%ions_dyn%is_active()) then
1453 call ion_dynamics_load(td%ions_dyn, restart, err)
1454 if (err /= 0) ierr = ierr + 64
1455 end if
1456
1457 message(1) = "Debug: Reading td restart done."
1458 call messages_info(1, namespace=namespace, debug_only=.true.)
1459
1460 pop_sub(td_load)
1461 end subroutine td_load
1462 ! ---------------------------------------------------------
1463 subroutine td_load_frozen(namespace, restart, space, mesh, st, hm, ierr)
1464 type(namespace_t), intent(in) :: namespace
1465 type(restart_t), intent(in) :: restart
1466 class(space_t), intent(in) :: space
1467 class(mesh_t), intent(in) :: mesh
1468 type(states_elec_t), intent(inout) :: st
1469 type(hamiltonian_elec_t), intent(inout) :: hm
1470 integer, intent(out) :: ierr
1471
1472 push_sub(td_load_frozen)
1473
1474 ierr = 0
1475
1476 if (restart%skip()) then
1477 ierr = -1
1478 pop_sub(td_load_frozen)
1479 return
1480 end if
1481
1482 message(1) = "Debug: Reading td frozen restart."
1483 call messages_info(1, namespace=namespace, debug_only=.true.)
1484
1485 safe_allocate(st%frozen_rho(1:mesh%np, 1:st%d%nspin))
1486 if (family_is_mgga(hm%xc%family)) then
1487 safe_allocate(st%frozen_tau(1:mesh%np, 1:st%d%nspin))
1488 safe_allocate(st%frozen_gdens(1:mesh%np, 1:space%dim, 1:st%d%nspin))
1489 safe_allocate(st%frozen_ldens(1:mesh%np, 1:st%d%nspin))
1490 end if
1491
1492 call states_elec_load_frozen(restart, space, st, mesh, ierr)
1493
1494 message(1) = "Debug: Reading td frozen restart done."
1495 call messages_info(1, namespace=namespace, debug_only=.true.)
1496
1497 pop_sub(td_load_frozen)
1498 end subroutine td_load_frozen
1499
1500 ! ---------------------------------------------------------
1501 logical function td_get_from_scratch(td)
1502 type(td_t), intent(in) :: td
1503
1504 push_sub(td_get_from_scratch)
1505
1506 td_get_from_scratch = td%from_scratch
1507
1508 pop_sub(td_get_from_scratch)
1509 end function td_get_from_scratch
1510
1511 ! ---------------------------------------------------------
1512 subroutine td_set_from_scratch(td, from_scratch)
1513 type(td_t), intent(inout) :: td
1514 logical, intent(in) :: from_scratch
1515
1516 push_sub(td_set_from_scratch)
1517
1518 td%from_scratch = from_scratch
1519
1520 pop_sub(td_set_from_scratch)
1521 end subroutine td_set_from_scratch
1522end module td_oct_m
1523
1524!! Local Variables:
1525!! mode: f90
1526!! coding: utf-8
1527!! End:
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:853
subroutine, public states_elec_freeze_redistribute_states(st, namespace, mesh, mc, nn)
Definition: density.F90:825
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:694
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:653
subroutine, public dm_propagation_init_run(dmp, namespace, space, gr, ions, st, hm, mc, from_scratch)
Initialise the adiabatic states prior to running TD propagation.
subroutine, public dm_end_run(system_grp, dmp)
subroutine, public dm_propagation_run(dmp, namespace, space, gr, ions, st, mc, hm, ks, iter, dt, ext_partners, update_energy)
Density matrix propagation.
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:169
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:200
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public m_one
Definition: global.F90:201
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:467
subroutine, public io_skip_header(iunit)
Definition: io.F90:646
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
Definition: io.F90:318
subroutine, public io_debug_on_the_fly(namespace)
check if debug mode should be enabled or disabled on the fly
Definition: io.F90:535
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
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:756
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:743
subroutine, public laser_write_info(lasers, namespace, dt, max_iter, iunit)
Definition: lasers.F90:936
subroutine, public lda_u_dump(restart, namespace, this, st, mesh, ierr)
Definition: lda_u_io.F90:642
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:528
subroutine, public lda_u_load(restart, this, st, dftu_energy, ierr, occ_only, u_only)
Definition: lda_u_io.F90:723
subroutine, public lda_u_write_v(this, iunit, namespace)
Definition: lda_u_io.F90:576
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
subroutine, public lda_u_init(this, namespace, space, level, gr, ions, st, mc, kpoints)
Definition: lda_u.F90:287
subroutine, public lda_u_freeze_occ(this)
Definition: lda_u.F90:984
subroutine, public lda_u_freeze_u(this)
Definition: lda_u.F90:991
subroutine, public lda_u_end(this)
Definition: lda_u.F90:660
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:895
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:205
This module implements fully polymorphic linked lists, and some specializations thereof.
System information (time, memory, sysname)
Definition: loct.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1068
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
type(namespace_t), public global_namespace
Definition: namespace.F90:135
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:463
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:631
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:337
integer, parameter, public restart_gs
Definition: restart.F90:156
integer, parameter, public restart_type_dump
Definition: restart.F90:184
integer, parameter, public restart_td
Definition: restart.F90:156
integer, parameter, public restart_type_load
Definition: restart.F90:184
subroutine, public scf_set_lower_bound_is_known(scf, known_lower_bound)
Set the flag lower_bound_is_known.
Definition: scf.F90:1683
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:260
subroutine, public scf_end(scf)
Definition: scf.F90:558
subroutine, public 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_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_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, 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_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, verbose)
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_init_with_wavefunctions(td, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, dmp, from_scratch)
Definition: td.F90:929
subroutine, public td_end(td)
Definition: td.F90:650
subroutine, public td_load_restart_from_gs(td, namespace, space, mc, gr, ext_partners, st, ks, hm)
Definition: td.F90:1251
subroutine td_print_message(td, namespace, ions, hm, iter, scsteps, etime)
Definition: td.F90:897
subroutine td_run_zero_iter(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, mc, dmp_st)
Definition: td.F90:1294
subroutine, public td_allocate_wavefunctions(td, namespace, mc, gr, ions, st, hm, space)
Definition: td.F90:585
subroutine, public td_init(td, namespace, space, gr, ions, st, ks, hm, ext_partners, outp, dmp)
Definition: td.F90:240
logical function, public td_get_from_scratch(td)
Definition: td.F90:1597
subroutine, public td_init_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, dmp, from_scratch)
Definition: td.F90:530
subroutine td_load_frozen(namespace, restart, space, mesh, st, hm, ierr)
Definition: td.F90:1559
subroutine td_print_header(namespace)
Definition: td.F90:798
subroutine, public td_run(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, dmp, from_scratch)
Definition: td.F90:690
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:1469
subroutine, public td_set_from_scratch(td, from_scratch)
Definition: td.F90:1608
subroutine, public td_dump(td, namespace, space, gr, st, hm, ks, ext_partners, iter, ierr)
Definition: td.F90:1384
subroutine td_read_coordinates(td, namespace, ions)
reads the pos and vel from coordinates file
Definition: td.F90:1341
subroutine, public td_init_gaugefield(td, namespace, gr, st, ks, hm, ext_partners, space)
Definition: td.F90:619
subroutine td_update_elapsed_time(etime)
Definition: td.F90:918
subroutine, public td_check_point(td, namespace, mc, gr, ions, st, ks, hm, ext_partners, outp, space, iter, scsteps, etime, dmp, stopping, from_scratch)
Definition: td.F90:813
subroutine, public td_load_restart_from_td(td, namespace, space, mc, gr, ext_partners, st, ks, hm, from_scratch)
Definition: td.F90:1215
subroutine td_init_ions_and_forces(td, namespace, space, gr, ions, ext_partners, st, ks, hm, outp)
Definition: td.F90:1169
subroutine, public td_end_run(td, st, hm, dmp)
Definition: td.F90:665
subroutine, public td_write_output(namespace, space, gr, st, hm, ks, outp, ions, ext_partners, iter, dt)
Definition: td_write.F90:1286
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, dmp)
Initialize files to write when prograting in time.
Definition: td_write.F90:380
subroutine, public td_write_iter(writ, namespace, space, outp, gr, st, hm, ions, ext_partners, kick, ks, dt, iter, mc, recalculate_gs, dmp_st)
Definition: td_write.F90:1069
subroutine, public td_write_data(writ)
Definition: td_write.F90:1252
subroutine, public td_write_kick(outp, namespace, space, mesh, kick, ions, iter)
Definition: td_write.F90:342
subroutine, public td_write_end(writ)
Definition: td_write.F90:1022
integer, parameter, public out_separate_forces
Definition: td_write.F90:204
type(type_t), parameter, 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:1516
subroutine, public v_ks_calculate_current(this, calc_cur)
Definition: v_ks.F90:1527
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:755
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:120
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:703
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:722
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)