Octopus
curv_gygi.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
23
24module curv_gygi_oct_m
26 use debug_oct_m
27 use global_oct_m
28 use, intrinsic :: iso_fortran_env
30 use math_oct_m
33 use parser_oct_m
36 use unit_oct_m
38
39 implicit none
40
41 private
42 public :: &
48
49 type, extends(coordinate_system_t) :: curv_gygi_t
50 private
51 real(real64), public :: A
52 real(real64), public :: alpha
53 real(real64), public :: beta
54 real(real64), allocatable :: pos(:, :)
55 integer :: npos
56 type(root_solver_t) :: rs
57 contains
58 procedure :: to_cartesian => curv_gygi_to_cartesian
59 procedure :: from_cartesian => curv_gygi_from_cartesian
60 procedure :: dvector_from_cartesian => dcurv_gygi_vector_from_cartesian
61 procedure :: zvector_from_cartesian => zcurv_gygi_vector_from_cartesian
62 procedure :: dcovector_to_cartesian => dcurv_gygi_covector_to_cartesian
63 procedure :: zcovector_to_cartesian => zcurv_gygi_covector_to_cartesian
64 procedure :: det_jac => curv_gygi_det_jac
65 procedure :: write_info => curv_gygi_write_info
66 procedure :: surface_element => curv_gygi_surface_element
67 final :: curv_gygi_finalize
68 end type curv_gygi_t
69
70 interface curv_gygi_t
71 procedure curv_gygi_constructor
72 end interface curv_gygi_t
73
74 ! Auxiliary variables for the root solver.
75 class(curv_gygi_t), pointer :: gygi_p
76 integer :: i_p
77 real(real64), allocatable :: chi_p(:)
78
79contains
80
81 ! ---------------------------------------------------------
82 function curv_gygi_constructor(namespace, dim, npos, pos) result(gygi)
83 type(namespace_t), intent(in) :: namespace
84 integer, intent(in) :: dim
85 integer, intent(in) :: npos
86 real(real64), intent(in) :: pos(1:dim,1:npos)
87 class(curv_gygi_t), pointer :: gygi
88
89 push_sub(curv_gygi_constructor)
90
91 safe_allocate(gygi)
92
93 gygi%dim = dim
94 gygi%local_basis = .true.
95 gygi%orthogonal = .true. ! This needs to be checked.
96
97 gygi%npos = npos
98 safe_allocate(gygi%pos(1:dim, 1:gygi%npos))
99 gygi%pos(:, :) = pos(:, :)
100
101 !%Variable CurvGygiA
102 !%Type float
103 !%Default 0.5
104 !%Section Mesh::Curvilinear::Gygi
105 !%Description
106 !% The grid spacing is reduced locally around each atom, and the reduction is
107 !% given by 1/(1+<i>A</i>), where <i>A</i> is specified by this variable. So, if
108 !% <i>A</i>=1/2 (the default), the grid spacing is reduced to two thirds = 1/(1+1/2).
109 !% [This is the <math>A_{\alpha}</math> variable in Eq. 2 of F. Gygi and G. Galli, <i>Phys.
110 !% Rev. B</i> <b>52</b>, R2229 (1995)]. It must be larger than zero.
111 !%End
112 call parse_variable(namespace, 'CurvGygiA', m_half, gygi%A)
113
114 !%Variable CurvGygiAlpha
115 !%Type float
116 !%Default 2.0 a.u.
117 !%Section Mesh::Curvilinear::Gygi
118 !%Description
119 !% This number determines the region over which the grid is enhanced (range of
120 !% enhancement of the resolution). That is, the grid is enhanced on a sphere
121 !% around each atom, whose radius is given by this variable. [This is the <math>a_{\alpha}</math>
122 !% variable in Eq. 2 of F. Gygi and G. Galli, <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)].
123 !% It must be larger than zero.
124 !%End
125
126 call parse_variable(namespace, 'CurvGygiAlpha', m_two, gygi%alpha, units_inp%length)
127 !%Variable CurvGygiBeta
128 !%Type float
129 !%Default 4.0 a.u.
130 !%Section Mesh::Curvilinear::Gygi
131 !%Description
132 !% This number determines the distance over which Euclidean coordinates are
133 !% recovered. [This is the <math>b_{\alpha}</math> variable in Eq. 2 of F. Gygi and G. Galli,
134 !% <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)]. It must be larger than zero.
135 !%End
136 call parse_variable(namespace, 'CurvGygiBeta', m_four, gygi%beta, units_inp%length)
137
138 if (gygi%a <= m_zero) call messages_input_error(namespace, 'CurvGygiA')
139 if (gygi%alpha <= m_zero) call messages_input_error(namespace, 'CurvGygiAlpha')
140 if (gygi%beta <= m_zero) call messages_input_error(namespace, 'CurvGygiBeta')
141
142 gygi%min_mesh_scaling_product = (m_one / (m_one + gygi%A))**gygi%dim
143
144 ! initialize root solver
145 call root_solver_init(gygi%rs, namespace, dim, solver_type = root_newton, maxiter = 500, abs_tolerance = 1.0e-10_real64)
150 ! ---------------------------------------------------------
151 subroutine curv_gygi_copy(this_out, this_in)
152 type(curv_gygi_t), intent(inout) :: this_out
153 type(curv_gygi_t), intent(in) :: this_in
157 this_out%A = this_in%A
158 this_out%alpha = this_in%alpha
159 this_out%beta = this_in%beta
160 safe_allocate_source_a(this_out%pos, this_in%pos)
161 this_out%npos = this_in%npos
162
163 pop_sub(curv_gygi_copy)
164 end subroutine curv_gygi_copy
165
166 ! ---------------------------------------------------------
167 subroutine curv_gygi_finalize(this)
168 type(curv_gygi_t), intent(inout) :: this
171
172 safe_deallocate_a(this%pos)
173
174 pop_sub(curv_gygi_finalize)
175 end subroutine curv_gygi_finalize
176
177 ! ---------------------------------------------------------
178 real(real64) pure function gygi_f(this, r)
179 class(curv_gygi_t), intent(in) :: this
180 real(real64), intent(in) :: r
181
182 ! no PUSH_SUB, called too often
183
184 if (r < 1.0e-4_real64) then
185 ! below 1e-4, the relative deviation of this Taylor expansion from the function is less than 1e-15
186 gygi_f = this%A * (m_one + (-m_one/(m_three*this%alpha**2) - m_one/this%beta**2)*r**2)
187 else
188 gygi_f = this%A * this%alpha/r * tanh(r/this%alpha) * exp(-(r/this%beta)**2)
189 end if
190 end function gygi_f
191
192 ! ---------------------------------------------------------
193 ! absorb division by r into the function to get the correct expansion
194 real(real64) pure function gygi_dfdr_over_r(this, r)
195 class(curv_gygi_t), intent(in) :: this
196 real(real64), intent(in) :: r
197
198 ! no PUSH_SUB, called too often
199
200 if (r < 3.0e-3_real64) then
201 ! below 3e-3, the expansion is closer than 1e-13 to the function and does not suffer from
202 ! rounding errors
203 gygi_dfdr_over_r = this%A * ((-m_two/(m_three*this%alpha**2) - m_two/this%beta**2) + &
204 (8._real64/(15._real64*this%alpha**4) + m_four/(m_three*this%alpha**2*this%beta**2) &
205 + m_two/this%beta**4) * r**2)
206 else
207 gygi_dfdr_over_r = this%A * (this%beta**2*r/cosh(r/this%alpha)**2 - &
208 this%alpha*(this%beta**2+m_two*r**2) * tanh(r/this%alpha)) &
209 / (this%beta**2*r**3) * exp(-(r/this%beta)**2)
210 end if
211 end function gygi_dfdr_over_r
212
213 ! ---------------------------------------------------------
214 ! only a combination of the second derivative is needed, namely
215 ! (d2f/d2r - 1/r df/dr)/r**2
216 real(real64) pure function gygi_d2fdr2_combination(this, r)
217 class(curv_gygi_t), intent(in) :: this
218 real(real64), intent(in) :: r
219
220 ! no PUSH_SUB, called too often
221
222 if (r < 3.0e-3_real64) then
223 ! below 3e-3, the expansion is closer than 1e-5 to the function and does not suffer from
224 ! rounding errors
225 gygi_d2fdr2_combination = this%A * ( &
226 m_two*(8._real64/(15._real64*this%alpha**4) + m_four/(m_three*this%alpha**2*this%beta**2) &
227 + m_two/this%beta**4) &
228 - m_two*(34._real64/(105._real64*this%alpha**6) + m_four/(m_five*this%alpha**4*this%beta**2) &
229 + m_one/(this%alpha**2*this%beta**4) + m_one/this%beta**6) * r**2 &
230 )
231 else
232 gygi_d2fdr2_combination = this%A * exp(-(r/this%beta)**2)*( &
233 -m_two*tanh(r/this%alpha)/(cosh(r/this%alpha)**2*this%alpha*r) &
234 + tanh(r/this%alpha)*(m_three*this%alpha/r**3+m_four*this%alpha/(this%beta**2*r)+m_four*this%alpha*r/this%beta**4) &
235 + m_one/cosh(r/this%alpha)**2*(-m_three/r**2-m_four/this%beta**2))/r**2
236 end if
237 end function gygi_d2fdr2_combination
238
239
240
241 ! ---------------------------------------------------------
242 function curv_gygi_to_cartesian(this, chi) result(xx)
243 class(curv_gygi_t), target, intent(in) :: this
244 real(real64), intent(in) :: chi(:)
245 real(real64) :: xx(1:this%dim)
246
247 integer :: i
248 logical :: conv
249
250 ! no PUSH_SUB, called too often
251
252 gygi_p => this
253 i_p = this%npos
254 safe_allocate(chi_p(1:this%dim))
255 chi_p(:) = chi(:)
256
257 call droot_solver_run(this%rs, getf, xx, conv, startval = chi)
258
259 if (.not. conv) then
260 do i = 1, this%npos
261 conv = .false.
262 i_p = i
263 call droot_solver_run(this%rs, getf, xx, conv, startval = xx(1:this%dim))
264 end do
265 end if
266
267 nullify(gygi_p)
268 safe_deallocate_a(chi_p)
269
270 if (.not. conv) then
271 message(1) = "During the construction of the adaptive grid, the Newton-Raphson"
272 message(2) = "method did not converge for point:"
273 write(message(3),'(9f14.6)') xx(1:this%dim)
274 message(4) = "Try varying the Gygi parameters -- usually reducing CurvGygiA or"
275 message(5) = "CurvGygiAlpha (or both) solves the problem."
276 call messages_fatal(5)
277 end if
278
279 end function curv_gygi_to_cartesian
280
281 ! ---------------------------------------------------------
282 pure function curv_gygi_from_cartesian(this, xx) result(chi)
283 class(curv_gygi_t), target, intent(in) :: this
284 real(real64), intent(in) :: xx(:)
285 real(real64) :: chi(1:this%dim)
286
287 integer :: i, ia
288 real(real64) :: diff(this%dim), f
289
290 ! no PUSH_SUB, called too often
291
292 chi(1:this%dim) = xx(1:this%dim)
293 do ia = 1, this%npos
294 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
295 f = gygi_f(this, norm2(diff))
296 do i = 1, this%dim
297 chi(i) = chi(i) + diff(i) * f
298 end do
299 end do
300
301 end function curv_gygi_from_cartesian
302
303 ! ---------------------------------------------------------
304 real(real64) function curv_gygi_det_jac(this, xx, chi) result(jdet)
305 class(curv_gygi_t), intent(in) :: this
306 real(real64), intent(in) :: xx(:)
307 real(real64), intent(in) :: chi(:)
308
309 real(real64) :: dummy(this%dim)
310 real(real64) :: jac(1:this%dim, 1:this%dim)
311
312 ! no PUSH_SUB, called too often
313
314 call curv_gygi_jacobian(this, xx, dummy, jac)
315 jdet = m_one/lalg_determinant(this%dim, jac, preserve_mat = .false.)
316
317 end function curv_gygi_det_jac
318
319 ! ---------------------------------------------------------
320 subroutine curv_gygi_write_info(this, iunit, namespace)
321 class(curv_gygi_t), intent(in) :: this
322 integer, optional, intent(in) :: iunit
323 type(namespace_t), optional, intent(in) :: namespace
324
325 push_sub(curv_gygi_write_info)
326
327 write(message(1), '(a)') ' Curvilinear Method = gygi'
328 write(message(2), '(a)') ' Gygi Parameters:'
329 write(message(3), '(4x,a,f6.3)') 'A = ', this%a
330 write(message(4), '(4x,3a,f6.3)') 'alpha [', trim(units_abbrev(units_out%length)), '] = ', &
331 units_from_atomic(units_out%length, this%alpha)
332 write(message(5), '(4x,3a,f6.3)') 'beta [', trim(units_abbrev(units_out%length)), '] = ', &
333 units_from_atomic(units_out%length, this%beta)
334 call messages_info(5, iunit=iunit, namespace=namespace)
336 pop_sub(curv_gygi_write_info)
337 end subroutine curv_gygi_write_info
338
339 ! ---------------------------------------------------------
340 real(real64) function curv_gygi_surface_element(this, idir) result(ds)
341 class(curv_gygi_t), intent(in) :: this
342 integer, intent(in) :: idir
343
344 ds = m_zero
345 message(1) = 'Surface element with gygi curvilinear coordinates not implemented'
346 call messages_fatal(1)
347
348 end function curv_gygi_surface_element
349
350 ! ---------------------------------------------------------
351 pure subroutine curv_gygi_jacobian(this, xx, chi, jac, natoms)
352 class(curv_gygi_t), intent(in) :: this
353 real(real64), intent(in) :: xx(:)
354 real(real64), intent(out) :: chi(:)
355 real(real64), intent(out) :: jac(:, :)
356 integer, optional, intent(in) :: natoms
357
358 integer :: i, ix, iy, natoms_
359 real(real64) :: r, diff(this%dim), dfdr_over_r, f
360
361 ! no PUSH_SUB, called too often
362
363 jac(1:this%dim, 1:this%dim) = diagonal_matrix(this%dim, m_one)
364 chi(1:this%dim) = xx(1:this%dim)
365
366 natoms_ = this%npos
367 if (present(natoms)) natoms_ = natoms
368
369 do i = 1, natoms_
370 diff = xx(1:this%dim) - this%pos(1:this%dim, i)
371 r = norm2(diff)
372
373 f = gygi_f(this, r)
374 dfdr_over_r = gygi_dfdr_over_r(this, r)
376 do ix = 1, this%dim
377 chi(ix) = chi(ix) + diff(ix) * f
378
379 jac(ix, ix) = jac(ix, ix) + f
380 do iy = 1, this%dim
381 jac(ix, iy) = jac(ix, iy) + diff(ix)*diff(iy)*dfdr_over_r
382 end do
383 end do
384 end do
385
386 end subroutine curv_gygi_jacobian
387
388 ! ---------------------------------------------------------
389 pure subroutine curv_gygi_hessian(this, xx, hessian, natoms)
390 class(curv_gygi_t), intent(in) :: this
391 real(real64), intent(in) :: xx(:)
392 real(real64), intent(out) :: hessian(:, :, :)
393 integer, optional, intent(in) :: natoms
394
395 integer :: ia, i, j, k, natoms_
396 real(real64) :: r, diff(this%dim), dfdr_over_r, d2fdr2_combination
398 ! no PUSH_SUB, called too often
399
400 hessian(1:this%dim, 1:this%dim, 1:this%dim) = m_zero
401
402 natoms_ = this%npos
403 if (present(natoms)) natoms_ = natoms
404
405 do ia = 1, natoms_
406 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
407 r = norm2(diff)
408
409 dfdr_over_r = gygi_dfdr_over_r(this, r)
410 d2fdr2_combination = gygi_d2fdr2_combination(this, r)
411
412 do i = 1, this%dim
413 do j = 1, this%dim
414 do k = 1, this%dim
415 if (i == j) then
416 hessian(i, j, k) = hessian(i, j, k) + diff(k)*dfdr_over_r
417 end if
418 if (i == k) then
419 hessian(i, j, k) = hessian(i, j, k) + diff(j)*dfdr_over_r
420 end if
421 if (j == k) then
422 hessian(i, j, k) = hessian(i, j, k) + diff(i)*dfdr_over_r
423 end if
424 hessian(i, j, k) = hessian(i, j, k) + diff(i)*diff(j)*diff(k)*d2fdr2_combination
425 end do
426 end do
427 end do
428 end do
429
430 end subroutine curv_gygi_hessian
431
432 ! ---------------------------------------------------------
433 pure subroutine curv_gygi_hessian_trace(this, xx, hessian_trace, natoms)
434 class(curv_gygi_t), intent(in) :: this
435 real(real64), intent(in) :: xx(:)
436 real(real64), intent(out) :: hessian_trace(:)
437 integer, optional, intent(in) :: natoms
438
439 integer :: ia, i, k, natoms_
440 real(real64) :: r, diff(this%dim), dfdr_over_r, d2fdr2_combination
441
442 ! no PUSH_SUB, called too often
443
444 hessian_trace(1:this%dim) = m_zero
445
446 natoms_ = this%npos
447 if (present(natoms)) natoms_ = natoms
448
449 do ia = 1, natoms_
450 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
451 r = norm2(diff)
452
453 dfdr_over_r = gygi_dfdr_over_r(this, r)
454 d2fdr2_combination = gygi_d2fdr2_combination(this, r)
455
456 do i = 1, this%dim
457 do k = 1, this%dim
458 if (i == k) then
459 hessian_trace(i) = hessian_trace(i) + m_two*diff(k)*dfdr_over_r
460 end if
461 hessian_trace(i) = hessian_trace(i) + diff(i)*dfdr_over_r
462 hessian_trace(i) = hessian_trace(i) + diff(i)*diff(k)*diff(k)*d2fdr2_combination
463 end do
464 end do
465 end do
466
467 end subroutine curv_gygi_hessian_trace
468
469 ! ---------------------------------------------------------
470 pure subroutine getf(y, f, jf)
471 real(real64), intent(in) :: y(:)
472 real(real64), intent(out) :: f(:), jf(:, :)
473
474 ! no PUSH_SUB, called too often
475
476 call curv_gygi_jacobian(gygi_p, y, f, jf, i_p)
477 f(1:gygi_p%dim) = f(1:gygi_p%dim) - chi_p(1:gygi_p%dim)
478
479 end subroutine getf
480
481#include "undef.F90"
482#include "real.F90"
483#include "curv_gygi_inc.F90"
484
485#include "undef.F90"
486#include "complex.F90"
487#include "curv_gygi_inc.F90"
488
489end module curv_gygi_oct_m
490
491!! Local Variables:
492!! mode: f90
493!! coding: utf-8
494!! End:
double exp(double __x) __attribute__((__nothrow__
double tanh(double __x) __attribute__((__nothrow__
double cosh(double __x) __attribute__((__nothrow__
This module implements the curvilinear coordinates given in F. Gygi and G. Galli, PRB 52 R2229 (1996)...
Definition: curv_gygi.F90:117
pure subroutine, public curv_gygi_jacobian(this, xx, chi, jac, natoms)
Definition: curv_gygi.F90:445
pure subroutine zcurv_gygi_vector_from_cartesian(this, xx, vv, src)
Definition: curv_gygi.F90:755
subroutine curv_gygi_write_info(this, iunit, namespace)
Definition: curv_gygi.F90:414
real(real64), dimension(:), allocatable chi_p
Definition: curv_gygi.F90:170
real(real64) pure function gygi_f(this, r)
Definition: curv_gygi.F90:272
real(real64) pure function gygi_d2fdr2_combination(this, r)
Definition: curv_gygi.F90:310
pure subroutine zcurv_gygi_covector_to_cartesian(this, xx, cv, src)
Definition: curv_gygi.F90:777
real(real64) pure function gygi_dfdr_over_r(this, r)
Definition: curv_gygi.F90:288
pure subroutine getf(y, f, jf)
Definition: curv_gygi.F90:564
real(real64) function curv_gygi_det_jac(this, xx, chi)
Definition: curv_gygi.F90:398
pure subroutine, public curv_gygi_hessian(this, xx, hessian, natoms)
Definition: curv_gygi.F90:483
pure subroutine dcurv_gygi_vector_from_cartesian(this, xx, vv, src)
Definition: curv_gygi.F90:643
pure real(real64) function, dimension(1:this%dim) curv_gygi_from_cartesian(this, xx)
Definition: curv_gygi.F90:376
subroutine, public curv_gygi_copy(this_out, this_in)
Definition: curv_gygi.F90:245
real(real64) function, dimension(1:this%dim) curv_gygi_to_cartesian(this, chi)
Definition: curv_gygi.F90:336
real(real64) function curv_gygi_surface_element(this, idir)
Definition: curv_gygi.F90:434
class(curv_gygi_t) function, pointer curv_gygi_constructor(namespace, dim, npos, pos)
Definition: curv_gygi.F90:176
class(curv_gygi_t), pointer gygi_p
Definition: curv_gygi.F90:168
pure subroutine dcurv_gygi_covector_to_cartesian(this, xx, cv, src)
Definition: curv_gygi.F90:665
subroutine curv_gygi_finalize(this)
Definition: curv_gygi.F90:261
pure subroutine, public curv_gygi_hessian_trace(this, xx, hessian_trace, natoms)
Definition: curv_gygi.F90:527
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
integer, parameter, public root_newton
subroutine, public root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, rel_tolerance, abs_tolerance)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
static double f(double w, void *p)
abstract class to describe coordinate systems
int true(void)