Octopus
projector.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
21module projector_oct_m
22 use accel_oct_m
23 use batch_oct_m
26 use comm_oct_m
27 use debug_oct_m
28 use global_oct_m
29 use grid_oct_m
31 use ions_oct_m
34 use mesh_oct_m
36 use mpi_oct_m
39 use ps_oct_m
40 use pseudo_oct_m
47
48 implicit none
49
50 private
51 public :: &
69
70 integer, parameter :: MAX_NPROJECTIONS = 4
71 integer, parameter :: MAX_L = 5
72
83
84 type projector_t
85 private
86 integer, public :: type = proj_none
87 integer :: nprojections
88 integer, public :: lmax
89 integer, public :: lloc
90 integer :: nik
91 integer :: reltype
92
93 type(submesh_t), public :: sphere
94
95
98 type(hgh_projector_t), allocatable, public :: hgh_p(:, :)
99 type(kb_projector_t), allocatable, public :: kb_p(:, :)
100 type(rkb_projector_t), allocatable, public :: rkb_p(:, :)
101 complex(real64), allocatable, public :: phase(:, :, :)
102 end type projector_t
103
104contains
105
106 !---------------------------------------------------------
107 logical elemental function projector_is_null(p)
108 type(projector_t), intent(in) :: p
109
110 projector_is_null = (p%type == proj_none)
111 end function projector_is_null
112
113 !---------------------------------------------------------
114 logical elemental function projector_is(p, type)
115 type(projector_t), intent(in) :: p
116 integer, intent(in) :: type
117 projector_is = (p%type == type)
118 end function projector_is
119
120 !---------------------------------------------------------
121 subroutine projector_init(p, pseudo, namespace, dim, reltype)
122 type(projector_t), intent(inout) :: p
123 type(pseudopotential_t), target, intent(in) :: pseudo
124 type(namespace_t), intent(in) :: namespace
125 integer, intent(in) :: dim
126 integer, intent(in) :: reltype
127
128 type(ps_t), pointer :: ps
129
130 push_sub(projector_init)
131
132 ps => pseudo%ps
133
134 p%reltype = reltype
135 p%lmax = ps%lmax
136
137 if (ps%local) then
138 p%type = proj_none
139 pop_sub(projector_init)
140 return
141 end if
142
143 p%lloc = ps%llocal
144
145 p%type = ps%projector_type
146
147 if (p%type == proj_kb .and. reltype == 1) then
148 if (ps%relativistic_treatment == proj_j_dependent) then
149 p%type = proj_rkb
150 else
151 call messages_write("Spin-orbit coupling for species '"//trim(pseudo%get_label())//" is not available.")
152 call messages_warning(namespace=namespace)
153 end if
154 end if
155
156 select case (p%type)
157 case (proj_kb, proj_rkb)
158 p%nprojections = ps%kbc
159 case (proj_hgh)
160 p%nprojections = 3
161 end select
162
164 end subroutine projector_init
165
166 !---------------------------------------------
167
168 subroutine projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
169 type(projector_t), intent(inout) :: this
170 integer, intent(in) :: dim
171 type(states_elec_dim_t), intent(in) :: std
172 type(boundaries_t), intent(in) :: bnd
173 type(kpoints_t), intent(in) :: kpoints
174 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
175 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
176
177 integer :: ns, iq, is, ikpoint
178 real(real64) :: kr, kpoint(dim)
179 integer :: nphase, iphase
180 real(real64), allocatable :: diff(:,:)
184 ns = this%sphere%np
185 nphase = 1
186 if (bnd%spiralBC) nphase = 3
187
188 if (.not. allocated(this%phase) .and. ns > 0) then
189 safe_allocate(this%phase(1:ns, 1:nphase, std%kpt%start:std%kpt%end))
190 end if
192 ! Construct vectors which translate the submesh point back into the unit cell:
193 ! The positions this%sphere%x can lie outside the unit cell, while
194 ! this%sphere%mesh%x(this%sphere%map(is), 1:ndim) by construction is the periodic image inside the unit cell.
195 ! If a point of the submesh is inside the unit cell, diff(:,is) = 0.
196 safe_allocate(diff(1:dim, 1:ns))
197 !$omp parallel private(ikpoint, kpoint, iphase, is, kr)
198 !$omp do
199 do is = 1, ns
200 diff(:, is) = this%sphere%rel_x(:,is) + this%sphere%center - this%sphere%mesh%x(this%sphere%map(is), :)
201 end do
202
203 do iq = std%kpt%start, std%kpt%end
204 ikpoint = std%get_kpoint_index(iq)
205
206 ! if this fails, it probably means that sb is not compatible with std
207 assert(ikpoint <= kpoints_number(kpoints))
208
209 kpoint = m_zero
210 kpoint(1:dim) = kpoints%get_point(ikpoint)
211
212 do iphase = 1, nphase
213 !$omp do
214 do is = 1, ns
215 ! this is only the correction to the global phase, that can
216 ! appear if the sphere crossed the boundary of the cell. (diff=0 otherwise)
217
218 kr = sum(kpoint(1:dim)*diff(1:dim, is))
219
220 if (present(vec_pot)) then
221 if (allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*diff(1:dim, is))
222 end if
223
224 if (present(vec_pot_var)) then
225 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, this%sphere%map(is)) &
226 *(this%sphere%rel_x(:, is)+this%sphere%center))
227 end if
228
229 if (bnd%spiralBC .and. iphase > 1) then
230 kr = kr + (2*(iphase-1)-3)*sum(bnd%spiral_q(1:dim)*diff(1:dim, is))
231 end if
232
233 this%phase(is, iphase, iq) = exp(-m_zi*kr)
234 end do
235 !$omp end do nowait
236 end do
237 end do
238 !$omp end parallel
239
240 safe_deallocate_a(diff)
241
242 pop_sub(projector_init_phases)
243
244 end subroutine projector_init_phases
245
246 !---------------------------------------------------------
247 subroutine projector_build(p, ps, so_strength)
248 type(projector_t), intent(inout) :: p
249 class(pseudopotential_t), intent(in) :: ps
250 real(real64), intent(in) :: so_strength
251
252 integer :: ll, mm
253
254 push_sub(projector_build)
255
256 select case (p%type)
257
258 case (proj_hgh)
259 safe_allocate(p%hgh_p(0:p%lmax, -p%lmax:p%lmax))
260 do ll = 0, p%lmax
261 if (ll == p%lloc) cycle
262 do mm = -ll, ll
263 call hgh_projector_init(p%hgh_p(ll, mm), p%sphere, p%reltype, ps, ll, mm, so_strength)
264 end do
265 end do
266
267 case (proj_kb)
268 safe_allocate(p%kb_p(0:p%lmax, -p%lmax:p%lmax))
269 do ll = 0, p%lmax
270 if (ll == p%lloc) cycle
271 do mm = -ll, ll
272 call kb_projector_init(p%kb_p(ll, mm), p%sphere, ps, ll, mm)
273 end do
274 end do
275
276 case (proj_rkb)
277 safe_allocate(p%rkb_p(1:p%lmax, -p%lmax:p%lmax))
278 do ll = 1, p%lmax
279 if (ll == p%lloc) cycle
280 do mm = -ll, ll
281 call rkb_projector_init(p%rkb_p(ll, mm), p%sphere, ps, ll, mm, so_strength)
282 end do
283 end do
284 ! for rkb, l = 0 is a normal kb
285 if (p%lloc /= 0) then
286 safe_allocate(p%kb_p(1, 1))
287 call kb_projector_init(p%kb_p(1, 1), p%sphere, ps, 0, 0)
288 end if
289
290 end select
291
292 pop_sub(projector_build)
293 end subroutine projector_build
294
295 !---------------------------------------------------------
296 subroutine projector_end(p)
297 type(projector_t), intent(inout) :: p
298
299 integer :: ll, mm
300
301 push_sub(projector_end)
302
303 call submesh_end(p%sphere)
304
305 select case (p%type)
306 case (proj_hgh)
307 do ll = 0, p%lmax
308 if (ll == p%lloc) cycle
309 do mm = -ll, ll
310 call hgh_projector_end(p%hgh_p(ll, mm))
311 end do
312 end do
313 safe_deallocate_a(p%hgh_p)
314
315 case (proj_kb)
316 do ll = 0, p%lmax
317 if (ll == p%lloc) cycle
318 do mm = -ll, ll
319 call kb_projector_end(p%kb_p(ll, mm))
320 end do
321 end do
322 safe_deallocate_a(p%kb_p)
323
324 case (proj_rkb)
325 do ll = 1, p%lmax
326 if (ll == p%lloc) cycle
327 do mm = -ll, ll
328 call rkb_projector_end(p%rkb_p(ll, mm))
329 end do
330 end do
331 safe_deallocate_a(p%rkb_p)
332 if (p%lloc /= 0) then
333 call kb_projector_end(p%kb_p(1, 1))
334 safe_deallocate_a(p%kb_p)
335 end if
336
337 end select
338
339 p%type = proj_none
341 safe_deallocate_a(p%phase)
342
343 pop_sub(projector_end)
344 end subroutine projector_end
345
346#include "undef.F90"
347#include "real.F90"
348#include "projector_inc.F90"
349
350#include "undef.F90"
351#include "complex.F90"
352#include "projector_inc.F90"
353
354end module projector_oct_m
355
356
357
358!! Local Variables:
359!! mode: f90
360!! coding: utf-8
361!! End:
double exp(double __x) __attribute__((__nothrow__
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 implements the underlying real-space grid.
Definition: grid.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
real(real64) function, public dprojector_matrix_element(pj, bnd, dim, ik, psia, psib)
dprojector_matrix_element calculates <psia|projector|psib>
Definition: projector.F90:845
subroutine, public projector_build(p, ps, so_strength)
Definition: projector.F90:341
subroutine, public dprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:964
logical elemental function, public projector_is(p, type)
Definition: projector.F90:208
subroutine, public dproject_psi(mesh, bnd, pj, npj, dim, psi, ppsi, ik)
dproject_psi calculates the action of a projector on the psi wavefunction. The result is summed up to...
Definition: projector.F90:510
subroutine, public projector_init(p, pseudo, namespace, dim, reltype)
Definition: projector.F90:215
subroutine, public projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
Definition: projector.F90:262
subroutine, public dproject_psi_batch(mesh, bnd, pj, npj, dim, psib, ppsib)
To optimize the application of the non-local operator in parallel, the projectors are applied in step...
Definition: projector.F90:545
complex(real64) function, public zprojector_matrix_element(pj, bnd, dim, ik, psia, psib)
zprojector_matrix_element calculates <psia|projector|psib>
Definition: projector.F90:1556
subroutine, public zprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:1760
subroutine, public projector_end(p)
Definition: projector.F90:390
subroutine, public zproject_psi(mesh, bnd, pj, npj, dim, psi, ppsi, ik)
zproject_psi calculates the action of a projector on the psi wavefunction. The result is summed up to...
Definition: projector.F90:1221
subroutine, public zproject_psi_batch(mesh, bnd, pj, npj, dim, psib, ppsib)
To optimize the application of the non-local operator in parallel, the projectors are applied in step...
Definition: projector.F90:1256
subroutine, public dprojector_commute_r_allatoms_alldir(pj, ions, mesh, dim, bnd, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:1049
subroutine, public zprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:1675
logical elemental function, public projector_is_null(p)
Definition: projector.F90:201
Definition: ps.F90:114
integer, parameter, public proj_none
Definition: ps.F90:167
This module handles spin dimensions of the states and the k-point distribution.
The projector data type is intended to hold the local and non-local parts of the pseudopotentials....
Definition: projector.F90:177