Octopus
perturbation_magnetic.F90
Go to the documentation of this file.
1!! Copyright (C) 2007 X. Andrade
2!! Copyright (C) 2021 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use batch_oct_m
26 use comm_oct_m
27 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
32 use ions_oct_m
35 use lda_u_oct_m
36 use math_oct_m
37 use mesh_oct_m
41 use mpi_oct_m
43 use parser_oct_m
48 use space_oct_m
55 use xc_oct_m
56
57 implicit none
58
59 private
60 public :: &
62
63 integer, public, parameter :: &
64 GAUGE_GIPAW = 1, &
65 gauge_icl = 2
66
68 private
69 integer :: gauge
70 type(ions_t), pointer :: ions => null()
71 contains
72 procedure :: copy_to => perturbation_magnetic_copy
73 generic :: assignment(=) => copy_to
74 procedure :: info => perturbation_magnetic_info
75 procedure :: dapply => dperturbation_magnetic_apply
76 procedure :: zapply => zperturbation_magnetic_apply
77 procedure :: dapply_order_2 => dperturbation_magnetic_apply_order_2
78 procedure :: zapply_order_2 => zperturbation_magnetic_apply_order_2
81
83 procedure perturbation_magnetic_constructor
84 end interface perturbation_magnetic_t
85
86
87contains
88 ! ---------------------------------------------------------
93 function perturbation_magnetic_constructor(namespace, ions) result(pert)
94 class(perturbation_magnetic_t), pointer :: pert
95 type(namespace_t), intent(in) :: namespace
96 type(ions_t), target, intent(in) :: ions
97
99
100 safe_allocate(pert)
101
102 call perturbation_magnetic_init(pert, namespace, ions)
103
106
107
108 ! --------------------------------------------------------------------
109 subroutine perturbation_magnetic_init(this, namespace, ions)
110 type(perturbation_magnetic_t), intent(out) :: this
111 type(namespace_t), intent(in) :: namespace
112 type(ions_t), target, intent(in) :: ions
113
116 this%dir = -1
117 this%dir2 = -1
118
119 this%ions => ions
120
121 !%Variable MagneticGaugeCorrection
122 !%Type integer
123 !%Default gipaw
124 !%Section Linear Response
125 !%Description
126 !% For magnetic linear response: how to handle gauge-invariance in the description
127 !% of the coupling of electrons to the magnetic field.
128 !%Option none 0
129 !% No correction.
130 !%Option gipaw 1
131 !% GIPAW correction: C Pickard and F Mauri, <i>Phys. Rev. Lett.</i> <b>91</b>, 196401 (2003).
132 !%Option icl 2
133 !% ICL correction: S Ismail-Beigi, EK Chang, and SG Louie, <i>Phys. Rev. Lett.</i> <b>87</b>, 087402 (2001).
134 !%End
135
136 call parse_variable(namespace, 'MagneticGaugeCorrection', gauge_gipaw, this%gauge)
137 if (.not. varinfo_valid_option('MagneticGaugeCorrection', this%gauge)) then
138 call messages_input_error(namespace, 'MagneticGaugeCorrection')
139 end if
140
142 end subroutine perturbation_magnetic_init
143
144 ! --------------------------------------------------------------------
145 subroutine perturbation_magnetic_finalize(this)
146 type(perturbation_magnetic_t), intent(inout) :: this
147
149
151 end subroutine perturbation_magnetic_finalize
152
153 ! --------------------------------------------------------------------
154 subroutine perturbation_magnetic_copy(this, source)
155 class(perturbation_magnetic_t), intent(out) :: this
156 class(perturbation_magnetic_t), intent(in) :: source
157
159
160 call perturbation_copy(this, source)
161 this%ions => source%ions
163 this%gauge = source%gauge
164
168 ! --------------------------------------------------------------------
170 class(perturbation_magnetic_t), intent(in) :: this
173
175 end subroutine perturbation_magnetic_info
176
177
178#include "undef.F90"
179#include "real.F90"
180#include "perturbation_magnetic_inc.F90"
181
182#include "undef.F90"
183#include "complex.F90"
184#include "perturbation_magnetic_inc.F90"
185
187
188!! Local Variables:
189!! mode: f90
190!! coding: utf-8
191!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
This module implements the underlying real-space grid.
Definition: grid.F90:117
A module to handle KS potential, without the external potential.
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
subroutine zperturbation_magnetic_apply(this, namespace, space, gr, hm, ik, f_in, f_out, set_bc)
Returns f_out = H' f_in, where H' is perturbation Hamiltonian Note that e^ikr phase is applied to f_i...
class(perturbation_magnetic_t) function, pointer perturbation_magnetic_constructor(namespace, ions)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
subroutine zperturbation_magnetic_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
subroutine dperturbation_magnetic_apply(this, namespace, space, gr, hm, ik, f_in, f_out, set_bc)
Returns f_out = H' f_in, where H' is perturbation Hamiltonian Note that e^ikr phase is applied to f_i...
integer, parameter, public gauge_icl
subroutine dperturbation_magnetic_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
subroutine perturbation_magnetic_copy(this, source)
subroutine perturbation_magnetic_finalize(this)
subroutine perturbation_magnetic_init(this, namespace, ions)
subroutine, public perturbation_copy(this, source)
This module handles spin dimensions of the states and the k-point distribution.