39 use,
intrinsic :: iso_fortran_env
70 real(real64),
allocatable :: vmagnus(:, :, :)
76 procedure magnus4_operator_constructor
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
91 procedure magnus3_linear_operator_constructor
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
107 procedure magnus3_operator_constructor
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
126 real(real64) :: atime(2)
127 real(real64),
allocatable :: vaux(:, :, :), vmagnus(:, :, :), pot(:)
128 type(lasers_t),
pointer :: lasers
129 class(operator_t),
pointer :: op
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))
152 if(
associated(lasers))
then
153 do i = 1, lasers%no_lasers
156 safe_allocate(pot(1:gr%np))
159 do is = 1, st%d%nspin
160 vaux(:, is, j) = vaux(:, is, j) + pot(:)
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.'
172 vmagnus(:, :, 2) =
m_half*(vaux(:, :, 1) + vaux(:, :, 2))
173 vmagnus(:, :, 1) = (
sqrt(
m_three)/12.0_real64)*dt*(vaux(:, :, 2) - vaux(:, :, 1))
177 safe_deallocate_p(op)
179 safe_deallocate_a(vaux)
180 safe_deallocate_a(vmagnus)
186 class(
mesh_t),
target,
intent(in) :: mesh
188 real(real64),
intent(in) :: vmagnus(:, :, :)
194 call this%init(namespace, mesh, hm)
198 assert(
size(vmagnus, 1) >= mesh%np)
199 assert(
size(vmagnus, 2) == hm%d%nspin)
200 assert(
size(vmagnus, 3) == 2)
202 write(
message(1),
'(a)')
'Magnus operators only implemented for electrons at the moment'
205 safe_allocate_source(this%vmagnus, vmagnus)
212 class(
batch_t),
intent(inout) :: psib
213 class(
batch_t),
intent(inout) :: hpsib
220 select type(hm => this%hm)
229 assert(psib%nst == hpsib%nst)
231 ispin = hm%d%get_spin_index(psib%ik)
233 call hpsib%copy_to(auxpsib, copy_data=.false.)
234 call hpsib%copy_to(aux2psib, copy_data=.false.)
241 call auxpsib%copy_data_to(this%mesh%np, hpsib)
244 call batch_mul(this%mesh%np, hm%abs_boundaries%mf(1:this%mesh%np), psib, aux2psib)
247 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 2), psib, aux2psib)
251 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 1), auxpsib, aux2psib)
255 call batch_mul(this%mesh%np, this%vmagnus(1:this%mesh%np, ispin, 1), psib, auxpsib)
263 message(1) =
"Internal error: imcompatible batch_t in magnus4_operator_apply for argument hpsib."
267 message(1) =
"Internal error: imcompatible batch_t in magnus4_operator_apply for argument psib."
279 class(
mesh_t),
target,
intent(in) :: mesh
282 real(real64),
intent(in) :: c1
283 real(real64),
intent(in) :: c2
289 call this%init(namespace, mesh, hm1)
302 class(
mesh_t),
intent(in) :: mesh
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)
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)
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)
331 class(
batch_t),
intent(inout) :: psib
332 class(
batch_t),
intent(inout) :: hpsib
338 select type (hm => this%hm)
346 assert(psib%nst == hpsib%nst)
348 call hpsib%copy_to(auxpsib, copy_data=.false.)
354 call batch_axpy(this%mesh%np, this%c2, auxpsib, hpsib)
358 message(1) =
"Internal error: imcompatible batch_t in magnus3_linear_operator_apply for argument hpsib."
362 message(1) =
"Internal error: imcompatible batch_t in magnus3_linear_operator_apply for argument psib."
374 class(
mesh_t),
target,
intent(in) :: mesh
379 real(real64),
intent(in) :: dt
385 call this%init(namespace, mesh, hm1)
402 class(
batch_t),
intent(inout) :: psib
403 class(
batch_t),
intent(inout) :: hpsib
405 type(
wfs_elec_t) :: k1psib, k2psib, k3psib, k4psib
406 type(
wfs_elec_t) :: work1psib, work2psib, sum12psib, commpsib
410 select type (hm => this%hm)
418 assert(psib%nst == hpsib%nst)
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.)
438 call k1psib%copy_data_to(this%mesh%np, hpsib)
452 call k1psib%copy_data_to(this%mesh%np, sum12psib)
463 call batch_axpy(this%mesh%np, (
m_zi*this%dt)/12.0_real64, work1psib, hpsib)
474 message(1) =
"Internal error: imcompatible batch_t in magnus3_operator_apply for argument hpsib."
478 message(1) =
"Internal error: imcompatible batch_t in magnus3_operator_apply for argument psib."
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
496 type(
grid_t),
target,
intent(inout) :: gr
499 real(real64),
intent(in) :: time
500 real(real64),
intent(in) :: dt
502 type(
ions_t),
intent(inout) :: ions
527 if (
associated(gfield))
then
532 call v_ks_calc(ks, namespace, space, hm, st2, ions, ext_partners, &
533 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
543 safe_deallocate_p(op)
546 call v_ks_calc(ks, namespace, space, hm, st3, ions, ext_partners, &
547 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
558 ext_partners, mc, time,
m_half*dt)
559 if (
associated(gfield))
then
564 call v_ks_calc(ks, namespace, space, hm, st4, ions, ext_partners, &
565 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
575 safe_deallocate_p(op)
579 if (
associated(gfield))
then
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
618 type(
grid_t),
target,
intent(inout) :: gr
621 real(real64),
intent(in) :: time
622 real(real64),
intent(in) :: dt
624 type(
ions_t),
intent(inout) :: ions
627 integer,
intent(in) :: iter
629 real(real64) :: alpha1, alpha2, c1, c2, t1, t2, dt1, dt2
640 t1 = time - dt + c1*dt
641 t2 = time - dt + c2*dt
648 ext_partners, mc, t1, dt1, save_pos=.
true.)
649 if (
associated(gfield))
then
653 call hm%ks_pot%interpolate_potentials(tr%vks_old, 4, time, dt, t1)
658 ext_partners, mc, t2, dt2)
659 if (
associated(gfield))
then
663 call hm%ks_pot%interpolate_potentials(tr%vks_old, 4, time, dt, t2)
668 if (
associated(gfield))
then
674 safe_deallocate_p(op)
678 safe_deallocate_p(op)
batchified version of the BLAS axpy routine:
scale a batch by a constant or vector
double sqrt(double __x) __attribute__((__nothrow__
integer, parameter, public imaginary_absorbing
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
This module implements a calculator for the density and defines related functions.
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
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
integer, parameter, public independent_particles
Theory level.
complex(real64), parameter, public m_zi
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module implements the underlying real-space grid.
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
integer, parameter, public e_field_vector_potential
subroutine, public laser_potential(laser, mesh, pot, time)
integer pure elemental function, public laser_kind(laser)
integer, parameter, public e_field_magnetic
This module defines the meshes, which are used in Octopus.
subroutine, public messages_not_implemented(feature, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
This module handles the communicators for the various parallelization strategies.
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)
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Class defining batches of mesh functions.
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
The abstract Hamiltonian class defines a skeleton for specific implementations.
Describes mesh distribution to nodes.
Stores all communicators and groups.
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