Octopus
propagator_magnus.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
23 use batch_oct_m
25 use debug_oct_m
31 use global_oct_m
32 use grid_oct_m
38 use ions_oct_m
39 use, intrinsic :: iso_fortran_env
41 use lasers_oct_m
43 use mesh_oct_m
46 use parser_oct_m
52 use space_oct_m
54 use v_ks_oct_m
56 use xc_oct_m
57
58 implicit none
59
60 private
61
62 public :: &
63 td_magnus, &
66
67
69 private
70 real(real64), allocatable :: vmagnus(:, :, :)
71 contains
72 procedure :: apply => magnus4_operator_apply
73 end type magnus4_operator_t
74
75 interface magnus4_operator_t
76 procedure magnus4_operator_constructor
77 end interface magnus4_operator_t
78
79
81 private
82 type(hamiltonian_elec_t), pointer :: hm1 => null()
83 type(hamiltonian_elec_t), pointer :: hm2 => null()
84 real(real64) :: c1 = m_zero
85 real(real64) :: c2 = m_zero
86 contains
87 procedure :: apply => magnus3_linear_operator_apply
89
91 procedure magnus3_linear_operator_constructor
92 end interface magnus3_linear_operator_t
93
94
96 private
97 type(hamiltonian_elec_t), pointer :: hm1 => null()
98 type(hamiltonian_elec_t), pointer :: hm2 => null()
99 type(hamiltonian_elec_t), pointer :: hm3 => null()
100 type(hamiltonian_elec_t), pointer :: hm4 => null()
101 real(real64) :: dt = m_zero
102 contains
103 procedure :: apply => magnus3_operator_apply
104 end type magnus3_operator_t
105
106 interface magnus3_operator_t
107 procedure magnus3_operator_constructor
108 end interface magnus3_operator_t
109
110
111contains
112
113 ! ---------------------------------------------------------
115 subroutine td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
116 type(hamiltonian_elec_t), intent(inout) :: hm
117 type(partner_list_t), intent(in) :: ext_partners
118 type(grid_t), intent(in) :: gr
119 type(states_elec_t), intent(inout) :: st
120 type(propagator_base_t), intent(inout) :: tr
121 type(namespace_t), intent(in) :: namespace
122 real(real64), intent(in) :: time
123 real(real64), intent(in) :: dt
124
125 integer :: j, is, i
126 real(real64) :: atime(2)
127 real(real64), allocatable :: vaux(:, :, :), vmagnus(:, :, :), pot(:)
128 type(lasers_t), pointer :: lasers
129 class(operator_t), pointer :: op
130
131 push_sub(propagator_dt.td_magnus)
132
133 assert(.not. family_is_mgga_with_exc(hm%xc))
134
135 safe_allocate(vaux(1:gr%np, 1:st%d%nspin, 1:2))
136 safe_allocate(vmagnus(1:gr%np, 1:st%d%nspin, 1:2))
137
138 atime(1) = (m_half-sqrt(m_three)/6.0_real64)*dt
139 atime(2) = (m_half+sqrt(m_three)/6.0_real64)*dt
140
141 if (hm%theory_level /= independent_particles) then
142 do j = 1, 2
143 call potential_interpolation_interpolate(tr%vks_old, 3, time, dt, time - dt + atime(j), vaux(:, :, j))
144 end do
145 else
146 vaux = m_zero
147 end if
148
149 do j = 1, 2
150 ! WARNING: This should be carefully tested, and extended to allow for velocity-gauge laser fields.
151 lasers => list_get_lasers(ext_partners)
152 if(associated(lasers)) then
153 do i = 1, lasers%no_lasers
154 select case (laser_kind(lasers%lasers(i)))
155 case (e_field_electric)
156 safe_allocate(pot(1:gr%np))
157 pot = m_zero
158 call laser_potential(lasers%lasers(i), gr, pot, time - dt + atime(j))
159 do is = 1, st%d%nspin
160 vaux(:, is, j) = vaux(:, is, j) + pot(:)
161 end do
162 safe_deallocate_a(pot)
164 write(message(1),'(a)') 'The Magnus propagator cannot be used with magnetic fields, or'
165 write(message(2),'(a)') 'with an electric field described in the velocity gauge.'
166 call messages_fatal(2, namespace=namespace)
167 end select
168 end do
169 end if
170 end do
171
172 vmagnus(:, :, 2) = m_half*(vaux(:, :, 1) + vaux(:, :, 2))
173 vmagnus(:, :, 1) = (sqrt(m_three)/12.0_real64)*dt*(vaux(:, :, 2) - vaux(:, :, 1))
174
175 op => magnus4_operator_t(namespace, gr, hm, vmagnus)
176 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, dt, op=op)
177 safe_deallocate_p(op)
179 safe_deallocate_a(vaux)
180 safe_deallocate_a(vmagnus)
181 pop_sub(propagator_dt.td_magnus)
182 end subroutine td_magnus
183
184 function magnus4_operator_constructor(namespace, mesh, hm, vmagnus) result(this)
185 type(namespace_t), target, intent(in) :: namespace
186 class(mesh_t), target, intent(in) :: mesh
187 class(hamiltonian_abst_t), target, intent(in) :: hm
188 real(real64), intent(in) :: vmagnus(:, :, :)
189 type(magnus4_operator_t), pointer :: this
193 allocate(this)
194 call this%init(namespace, mesh, hm)
196 select type(hm)
197 class is (hamiltonian_elec_t)
198 assert(size(vmagnus, 1) >= mesh%np)
199 assert(size(vmagnus, 2) == hm%d%nspin)
200 assert(size(vmagnus, 3) == 2)
201 class default
202 write(message(1), '(a)') 'Magnus operators only implemented for electrons at the moment'
203 call messages_fatal(1, namespace=namespace)
204 end select
205 safe_allocate_source(this%vmagnus, vmagnus)
206
209
210 subroutine magnus4_operator_apply(this, psib, hpsib)
211 class(magnus4_operator_t), intent(in) :: this
212 class(batch_t), intent(inout) :: psib
213 class(batch_t), intent(inout) :: hpsib
214
215 integer :: ispin
216 type(wfs_elec_t) :: auxpsib, aux2psib
217
218 push_sub(magnus4_operator_apply)
219
220 select type(hm => this%hm)
221 class is (hamiltonian_elec_t)
222 select type (psib)
223 class is (wfs_elec_t)
224 select type (hpsib)
225 class is (wfs_elec_t)
226 ! We will assume, for the moment, no spinors.
227 if (hm%d%dim /= 1) call messages_not_implemented("Magnus with spinors", namespace=this%namespace)
228
229 assert(psib%nst == hpsib%nst)
230
231 ispin = hm%d%get_spin_index(psib%ik)
232
233 call hpsib%copy_to(auxpsib, copy_data=.false.)
234 call hpsib%copy_to(aux2psib, copy_data=.false.)
235
236 ! Compute (T + Vnl)|psi> and store it
237 call zhamiltonian_elec_apply_batch(hm, this%namespace, this%mesh, psib, auxpsib, &
239
240 ! H|psi> = (T + Vnl)|psi> + Vpsl|psi> + Vmagnus(t2)|psi> + Vborders|psi>
241 call auxpsib%copy_data_to(this%mesh%np, hpsib)
242 call zhamiltonian_elec_external(hm, this%mesh, psib, hpsib)
243 if (hm%abs_boundaries%abtype == imaginary_absorbing) then
244 call batch_mul(this%mesh%np, hm%abs_boundaries%mf(1:this%mesh%np), psib, aux2psib)
245 call batch_axpy(this%mesh%np, m_zi, aux2psib, hpsib)
246 end if
247 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 2), psib, aux2psib)
248 call batch_axpy(this%mesh%np, m_one, aux2psib, hpsib)
249
250 ! Add first term of the commutator: - i Vmagnus(t1) (T + Vnl) |psi>
251 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 1), auxpsib, aux2psib)
252 call batch_axpy(this%mesh%np, -m_zi, aux2psib, hpsib)
253
254 ! Add second term of commutator: i (T + Vnl) Vmagnus(t1) |psi>
255 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 1), psib, auxpsib)
256 call zhamiltonian_elec_apply_batch(hm, this%namespace, this%mesh, auxpsib, aux2psib, &
258 call batch_axpy(this%mesh%np, m_zi, aux2psib, hpsib)
259
260 call auxpsib%end()
261 call aux2psib%end()
262 class default
263 message(1) = "Internal error: imcompatible batch_t in magnus4_operator_apply for argument hpsib."
264 call messages_fatal(1, namespace=this%namespace)
265 end select
266 class default
267 message(1) = "Internal error: imcompatible batch_t in magnus4_operator_apply for argument psib."
268 call messages_fatal(1, namespace=this%namespace)
269 end select
270 end select
271
273 end subroutine magnus4_operator_apply
274
275 ! ---------------------------------------------------------
277 function magnus3_linear_operator_constructor(namespace, mesh, hm1, hm2, c1, c2) result(this)
278 type(namespace_t), target, intent(in) :: namespace
279 class(mesh_t), target, intent(in) :: mesh
280 type(hamiltonian_elec_t), target, intent(in) :: hm1
281 type(hamiltonian_elec_t), target, intent(in) :: hm2
282 real(real64), intent(in) :: c1
283 real(real64), intent(in) :: c2
284 type(magnus3_linear_operator_t), pointer :: this
285
287
288 allocate(this)
289 call this%init(namespace, mesh, hm1)
290 this%hm1 => hm1
291 this%hm2 => hm2
292 this%c1 = c1
293 this%c2 = c2
294
297
298 ! ---------------------------------------------------------
300 subroutine magnus3_apply_stage_in_ref_phase(namespace, mesh, hm_ref, hm_stage, psib, hpsib)
301 type(namespace_t), intent(in) :: namespace
302 class(mesh_t), intent(in) :: mesh
303 type(hamiltonian_elec_t), intent(in) :: hm_ref
304 type(hamiltonian_elec_t), intent(in) :: hm_stage
305 type(wfs_elec_t), intent(inout) :: psib
306 type(wfs_elec_t), intent(inout) :: hpsib
307
309
310 if (hm_ref%phase%is_allocated()) call hm_ref%phase%unset_phase_corr(mesh, psib)
311 if (hm_stage%phase%is_allocated()) call hm_stage%phase%set_phase_corr(mesh, psib)
312
313 call zhamiltonian_elec_apply_batch(hm_stage, namespace, mesh, psib, hpsib)
314
315 if (hm_stage%phase%is_allocated()) then
316 call hm_stage%phase%unset_phase_corr(mesh, hpsib)
317 call hm_stage%phase%unset_phase_corr(mesh, psib)
318 end if
319 if (hm_ref%phase%is_allocated()) then
320 call hm_ref%phase%set_phase_corr(mesh, hpsib)
321 call hm_ref%phase%set_phase_corr(mesh, psib)
322 end if
323
326
327 ! ---------------------------------------------------------
329 subroutine magnus3_linear_operator_apply(this, psib, hpsib)
330 class(magnus3_linear_operator_t), intent(in) :: this
331 class(batch_t), intent(inout) :: psib
332 class(batch_t), intent(inout) :: hpsib
333
334 type(wfs_elec_t) :: auxpsib
335
337
338 select type (hm => this%hm)
339 class is (hamiltonian_elec_t)
340 select type (psib)
341 class is (wfs_elec_t)
342 select type (hpsib)
343 class is (wfs_elec_t)
344 if (hm%d%dim /= 1) call messages_not_implemented("Explicit Magnus 3 with spinors", namespace=this%namespace)
345
346 assert(psib%nst == hpsib%nst)
347
348 call hpsib%copy_to(auxpsib, copy_data=.false.)
349
350 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm1, psib, hpsib)
351 call batch_scal(this%mesh%np, this%c1, hpsib)
352
353 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm2, psib, auxpsib)
354 call batch_axpy(this%mesh%np, this%c2, auxpsib, hpsib)
355
356 call auxpsib%end()
357 class default
358 message(1) = "Internal error: imcompatible batch_t in magnus3_linear_operator_apply for argument hpsib."
359 call messages_fatal(1, namespace=this%namespace)
360 end select
361 class default
362 message(1) = "Internal error: imcompatible batch_t in magnus3_linear_operator_apply for argument psib."
363 call messages_fatal(1, namespace=this%namespace)
364 end select
365 end select
366
368 end subroutine magnus3_linear_operator_apply
369
370 ! ---------------------------------------------------------
372 function magnus3_operator_constructor(namespace, mesh, hm1, hm2, hm3, hm4, dt) result(this)
373 type(namespace_t), target, intent(in) :: namespace
374 class(mesh_t), target, intent(in) :: mesh
375 type(hamiltonian_elec_t), target, intent(in) :: hm1
376 type(hamiltonian_elec_t), target, intent(in) :: hm2
377 type(hamiltonian_elec_t), target, intent(in) :: hm3
378 type(hamiltonian_elec_t), target, intent(in) :: hm4
379 real(real64), intent(in) :: dt
380 type(magnus3_operator_t), pointer :: this
381
383
384 allocate(this)
385 call this%init(namespace, mesh, hm1)
386 this%hm1 => hm1
387 this%hm2 => hm2
388 this%hm3 => hm3
389 this%hm4 => hm4
390 this%dt = dt
391
394
395 ! ---------------------------------------------------------
400 subroutine magnus3_operator_apply(this, psib, hpsib)
401 class(magnus3_operator_t), intent(in) :: this
402 class(batch_t), intent(inout) :: psib
403 class(batch_t), intent(inout) :: hpsib
404
405 type(wfs_elec_t) :: k1psib, k2psib, k3psib, k4psib
406 type(wfs_elec_t) :: work1psib, work2psib, sum12psib, commpsib
407
408 push_sub(magnus3_operator_apply)
409
410 select type (hm => this%hm)
411 class is (hamiltonian_elec_t)
412 select type (psib)
413 class is (wfs_elec_t)
414 select type (hpsib)
415 class is (wfs_elec_t)
416 if (hm%d%dim /= 1) call messages_not_implemented("Explicit Magnus 3 with spinors", namespace=this%namespace)
417
418 assert(psib%nst == hpsib%nst)
419
420 call hpsib%copy_to(k1psib, copy_data=.false.)
421 call hpsib%copy_to(k2psib, copy_data=.false.)
422 call hpsib%copy_to(k3psib, copy_data=.false.)
423 call hpsib%copy_to(k4psib, copy_data=.false.)
424 call hpsib%copy_to(work1psib, copy_data=.false.)
425 call hpsib%copy_to(work2psib, copy_data=.false.)
426 call hpsib%copy_to(sum12psib, copy_data=.false.)
427 call hpsib%copy_to(commpsib, copy_data=.false.)
428
429 ! Precompute stage derivatives once:
430 ! k1 = H1 psi, k2 = H2 psi, k3 = H3 psi, k4 = H4 psi.
431 ! no phase treatment needed for first stage
432 call zhamiltonian_elec_apply_batch(this%hm1, this%namespace, this%mesh, psib, k1psib)
433 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm2, psib, k2psib)
434 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm3, psib, k3psib)
435 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm4, psib, k4psib)
436
437 ! Weighted sum in eq. 22: (k1 + 4 k3 + k4)/6.
438 call k1psib%copy_data_to(this%mesh%np, hpsib)
439 call batch_scal(this%mesh%np, m_one/6.0_real64, hpsib)
440 call batch_axpy(this%mesh%np, m_four/6.0_real64, k3psib, hpsib)
441 call batch_axpy(this%mesh%np, m_one/6.0_real64, k4psib, hpsib)
442
443 ! Commutator term from eq. 22:
444 ! -1/3 [u3, k3] - 1/12 [u4, k4]
445 ! where u3 = (k1 + k2)/4 and u4 = k2.
446
447 ! [H1 + H2, H3] psi
448 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm1, k3psib, work1psib)
449 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm2, k3psib, work2psib)
450 call batch_axpy(this%mesh%np, m_one, work2psib, work1psib)
451
452 call k1psib%copy_data_to(this%mesh%np, sum12psib)
453 call batch_axpy(this%mesh%np, m_one, k2psib, sum12psib)
454 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm3, sum12psib, commpsib)
455 call batch_axpy(this%mesh%np, -m_one, commpsib, work1psib)
456
457 ! [H2, H4] psi
458 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm2, k4psib, work2psib)
459 call magnus3_apply_stage_in_ref_phase(this%namespace, this%mesh, this%hm1, this%hm4, k2psib, commpsib)
460 call batch_axpy(this%mesh%np, -m_one, commpsib, work2psib)
461
462 call batch_axpy(this%mesh%np, m_one, work2psib, work1psib)
463 call batch_axpy(this%mesh%np, (m_zi*this%dt)/12.0_real64, work1psib, hpsib)
464
465 call k1psib%end()
466 call k2psib%end()
467 call k3psib%end()
468 call k4psib%end()
469 call work1psib%end()
470 call work2psib%end()
471 call sum12psib%end()
472 call commpsib%end()
473 class default
474 message(1) = "Internal error: imcompatible batch_t in magnus3_operator_apply for argument hpsib."
475 call messages_fatal(1, namespace=this%namespace)
476 end select
477 class default
478 message(1) = "Internal error: imcompatible batch_t in magnus3_operator_apply for argument psib."
479 call messages_fatal(1, namespace=this%namespace)
480 end select
481 end select
482
484 end subroutine magnus3_operator_apply
485
486 ! ---------------------------------------------------------
491 subroutine td_explicit_magnus3(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
492 type(v_ks_t), intent(inout) :: ks
493 type(namespace_t), intent(in) :: namespace
494 type(electron_space_t), intent(in) :: space
495 type(hamiltonian_elec_t), target, intent(inout) :: hm
496 type(grid_t), target, intent(inout) :: gr
497 type(states_elec_t), target, intent(inout) :: st
498 type(propagator_base_t), target, intent(inout) :: tr
499 real(real64), intent(in) :: time
500 real(real64), intent(in) :: dt
501 type(ion_dynamics_t), intent(inout) :: ions_dyn
502 type(ions_t), intent(inout) :: ions
503 type(partner_list_t), intent(in) :: ext_partners
504 type(multicomm_t), intent(inout) :: mc
505
506 type(gauge_field_t), pointer :: gfield
507 type(states_elec_t) :: st0, st2, st3, st4
508 type(hamiltonian_elec_t), target :: hm1, hm2, hm3, hm4
509 class(operator_t), pointer :: op
510
511 push_sub(propagator_dt.td_explicit_magnus3)
512
513 ! Stage u1 / k1 (eq. 22):
514 ! Keep Y_n and snapshot H1 = H(t_n, Y_n).
515 call states_elec_copy(st0, st)
516 call hamiltonian_elec_copy(hm1, hm)
517
518 ! Stage u2 = 1/2 k1 -> k2 (eq. 22):
519 ! Predict Y at t_n + dt/2 using H1, then build hm2 from that Y at the midpoint
520 call states_elec_copy(st2, st0)
521 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st2, gr, hm1, m_half*dt)
522
523 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
524 ext_partners, mc, time - m_half*dt, m_half*dt, save_pos = .true.)
525
526 gfield => list_get_gauge_field(ext_partners)
527 if (associated(gfield)) then
528 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, m_half*dt, time, save_gf=.true.)
529 end if
530
531 if (hm%theory_level /= independent_particles) then
532 call v_ks_calc(ks, namespace, space, hm, st2, ions, ext_partners, &
533 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
534 end if
535 call propagation_ops_elec_update_hamiltonian(namespace, space, st2, gr, hm, ext_partners, time - m_half*dt)
536 call hamiltonian_elec_copy(hm2, hm)
537
538 ! Stage u3 = 1/4 (k1 + k2) -> k3 (eq. 22):
539 ! Propagate with 1/4 H1 + 1/4 H2 and build hm3 at the midpoint from st3
540 call states_elec_copy(st3, st0)
541 op => magnus3_linear_operator_t(namespace, gr, hm1, hm2, m_one/4.0_real64, m_one/4.0_real64)
542 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st3, gr, hm1, dt, op=op)
543 safe_deallocate_p(op)
544
545 if (hm%theory_level /= independent_particles) then
546 call v_ks_calc(ks, namespace, space, hm, st3, ions, ext_partners, &
547 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
548 end if
549 call propagation_ops_elec_update_hamiltonian(namespace, space, st3, gr, hm, ext_partners, time - m_half*dt)
550 call hamiltonian_elec_copy(hm3, hm)
551
552 ! Stage u4 = k2 -> k4 (eq. 22):
553 ! Propagate Y_n with H2 to t_n + dt, then rebuild H4 at the end of the step.
554 call states_elec_copy(st4, st0)
555 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st4, gr, hm2, dt)
556
557 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
558 ext_partners, mc, time, m_half*dt)
559 if (associated(gfield)) then
560 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, m_half*dt, time)
561 end if
562
563 if (hm%theory_level /= independent_particles) then
564 call v_ks_calc(ks, namespace, space, hm, st4, ions, ext_partners, &
565 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
566 end if
567 call propagation_ops_elec_update_hamiltonian(namespace, space, st4, gr, hm, ext_partners, time)
568 call hamiltonian_elec_copy(hm4, hm)
569
570 ! Final stage v3 (eq. 22):
571 ! Combine k1, k3, k4 and the commutator correction
572 ! -1/3 [u3, k3] - 1/12 [u4, k4].
573 op => magnus3_operator_t(namespace, gr, hm1, hm2, hm3, hm4, dt)
574 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm1, dt, op=op)
575 safe_deallocate_p(op)
576
577 ! Restore ions and gauge field to their saved values at t_n.
578 call propagation_ops_elec_restore_ions(tr%propagation_ops_elec, ions_dyn, ions)
579 if (associated(gfield)) then
580 call propagation_ops_elec_restore_gauge_field(tr%propagation_ops_elec, namespace, space, hm, gr, ext_partners)
581 end if
582
583 call states_elec_end(st0)
584 call states_elec_end(st2)
585 call states_elec_end(st3)
587 call hamiltonian_elec_end(hm1)
588 call hamiltonian_elec_end(hm2)
589 call hamiltonian_elec_end(hm3)
590 call hamiltonian_elec_end(hm4)
591
592 pop_sub(propagator_dt.td_explicit_magnus3)
593 end subroutine td_explicit_magnus3
594
595 ! ---------------------------------------------------------
613 subroutine td_cfmagnus4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc, iter)
614 type(v_ks_t), target, intent(inout) :: ks
615 type(namespace_t), intent(in) :: namespace
616 type(electron_space_t), intent(in) :: space
617 type(hamiltonian_elec_t), target, intent(inout) :: hm
618 type(grid_t), target, intent(inout) :: gr
619 type(states_elec_t), target, intent(inout) :: st
620 type(propagator_base_t), target, intent(inout) :: tr
621 real(real64), intent(in) :: time
622 real(real64), intent(in) :: dt
623 type(ion_dynamics_t), intent(inout) :: ions_dyn
624 type(ions_t), intent(inout) :: ions
625 type(partner_list_t), intent(in) :: ext_partners
626 type(multicomm_t), intent(inout) :: mc
627 integer, intent(in) :: iter
628
629 real(real64) :: alpha1, alpha2, c1, c2, t1, t2, dt1, dt2
630 type(hamiltonian_elec_t), target :: hm1, hm2
631 class(operator_t), pointer :: op
632 type(gauge_field_t), pointer :: gfield
633
634 push_sub(propagator_dt.td_cfmagnus4)
635
636 alpha1 = (m_three - m_two * sqrt(m_three))/12.0_real64
637 alpha2 = (m_three + m_two * sqrt(m_three))/12.0_real64
638 c1 = m_half - sqrt(m_three)/6.0_real64
639 c2 = m_half + sqrt(m_three)/6.0_real64
640 t1 = time - dt + c1*dt
641 t2 = time - dt + c2*dt
642
643 gfield => list_get_gauge_field(ext_partners)
644 dt1 = c1*dt
645 dt2 = (c2 - c1)*dt
646
647 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
648 ext_partners, mc, t1, dt1, save_pos=.true.)
649 if (associated(gfield)) then
650 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt1, time, save_gf=.true.)
651 end if
652
653 call hm%ks_pot%interpolate_potentials(tr%vks_old, 4, time, dt, t1)
654 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, t1)
655 call hamiltonian_elec_copy(hm1, hm)
656
657 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
658 ext_partners, mc, t2, dt2)
659 if (associated(gfield)) then
660 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt2, time)
661 end if
662
663 call hm%ks_pot%interpolate_potentials(tr%vks_old, 4, time, dt, t2)
664 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, t2)
665 call hamiltonian_elec_copy(hm2, hm)
666
667 call propagation_ops_elec_restore_ions(tr%propagation_ops_elec, ions_dyn, ions)
668 if (associated(gfield)) then
669 call propagation_ops_elec_restore_gauge_field(tr%propagation_ops_elec, namespace, space, hm, gr, ext_partners)
670 end if
671
672 op => magnus3_linear_operator_t(namespace, gr, hm1, hm2, m_two*alpha2, m_two*alpha1)
673 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm1, m_half*dt, op=op)
674 safe_deallocate_p(op)
675
676 op => magnus3_linear_operator_t(namespace, gr, hm1, hm2, m_two*alpha1, m_two*alpha2)
677 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm1, m_half*dt, op=op)
678 safe_deallocate_p(op)
679
680 call hamiltonian_elec_end(hm1)
681 call hamiltonian_elec_end(hm2)
682
683 pop_sub(propagator_dt.td_cfmagnus4)
684 end subroutine td_cfmagnus4
685
687
688!! Local Variables:
689!! mode: f90
690!! coding: utf-8
691!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:156
scale a batch by a constant or vector
Definition: batch_ops.F90:164
double sqrt(double __x) __attribute__((__nothrow__
integer, parameter, public imaginary_absorbing
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:241
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines an abstract class for Hamiltonians.
integer, parameter, public term_non_local_potential
integer, parameter, public term_kinetic
subroutine, public zhamiltonian_elec_apply_batch(hm, namespace, mesh, psib, hpsib, terms, set_bc)
subroutine, public hamiltonian_elec_end(hm)
subroutine, public hamiltonian_elec_copy(hm_out, hm_in)
Deep-copy a hamiltonian_elec_t snapshot.
subroutine, public zhamiltonian_elec_external(this, mesh, psib, vpsib)
This module defines classes and functions for interaction partners.
A module to handle KS potential, without the external potential.
integer, parameter, public e_field_electric
Definition: lasers.F90:179
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:179
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1045
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:717
integer, parameter, public e_field_magnetic
Definition: lasers.F90:179
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
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
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public propagation_ops_elec_restore_ions(wo, ions_dyn, ions)
subroutine, public propagation_ops_elec_propagate_gauge_field(wo, gfield, dt, time, save_gf)
subroutine, public propagation_ops_elec_update_hamiltonian(namespace, space, st, mesh, hm, ext_partners, time)
subroutine, public propagation_ops_elec_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, op)
subroutine, public propagation_ops_elec_restore_gauge_field(wo, namespace, space, hm, mesh, ext_partners)
subroutine, public propagation_ops_elec_move_ions(wo, gr, hm, st, namespace, space, ions_dyn, ions, ext_partners, mc, time, dt, save_pos)
subroutine, public propagation_ops_elec_exp_apply(te, namespace, st, mesh, hm, dt, op)
type(magnus3_linear_operator_t) function, pointer magnus3_linear_operator_constructor(namespace, mesh, hm1, hm2, c1, c2)
Build a linear Hamiltonian operator c1*H1 + c2*H2.
subroutine magnus4_operator_apply(this, psib, hpsib)
subroutine magnus3_operator_apply(this, psib, hpsib)
Apply the explicit Magnus-3 operator (weighted stage sum plus commutator term). The equation refers t...
subroutine magnus3_linear_operator_apply(this, psib, hpsib)
Apply the linear Hamiltonian operator for the explicit Magnus-3 stages.
subroutine, public td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
Magnus propagator.
type(magnus4_operator_t) function, pointer magnus4_operator_constructor(namespace, mesh, hm, vmagnus)
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....
type(magnus3_operator_t) function, pointer magnus3_operator_constructor(namespace, mesh, hm1, hm2, hm3, hm4, dt)
Build the explicit Magnus-3 final-stage operator from stage Hamiltonians.
subroutine magnus3_apply_stage_in_ref_phase(namespace, mesh, hm_ref, hm_stage, psib, hpsib)
Apply a stage Hamiltonian in the reference Hamiltonian phase convention.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
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, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:599
Class defining batches of mesh functions.
Definition: batch.F90:161
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
The abstract Hamiltonian class defines a skeleton for specific implementations.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
Class to encapsulate the operation for the explicit Magnus 3 propagator.
Class to encapsulate linear combinations of two Hamiltonians.
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)