Octopus
perturbation.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 epot_oct_m
30 use global_oct_m
31 use grid_oct_m
34 use mpi_oct_m
35 use mesh_oct_m
40 use parser_oct_m
42 use space_oct_m
45 use types_oct_m
47
48 implicit none
49
50 private
51 public :: &
54
55 type, abstract :: perturbation_t
56 integer :: dir
57 integer :: dir2
58 contains
59 procedure :: setup_dir => perturbation_setup_dir
60 procedure :: apply_batch => perturbation_apply_batch
61 procedure :: dexpectation_value => dperturbation_expectation_value
62 procedure :: zexpectation_value => zperturbation_expectation_value
63 procedure :: dstates_elec_expectation_value => dperturbation_states_elec_expectation_value
64 procedure :: zstates_elec_expectation_value => zperturbation_states_elec_expectation_value
65 procedure :: dexpectation_density => dperturbation_expectation_density
66 procedure :: zexpectation_density => zperturbation_expectation_density
67 procedure(perturbation_info), deferred :: info
68 procedure(dperturbation_apply), deferred :: dapply
69 procedure(zperturbation_apply), deferred :: zapply
70 procedure(dperturbation_apply_order_2), deferred :: dapply_order_2
71 procedure(zperturbation_apply_order_2), deferred :: zapply_order_2
72 end type perturbation_t
73
74 abstract interface
75 ! ---------------------------------------------------------
76 subroutine dperturbation_apply(this, namespace, space, gr, hm, ik, f_in, f_out, set_bc)
77 import perturbation_t
78 import namespace_t
79 import space_t
80 import grid_t
82 import real64
83 class(perturbation_t), intent(in) :: this
84 type(namespace_t), intent(in) :: namespace
85 class(space_t), intent(in) :: space
86 type(grid_t), intent(in) :: gr
87 type(hamiltonian_elec_t), intent(in) :: hm
88 integer, intent(in) :: ik
89 real(real64), contiguous, intent(in) :: f_in(:, :)
90 real(real64), contiguous, intent(out) :: f_out(:, :)
91 logical, optional, intent(in) :: set_bc
92 end subroutine dperturbation_apply
93
94 ! ---------------------------------------------------------
95 subroutine zperturbation_apply(this, namespace, space, gr, hm, ik, f_in, f_out, set_bc)
96 import perturbation_t
97 import namespace_t
98 import space_t
99 import grid_t
100 import hamiltonian_elec_t
101 import real64
102 class(perturbation_t), intent(in) :: this
103 type(namespace_t), intent(in) :: namespace
104 class(space_t), intent(in) :: space
105 type(grid_t), intent(in) :: gr
106 type(hamiltonian_elec_t), intent(in) :: hm
107 integer, intent(in) :: ik
108 complex(real64), contiguous, intent(in) :: f_in(:, :)
109 complex(real64), contiguous, intent(out) :: f_out(:, :)
110 logical, optional, intent(in) :: set_bc
111 end subroutine zperturbation_apply
112
113 ! ---------------------------------------------------------
114 subroutine dperturbation_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
116 import namespace_t
117 import space_t
118 import grid_t
119 import hamiltonian_elec_t
120 import real64
121 class(perturbation_t), intent(in) :: this
122 type(namespace_t), intent(in) :: namespace
123 class(space_t), intent(in) :: space
124 type(grid_t), intent(in) :: gr
125 type(hamiltonian_elec_t), intent(in) :: hm
126 integer, intent(in) :: ik
127 real(real64), contiguous, intent(in) :: f_in(:, :)
128 real(real64), contiguous, intent(out) :: f_out(:, :)
129 end subroutine dperturbation_apply_order_2
130
131 ! ---------------------------------------------------------
132 subroutine zperturbation_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
133 import perturbation_t
134 import namespace_t
135 import space_t
136 import grid_t
137 import hamiltonian_elec_t
138 import real64
139 class(perturbation_t), intent(in) :: this
140 type(namespace_t), intent(in) :: namespace
141 class(space_t), intent(in) :: space
142 type(grid_t), intent(in) :: gr
143 type(hamiltonian_elec_t), intent(in) :: hm
144 integer, intent(in) :: ik
145 complex(real64), contiguous, intent(in) :: f_in(:, :)
146 complex(real64), contiguous, intent(out) :: f_out(:, :)
147 end subroutine zperturbation_apply_order_2
150 ! ---------------------------------------------------------
151 subroutine perturbation_info(this)
153 class(perturbation_t), intent(in) :: this
154 end subroutine perturbation_info
155 end interface
157contains
159 ! --------------------------------------------------------------------
160 subroutine perturbation_copy(this, source)
161 class(perturbation_t), intent(out) :: this
162 class(perturbation_t), intent(in) :: source
165
166 this%dir = source%dir
167 this%dir2 = source%dir2
168
170 end subroutine perturbation_copy
171
172 ! --------------------------------------------------------------------
173 subroutine perturbation_setup_dir(this, dir, dir2)
174 class(perturbation_t), intent(inout) :: this
175 integer, intent(in) :: dir
176 integer, optional, intent(in) :: dir2
177
178 push_sub(perturbation_setup_dir)
179
180 this%dir = dir
181 if (present(dir2)) this%dir2 = dir2
182
184 end subroutine perturbation_setup_dir
185
186 ! --------------------------------------------------------------------------
187 subroutine perturbation_apply_batch(this, namespace, space, gr, hm, f_in, f_out)
188 class(perturbation_t), intent(in) :: this
189 type(namespace_t), intent(in) :: namespace
190 class(space_t), intent(in) :: space
191 type(grid_t), intent(in) :: gr
192 type(hamiltonian_elec_t), intent(in) :: hm
193 type(wfs_elec_t), intent(in) :: f_in
194 type(wfs_elec_t), intent(inout) :: f_out
195
196 integer :: ist
197 real(real64), allocatable :: dfi(:, :), dfo(:, :)
198 complex(real64), allocatable :: zfi(:, :), zfo(:, :)
199
201
202 assert(f_in%status() == f_out%status())
203
204 if (f_in%type() == type_float) then
205 safe_allocate(dfi(1:gr%np, 1:hm%d%dim))
206 safe_allocate(dfo(1:gr%np, 1:hm%d%dim))
208 do ist = 1, f_in%nst
209 call batch_get_state(f_in, ist, gr%np, dfi)
210 call this%dapply(namespace, space, gr, hm, f_in%ik, dfi, dfo)
211 call batch_set_state(f_out, ist, gr%np, dfo)
212 end do
213
214 safe_deallocate_a(dfi)
215 safe_deallocate_a(dfo)
216
217 else
218
219 safe_allocate(zfi(1:gr%np, 1:hm%d%dim))
220 safe_allocate(zfo(1:gr%np, 1:hm%d%dim))
221
222 do ist = 1, f_in%nst
223 call batch_get_state(f_in, ist, gr%np, zfi)
224 call this%zapply(namespace, space, gr, hm, f_in%ik, zfi, zfo)
225 call batch_set_state(f_out, ist, gr%np, zfo)
226 end do
227
228 safe_deallocate_a(zfi)
229 safe_deallocate_a(zfo)
230
231 end if
232
234 end subroutine perturbation_apply_batch
235
236
237#include "undef.F90"
238#include "real.F90"
239#include "perturbation_inc.F90"
240
241#include "undef.F90"
242#include "complex.F90"
243#include "perturbation_inc.F90"
245end module perturbation_oct_m
246
247!! Local Variables:
248!! mode: f90
249!! coding: utf-8
250!! End:
There are several ways how to call batch_set_state and batch_get_state:
Definition: batch_ops.F90:199
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
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
real(real64) function dperturbation_expectation_value(this, namespace, space, gr, hm, st, psia, psib, perturbation_order)
subroutine zperturbation_expectation_density(this, namespace, space, gr, hm, st, psia, psib, density, perturbation_order)
This routine includes occupations for psib if perturbation_order == 2, correct if used as ....
subroutine perturbation_apply_batch(this, namespace, space, gr, hm, f_in, f_out)
complex(real64) function zperturbation_states_elec_expectation_value(this, namespace, space, gr, hm, st, perturbation_order)
subroutine, public perturbation_copy(this, source)
complex(real64) function zperturbation_expectation_value(this, namespace, space, gr, hm, st, psia, psib, perturbation_order)
subroutine perturbation_setup_dir(this, dir, dir2)
real(real64) function dperturbation_states_elec_expectation_value(this, namespace, space, gr, hm, st, perturbation_order)
subroutine dperturbation_expectation_density(this, namespace, space, gr, hm, st, psia, psib, density, perturbation_order)
This routine includes occupations for psib if perturbation_order == 2, correct if used as ....
This module handles spin dimensions of the states and the k-point distribution.
type(type_t), public type_float
Definition: types.F90:133
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
batches of electronic states
Definition: wfs_elec.F90:138