Octopus
logrid.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 logrid_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use, intrinsic :: iso_fortran_env
28
29 implicit none
30
31 private
32 public :: &
33 logrid_t, &
35 logrid_end, &
41
42 integer, parameter, public :: &
43 LOGRID_PSF = 1, & !< log grid used in Troullier-Martins code
44 logrid_cpi = 2
45
46 type logrid_t
47 ! Components are public by default
48 integer :: flavor
49
50 real(real64) :: a, b
51 integer :: nrval
52
53 real(real64), allocatable :: rofi(:)
54 real(real64), allocatable :: r2ofi(:)
55 real(real64), allocatable :: drdi(:)
56 real(real64), allocatable :: s(:)
57 end type logrid_t
58
59contains
60
61 ! ---------------------------------------------------------
62 subroutine logrid_init(grid, flavor, aa, bb, nrval)
63 type(logrid_t), intent(out) :: grid
64 integer, intent(in) :: flavor
65 real(real64), intent(in) :: aa, bb
66 integer, intent(in) :: nrval
67
68 real(real64) :: rpb, ea
69 integer :: ir
70
71 push_sub(logrid_init)
72
73 assert(flavor == logrid_psf .or. flavor == logrid_cpi)
74 assert(nrval >= 2)
75 assert(aa > m_zero)
76 assert(bb > m_zero)
77
78 grid%flavor = flavor
79 grid%a = aa
80 grid%b = bb
81 grid%nrval = nrval
82
83 safe_allocate(grid%rofi(1:nrval))
84 safe_allocate(grid%r2ofi(1:nrval))
85 safe_allocate(grid%drdi(1:nrval))
86 safe_allocate(grid%s(1:nrval))
87
88 select case (grid%flavor)
89 case (logrid_psf)
90 rpb = bb
91 ea = exp(aa)
92 do ir = 1, nrval
93 grid%drdi(ir) = aa*rpb
94 rpb = rpb*ea
95 grid%rofi(ir) = bb*(exp(aa*(ir-1)) - m_one)
96 end do
97
98 case (logrid_cpi)
99 grid%rofi(1) = m_zero
100 grid%drdi(1) = m_zero
101
102 rpb = log(aa)
103 grid%rofi(2) = bb
104 grid%drdi(2) = bb*rpb
105 do ir = 3, grid%nrval
106 grid%rofi(ir) = grid%rofi(ir-1)*aa
107 grid%drdi(ir) = grid%rofi(ir)*rpb
108 end do
109 end select
110
111 ! calculate sqrt(drdi)
112 do ir = 1, grid%nrval
113 grid%s(ir) = sqrt(grid%drdi(ir))
114 grid%r2ofi(ir) = grid%rofi(ir)**2
115 end do
116
117 pop_sub(logrid_init)
118 end subroutine logrid_init
119
120
121 ! ---------------------------------------------------------
122 subroutine logrid_end(grid)
123 type(logrid_t), intent(inout) :: grid
124
125 push_sub(logrid_end)
126
127 safe_deallocate_a(grid%rofi)
128 safe_deallocate_a(grid%r2ofi)
129 safe_deallocate_a(grid%drdi)
130 safe_deallocate_a(grid%s)
131
132 pop_sub(logrid_end)
133 end subroutine logrid_end
134
136 ! ---------------------------------------------------------
137 subroutine logrid_copy(grid_in, grid_out)
138 type(logrid_t), intent(in) :: grid_in
139 type(logrid_t), intent(inout) :: grid_out
140
141 push_sub(logrid_copy)
142
143 call logrid_end(grid_out)
145 grid_out%flavor = grid_in%flavor
146 grid_out%a = grid_in%a
147 grid_out%b = grid_in%b
148 grid_out%nrval = grid_in%nrval
150 safe_allocate(grid_out%rofi (1:grid_out%nrval))
151 safe_allocate(grid_out%r2ofi(1:grid_out%nrval))
152 safe_allocate(grid_out%drdi (1:grid_out%nrval))
153 safe_allocate(grid_out%s (1:grid_out%nrval))
154
155 grid_out%rofi(:) = grid_in%rofi(:)
156 grid_out%r2ofi(:) = grid_in%r2ofi(:)
157 grid_out%drdi(:) = grid_in%drdi(:)
158 grid_out%s(:) = grid_in%s(:)
159
160 pop_sub(logrid_copy)
161 end subroutine logrid_copy
162
163
164 ! ---------------------------------------------------------
165 integer function logrid_index(grid, rofi) result(ii)
166 type(logrid_t), intent(in) :: grid
167 real(real64), intent(in) :: rofi
168
169 integer :: ir
170
171 push_sub(logrid_index)
172
173 ii = 0
174 do ir = 1, grid%nrval-1
175
176 if (rofi >= grid%rofi(ir).and.rofi < grid%rofi(ir+1)) then
177 if (abs(rofi-grid%rofi(ir)) < abs(rofi-grid%rofi(ir+1))) then
178 ii = ir
179 else
180 ii = ir + 1
181 end if
182 exit
183 end if
184
185 end do
186
187 pop_sub(logrid_index)
188 end function logrid_index
189
190
191 ! ---------------------------------------------------------
192 subroutine derivative_in_log_grid(grid, ff, dfdr)
193 type(logrid_t), intent(in) :: grid
194 real(real64), intent(in) :: ff(:)
195 real(real64), intent(out) :: dfdr(:)
196
197 integer :: ii
198
199 push_sub(derivative_in_log_grid)
200
201 dfdr(1) = (ff(2) - ff(1))/(grid%rofi(2) - grid%rofi(1))
202 do ii = 2, grid%nrval-1
203 dfdr(ii) = (ff(ii+1) - ff(ii-1))/(grid%rofi(ii+1) - grid%rofi(ii-1))
204 end do
205 dfdr(grid%nrval) = (ff(grid%nrval) - ff(grid%nrval-1))/(grid%rofi(grid%nrval) - grid%rofi(grid%nrval-1))
206
208 end subroutine derivative_in_log_grid
209
210 ! ----------------------------------------------------------
211 real(real64) pure function logrid_radius(grid) result(radius)
212 type(logrid_t), intent(in) :: grid
213
214 radius = grid%rofi(grid%nrval)
215 end function logrid_radius
216
217
218 subroutine logrid_find_parameters(namespace, zz, aa, bb, np)
219 type(namespace_t), intent(in) :: namespace
220 integer, intent(in) :: zz
221 real(real64), intent(out) :: aa, bb
222 integer, intent(out) :: np
223
224 real(real64) :: xmin, xmax, a1, a2, f1, fm
225
226 push_sub(logrid_find_parameters)
227
228 ! Initializes the logarithmic grid.
229 ! Parameters are obtained using the default values for the first non-zero point xmin,
230 ! the last point xmax, and the number of points np
231 ! These values have a default value obtained from the atomic number
232 ! Adapted from APE
233 xmin = sqrt(real(zz, real64) )*1e-5_real64
234 xmax = sqrt(real(zz, real64) )*30.0_real64
235 np = floor(sqrt(real(zz, real64) )*200_real64)
236 ! The code wants np to be an odd number
237 np = floor(np/m_two)*2+1
238
239 a1 = 1e-8_real64
240 f1 = func(xmin, xmax, real(np, real64) , a1)
241 a2 = m_one
242 do
243 aa = (a2 + a1)*m_half
244 fm = func(xmin, xmax, real(np, real64) , aa)
245 if (m_half*abs(a1 - a2) < 1.0e-16_real64) exit
246 if (fm*f1 > m_zero) then
247 a1 = aa
248 f1 = fm
249 else
250 a2 = aa
251 end if
252 end do
253
254 bb = xmin/(exp(aa)-m_one)
255
256 write(message(1), '(a,es13.6,a,es13.6,a,i4)') 'Debug: Log grid parameters: a = ', aa, &
257 ' b = ', bb, ' np = ', np
258 call messages_info(1, namespace=namespace, debug_only=.true.)
259
261 contains
262 real(real64) function func(r1, rn, n, a)
263 real(real64), intent(in) :: r1, rn, a, n
264 if((n-m_one)*a < m_max_exp_arg) then ! To avoid FPE
265 func = exp((n-m_one)*a)*r1 - m_one*r1 - rn*exp(a) + rn*m_one
266 else
267 func = m_huge*r1 - m_one*r1 - rn*exp(a) + rn*m_one
268 end if
269 end function func
270 end subroutine logrid_find_parameters
271
272end module logrid_oct_m
273
274!! Local Variables:
275!! mode: f90
276!! coding: utf-8
277!! End:
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double floor(double __x) __attribute__((__nothrow__
real(real64) function func(r1, rn, n, a)
Definition: logrid.F90:356
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_one
Definition: global.F90:189
subroutine, public logrid_find_parameters(namespace, zz, aa, bb, np)
Definition: logrid.F90:312
real(real64) pure function, public logrid_radius(grid)
Definition: logrid.F90:305
subroutine, public logrid_copy(grid_in, grid_out)
Definition: logrid.F90:231
integer, parameter, public logrid_cpi
log grid used in FHI code
Definition: logrid.F90:135
integer function, public logrid_index(grid, rofi)
Definition: logrid.F90:259
subroutine, public logrid_end(grid)
Definition: logrid.F90:216
subroutine, public derivative_in_log_grid(grid, ff, dfdr)
Definition: logrid.F90:286
subroutine, public logrid_init(grid, flavor, aa, bb, nrval)
Definition: logrid.F90:156
int true(void)