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