Octopus
rkb_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
22 use debug_oct_m
23 use global_oct_m
24 use mesh_oct_m
26 use mpi_oct_m
29 use ps_oct_m
32
33 implicit none
34
35 private
36 public :: &
43
48 private
49 integer :: n_s
50 integer, public :: n_c
51 complex(real64), allocatable :: bra(:, :)
52 complex(real64), public, allocatable :: ket(:, :, :, :)
53 real(real64), public, allocatable :: f(:, :, :)
54 end type rkb_projector_t
55
56
57contains
58
59 ! ---------------------------------------------------------
60 subroutine rkb_projector_init(rkb_p, sm, spec, l, lm, so_strength)
61 type(rkb_projector_t), intent(inout) :: rkb_p
62 type(submesh_t), intent(in) :: sm
63 class(pseudopotential_t), intent(in) :: spec
64 integer, intent(in) :: l, lm
65 real(real64), intent(in) :: so_strength
66
67 integer :: i
68
69 push_sub(rkb_projector_init)
70
71 assert(l >= 0)
72 assert(lm >= -l .and. lm <= l)
73
74 rkb_p%n_s = sm%np
75 rkb_p%n_c = spec%ps%projectors_per_l(l+1)
76 assert(mod(rkb_p%n_c, 2) == 0)
77 assert(spec%ps%relativistic_treatment == proj_j_dependent)
78
79 !Allocate memory
80 safe_allocate(rkb_p%bra(1:rkb_p%n_s, 1:rkb_p%n_c))
81 safe_allocate(rkb_p%ket(1:rkb_p%n_s, 1:rkb_p%n_c, 1:2, 1:2))
82 safe_allocate(rkb_p%f(1:rkb_p%n_c, 1:2, 1:2))
83
84 ! Build projectors
85 ! He only need P^{m,m}, P^{m+1,m}, and P^{m-1,m}, for l+1/2 and for l-1/2.
86 ! If we use Vanderbilt-type multi-projector per l, we repeat this for each pair of (l+1/2,l-1/2) projectors.
87 ! For instead ONCVPSP uses 2 projector per l, so rkb_p%n_c is 4
88 do i = 1, rkb_p%n_c
89 call pseudopotential_nl_projector(spec, rkb_p%n_s, sm%rel_x, sm%r, l, lm, i, rkb_p%ket(:, i, 1, 1))
90 rkb_p%bra(:, i) = conjg(rkb_p%ket(:, i, 1, 1))
91 rkb_p%ket(:, i, 2, 2) = rkb_p%ket(:, i, 1, 1)
92
93 if (lm /= l) then
94 call pseudopotential_nl_projector(spec, rkb_p%n_s, sm%rel_x, sm%r, l, lm+1, i, rkb_p%ket(:, i, 2, 1))
95 else
96 rkb_p%ket(:, i, 2, 1) = m_z0
97 end if
98 if (lm /= -l) then
99 call pseudopotential_nl_projector(spec, rkb_p%n_s, sm%rel_x, sm%r, l, lm-1, i, rkb_p%ket(:, i, 1, 2))
100 else
101 rkb_p%ket(:, i, 1, 2) = m_z0
102 end if
103 end do
104
105 ! The l- and m-dependent prefactors are included in the KB energies
106 ! See Eq. B.7 in the thesis of M. Oliveira
107 do i = 0, rkb_p%n_c/2-1
108 rkb_p%f(i*2+1, 1, 1) = real(l + so_strength*lm + 1, real64) ! S_{\up\up}
109 rkb_p%f(i*2+1, 2, 1) = so_strength*sqrt(real((l + lm + 1)*(l - lm), real64)) ! S_{\down\up}
110 rkb_p%f(i*2+1, 1, 2) = so_strength*sqrt(real((l - lm + 1)*(l + lm), real64)) ! S_{\up\down}
111 rkb_p%f(i*2+1, 2, 2) = real(l - so_strength*lm + 1, real64) ! S_{\down\down}
112 rkb_p%f(i*2+2, 1, 1) = real(l - so_strength*lm, real64) ! S_{\up\up}
113 rkb_p%f(i*2+2, 2, 1) = -so_strength*sqrt(real((l + lm + 1)*(l - lm), real64)) ! S_{\down\up}
114 rkb_p%f(i*2+2, 1, 2) = -so_strength*sqrt(real((l - lm + 1)*(l + lm), real64)) ! S_{\up\down}
115 rkb_p%f(i*2+2, 2, 2) = real(l + so_strength*lm, real64) ! S_{\down\down}
116 end do
117 rkb_p%f = rkb_p%f/real(2*l + 1, real64)
118
119 ! Multiply by the KB energies
120 do i = 1, rkb_p%n_c
121 rkb_p%f(i, :, :) = rkb_p%f(i, :, :) * spec%ps%h(l, i, i)
122 end do
123
124 pop_sub(rkb_projector_init)
125 end subroutine rkb_projector_init
126
127 ! ---------------------------------------------------------
128 subroutine rkb_projector_end(rkb_p)
129 type(rkb_projector_t), intent(inout) :: rkb_p
130
131 push_sub(rkb_projector_end)
132
133 safe_deallocate_a(rkb_p%bra)
134 safe_deallocate_a(rkb_p%ket)
135 safe_deallocate_a(rkb_p%f)
136
137 pop_sub(rkb_projector_end)
138 end subroutine rkb_projector_end
139
140 ! ---------------------------------------------------------
141 subroutine rkb_project(mesh, sm, rkb_p, psi, ppsi)
142 type(mesh_t), intent(in) :: mesh
143 type(submesh_t), intent(in) :: sm
144 type(rkb_projector_t), intent(in) :: rkb_p
145 complex(real64), intent(in) :: psi(:, :)
146 complex(real64), intent(inout) :: ppsi(:, :)
147
148 complex(real64) :: uvpsi(1:2, 1:rkb_p%n_c)
149
150 push_sub(rkb_project)
151
152 call rkb_project_bra(mesh, sm, rkb_p, psi, uvpsi)
153 call mesh%allreduce(uvpsi)
154
155 call rkb_project_ket(rkb_p, uvpsi, ppsi)
156
157 pop_sub(rkb_project)
158 end subroutine rkb_project
159
160 ! ---------------------------------------------------------
162 subroutine rkb_project_bra(mesh, sm, rkb_p, psi, uvpsi)
163 type(mesh_t), intent(in) :: mesh
164 type(submesh_t), intent(in) :: sm
165 type(rkb_projector_t), intent(in) :: rkb_p
166 complex(real64), intent(in) :: psi(:, :)
167 complex(real64), intent(out) :: uvpsi(:,:)
168
169 integer :: idim, n_s, is, ic
170
171 complex(real64), allocatable :: bra(:, :)
172
173#ifndef HAVE_OPENMP
174 push_sub(rkb_project_bra)
175#endif
176
177 call profiling_in("RKB_PROJECT_BRA")
178
179 uvpsi(1:2, 1:rkb_p%n_c) = m_zero
180
181 n_s = rkb_p%n_s
182
183 if (mesh%use_curvilinear) then
184 safe_allocate(bra(1:n_s, 1:rkb_p%n_c))
185 do ic = 1, rkb_p%n_c
186 bra(1:n_s, ic) = rkb_p%bra(1:n_s, ic)*mesh%vol_pp(sm%map(1:n_s))
187 do idim = 1, 2
188 do is = 1, n_s
189 uvpsi(idim, ic) = uvpsi(idim, ic) + psi(is, idim)*bra(is, ic)
190 end do
191 end do
192 end do
193 safe_deallocate_a(bra)
194 else
195 do ic = 1, rkb_p%n_c
196 do idim = 1, 2
197 do is = 1, n_s
198 uvpsi(idim, ic) = uvpsi(idim, ic) + psi(is, idim)*rkb_p%bra(is, ic)
199 end do
200 end do
201 end do
202 end if
203
204 uvpsi(1:2, 1:rkb_p%n_c) = uvpsi(1:2, 1:rkb_p%n_c)*mesh%volume_element
205
206 call profiling_out("RKB_PROJECT_BRA")
207
208#ifndef HAVE_OPENMP
209 pop_sub(rkb_project_bra)
210#endif
211 end subroutine rkb_project_bra
212
213 ! ---------------------------------------------------------
215 subroutine rkb_project_ket(rkb_p, uvpsi, psi)
216 type(rkb_projector_t), intent(in) :: rkb_p
217 complex(real64), intent(in) :: uvpsi(:, :)
218 complex(real64), intent(inout) :: psi(:, :)
219
220 integer :: idim, jdim, n_s, is, ic
221 complex(real64) :: aa
222 complex(real64) :: weight(2,2,rkb_p%n_c)
223
224#ifndef HAVE_OPENMP
225 push_sub(rkb_project_ket)
226#endif
227
228 call profiling_in("RKB_PROJECT_KET")
229
230 n_s = rkb_p%n_s
231
232 ! Weight the projectors
233 do jdim = 1, 2
234 do ic = 1, rkb_p%n_c
235 do idim = 1, 2
236 weight(jdim, idim, ic) = rkb_p%f(ic, jdim, idim)*uvpsi(idim, ic)
237 end do
238 end do
239 end do
240
241 do ic = 1, rkb_p%n_c
242 do jdim = 1, 2
243 do is = 1, n_s
244 aa = m_z0
245 do idim = 1, 2
246 aa = aa + weight(jdim, idim, ic)*rkb_p%ket(is, ic, jdim, idim)
247 end do
248 psi(is, jdim) = psi(is, jdim) + aa
249 end do
250 end do
251 end do
252
253 call profiling_out("RKB_PROJECT_KET")
254
255#ifndef HAVE_OPENMP
256 pop_sub(rkb_project_ket)
257#endif
258 end subroutine rkb_project_ket
259
260end module rkb_projector_oct_m
261
262!! Local Variables:
263!! mode: f90
264!! coding: utf-8
265!! End:
double sqrt(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
Definition: ps.F90:114
integer, parameter, public proj_j_dependent
Fully-relativistic j-dependent pseudopotentials.
Definition: ps.F90:173
subroutine, public pseudopotential_nl_projector(spec, np, x, r, l, lm, i, uV)
This routine returns the non-local projector, built using spherical harmonics.
subroutine, public rkb_project_bra(mesh, sm, rkb_p, psi, uvpsi)
THREADSAFE.
subroutine, public rkb_projector_init(rkb_p, sm, spec, l, lm, so_strength)
subroutine, public rkb_project_ket(rkb_p, uvpsi, psi)
THREADSAFE.
subroutine, public rkb_projector_end(rkb_p)
subroutine, public rkb_project(mesh, sm, rkb_p, psi, ppsi)
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The rkb_projector data type holds the KB projectors build with total angular momentum eigenfunctions....
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:175