Octopus
perturbation_kdotp.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
24 use comm_oct_m
25 use debug_oct_m
27 use global_oct_m
28 use grid_oct_m
30 use ions_oct_m
32 use math_oct_m
33 use mesh_oct_m
38 use parser_oct_m
42 use space_oct_m
49
50 implicit none
51
52 private
53 public :: &
55
56 type, extends(perturbation_t) :: perturbation_kdotp_t
57 private
58 integer :: vel_method
59 logical :: use_nonlocalpps
60 type(ions_t), pointer :: ions => null()
61 contains
62 procedure :: copy_to => perturbation_kdotp_copy
63 generic :: assignment(=) => copy_to
64 procedure :: info => perturbation_kdotp_info
65 procedure :: dapply => dperturbation_kdotp_apply
66 procedure :: zapply => zperturbation_kdotp_apply
67 procedure :: dapply_order_2 => dperturbation_kdotp_apply_order_2
68 procedure :: zapply_order_2 => zperturbation_kdotp_apply_order_2
71
72 interface perturbation_kdotp_t
73 procedure perturbation_kdotp_constructor
74 end interface perturbation_kdotp_t
75
76
77contains
78
79 ! ---------------------------------------------------------
84 function perturbation_kdotp_constructor(namespace, ions) result(pert)
85 class(perturbation_kdotp_t), pointer :: pert
86 type(namespace_t), intent(in) :: namespace
87 type(ions_t), target, intent(in) :: ions
88
90
91 safe_allocate(pert)
92
93 call perturbation_kdotp_init(pert, namespace, ions)
94
97
98
99 ! --------------------------------------------------------------------
100 subroutine perturbation_kdotp_init(this, namespace, ions)
101 type(perturbation_kdotp_t), intent(out) :: this
102 type(namespace_t), intent(in) :: namespace
103 type(ions_t), target, intent(in) :: ions
104
106
107 this%dir = -1
108 this%dir2 = -1
109
110 this%ions => ions
111
112 !%Variable KdotPUseNonLocalPseudopotential
113 !%Type logical
114 !%Default true
115 !%Section Linear Response::KdotP
116 !%Description
117 !% For testing purposes, set to false to ignore the term <math>-i \left[\vec{r}, V\right]</math> in
118 !% the <math>\vec{k} \cdot \vec{p}</math> perturbation, which is due to non-local pseudopotentials.
119 !%End
120 call messages_obsolete_variable(namespace, 'KdotP_UseNonLocalPseudopotential', 'KdotPUseNonLocalPseudopotential')
121 call parse_variable(namespace, 'KdotPUseNonLocalPseudopotential', .true., this%use_nonlocalpps)
122
123 !%Variable KdotPVelMethod
124 !%Type integer
125 !%Default grad_vel
126 !%Section Linear Response::KdotP
127 !%Description
128 !% Method of velocity calculation.
129 !%Option grad_vel 0
130 !% <math>-i \left(\nabla + \left[r, V_{\rm nl} \right] \right)</math>
131 !%Option hcom_vel 1
132 !% As a commutator of the position operator and Hamiltonian, <math>-i \left[ r, H \right]</math>.
133 !%End
134 call parse_variable(namespace, 'KdotPVelMethod', option__kdotpvelmethod__grad_vel, this%vel_method)
135
137 end subroutine perturbation_kdotp_init
138
139 ! --------------------------------------------------------------------
140 subroutine perturbation_kdotp_finalize(this)
141 type(perturbation_kdotp_t), intent(inout) :: this
142
144
145 nullify(this%ions)
146
148 end subroutine perturbation_kdotp_finalize
150 ! --------------------------------------------------------------------
151 subroutine perturbation_kdotp_copy(this, source)
152 class(perturbation_kdotp_t), intent(out) :: this
153 class(perturbation_kdotp_t), intent(in) :: source
154
157 call perturbation_copy(this, source)
158 this%ions => source%ions
160 this%use_nonlocalpps = source%use_nonlocalpps
161 this%vel_method = source%vel_method
164 end subroutine perturbation_kdotp_copy
165
166 ! --------------------------------------------------------------------
167 subroutine perturbation_kdotp_info(this)
168 class(perturbation_kdotp_t), intent(in) :: this
169
171
172 if (.not. this%use_nonlocalpps) then
173 write(message(1), '(a)') 'Ignoring non-local pseudopotential term.'
174 call messages_info(1)
175 end if
176
178 end subroutine perturbation_kdotp_info
179
180#include "undef.F90"
181#include "real.F90"
182#include "perturbation_kdotp_inc.F90"
183
184#include "undef.F90"
185#include "complex.F90"
186#include "perturbation_kdotp_inc.F90"
187
189
190!! Local Variables:
191!! mode: f90
192!! coding: utf-8
193!! End:
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 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_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
subroutine zperturbation_kdotp_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
d^2/dki dkj (-(1/2) ki kj [ri,[rj,H]]) = for i = j : 1 - [ri,[rj,Vnl]] for i != j : -(1/2) [ri,...
subroutine perturbation_kdotp_copy(this, source)
subroutine zperturbation_kdotp_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...
subroutine dperturbation_kdotp_apply_order_2(this, namespace, space, gr, hm, ik, f_in, f_out)
d^2/dki dkj (-(1/2) ki kj [ri,[rj,H]]) = for i = j : 1 - [ri,[rj,Vnl]] for i != j : -(1/2) [ri,...
class(perturbation_kdotp_t) function, pointer perturbation_kdotp_constructor(namespace, ions)
The factory routine (or constructor) allocates a pointer of the corresponding type and then calls the...
subroutine perturbation_kdotp_init(this, namespace, ions)
subroutine perturbation_kdotp_finalize(this)
subroutine perturbation_kdotp_info(this)
subroutine dperturbation_kdotp_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...
subroutine, public perturbation_copy(this, source)
This module handles spin dimensions of the states and the k-point distribution.
int true(void)