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 :: &
71
72 integer, parameter :: MAX_NPROJECTIONS = 4
73 integer, parameter :: MAX_L = 5
74
85
86 type projector_t
87 private
88 integer, public :: type = proj_none
89 integer :: nprojections
90 integer, public :: lmax
91 integer, public :: lloc
92 integer :: nik
93 integer :: reltype
94
95 type(submesh_t), public :: sphere
96
97
100 type(hgh_projector_t), allocatable, public :: hgh_p(:, :)
101 type(kb_projector_t), allocatable, public :: kb_p(:, :)
102 type(rkb_projector_t), allocatable, public :: rkb_p(:, :)
103 complex(real64), allocatable, public :: phase(:, :, :)
104 end type projector_t
105
106contains
107
108 !---------------------------------------------------------
109 logical elemental function projector_is_null(p)
110 type(projector_t), intent(in) :: p
111
112 projector_is_null = (p%type == proj_none)
113 end function projector_is_null
115 !---------------------------------------------------------
116 logical elemental function projector_is(p, type)
117 type(projector_t), intent(in) :: p
118 integer, intent(in) :: type
119 projector_is = (p%type == type)
120 end function projector_is
121
122 !---------------------------------------------------------
123 subroutine projector_init(p, pseudo, namespace, dim, reltype)
124 type(projector_t), intent(inout) :: p
125 type(pseudopotential_t), target, intent(in) :: pseudo
126 type(namespace_t), intent(in) :: namespace
127 integer, intent(in) :: dim
128 integer, intent(in) :: reltype
129
130 type(ps_t), pointer :: ps
131
132 push_sub(projector_init)
133
134 ps => pseudo%ps
135
136 p%reltype = reltype
137 p%lmax = ps%lmax
138
139 if (ps%local) then
140 p%type = proj_none
141 pop_sub(projector_init)
142 return
143 end if
144
145 p%lloc = ps%llocal
146
147 p%type = ps%projector_type
148
149 if (p%type == proj_kb .and. reltype == 1) then
150 if (ps%relativistic_treatment == proj_j_dependent) then
151 p%type = proj_rkb
152 else
153 call messages_write("Spin-orbit coupling for species '"//trim(pseudo%get_label())//" is not available.")
154 call messages_warning(namespace=namespace)
155 end if
156 end if
157
158 select case (p%type)
159 case (proj_kb, proj_rkb)
160 p%nprojections = ps%kbc
161 case (proj_hgh)
162 p%nprojections = 3
163 end select
164
166 end subroutine projector_init
167
168 !---------------------------------------------
169
170 subroutine projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
171 type(projector_t), intent(inout) :: this
172 integer, intent(in) :: dim
173 type(states_elec_dim_t), intent(in) :: std
174 type(boundaries_t), intent(in) :: bnd
175 type(kpoints_t), intent(in) :: kpoints
176 real(real64), optional, allocatable, intent(in) :: vec_pot(:)
177 real(real64), optional, allocatable, intent(in) :: vec_pot_var(:, :)
178
179 integer :: ns, iq, is, ikpoint
180 real(real64) :: kr, kpoint(dim)
181 integer :: nphase, iphase
182 real(real64), allocatable :: diff(:,:)
186 ns = this%sphere%np
187 nphase = 1
188 if (bnd%spiralBC) nphase = 3
189
190 if (.not. allocated(this%phase) .and. ns > 0) then
191 safe_allocate(this%phase(1:ns, 1:nphase, std%kpt%start:std%kpt%end))
192 end if
194 ! Conctruct vectors which translate the submesh point back into the unit cell:
195 ! The positions this%sphere%x can lie outside the unit cell, while
196 ! this%sphere%mesh%x(this%sphere%map(is), 1:ndim) by construction is the periodic image inside the unit cell.
197 ! If a point of the submesh is inside the unit cell, diff(:,is) = 0.
198 safe_allocate(diff(1:dim, 1:ns))
199 !$omp parallel private(ikpoint, kpoint, iphase, is, kr)
200 !$omp do
201 do is = 1, ns
202 diff(:, is) = this%sphere%rel_x(:,is) + this%sphere%center - this%sphere%mesh%x(this%sphere%map(is), :)
203 end do
204
205 do iq = std%kpt%start, std%kpt%end
206 ikpoint = std%get_kpoint_index(iq)
207
208 ! if this fails, it probably means that sb is not compatible with std
209 assert(ikpoint <= kpoints_number(kpoints))
210
211 kpoint = m_zero
212 kpoint(1:dim) = kpoints%get_point(ikpoint)
213
214 do iphase = 1, nphase
215 !$omp do
216 do is = 1, ns
217 ! this is only the correction to the global phase, that can
218 ! appear if the sphere crossed the boundary of the cell. (diff=0 otherwise)
219
220 kr = sum(kpoint(1:dim)*diff(1:dim, is))
221
222 if (present(vec_pot)) then
223 if (allocated(vec_pot)) kr = kr + sum(vec_pot(1:dim)*diff(1:dim, is))
224 end if
225
226 if (present(vec_pot_var)) then
227 if (allocated(vec_pot_var)) kr = kr + sum(vec_pot_var(1:dim, this%sphere%map(is)) &
228 *(this%sphere%rel_x(:, is)+this%sphere%center))
229 end if
230
231 if (bnd%spiralBC .and. iphase > 1) then
232 kr = kr + (2*(iphase-1)-3)*sum(bnd%spiral_q(1:dim)*diff(1:dim, is))
233 end if
234
235 this%phase(is, iphase, iq) = exp(-m_zi*kr)
236 end do
237 !$omp end do nowait
238 end do
239 end do
240 !$omp end parallel
241
242 safe_deallocate_a(diff)
243
244 pop_sub(projector_init_phases)
245
246 end subroutine projector_init_phases
247
248 !---------------------------------------------------------
249 subroutine projector_build(p, ps, so_strength)
250 type(projector_t), intent(inout) :: p
251 class(pseudopotential_t), intent(in) :: ps
252 real(real64), intent(in) :: so_strength
253
254 integer :: ll, mm
255
256 push_sub(projector_build)
257
258 select case (p%type)
259
260 case (proj_hgh)
261 safe_allocate(p%hgh_p(0:p%lmax, -p%lmax:p%lmax))
262 do ll = 0, p%lmax
263 if (ll == p%lloc) cycle
264 do mm = -ll, ll
265 call hgh_projector_init(p%hgh_p(ll, mm), p%sphere, p%reltype, ps, ll, mm, so_strength)
266 end do
267 end do
268
269 case (proj_kb)
270 safe_allocate(p%kb_p(0:p%lmax, -p%lmax:p%lmax))
271 do ll = 0, p%lmax
272 if (ll == p%lloc) cycle
273 do mm = -ll, ll
274 call kb_projector_init(p%kb_p(ll, mm), p%sphere, ps, ll, mm)
275 end do
276 end do
277
278 case (proj_rkb)
279 safe_allocate(p%rkb_p(1:p%lmax, -p%lmax:p%lmax))
280 do ll = 1, p%lmax
281 if (ll == p%lloc) cycle
282 do mm = -ll, ll
283 call rkb_projector_init(p%rkb_p(ll, mm), p%sphere, ps, ll, mm, so_strength)
284 end do
285 end do
286 ! for rkb, l = 0 is a normal kb
287 if (p%lloc /= 0) then
288 safe_allocate(p%kb_p(1, 1))
289 call kb_projector_init(p%kb_p(1, 1), p%sphere, ps, 0, 0)
290 end if
291
292 end select
293
294 pop_sub(projector_build)
295 end subroutine projector_build
296
297 !---------------------------------------------------------
298 subroutine projector_end(p)
299 type(projector_t), intent(inout) :: p
300
301 integer :: ll, mm
302
303 push_sub(projector_end)
304
305 call submesh_end(p%sphere)
306
307 select case (p%type)
308 case (proj_hgh)
309 do ll = 0, p%lmax
310 if (ll == p%lloc) cycle
311 do mm = -ll, ll
312 call hgh_projector_end(p%hgh_p(ll, mm))
313 end do
314 end do
315 safe_deallocate_a(p%hgh_p)
316
317 case (proj_kb)
318 do ll = 0, p%lmax
319 if (ll == p%lloc) cycle
320 do mm = -ll, ll
321 call kb_projector_end(p%kb_p(ll, mm))
322 end do
323 end do
324 safe_deallocate_a(p%kb_p)
325
326 case (proj_rkb)
327 do ll = 1, p%lmax
328 if (ll == p%lloc) cycle
329 do mm = -ll, ll
330 call rkb_projector_end(p%rkb_p(ll, mm))
331 end do
332 end do
333 safe_deallocate_a(p%rkb_p)
334 if (p%lloc /= 0) then
335 call kb_projector_end(p%kb_p(1, 1))
336 safe_deallocate_a(p%kb_p)
337 end if
338
339 end select
340
341 p%type = proj_none
343 safe_deallocate_a(p%phase)
344
345 pop_sub(projector_end)
346 end subroutine projector_end
347
348#include "undef.F90"
349#include "real.F90"
350#include "projector_inc.F90"
351
352#include "undef.F90"
353#include "complex.F90"
354#include "projector_inc.F90"
355
356end module projector_oct_m
357
358
359
360!! Local Variables:
361!! mode: f90
362!! coding: utf-8
363!! 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:847
subroutine, public projector_build(p, ps, so_strength)
Definition: projector.F90:343
subroutine, public dprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:966
logical elemental function, public projector_is(p, type)
Definition: projector.F90:210
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:512
subroutine, public projector_init(p, pseudo, namespace, dim, reltype)
Definition: projector.F90:217
subroutine, public projector_init_phases(this, dim, std, bnd, kpoints, vec_pot, vec_pot_var)
Definition: projector.F90:264
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:547
complex(real64) function, public zprojector_matrix_element(pj, bnd, dim, ik, psia, psib)
zprojector_matrix_element calculates <psia|projector|psib>
Definition: projector.F90:1558
subroutine, public zproject_sphere(mesh, pj, dim, psi, ppsi)
Definition: projector.F90:1635
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:1762
subroutine, public projector_end(p)
Definition: projector.F90:392
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:1223
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:1258
subroutine, public dproject_sphere(mesh, pj, dim, psi, ppsi)
Definition: projector.F90:924
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:1051
subroutine, public zprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
Definition: projector.F90:1677
logical elemental function, public projector_is_null(p)
Definition: projector.F90:203
Definition: ps.F90:114
integer, parameter, public proj_none
Definition: ps.F90:165
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:179