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
22 use batch_oct_m
23 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
34 use ions_oct_m
35 use, intrinsic :: iso_fortran_env
37 use lasers_oct_m
40 use parser_oct_m
46 use space_oct_m
48 use v_ks_oct_m
49 use xc_oct_m
50
51 implicit none
52
53 private
54
55 public :: &
56 td_magnus, &
58
59contains
60
61 ! ---------------------------------------------------------
63 subroutine td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
64 type(hamiltonian_elec_t), intent(inout) :: hm
65 type(partner_list_t), intent(in) :: ext_partners
66 type(grid_t), intent(in) :: gr
67 type(states_elec_t), intent(inout) :: st
68 type(propagator_base_t), intent(inout) :: tr
69 type(namespace_t), intent(in) :: namespace
70 real(real64), intent(in) :: time
71 real(real64), intent(in) :: dt
72
73 integer :: j, is, i
74 real(real64) :: atime(2)
75 real(real64), allocatable :: vaux(:, :, :), pot(:)
76 type(lasers_t), pointer :: lasers
77
78 push_sub(propagator_dt.td_magnus)
79
80 assert(.not. family_is_mgga_with_exc(hm%xc))
81
82 safe_allocate(vaux(1:gr%np, 1:st%d%nspin, 1:2))
83
84 atime(1) = (m_half-sqrt(m_three)/6.0_real64)*dt
85 atime(2) = (m_half+sqrt(m_three)/6.0_real64)*dt
86
87 if (hm%theory_level /= independent_particles) then
88 do j = 1, 2
89 call potential_interpolation_interpolate(tr%vks_old, 3, time, dt, time - dt + atime(j), vaux(:, :, j))
90 end do
91 else
92 vaux = m_zero
93 end if
94
95 do j = 1, 2
96 ! WARNING: This should be carefully tested, and extended to allow for velocity-gauge laser fields.
97 lasers => list_get_lasers(ext_partners)
98 if(associated(lasers)) then
99 do i = 1, lasers%no_lasers
100 select case (laser_kind(lasers%lasers(i)))
101 case (e_field_electric)
102 safe_allocate(pot(1:gr%np))
103 pot = m_zero
104 call laser_potential(lasers%lasers(i), gr, pot, time - dt + atime(j))
105 do is = 1, st%d%nspin
106 vaux(:, is, j) = vaux(:, is, j) + pot(:)
107 end do
108 safe_deallocate_a(pot)
110 write(message(1),'(a)') 'The Magnus propagator cannot be used with magnetic fields, or'
111 write(message(2),'(a)') 'with an electric field described in the velocity gauge.'
112 call messages_fatal(2, namespace=namespace)
113 end select
114 end do
115 end if
116 end do
117
118 tr%vmagnus(:, :, 2) = m_half*(vaux(:, :, 1) + vaux(:, :, 2))
119 tr%vmagnus(:, :, 1) = (sqrt(m_three)/12.0_real64)*dt*(vaux(:, :, 2) - vaux(:, :, 1))
120
121 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, dt, vmagnus = tr%vmagnus)
122
123 safe_deallocate_a(vaux)
124 pop_sub(propagator_dt.td_magnus)
125 end subroutine td_magnus
126
127
128 ! ---------------------------------------------------------
130 subroutine td_cfmagnus4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, iter)
131 type(v_ks_t), target, intent(inout) :: ks
132 type(namespace_t), intent(in) :: namespace
133 type(electron_space_t), intent(in) :: space
134 type(hamiltonian_elec_t), target, intent(inout) :: hm
135 type(grid_t), target, intent(in) :: gr
136 type(states_elec_t), target, intent(inout) :: st
137 type(propagator_base_t), target, intent(inout) :: tr
138 real(real64), intent(in) :: time
139 real(real64), intent(in) :: dt
140 type(ion_dynamics_t), intent(inout) :: ions_dyn
141 type(ions_t), intent(inout) :: ions
142 type(partner_list_t), intent(in) :: ext_partners
143 integer, intent(in) :: iter
144
145 real(real64) :: alpha1, alpha2, c1, c2, t1, t2
146 real(real64), allocatable :: vhxc1(:, :), vhxc2(:, :)
147
148 if (ion_dynamics_ions_move(ions_dyn) .or. list_has_gauge_field(ext_partners)) then
149 message(1) = "The commutator-free Magnus expansion cannot be used with moving ions or gauge fields"
150 call messages_fatal(1, namespace=namespace)
151 end if
152
153 push_sub(propagator_dt.td_cfmagnus4)
154
155 if (iter < 4) then
156 call td_explicit_runge_kutta4(ks, namespace, space, hm, gr, st, time, dt, ions_dyn, ions, ext_partners)
157 pop_sub(propagator_dt.td_cfmagnus4)
158 return
159 end if
160
161
162 alpha1 = (m_three - m_two * sqrt(m_three))/12.0_real64
163 alpha2 = (m_three + m_two * sqrt(m_three))/12.0_real64
164 c1 = m_half - sqrt(m_three)/6.0_real64
165 c2 = m_half + sqrt(m_three)/6.0_real64
166 t1 = time - dt + c1*dt
167 t2 = time - dt + c2*dt
168
169 safe_allocate(vhxc1(1:gr%np, 1:st%d%nspin))
170 safe_allocate(vhxc2(1:gr%np, 1:st%d%nspin))
171
172 call potential_interpolation_interpolate(tr%vks_old, 4, time, dt, t1, vhxc1)
173 call potential_interpolation_interpolate(tr%vks_old, 4, time, dt, t2, vhxc2)
174
175 hm%ks_pot%vhxc = m_two * (alpha2 * vhxc1 + alpha1 * vhxc2)
176 call hamiltonian_elec_update_with_ext_pot(hm, gr, space, ext_partners, (/ t1, t2 /), (/ m_two * alpha2, m_two * alpha1 /))
177 ! propagate by dt/2
178 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
179
180 hm%ks_pot%vhxc = m_two * (alpha1 * vhxc1 + alpha2 * vhxc2)
181 call hamiltonian_elec_update_with_ext_pot(hm, gr, space, ext_partners, (/ t1, t2 /), (/ m_two * alpha1, m_two * alpha2 /))
182 ! propagate by dt/2
183 !TODO: fuse this with density calc
184 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half * dt)
185
186 call density_calc(st, gr, st%rho)
187
188 safe_deallocate_a(vhxc1)
189 safe_deallocate_a(vhxc2)
190 pop_sub(propagator_dt.td_cfmagnus4)
191 end subroutine td_cfmagnus4
192
194
195!! Local Variables:
196!! mode: f90
197!! coding: utf-8
198!! End:
double sqrt(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
logical function, public list_has_gauge_field(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_three
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public hamiltonian_elec_update_with_ext_pot(this, mesh, space, ext_partners, time, mu)
This is an extension of "hamiltonian_elec_update_pot" to be used by the CFM4 propagator....
This module defines classes and functions for interaction partners.
logical pure function, public ion_dynamics_ions_move(this)
A module to handle KS potential, without the external potential.
integer, parameter, public independent_particles
integer, parameter, public e_field_electric
Definition: lasers.F90:177
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:177
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1043
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:715
integer, parameter, public e_field_magnetic
Definition: lasers.F90:177
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 potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public propagation_ops_elec_exp_apply(te, namespace, st, mesh, hm, dt)
subroutine, public propagation_ops_elec_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, vmagnus)
subroutine, public td_cfmagnus4(ks, namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, iter)
Commutator-free Magnus propagator of order 4.
subroutine, public td_magnus(hm, ext_partners, gr, st, tr, namespace, time, dt)
Magnus propagator.
subroutine, public td_explicit_runge_kutta4(ks, namespace, space, hm, gr, st, time, dt, ions_dyn, ions, ext_partners, qcchi)