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