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