Octopus
propagator_elec.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
22 use debug_oct_m
27 use forces_oct_m
29 use grid_oct_m
30 use global_oct_m
34 use ions_oct_m
35 use, intrinsic :: iso_fortran_env
36 use lasers_oct_m
37 use lda_u_oct_m
39 use parser_oct_m
57 use scf_oct_m
59 use space_oct_m
61 use stress_oct_m
63 use v_ks_oct_m
65 use xc_oct_m
66
67 implicit none
68
69 private
70 public :: &
79
80contains
81
82 ! ---------------------------------------------------------
83 subroutine propagator_elec_copy(tro, tri)
84 type(propagator_base_t), intent(inout) :: tro
85 type(propagator_base_t), intent(in) :: tri
86
87 push_sub(propagator_elec_copy)
88
89 tro%method = tri%method
90
91 select case (tro%method)
93 tro%tdsk_size = tri%tdsk_size
94 call sparskit_solver_copy(tro%tdsk, tri%tdsk)
95
97 tro%tdsk_size = tri%tdsk_size
98 call sparskit_solver_copy(tro%tdsk, tri%tdsk)
99
100 case (prop_runge_kutta2)
101 tro%tdsk_size = tri%tdsk_size
102 call sparskit_solver_copy(tro%tdsk, tri%tdsk)
103
104 end select
105
106 call potential_interpolation_copy(tro%vks_old, tri%vks_old)
107
108 call exponential_copy(tro%te, tri%te)
109 tro%scf_propagation_steps = tri%scf_propagation_steps
110
111 tro%scf_threshold = tri%scf_threshold
112 pop_sub(propagator_elec_copy)
113 end subroutine propagator_elec_copy
114 ! ---------------------------------------------------------
115
117 ! ---------------------------------------------------------
118 subroutine propagator_elec_init(gr, namespace, st, tr, ks_pot, have_fields, family_is_mgga_with_exc, relax_cell)
119 type(grid_t), intent(in) :: gr
120 type(namespace_t), intent(in) :: namespace
121 type(states_elec_t), intent(in) :: st
122 type(propagator_base_t), intent(inout) :: tr
123 type(ks_potential_t), intent(in) :: ks_pot
127 logical, intent(in) :: have_fields
128 logical, intent(in) :: family_is_mgga_with_exc
129 logical, intent(in) :: relax_cell
130
131 push_sub(propagator_elec_init)
132
133 !%Variable TDPropagator
134 !%Type integer
135 !%Default etrs
136 !%Section Time-Dependent::Propagation
137 !%Description
138 !% This variable determines which algorithm will be used to approximate
139 !% the evolution operator <math>U(t+\delta t, t)</math>. That is, given
140 !% <math>\psi(\tau)</math> and <math>H(\tau)</math> for <math>\tau \le t</math>,
141 !% calculate <math>t+\delta t</math>. Note that in general the Hamiltonian
142 !% is not known at times in the interior of the interval <math>[t,t+\delta t]</math>.
143 !% This is due to the self-consistent nature of the time-dependent Kohn-Sham problem:
144 !% the Hamiltonian at a given time <math>\tau</math> is built from the
145 !% "solution" wavefunctions at that time.
146 !%
147 !% Some methods, however, do require the knowledge of the Hamiltonian at some
148 !% point of the interval <math>[t,t+\delta t]</math>. This problem is solved by making
149 !% use of extrapolation: given a number <math>l</math> of time steps previous to time
150 !% <math>t</math>, this information is used to build the Hamiltonian at arbitrary times
151 !% within <math>[t,t+\delta t]</math>. To be fully precise, one should then proceed
152 !% <i>self-consistently</i>: the obtained Hamiltonian at time <math>t+\delta t</math>
153 !% may then be used to interpolate the Hamiltonian, and repeat the evolution
154 !% algorithm with this new information. Whenever iterating the procedure does
155 !% not change the solution wavefunctions, the cycle is stopped. In practice,
156 !% in <tt>Octopus</tt> we perform a second-order extrapolation without a
157 !% self-consistency check, except for the first two iterations, where obviously
158 !% the extrapolation is not reliable.
159 !%
160 !% The proliferation of methods is certainly excessive. The reason for it is that
161 !% the propagation algorithm is currently a topic of active development. We
162 !% hope that in the future the optimal schemes are clearly identified. In the
163 !% mean time, if you do not feel like testing, use the default choices and
164 !% make sure the time step is small enough.
165 !%
166 !% More details about the propagation methods could be found in
167 !% A. Castro, M.A.L. Marques, and A. Rubio, J. Chem. Phys 121 3425-3433 (2004)
168 !% and
169 !% A. Pueyo et al., J. Chem. Theory Comput. 14, 6, 3040 (2018)
170 !%
171 !%Option etrs 2
172 !% The idea is to make use of time-reversal symmetry from the beginning:
173 !%
174 !% <math>
175 !% \exp \left(-i\delta t H_{n} / 2 \right)\psi_n = \exp \left(i\delta t H_{n+1} / 2 \right)\psi_{n+1},
176 !% </math>
177 !%
178 !% and then invert to obtain:
179 !%
180 !% <math>
181 !% \psi_{n+1} = \exp \left(-i\delta t H_{n+1} / 2 \right) \exp \left(-i\delta t H_{n} / 2 \right)\psi_{n}.
182 !% </math>
183 !%
184 !% But we need to know <math>H_{n+1}</math>, which can only be known exactly through the solution
185 !% <math>\psi_{n+1}</math>. What we do is to estimate it by performing a single exponential:
186 !% <math>\psi\^{\*}\_{n+1}=\exp \left( -i\delta t H_{n} \right) \psi_n</math>, and then
187 !% <math>H_{n+1} = H[\psi\^{\*}_{n+1}]</math>. Thus no extrapolation is performed in this case.
188 !%Option aetrs 3
189 !% Approximated Enforced Time-Reversal Symmetry (AETRS).
190 !% A modification of previous method to make it faster.
191 !% It is based on extrapolation of the time-dependent potentials. It is faster
192 !% by about 40%.
193 !% The only difference is the procedure to estimate <math>H_{n+1}</math>: in this case
194 !% it is extrapolated via a second-order polynomial by making use of the
195 !% Hamiltonian at time <math>t-2\delta t</math>, <math>t-\delta t</math> and <math>t</math>.
196 !%Option caetrs 12
197 !% (experimental) Corrected Approximated Enforced Time-Reversal
198 !% Symmetry (AETRS), this is the previous propagator but including
199 !% a correction step to the exponential.
200 !%Option exp_mid 4
201 !% Exponential Midpoint Rule (EM).
202 !% This is maybe the simplest method, but it is very well grounded theoretically:
203 !% it is unitary (if the exponential is performed correctly) and preserves
204 !% time-reversal symmetry (if the self-consistency problem is dealt with correctly).
205 !% It is defined as:
206 !% <math>
207 !% U_{\rm EM}(t+\delta t, t) = \exp \left( -i\delta t H_{t+\delta t/2}\right)\,.
208 !% </math>
209 !%Option crank_nicholson 5
210 !%Option crank_nicolson 5
211 !% Classical Crank-Nicolson propagator.
212 !% <math>
213 !% (1 + i\delta t H_{n+1/2} / 2) \psi_{n+1} = (1 - i\delta t H_{n+1/2} / 2) \psi_{n}
214 !% </math>
215 !%Option crank_nicholson_sparskit 6
216 !%Option crank_nicolson_sparskit 6
217 !% Classical Crank-Nicolson propagator. Requires the SPARSKIT library.
218 !% <math>
219 !% (1 + i\delta t H_{n+1/2} / 2) \psi_{n+1} = (1 - i\delta t H_{n+1/2} / 2) \psi_{n}
220 !% </math>
221 !%Option magnus 7
222 !% Magnus Expansion (M4).
223 !% This is the most sophisticated approach. It is a fourth-order scheme (a feature
224 !% which it shares with the ST scheme; the other schemes are in principle second-order).
225 !% It is tailored for making use of very large time steps, or equivalently,
226 !% dealing with problem with very high-frequency time-dependence.
227 !% It is still in a experimental state; we are not yet sure of when it is
228 !% advantageous.
229 !%Option qoct_tddft_propagator 10
230 !% WARNING: EXPERIMENTAL
231 !%Option runge_kutta4 13
232 !% WARNING: EXPERIMENTAL. Implicit Gauss-Legendre 4th order Runge-Kutta.
233 !%Option runge_kutta2 14
234 !% WARNING: EXPERIMENTAL. Implicit 2nd order Runge-Kutta (trapezoidal rule).
235 !% Similar, but not identical, to Crank-Nicolson method.
236 !%Option expl_runge_kutta4 15
237 !% WARNING: EXPERIMENTAL. Explicit RK4 method.
238 !%Option cfmagnus4 16
239 !% Commutator-free magnus expansion. This method has been developed for linear
240 !% systems and involves products of exponentials to avoid commutators. To apply it
241 !% to our non-linear equations, we extrapolate the potential from previous timesteps
242 !% to the times needed.The fourth-ordermethod implemented here follows equation 54
243 !% in Gómez Pueyo, A., Marques, M. A. L., Rubio, A. & Castro, A. Propagators for the
244 !% Time-Dependent Kohn–Sham Equations: Multistep, Runge–Kutta, Exponential Runge–Kutta,
245 !% and Commutator Free Magnus Methods. J. Chem. Theory Comput. 14, 3040–3052 (2018)
246 !% which corresponds to equation 43 in Blanes, S. & Moan, P. C. Fourth- and sixth-order
247 !% commutator-free Magnus integrators for linear and non-linear dynamical systems.
248 !% Applied Numerical Mathematics 56, 1519–1537 (2006):
249 !% <math>
250 !% \psi(t+\Delta t) = \exp(-i \Delta t ( \alpha_1 H(t_1) + \alpha_2 H(t_2) ) )
251 !% \exp(-i \Delta t ( \alpha_2 H(t_1) + \alpha_1 H(t_2) ) ) \psi(t),
252 !% </math>
253 !% with $\alpha_{1,2} = \frac{3\mp 2\sqrt{3}}{12}$ and $t_i = t_0 + c_i \Delta t$,
254 !% $c_{1,2} = \frac{1}{2} \mp \frac{\sqrt{3}}{6}$.
255 !% Still in experimental state.
256 !%Option exp_mid_predictor 17
257 !% Exponential Midpoint Rule (EM).
258 !% This is maybe the simplest method, but it is very well grounded theoretically:
259 !% it is unitary (if the exponential is performed correctly) and preserves
260 !% time-reversal symmetry (if the self-consistency problem is dealt with correctly).
261 !% It is defined as:
262 !% <math>
263 !% U_{\rm EM}(t+\delta t, t) = \exp \left( -i\delta t H_{t+\delta t/2}\right)\,.
264 !% </math>
265 !% As opposed to the exp_mid propagator, $H_{t+\delta t/2}$ is not obtained by
266 !% extrapolation from previous timesteps, but by predicting by applying the exponential with
267 !% $\delta t/2$ and computing the corresponding density and potential.
268 !% It is equivalent to an explicit Magnus integration of second order, see equation 20 in
269 !% Casas, F. & Iserles, A. Explicit Magnus expansions for nonlinear equations. J. Phys. A: Math. Gen. 39,
270 !% 5445 (2006) [https://dx.doi.org/10.1088/0305-4470/39/19/S07]
271 !%Option explicit_magnus3 18
272 !% Third-order explicit Magnus expansion.
273 !% This method follows equations 22 in
274 !% Casas, F. & Iserles, A. Explicit Magnus expansions for nonlinear equations.
275 !% J. Phys. A: Math. Gen. 39, 5445 (2006) [https://dx.doi.org/10.1088/0305-4470/39/19/S07]
276 !% Applying eq. 22 to our notation for the Schrödinger equation yields:
277 !% <math>
278 !% \psi_0 = \psi(t) \\
279 !% H_1 = H(t, \psi_0) \\
280 !% \psi_2 = \exp(-i \frac{\Delta t}{2} H_1) \psi_0 \\
281 !% H_2 = H(t+\frac{\Delta t}{2}, \psi_2) \\
282 !% \psi_3 = \exp(-i \frac{\Delta t}{4} (H_1 + H_2)) \psi_0 \\
283 !% H_3 = H(t+\frac{\Delta t}{2}, \psi_3) \\
284 !% \psi_4 = \exp(-i \Delta t H_2) \psi_0 \\
285 !% H_4 = H(t+\Delta t, \psi_4) \\
286 !% \psi(t+\Delta t) = \exp(-i \Delta t ( \frac{1}{6}(H_1 + 4 H_3 + H_4)
287 !% + \frac{i \Delta t}{12} ( [H_1 + H_2, H_3] + [H_2, H_4]))) \psi_0
288 !% </math>
289 !% Still in experimental state.
290 !%End
291 call messages_obsolete_variable(namespace, 'TDEvolutionMethod', 'TDPropagator')
292
293 call parse_variable(namespace, 'TDPropagator', prop_etrs, tr%method)
294 if (.not. varinfo_valid_option('TDPropagator', tr%method)) call messages_input_error(namespace, 'TDPropagator')
295
296 select case (tr%method)
297 case (prop_etrs)
298 case (prop_aetrs)
299 case (prop_caetrs)
300 call messages_experimental("CAETRS propagator", namespace=namespace)
304 ! set up pointer for zmf_dotu_aux, zmf_nrm2_aux
305 call mesh_init_mesh_aux(gr)
306 case (prop_runge_kutta4)
307 ! set up pointer for zmf_dotu_aux, zmf_nrm2_aux
308 call mesh_init_mesh_aux(gr)
310 tr%tdsk_size = 2 * st%d%dim * gr%np * (st%st_end - st%st_start + 1) * (st%d%kpt%end - st%d%kpt%start + 1)
311 call sparskit_solver_init(namespace, tr%tdsk_size, tr%tdsk, .true.)
312
313#ifndef HAVE_SPARSKIT
314 message(1) = 'Octopus was not compiled with support for the SPARSKIT library. This'
315 message(2) = 'library is required if the "runge_kutta4" propagator is selected.'
316 message(3) = 'Try using a different propagation scheme or recompile with SPARSKIT support.'
317 call messages_fatal(3, namespace=namespace)
318#endif
319
320 call messages_experimental("Runge-Kutta 4 propagator", namespace=namespace)
321 case (prop_runge_kutta2)
322
323 ! set up pointer for zmf_dotu_aux, zmf_nrm2_aux
324 call mesh_init_mesh_aux(gr)
326 tr%tdsk_size = st%d%dim * gr%np * (st%st_end - st%st_start + 1) * (st%d%kpt%end - st%d%kpt%start + 1)
327 call sparskit_solver_init(namespace, tr%tdsk_size, tr%tdsk, .true.)
328
329#ifndef HAVE_SPARSKIT
330 message(1) = 'Octopus was not compiled with support for the SPARSKIT library. This'
331 message(2) = 'library is required if the "runge_kutta2" propagator is selected.'
332 message(3) = 'Try using a different propagation scheme or recompile with SPARSKIT support.'
333 call messages_fatal(3, namespace=namespace)
334#endif
335
336 call messages_experimental("Runge-Kutta 2 propagator", namespace=namespace)
338 ! set up pointer for zmf_dotu_aux
339 call mesh_init_mesh_aux(gr)
341 tr%tdsk_size = st%d%dim*gr%np
342 call sparskit_solver_init(namespace, st%d%dim*gr%np, tr%tdsk, .true.)
343
344#ifndef HAVE_SPARSKIT
345 message(1) = 'Octopus was not compiled with support for the SPARSKIT library. This'
346 message(2) = 'library is required if the "crank_nicolson_sparskit" propagator is selected.'
347 message(3) = 'Try using a different propagation scheme or recompile with SPARSKIT support.'
348 call messages_fatal(3, namespace=namespace)
349#endif
350 case (prop_magnus)
351 call messages_experimental("Magnus propagator", namespace=namespace)
352 if (family_is_mgga_with_exc) then
353 message(1) = "Magnus propagator with MGGA"
354 call messages_fatal(1, namespace=namespace)
355 end if
357 call messages_experimental("QOCT+TDDFT propagator", namespace=namespace)
359 call messages_experimental("explicit Runge-Kutta 4 propagator", namespace=namespace)
360 case (prop_cfmagnus4)
361 call messages_experimental("Commutator-Free Magnus propagator", namespace=namespace)
363 call messages_experimental("Explicit Magnus 3 propagator", namespace=namespace)
364 case default
365 call messages_input_error(namespace, 'TDPropagator')
366 end select
367 call messages_print_var_option('TDPropagator', tr%method, namespace=namespace)
368
369 if (have_fields) then
370 if (tr%method /= prop_etrs .and. &
371 tr%method /= prop_aetrs .and. &
372 tr%method /= prop_exponential_midpoint .and. &
373 tr%method /= prop_exponential_midpoint_predictor .and. &
374 tr%method /= prop_cfmagnus4 .and. &
375 tr%method /= prop_explicit_magnus3 .and. &
376 tr%method /= prop_qoct_tddft_propagator .and. &
377 tr%method /= prop_crank_nicolson .and. &
378 tr%method /= prop_runge_kutta4 .and. &
379 tr%method /= prop_explicit_runge_kutta4 .and. &
380 tr%method /= prop_runge_kutta2 .and. &
381 tr%method /= prop_crank_nicolson_sparskit) then
382 message(1) = "To move the ions or put in a gauge field, use one of the following propagators:"
383 message(2) = "etrs, aetrs, cfmagnus4, crank-nicolson, rk2, rk4, or exp_mid."
384 call messages_fatal(2, namespace=namespace)
385 end if
386 end if
387
388 if (relax_cell) then
389 if (tr%method /= prop_etrs .and. &
390 tr%method /= prop_aetrs .and. &
391 tr%method /= prop_exponential_midpoint .and. &
392 tr%method /= prop_exponential_midpoint_predictor .and. &
393 tr%method /= prop_cfmagnus4 .and. &
394 tr%method /= prop_explicit_magnus3 .and. &
395 tr%method /= prop_crank_nicolson .and. &
396 tr%method /= prop_crank_nicolson_sparskit) then
397 message(1) = "To relax the cell, use the one of the following propagators:"
398 message(2) = "etrs, aetrs, cfmagnus4, crank-nicolson, or exp_mid."
399 call messages_fatal(2, namespace=namespace)
400 end if
401 end if
402
403 select case (tr%method)
404 case (prop_cfmagnus4)
405 call ks_pot%init_interpolation(tr%vks_old, order = 4)
406 case default
407 call ks_pot%init_interpolation(tr%vks_old)
408 end select
409
410 call exponential_init(tr%te, namespace) ! initialize propagator
411
412 call messages_obsolete_variable(namespace, 'TDSelfConsistentSteps', 'TDStepsWithSelfConsistency')
413
414 !%Variable TDStepsWithSelfConsistency
415 !%Type integer
416 !%Default 0
417 !%Section Time-Dependent::Propagation
418 !%Description
419 !% Since the KS propagator is non-linear, each propagation step
420 !% should be performed self-consistently. In practice, for most
421 !% purposes this is not necessary, except perhaps in the first
422 !% iterations. This variable holds the number of propagation steps
423 !% for which the propagation is done self-consistently.
424 !%
425 !% The special value <tt>all_steps</tt> forces self-consistency to
426 !% be imposed on all propagation steps. A value of 0 means that
427 !% self-consistency will not be imposed. The default is 0.
428 !%Option all_steps -1
429 !% Self-consistency is imposed for all propagation steps.
430 !%End
431
432 call parse_variable(namespace, 'TDStepsWithSelfConsistency', 0, tr%scf_propagation_steps)
433
434 if (tr%scf_propagation_steps == -1) tr%scf_propagation_steps = huge(tr%scf_propagation_steps)
435 if (tr%scf_propagation_steps < 0) call messages_input_error(namespace, 'TDStepsWithSelfConsistency', 'Cannot be negative')
436
437 if (tr%scf_propagation_steps /= 0) then
438 call messages_experimental('TDStepsWithSelfConsistency', namespace=namespace)
439
440 if (tr%method /= prop_etrs) then
441 call messages_write('TDStepsWithSelfConsistency only works with the ETRS propagator')
442 call messages_fatal(namespace=namespace)
443 end if
444 end if
445
446 !%Variable TDSCFThreshold
447 !%Type float
448 !%Default 1.0e-6
449 !%Section Time-Dependent::Propagation
450 !%Description
451 !% Since the KS propagator is non-linear, each propagation step
452 !% should be performed self-consistently. In practice, for most
453 !% purposes this is not necessary, except perhaps in the first
454 !% iterations. This variable holds the number of propagation steps
455 !% for which the propagation is done self-consistently.
456 !%
457 !% The self consistency has to be measured against some accuracy
458 !% threshold. This variable controls the value of that threshold.
459 !%End
460 call parse_variable(namespace, 'TDSCFThreshold', 1.0e-6_real64, tr%scf_threshold)
461
462 pop_sub(propagator_elec_init)
463 end subroutine propagator_elec_init
464 ! ---------------------------------------------------------
465
466
467 ! ---------------------------------------------------------
468 subroutine propagator_elec_set_scf_prop(tr, threshold)
469 type(propagator_base_t), intent(inout) :: tr
470 real(real64), optional, intent(in) :: threshold
471
473
474 tr%scf_propagation_steps = huge(1)
475 if (present(threshold)) then
476 tr%scf_threshold = threshold
477 end if
478
480 end subroutine propagator_elec_set_scf_prop
481 ! ---------------------------------------------------------
482
483
484 ! ---------------------------------------------------------
486 type(propagator_base_t), intent(inout) :: tr
487
489
490 tr%scf_propagation_steps = -1
491
494 ! ---------------------------------------------------------
495
496
497 ! ---------------------------------------------------------
498 subroutine propagator_elec_end(tr)
499 type(propagator_base_t), intent(inout) :: tr
500
501 push_sub(propagator_elec_end)
502
503 call potential_interpolation_end(tr%vks_old)
504
505 select case (tr%method)
507 call propagator_rk_end()
508 call sparskit_solver_end(tr%tdsk)
509
511 call sparskit_solver_end(tr%tdsk)
512
513 end select
514
515 pop_sub(propagator_elec_end)
516 end subroutine propagator_elec_end
517 ! ---------------------------------------------------------
518
519
520 ! ---------------------------------------------------------
523 ! ---------------------------------------------------------
524 subroutine propagator_elec_dt(ks, namespace, space, hm, gr, st, tr, time, dt, nt, &
525 ions_dyn, ions, ext_partners, mc, outp, write_handler, scsteps, update_energy, qcchi)
526 type(v_ks_t), target, intent(inout) :: ks
527 type(namespace_t), intent(in) :: namespace
528 type(electron_space_t), intent(in) :: space
529 type(hamiltonian_elec_t), target, intent(inout) :: hm
530 type(grid_t), target, intent(inout) :: gr
531 type(states_elec_t), target, intent(inout) :: st
532 type(propagator_base_t), target, intent(inout) :: tr
533 real(real64), intent(in) :: time
534 real(real64), intent(in) :: dt
535 integer, intent(in) :: nt
536 type(ion_dynamics_t), intent(inout) :: ions_dyn
537 type(ions_t), intent(inout) :: ions
538 type(partner_list_t), intent(in) :: ext_partners
539 type(multicomm_t), intent(inout) :: mc
540 type(output_t), intent(in) :: outp
541 type(td_write_t), intent(in) :: write_handler
542 integer, optional, intent(out) :: scsteps
543 logical, optional, intent(in) :: update_energy
544 type(opt_control_state_t), optional, target, intent(inout) :: qcchi
545 type(gauge_field_t), pointer :: gfield
546 type(lasers_t), pointer :: lasers
547 logical :: generate, update_energy_
548 real(real64) :: am(space%dim)
549
550 call profiling_in("TD_PROPAGATOR")
551 push_sub(propagator_elec_dt)
552
553 update_energy_ = optional_default(update_energy, .true.)
554
555 call hm%ks_pot%interpolation_new(tr%vks_old, time, dt)
556
557 if (present(scsteps)) scsteps = 1
558
559 select case (tr%method)
560 case (prop_etrs)
561 if (self_consistent_step()) then
562 call td_etrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
563 ions_dyn, ions, mc, tr%scf_threshold, scsteps)
564 else
565 call td_etrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
566 ions_dyn, ions, mc)
567 end if
568 case (prop_aetrs)
569 call td_aetrs(namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
570 case (prop_caetrs)
571 call td_caetrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc)
573 call exponential_midpoint(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
575 call exponential_midpoint_predictor(hm, ks, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, &
576 ext_partners, mc)
578 call td_crank_nicolson(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, .false.)
579 case (prop_runge_kutta4)
580 call td_runge_kutta4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
581 case (prop_runge_kutta2)
582 call td_runge_kutta2(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
584 call td_crank_nicolson(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, .true.)
585 case (prop_magnus)
586 call td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
588 call td_qoct_tddft_propagator(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
590 if (present(qcchi)) then
591 call td_explicit_runge_kutta4(ks, namespace, space, hm, gr, st, time, dt, ions_dyn, ions, ext_partners, qcchi)
592 else
593 call td_explicit_runge_kutta4(ks, namespace, space, hm, gr, st, time, dt, ions_dyn, ions, ext_partners)
594 end if
595 case (prop_cfmagnus4)
596 call td_cfmagnus4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, nt)
598 call td_explicit_magnus3(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
599 end select
600
601 generate = .false.
602 if (ions_dyn%is_active()) then
603 if (.not. propagator_elec_ions_are_propagated(tr)) then
604 call propagation_ops_elec_propagate_ions_and_cell(gr, hm, st, namespace, space, &
605 ions_dyn, ions, mc, abs(nt*dt), ions_dyn%ionic_scale*dt)
606 generate = .true.
607 end if
608 end if
609
610 gfield => list_get_gauge_field(ext_partners)
611 if(associated(gfield)) then
612 if (gauge_field_is_propagated(gfield) .and. .not. propagator_elec_ions_are_propagated(tr)) then
614 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
615 end if
616 end if
617
618 if (generate .or. ions%has_time_dependent_species()) then
619 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = abs(nt*dt))
620 end if
621
622 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
623 calc_eigenval = update_energy_, time = abs(nt*dt), calc_energy = update_energy_)
624
625 ! compute time integral of the paramagnetic current
626 if (ks%xc_photon /= 0) then
627 call ks%xc_photons%add_mean_field(ks%gr, space, hm%kpoints, st, abs(nt*dt), dt)
628 endif
629
630 if (update_energy_) call energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit = -1)
631
632 ! Recalculate forces, update velocities...
633 if ((ions_dyn%ions_move() .and. tr%method .ne. prop_explicit_runge_kutta4) .or. &
634 outp%what(option__output__forces) .or. write_handler%out(out_separate_forces)%write) then
635 call forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, t = abs(nt*dt), dt = dt)
636 end if
637
638 if (ions_dyn%cell_relax() .or. outp%what(option__output__stress)) then
639 call stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
640 if(ions_dyn%cell_relax()) then
641 call ions_dyn%update_stress(ions%space, st%stress_tensors%total, ions%latt%rlattice, ions%latt%rcell_volume)
642 end if
643 end if
644
645 if( ions_dyn%is_active()) then
646 call ion_dynamics_propagate_vel(ions_dyn, ions, atoms_moved = generate)
647 call ions%update_kinetic_energy()
648 end if
649
650 if(associated(gfield)) then
651 if(gauge_field_is_propagated(gfield)) then
652 call gauge_field_get_force(gfield, gr, st%d%spin_channels, st%current, ks%xc%lrc)
654 end if
655 end if
656
657 !We update the occupation matrices
658 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
659
660 !Updates Nondipole Vector potential
661 lasers => list_get_lasers(ext_partners)
662 if(associated(lasers)) then
663 if (lasers_with_nondipole_field(lasers)) then
664 call lasers_nondipole_laser_field_step(lasers, am, time)
665 call lasers_set_nondipole_parameters(lasers,am,time)
666 end if
667 end if
668 pop_sub(propagator_elec_dt)
669 call profiling_out("TD_PROPAGATOR")
670
671 contains
672
673 ! ---------------------------------------------------------
674 logical pure function self_consistent_step() result(scs)
675 scs = (hm%theory_level /= independent_particles) .and. (time <= tr%scf_propagation_steps*abs(dt) + m_epsilon)
676 end function self_consistent_step
677 ! ---------------------------------------------------------
678
679 end subroutine propagator_elec_dt
680 ! ---------------------------------------------------------
681
682
683
684
685 ! ---------------------------------------------------------
686 logical pure function propagator_elec_ions_are_propagated(tr) result(propagated)
687 type(propagator_base_t), intent(in) :: tr
688
689 select case (tr%method)
690 case (prop_etrs, prop_aetrs, prop_caetrs, prop_explicit_runge_kutta4)
691 propagated = .true.
692 case default
693 propagated = .false.
694 end select
695
697 ! ---------------------------------------------------------
698
699
700 ! ---------------------------------------------------------
701
702 subroutine propagator_elec_dt_bo(scf, namespace, space, gr, ks, st, hm, ions, ext_partners, mc, iter, dt, ions_dyn, scsteps)
703 type(scf_t), intent(inout) :: scf
704 type(namespace_t), intent(in) :: namespace
705 type(electron_space_t), intent(in) :: space
706 type(grid_t), intent(inout) :: gr
707 type(v_ks_t), intent(inout) :: ks
708 type(states_elec_t), intent(inout) :: st
709 type(hamiltonian_elec_t), intent(inout) :: hm
710 type(ions_t), intent(inout) :: ions
711 type(partner_list_t), intent(in) :: ext_partners
712 type(multicomm_t), intent(inout) :: mc
713 integer, intent(in) :: iter
714 real(real64), intent(in) :: dt
715 type(ion_dynamics_t), intent(inout) :: ions_dyn
716 integer, intent(inout) :: scsteps
717
718 type(gauge_field_t), pointer :: gfield
719
720 push_sub(propagator_elec_dt_bo)
721
722 ! move the hamiltonian to time t
723 call propagation_ops_elec_propagate_ions_and_cell(gr, hm, st, namespace, space, &
724 ions_dyn, ions, mc, iter*dt, dt)
725
726 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = iter*dt)
727 ! now calculate the eigenfunctions
728 call scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, &
729 verbosity = verb_compact, iters_done = scsteps)
730
731 ! Stress tensor calculation and uptade
732 if (ions_dyn%cell_relax()) then
733 if (.not. scf%calc_stress) then
734 call stress_calculate(namespace, gr, hm, st, ions, ks, ext_partners)
735 end if
736 call ions_dyn%update_stress(ions%space, st%stress_tensors%total, ions%latt%rlattice, ions%latt%rcell_volume)
737 end if
738
739 gfield => list_get_gauge_field(ext_partners)
740 if(associated(gfield)) then
741 if (gauge_field_is_propagated(gfield)) then
742 call gauge_field_do_algorithmic_operation(gfield, op_verlet_compute_acc, dt, iter*dt)
743 end if
744 end if
745
746 !TODO: we should update the occupation matrices here
747 if (hm%lda_u_level /= dft_u_none) then
748 call messages_not_implemented("DFT+U with propagator_elec_dt_bo", namespace=namespace)
749 end if
750
751 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = iter*dt)
752
753 ! update Hamiltonian and eigenvalues (fermi is *not* called)
754 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
755 calc_eigenval = .true., time = iter*dt, calc_energy = .true.)
756
757 ! Get the energies.
758 call energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit = -1)
759
760 call ion_dynamics_propagate_vel(ions_dyn, ions)
761 call hamiltonian_elec_epot_generate(hm, namespace, space, gr, ions, ext_partners, st, time = iter*dt)
762 call ions%update_kinetic_energy()
763
764 if(associated(gfield)) then
765 if (gauge_field_is_propagated(gfield)) then
766 call gauge_field_do_algorithmic_operation(gfield, op_verlet_compute_vel, dt, iter*dt)
767 end if
768 end if
770 pop_sub(propagator_elec_dt_bo)
771 end subroutine propagator_elec_dt_bo
772
773end module propagator_elec_oct_m
774
775
776!! Local Variables:
777!! mode: f90
778!! coding: utf-8
779!! End:
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,...
subroutine, public exponential_copy(teo, tei)
subroutine, public exponential_init(te, namespace, full_batch)
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_get_force(this, gr, spin_channels, current, lrc)
subroutine, public gauge_field_do_algorithmic_operation(this, operation, dt, time)
logical pure function, public gauge_field_is_propagated(this)
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:241
real(real64), parameter, public m_epsilon
Definition: global.F90:207
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
This module defines classes and functions for interaction partners.
subroutine, public ion_dynamics_propagate_vel(this, ions, atoms_moved)
A module to handle KS potential, without the external potential.
subroutine, public lasers_nondipole_laser_field_step(this, field, time)
Retrieves the NDSFA vector_potential correction. The nondipole field is obtained for consecutive time...
Definition: lasers.F90:1152
subroutine, public lasers_set_nondipole_parameters(this, ndfield, nd_integration_time)
Set parameters for nondipole SFA calculation.
Definition: lasers.F90:754
logical function, public lasers_with_nondipole_field(lasers)
Check if a nondipole SFA correction should be computed for the given laser.
Definition: lasers.F90:741
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:892
This module defines various routines, operating on mesh functions.
integer, public sp_distdot_mode
subroutine, public mesh_init_mesh_aux(mesh)
Initialise a pointer to the grid/mesh, that is globally exposed, such that low level mesh operations ...
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
This module holds the "opt_control_state_t" datatype, which contains a quantum-classical state.
this module contains the low-level part of the output system
Definition: output_low.F90:117
subroutine, public potential_interpolation_copy(vkso, vksi)
subroutine, public potential_interpolation_end(potential_interpolation)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public propagation_ops_elec_propagate_ions_and_cell(gr, hm, st, namespace, space, ions_dyn, ions, mc, time, dt_ions)
subroutine, public propagation_ops_elec_update_hamiltonian(namespace, space, st, mesh, hm, ext_partners, time)
integer, parameter, public prop_aetrs
integer, parameter, public prop_crank_nicolson_sparskit
integer, parameter, public prop_crank_nicolson
integer, parameter, public prop_runge_kutta2
integer, parameter, public prop_magnus
integer, parameter, public prop_cfmagnus4
integer, parameter, public prop_caetrs
integer, parameter, public prop_qoct_tddft_propagator
integer, parameter, public prop_exponential_midpoint
integer, parameter, public prop_explicit_magnus3
integer, parameter, public prop_runge_kutta4
integer, parameter, public prop_explicit_runge_kutta4
integer, parameter, public prop_etrs
integer, parameter, public prop_exponential_midpoint_predictor
subroutine, public td_crank_nicolson(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, use_sparskit)
Crank-Nicolson propagator.
subroutine, public propagator_elec_copy(tro, tri)
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_set_scf_prop(tr, threshold)
subroutine, public propagator_elec_init(gr, namespace, st, tr, ks_pot, have_fields, family_is_mgga_with_exc, relax_cell)
subroutine, public propagator_elec_remove_scf_prop(tr)
logical pure function, public propagator_elec_ions_are_propagated(tr)
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)
subroutine, public td_aetrs(namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Propagator with approximate enforced time-reversal symmetry.
subroutine, public td_etrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc, sctol, scsteps)
Propagator with enforced time-reversal symmetry and self-consistency.
subroutine, public td_etrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc)
Propagator with enforced time-reversal symmetry.
subroutine, public td_caetrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc)
Propagator with approximate enforced time-reversal symmetry.
subroutine, public exponential_midpoint_predictor(hm, ks, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Exponential midpoint with prediction of the half time step This is equivalent to an explicit Magnus i...
subroutine, public exponential_midpoint(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Exponential midpoint.
subroutine, public td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
Magnus propagator.
subroutine, public td_cfmagnus4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, iter)
Commutator-free Magnus propagator of order 4. This method has been developed for linear systems and i...
subroutine, public td_explicit_magnus3(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Third-order explicit Magnus propagator. This implements equation (22) from Casas, F....
subroutine, public td_qoct_tddft_propagator(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Propagator specifically designed for the QOCT+TDDFT problem.
subroutine, public td_runge_kutta2(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
subroutine, public td_runge_kutta4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
subroutine, public td_explicit_runge_kutta4(ks, namespace, space, hm, gr, st, time, dt, ions_dyn, ions, ext_partners, qcchi)
subroutine, public propagator_rk_end()
type(algorithmic_operation_t), parameter, public op_verlet_compute_acc
type(algorithmic_operation_t), parameter, public op_verlet_compute_vel
subroutine, public sparskit_solver_init(namespace, n, sk, is_complex)
Definition: sparskit.F90:240
subroutine, public sparskit_solver_end(sk)
Definition: sparskit.F90:437
subroutine, public sparskit_solver_copy(sko, ski)
Definition: sparskit.F90:450
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
integer, parameter, public out_separate_forces
Definition: td_write.F90:203
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:754
Definition: xc.F90:116
logical pure function self_consistent_step()
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
Stores all communicators and groups.
Definition: multicomm.F90:208
This is the datatype that contains the objects that are propagated: in principle this could be both t...
output handler class
Definition: output_low.F90:166
The states_elec_t class contains all electronic wave functions.
Time-dependent Write Properties.
Definition: td_write.F90:321
int true(void)