Octopus
ps_in_grid.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 atomic_oct_m
23 use debug_oct_m
24 use global_oct_m
25 use, intrinsic :: iso_fortran_env
26 use logrid_oct_m
30
31 implicit none
32
33 private
34 public :: &
44
45 type ps_in_grid_t
46 ! Components are public by default
47 type(logrid_t) :: g
48
49 real(real64) :: zval
50
51 integer :: no_l_channels
52 real(real64), allocatable :: vps(:, :)
53 real(real64), allocatable :: KB(:,:)
54 real(real64), allocatable :: dkbcos(:)
55 real(real64), allocatable, private :: dknorm(:)
56 real(real64), allocatable :: kb_radius(:)
57 integer, allocatable :: kb_nrval(:)
58
59 integer :: so_no_l_channels
60 real(real64), allocatable :: so_vps(:,:)
61 real(real64), allocatable, private :: so_KB(:,:)
62 real(real64), allocatable, private :: so_dkbcos(:)
63 real(real64), allocatable, private :: so_dknorm(:)
64 real(real64), allocatable, private :: so_kb_radius(:)
65 integer, allocatable :: so_kb_nrval(:)
66
67 real(real64), allocatable :: vlocal(:)
68 real(real64), allocatable :: rphi(:, :,:)
69
70 logical :: core_corrections
71 real(real64), allocatable :: chcore(:)
72
73 end type ps_in_grid_t
74
75contains
76
77 subroutine ps_in_grid_init(ps, flavor, a, b, nrval, no_l, so_no_l)
78 type(ps_in_grid_t), intent(out) :: ps
79 integer, intent(in) :: flavor, nrval
80 real(real64), intent(in) :: a, b
81 integer, intent(in) :: no_l, so_no_l
82
83 push_sub(ps_in_grid_init)
84
85 ! initialize logarithmic grid
86 call logrid_init(ps%g, flavor, a, b, nrval)
87
88 ! copy data
89 ps% no_l_channels = no_l
90 ps%so_no_l_channels = so_no_l
91
92 ! Allocate some stuff
93 safe_allocate(ps%vps (1:ps%g%nrval, 1:no_l))
94 safe_allocate(ps%chcore (1:ps%g%nrval))
95 safe_allocate(ps%vlocal (1:ps%g%nrval))
96
97 safe_allocate(ps%rphi (1:ps%g%nrval, 1:no_l, 1:3))
98 safe_allocate(ps%KB (1:ps%g%nrval, 1:no_l))
99 safe_allocate(ps%dkbcos(1:no_l))
100 safe_allocate(ps%dknorm(1:no_l))
101 safe_allocate(ps%kb_radius(1:no_l+1))
102 safe_allocate(ps%kb_nrval(1:no_l+1))
103
104 if (so_no_l > 0) then
105 safe_allocate(ps%so_vps(1:ps%g%nrval, 1:so_no_l))
106 safe_allocate(ps%so_KB (1:ps%g%nrval, 1:so_no_l))
107 safe_allocate(ps%so_dkbcos(1:so_no_l))
108 safe_allocate(ps%so_dknorm(1:so_no_l))
109 safe_allocate(ps%so_kb_radius(1:so_no_l))
110 safe_allocate(ps%so_kb_nrval(1:so_no_l))
111 end if
112
113 pop_sub(ps_in_grid_init)
114 end subroutine ps_in_grid_init
115
117 ! ---------------------------------------------------------
118 subroutine ps_in_grid_end(ps)
119 type(ps_in_grid_t), intent(inout) :: ps
120
121 push_sub(ps_in_grid_end)
122
123 safe_deallocate_a(ps%vps)
124 safe_deallocate_a(ps%chcore)
125 safe_deallocate_a(ps%vlocal)
126
127 safe_deallocate_a(ps%rphi)
128 safe_deallocate_a(ps%KB)
129 safe_deallocate_a(ps%dkbcos)
130 safe_deallocate_a(ps%dknorm)
131 safe_deallocate_a(ps%kb_radius)
132 safe_deallocate_a(ps%kb_nrval)
133
134 if (ps%so_no_l_channels > 0) then
135 safe_deallocate_a(ps%so_vps)
136 safe_deallocate_a(ps%so_KB)
137 safe_deallocate_a(ps%so_dkbcos)
138 safe_deallocate_a(ps%so_dknorm)
139 safe_deallocate_a(ps%so_kb_radius)
140 safe_deallocate_a(ps%so_kb_nrval)
141 end if
143 call logrid_end(ps%g)
145 pop_sub(ps_in_grid_end)
146 end subroutine ps_in_grid_end
149 ! ---------------------------------------------------------
150 subroutine ps_in_grid_vlocal(ps, l_loc, rcore, namespace)
151 type(ps_in_grid_t), intent(inout) :: ps
152 integer, intent(in) :: l_loc
153 real(real64), intent(in) :: rcore
154 type(namespace_t), intent(in) :: namespace
156 integer :: ir
157 real(real64) :: a, b, qtot
158 real(real64), allocatable :: rho(:)
161
162 if (l_loc >= 0) then
163 ps%vlocal(:) = ps%vps(:, l_loc+1)
164
165 else if (l_loc == -1) then
166 if (ps%g%flavor /= logrid_psf) then
167 message(1) = "For the moment, Vanderbilt local potentials are only possible with tm grids."
168 call messages_fatal(1, namespace=namespace)
169 end if
170
171 a = 1.82_real64 / rcore
173 safe_allocate(rho(1:ps%g%nrval))
174
175 do ir = 1, ps%g%nrval
176 rho(ir) = exp(-(sinh(a*b*ps%g%rofi(ir)) / sinh(b))**2)
177 rho(ir) = m_four * m_pi * rho(ir) * ps%g%r2ofi(ir)
178 end do
179! do ir =1,2
180 qtot = sum(rho(:)*ps%g%drdi(:))
181 rho(:) = - rho(:)*(ps%zval/qtot)
182
183 call vhrtre(rho, ps%vlocal, ps%g%rofi, ps%g%drdi, ps%g%s, ps%g%nrval, ps%g%a)
184 ps%vlocal(1) = ps%vlocal(2)
185
186 safe_deallocate_a(rho)
187 end if
188
189 pop_sub(ps_in_grid_vlocal)
190 end subroutine ps_in_grid_vlocal
191
192
193 ! ---------------------------------------------------------
196 subroutine ps_in_grid_kb_projectors(ps)
197 type(ps_in_grid_t), intent(inout) :: ps
198
199 integer :: l
200
202
203 do l = 1, ps%no_l_channels
204 ps%KB(2:, l) = (ps%vps(2:, l) - ps%vlocal(2:))*(ps%rphi(2:, l, 1)/ps%g%rofi(2:))*ps%dknorm(l)
205
206 ps%KB(1, l) = first_point_extrapolate(ps%g%rofi, ps%KB(:, l))
207 end do
208
209 do l = 1, ps%so_no_l_channels
210 ps%so_KB(2:, l) = ps%so_vps(2:, l)*(ps%rphi(2:, l, 1)/ps%g%rofi(2:))*ps%dknorm(l)
211
212 ps%so_KB(1, l) = first_point_extrapolate(ps%g%rofi, ps%so_KB(:, l))
213 end do
214
216 end subroutine ps_in_grid_kb_projectors
217
218
219 ! ---------------------------------------------------------
225 subroutine ps_in_grid_kb_cosines(ps, lloc)
226 type(ps_in_grid_t), intent(inout) :: ps
227 integer, intent(in) :: lloc
228
229 integer :: ir, l
230 real(real64) :: dnrm, avgv, vphi
231 real(real64), parameter :: tol = 1.0e-20_real64
232
233 push_sub(ps_in_grid_kb_cosines)
234
235 do l = 1, ps%no_l_channels
236 if (l-1 == lloc) then
237 ps%dkbcos(l) = m_zero
238 ps%dknorm(l) = m_zero
239 cycle
240 end if
241
242 dnrm = m_zero
243 avgv = m_zero
244 do ir = 1, ps%g%nrval
245 vphi = (ps%vps(ir, l) - ps%vlocal(ir))*ps%rphi(ir, l, 1)
246 dnrm = dnrm + vphi*vphi*ps%g%drdi(ir)
247 avgv = avgv + vphi*ps%rphi(ir, l, 1)*ps%g%drdi(ir)
248 end do
249 ps%dkbcos(l) = dnrm/safe_tol(avgv, tol)
250 ps%dknorm(l) = m_one/safe_tol(sqrt(dnrm), tol)
251 end do
252
253 do l = 1, ps%so_no_l_channels
254 dnrm = m_zero
255 avgv = m_zero
256 do ir = 1, ps%g%nrval
257 vphi = ps%so_vps(ir, l)*ps%rphi(ir, l, 1)
258 dnrm = dnrm + vphi*vphi*ps%g%drdi(ir)
259 avgv = avgv + vphi*ps%rphi(ir, l, 1)*ps%g%drdi(ir)
260 end do
261 ps%so_dkbcos(l) = dnrm/safe_tol(avgv, tol)
262 ps%so_dknorm(l) = m_one/safe_tol(sqrt(dnrm), tol)
263 end do
264
265 pop_sub(ps_in_grid_kb_cosines)
266 end subroutine ps_in_grid_kb_cosines
267
268
269 ! ---------------------------------------------------------
270 subroutine ps_in_grid_cutoff_radii(ps, lloc)
271 type(ps_in_grid_t), intent(inout) :: ps
272 integer, intent(in) :: lloc
273
274 integer :: l, ir
275 real(real64) :: dincv, phi
276 real(real64), parameter :: threshold = 1.0e-6_real64
277
279
280 ! local part ....
281 do ir = ps%g%nrval-1, 2, -1
282 dincv = abs(ps%vlocal(ir)*ps%g%rofi(ir) + m_two*ps%zval)
283 if (dincv > threshold) exit
284 end do
285 ps%kb_radius(ps%no_l_channels+1) = ps%g%rofi(ir + 1)
286 ps%kb_nrval(ps%no_l_channels+1) = ir + 1
287
288 ! non-local part....
289 do l = 1, ps%no_l_channels
290 if (l-1 == lloc) then
291 ps%kb_radius(l) = m_zero
292 ps%kb_nrval(l) = 0
293 cycle
294 end if
295
296 do ir = ps%g%nrval-1, 2, -1
297 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%dknorm(l)
298 dincv = abs((ps%vps(ir, l) - ps%vlocal(ir))*phi)
299
300 if (dincv > threshold) exit
301
302 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%dknorm(l)
303 end do
304
305 ps%kb_radius(l) = ps%g%rofi(ir + 1)
306 ps%kb_nrval(l) = ir + 1
307 end do
308
309 ! now the SO part
310 do l = 1, ps%so_no_l_channels
311 do ir = ps%g%nrval-1, 2, -1
312 phi = (ps%rphi(ir, l, 1)/ps%g%rofi(ir))*ps%so_dknorm(l)
313 dincv = abs((ps%so_vps(ir, l))*phi)
314
315 if (dincv > threshold) exit
316 end do
317
318 ps%so_kb_radius(l) = ps%g%rofi(ir + 1)
319 ps%so_kb_nrval(l) = ir + 1
320 end do
321
323 end subroutine ps_in_grid_cutoff_radii
324
325
326 ! ---------------------------------------------------------
328 subroutine ps_in_grid_check_rphi(ps, namespace)
329 type(ps_in_grid_t), intent(in) :: ps
330 type(namespace_t), intent(in) :: namespace
331
332 integer :: l
333 real(real64) :: nrm
334
335 push_sub(ps_in_grid_check_rphi)
336
337 ! checking normalization of the wavefunctions
338 do l = 1, ps%no_l_channels
339 nrm = sqrt(sum(ps%g%drdi(:)*ps%rphi(:, l, 1)**2))
340 nrm = abs(nrm - m_one)
341 if (nrm > 1.0e-5_real64) then
342 write(message(1), '(a,i2,a)') "Eigenstate for l = ", l-1, ' is not normalized.'
343 write(message(2), '(a, f12.6,a)') '(abs(1 - norm) = ', nrm, ')'
344 call messages_warning(2, namespace=namespace)
345 end if
346 end do
347
348 pop_sub(ps_in_grid_check_rphi)
349 end subroutine ps_in_grid_check_rphi
350
351 ! ---------------------------------------------------------
352
353 real(real64) function first_point_extrapolate(x, y, high_order) result(y0)
354 real(real64), intent(in) :: x(:)
355 real(real64), intent(in) :: y(:)
356 logical, optional, intent(in) :: high_order
357
358 real(real64) :: x1, x2, x3
359 real(real64) :: y1, y2, y3
360
361 x1 = x(2) - x(1)
362 x2 = x(3) - x(1)
363 x3 = x(4) - x(1)
364 y1 = y(2)
365 y2 = y(3)
366 y3 = y(4)
367
368 if (optional_default(high_order, .false.)) then
369
370 y0 = y1*x2*x3*(x2 - x3) + y2*x1*x3*(x3 - x1) + y3*x1*x2*(x1 - x2)
371 y0 = y0/((x1 - x2)*(x1 - x3)*(x2 - x3))
372
373 else
374
375 y0 = y1 - (y2 - y1)*x1/(x2 - x1)
376
377 end if
378
379 end function first_point_extrapolate
380
381end module ps_in_grid_oct_m
382
383!! Local Variables:
384!! mode: f90
385!! coding: utf-8
386!! End:
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double sinh(double __x) __attribute__((__nothrow__
subroutine, public vhrtre(rho, v, r, drdi, srdrdi, nr, a)
VHRTRE CONSTRUCTS THE ELECTROSTATIC POTENTIAL DUE TO A SUPPLIED ! ELECTRON DENSITY....
Definition: atomic.F90:603
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_one
Definition: global.F90:192
integer, parameter, public logrid_psf
log grid used in Troullier-Martins code
Definition: logrid.F90:135
subroutine, public logrid_end(grid)
Definition: logrid.F90:213
subroutine, public logrid_init(grid, flavor, aa, bb, nrval)
Definition: logrid.F90:156
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public ps_in_grid_kb_projectors(ps)
KB-projectors kb = (vps - vlocal) |phi> * dknorm.
Definition: ps_in_grid.F90:292
subroutine, public ps_in_grid_end(ps)
Definition: ps_in_grid.F90:214
subroutine, public ps_in_grid_cutoff_radii(ps, lloc)
Definition: ps_in_grid.F90:366
subroutine, public ps_in_grid_init(ps, flavor, a, b, nrval, no_l, so_no_l)
Definition: ps_in_grid.F90:173
subroutine, public ps_in_grid_kb_cosines(ps, lloc)
KB-cosines and KB-norms: dkbcos stores the KB "cosines:" || (v_l - v_local) phi_l ||^2 / < (v_l - v_l...
Definition: ps_in_grid.F90:321
subroutine, public ps_in_grid_vlocal(ps, l_loc, rcore, namespace)
Definition: ps_in_grid.F90:246
real(real64) function, public first_point_extrapolate(x, y, high_order)
Definition: ps_in_grid.F90:449
subroutine, public ps_in_grid_check_rphi(ps, namespace)
checks normalization of the pseudo wavefunctions
Definition: ps_in_grid.F90:424