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