Octopus
propagator_expmid.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
25 use grid_oct_m
27 use global_oct_m
31 use ions_oct_m
33 use parser_oct_m
37 use space_oct_m
39 use xc_oct_m
40
41 implicit none
42
43 private
44
45 public :: &
47
48contains
49
50 ! ---------------------------------------------------------
52 subroutine exponential_midpoint(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
53 type(hamiltonian_elec_t), target, intent(inout) :: hm
54 type(namespace_t), intent(in) :: namespace
55 type(electron_space_t), intent(in) :: space
56 type(grid_t), target, intent(in) :: gr
57 type(states_elec_t), target, intent(inout) :: st
58 type(propagator_base_t), target, intent(inout) :: tr
59 real(real64), intent(in) :: time
60 real(real64), intent(in) :: dt
61 type(ion_dynamics_t), intent(inout) :: ions_dyn
62 type(ions_t), intent(inout) :: ions
63 type(partner_list_t), intent(in) :: ext_partners
64
65 type(gauge_field_t), pointer :: gfield
66
67 push_sub(propagator_dt.exponential_midpoint)
68
69 ! the half step of this propagator screws with the gauge field kick
70
71 gfield => list_get_gauge_field(ext_partners)
72 if(associated(gfield)) then
73 assert(gauge_field_is_propagated(gfield) .eqv. .false.)
74 end if
75
76 if (hm%theory_level /= independent_particles) then
77 call potential_interpolation_interpolate(tr%vksold, 3, &
78 time, dt, time - dt/m_two, hm%vhxc, vtau = hm%vtau)
79 end if
80
81 !move the ions to time 'time - dt/2'
82 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
83 ext_partners, time - m_half*dt, m_half*dt, save_pos = .true.)
84
85 if(associated(gfield)) then
86 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, &
87 m_half*dt, time, save_gf = .true.)
88 end if
89
90 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time - dt*m_half)
91
92 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, dt)
93
94 !restore to time 'time - dt'
95 call propagation_ops_elec_restore_ions(tr%propagation_ops_elec, ions_dyn, ions)
96
97 call propagation_ops_elec_restore_gauge_field(tr%propagation_ops_elec, namespace, space, hm, gr, ext_partners)
98
99 pop_sub(propagator_dt.exponential_midpoint)
100 end subroutine exponential_midpoint
101! ---------------------------------------------------------
102
104
105!! Local Variables:
106!! mode: f90
107!! coding: utf-8
108!! End:
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
logical pure function, public gauge_field_is_propagated(this)
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_half
Definition: global.F90:193
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter, public independent_particles
This module defines classes and functions for interaction partners.
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_move_ions(wo, gr, hm, st, namespace, space, ions_dyn, ions, ext_partners, time, dt, save_pos)
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_restore_gauge_field(wo, namespace, space, hm, mesh, ext_partners)
subroutine, public propagation_ops_elec_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, vmagnus)
subroutine, public exponential_midpoint(hm, namespace, space, gr, st, tr, time, dt, ions_dyn, ions, ext_partners)
Exponential midpoint.
Definition: xc.F90:114
int true(void)