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!! Copyright (C) 2025 S. Ohlmann
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
24
25module curv_gygi_oct_m
27 use debug_oct_m
28 use global_oct_m
29 use, intrinsic :: iso_fortran_env
31 use math_oct_m
34 use parser_oct_m
37 use unit_oct_m
39
40 implicit none
41
42 private
43 public :: &
47
48 type, extends(coordinate_system_t) :: curv_gygi_t
49 private
50 real(real64), public :: A
51 real(real64), public :: alpha
52 real(real64), public :: beta
53 real(real64), allocatable :: pos(:, :)
54 integer :: npos
55 type(root_solver_t) :: rs
56 contains
57 procedure :: to_cartesian => curv_gygi_to_cartesian
58 procedure :: from_cartesian => curv_gygi_from_cartesian
59 procedure :: write_info => curv_gygi_write_info
60 procedure :: surface_element => curv_gygi_surface_element
61 procedure :: jacobian => curv_gygi_jacobian
62 procedure :: jacobian_inverse => curv_gygi_jacobian_inverse
63 procedure :: trace_hessian => curv_gygi_trace_hessian
64 final :: curv_gygi_finalize
65 end type curv_gygi_t
66
67 interface curv_gygi_t
68 procedure curv_gygi_constructor
69 end interface curv_gygi_t
70
71 ! Auxiliary variables for the root solver.
72 class(curv_gygi_t), pointer :: gygi_p
73 type(curv_gygi_t), target :: gygi_global
74 real(real64), allocatable :: chi_p(:)
75
76contains
77
78 ! ---------------------------------------------------------
79 function curv_gygi_constructor(namespace, dim, npos, pos) result(gygi)
80 type(namespace_t), intent(in) :: namespace
81 integer, intent(in) :: dim
82 integer, intent(in) :: npos
83 real(real64), intent(in) :: pos(1:dim,1:npos)
84 class(curv_gygi_t), pointer :: gygi
85
86 push_sub(curv_gygi_constructor)
87
88 safe_allocate(gygi)
89
90 gygi%dim = dim
91 gygi%local_basis = .true.
92 gygi%orthogonal = .true. ! This needs to be checked.
93
94 gygi%npos = npos
95 safe_allocate(gygi%pos(1:dim, 1:gygi%npos))
96 gygi%pos(:, :) = pos(:, :)
97
98 !%Variable CurvGygiA
99 !%Type float
100 !%Default 0.5
101 !%Section Mesh::Curvilinear::Gygi
102 !%Description
103 !% The grid spacing is reduced locally around each atom, and the reduction is
104 !% given by 1/(1+<i>A</i>), where <i>A</i> is specified by this variable. So, if
105 !% <i>A</i>=1/2 (the default), the grid spacing is reduced to two thirds = 1/(1+1/2).
106 !% [This is the <math>A_{\alpha}</math> variable in Eq. 2 of F. Gygi and G. Galli, <i>Phys.
107 !% Rev. B</i> <b>52</b>, R2229 (1995)]. It must be larger than zero.
108 !%End
109 call parse_variable(namespace, 'CurvGygiA', m_half, gygi%A)
110
111 !%Variable CurvGygiAlpha
112 !%Type float
113 !%Default 2.0 a.u.
114 !%Section Mesh::Curvilinear::Gygi
115 !%Description
116 !% This number determines the region over which the grid is enhanced (range of
117 !% enhancement of the resolution). That is, the grid is enhanced on a sphere
118 !% around each atom, whose radius is given by this variable. [This is the <math>a_{\alpha}</math>
119 !% variable in Eq. 2 of F. Gygi and G. Galli, <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)].
120 !% It must be larger than zero.
121 !%End
122
123 call parse_variable(namespace, 'CurvGygiAlpha', m_two, gygi%alpha, units_inp%length)
124 !%Variable CurvGygiBeta
125 !%Type float
126 !%Default 4.0 a.u.
127 !%Section Mesh::Curvilinear::Gygi
128 !%Description
129 !% This number determines the distance over which Euclidean coordinates are
130 !% recovered. [This is the <math>b_{\alpha}</math> variable in Eq. 2 of F. Gygi and G. Galli,
131 !% <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)]. It must be larger than zero.
132 !%End
133 call parse_variable(namespace, 'CurvGygiBeta', m_four, gygi%beta, units_inp%length)
134
135 if (gygi%a <= m_zero) call messages_input_error(namespace, 'CurvGygiA')
136 if (gygi%alpha <= m_zero) call messages_input_error(namespace, 'CurvGygiAlpha')
137 if (gygi%beta <= m_zero) call messages_input_error(namespace, 'CurvGygiBeta')
138
139 gygi%min_mesh_scaling_product = (m_one / (m_one + gygi%A))**gygi%dim
140
141 ! initialize root solver
142 call root_solver_init(gygi%rs, namespace, dim, solver_type = root_newton, maxiter = 500, abs_tolerance = 1.0e-10_real64)
147 ! ---------------------------------------------------------
148 subroutine curv_gygi_copy(this_out, this_in)
149 type(curv_gygi_t), intent(inout) :: this_out
150 type(curv_gygi_t), intent(in) :: this_in
154 this_out%A = this_in%A
155 this_out%alpha = this_in%alpha
156 this_out%beta = this_in%beta
157 safe_allocate_source_a(this_out%pos, this_in%pos)
158 this_out%pos = this_in%pos
159 this_out%npos = this_in%npos
160 this_out%dim = this_in%dim
161 this_out%local_basis = this_in%local_basis
162 this_out%orthogonal = this_in%orthogonal
163 call root_solver_init(this_out%rs, global_namespace, this_in%dim, solver_type = root_newton, &
164 maxiter = 500, abs_tolerance = 1.0e-10_real64)
167 end subroutine curv_gygi_copy
168
169 ! ---------------------------------------------------------
170 subroutine curv_gygi_finalize(this)
171 type(curv_gygi_t), intent(inout) :: this
173 push_sub(curv_gygi_finalize)
174
175 safe_deallocate_a(this%pos)
176
177 pop_sub(curv_gygi_finalize)
178 end subroutine curv_gygi_finalize
179
180 ! ---------------------------------------------------------
181 real(real64) pure function gygi_f(this, r)
182 class(curv_gygi_t), intent(in) :: this
183 real(real64), intent(in) :: r
184
185 ! no PUSH_SUB, called too often
186
187 if (r < 1.0e-4_real64) then
188 ! below 1e-4, the relative deviation of this Taylor expansion from the function is less than 1e-15
189 gygi_f = this%A * (m_one + (-m_one/(m_three*this%alpha**2) - m_one/this%beta**2)*r**2)
190 else if (-(r/this%beta)**2 <= m_min_exp_arg) then
191 gygi_f = m_zero
192 else
193 gygi_f = this%A * this%alpha/r * tanh(r/this%alpha) * exp(-(r/this%beta)**2)
194 end if
195 end function gygi_f
196
197 ! ---------------------------------------------------------
198 ! absorb division by r into the function to get the correct expansion
199 real(real64) pure function gygi_dfdr_over_r(this, r)
200 class(curv_gygi_t), intent(in) :: this
201 real(real64), intent(in) :: r
202
203 ! no PUSH_SUB, called too often
204
205 if (r < 3.0e-3_real64) then
206 ! below 3e-3, the expansion is closer than 1e-13 to the function and does not suffer from
207 ! rounding errors
208 gygi_dfdr_over_r = this%A * ((-m_two/(m_three*this%alpha**2) - m_two/this%beta**2) + &
209 (8._real64/(15._real64*this%alpha**4) + m_four/(m_three*this%alpha**2*this%beta**2) &
210 + m_two/this%beta**4) * r**2)
211 else if (-(r/this%beta)**2 <= m_min_exp_arg) then
212 gygi_dfdr_over_r = m_zero
213 else
214 gygi_dfdr_over_r = this%A * (this%beta**2*r/cosh(r/this%alpha)**2 - &
215 this%alpha*(this%beta**2+m_two*r**2) * tanh(r/this%alpha)) &
216 / (this%beta**2*r**3) * exp(-(r/this%beta)**2)
217 end if
218 end function gygi_dfdr_over_r
219
220 ! ---------------------------------------------------------
221 ! only a combination of the second derivative is needed, namely
222 ! (d2f/d2r - 1/r df/dr)/r**2
223 real(real64) pure function gygi_d2fdr2_combination(this, r)
224 class(curv_gygi_t), intent(in) :: this
225 real(real64), intent(in) :: r
226
227 ! no PUSH_SUB, called too often
228
229 if (r < 3.0e-3_real64) then
230 ! below 3e-3, the expansion is closer than 1e-5 to the function and does not suffer from
231 ! rounding errors
232 gygi_d2fdr2_combination = this%A * ( &
233 m_two*(8._real64/(15._real64*this%alpha**4) + m_four/(m_three*this%alpha**2*this%beta**2) &
234 + m_two/this%beta**4) &
235 - m_two*(34._real64/(105._real64*this%alpha**6) + m_four/(m_five*this%alpha**4*this%beta**2) &
236 + m_one/(this%alpha**2*this%beta**4) + m_one/this%beta**6) * r**2 &
237 )
238 else if (-(r/this%beta)**2 <= m_min_exp_arg) then
240 else
241 gygi_d2fdr2_combination = this%A * exp(-(r/this%beta)**2)*( &
242 -m_two*tanh(r/this%alpha)/(cosh(r/this%alpha)**2*this%alpha*r) &
243 + 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) &
244 + m_one/cosh(r/this%alpha)**2*(-m_three/r**2-m_four/this%beta**2))/r**2
245 end if
246 end function gygi_d2fdr2_combination
247
248
249 ! ---------------------------------------------------------
250 function curv_gygi_to_cartesian(this, chi) result(xx)
251 class(curv_gygi_t), target, intent(in) :: this
252 real(real64), intent(in) :: chi(:)
253 real(real64) :: xx(1:this%dim), xx_start(1:this%dim)
254
255 logical :: conv
256 integer :: i_conv, n_conv
257
258 ! no PUSH_SUB, called too often
259
260 gygi_p => this
261 safe_allocate(chi_p(1:this%dim))
262 chi_p(:) = chi(:)
264 call droot_solver_run(this%rs, getf, xx, conv, startval = chi)
265 nullify(gygi_p)
266
267 if (.not. conv) then
268 call curv_gygi_copy(gygi_global, this)
270 ! try to converge with decreasing A
271 n_conv = 8
272 do i_conv = 1, n_conv
273 gygi_p%A = gygi_p%A / m_two
274 call droot_solver_run(gygi_p%rs, getf, xx, conv, startval = chi)
275 if (conv) then
276 exit
277 end if
278 end do
279 if (conv) then
280 ! increase A again and take previous xx as starting point
281 n_conv = i_conv
282 do i_conv = n_conv, 1, -1
283 gygi_p%A = gygi_p%A * m_two
284 xx_start = xx
285 call droot_solver_run(gygi_p%rs, getf, xx, conv, startval = xx_start)
286 end do
287 end if
288 nullify(gygi_p)
290 end if
291
292 safe_deallocate_a(chi_p)
293
294 if (.not. conv) then
295 message(1) = "During the construction of the adaptive grid, the Newton-Raphson"
296 message(2) = "method did not converge for point:"
297 write(message(3),'(9f14.6)') xx(1:this%dim)
298 message(4) = "Try varying the Gygi parameters -- usually reducing CurvGygiA or"
299 message(5) = "CurvGygiAlpha (or both) solves the problem."
300 call messages_fatal(5)
301 end if
302
303 end function curv_gygi_to_cartesian
304
305 ! ---------------------------------------------------------
306 pure function curv_gygi_from_cartesian(this, xx) result(chi)
307 class(curv_gygi_t), target, intent(in) :: this
308 real(real64), intent(in) :: xx(:)
309 real(real64) :: chi(1:this%dim)
310
311 integer :: i, ia
312 real(real64) :: diff(this%dim), f
313
314 ! no PUSH_SUB, called too often
315
316 chi(1:this%dim) = xx(1:this%dim)
317 do ia = 1, this%npos
318 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
319 f = gygi_f(this, norm2(diff))
320 do i = 1, this%dim
321 chi(i) = chi(i) + diff(i) * f
322 end do
323 end do
324
325 end function curv_gygi_from_cartesian
326
327 ! ---------------------------------------------------------
328 subroutine curv_gygi_write_info(this, iunit, namespace)
329 class(curv_gygi_t), intent(in) :: this
330 integer, optional, intent(in) :: iunit
331 type(namespace_t), optional, intent(in) :: namespace
332
333 push_sub(curv_gygi_write_info)
334
335 write(message(1), '(a)') ' Curvilinear Method = gygi'
336 write(message(2), '(a)') ' Gygi Parameters:'
337 write(message(3), '(4x,a,f6.3)') 'A = ', this%a
338 write(message(4), '(4x,3a,f6.3)') 'alpha [', trim(units_abbrev(units_out%length)), '] = ', &
339 units_from_atomic(units_out%length, this%alpha)
340 write(message(5), '(4x,3a,f6.3)') 'beta [', trim(units_abbrev(units_out%length)), '] = ', &
341 units_from_atomic(units_out%length, this%beta)
342 call messages_info(5, iunit=iunit, namespace=namespace)
344 pop_sub(curv_gygi_write_info)
345 end subroutine curv_gygi_write_info
346
347 ! ---------------------------------------------------------
348 real(real64) function curv_gygi_surface_element(this, idir) result(ds)
349 class(curv_gygi_t), intent(in) :: this
350 integer, intent(in) :: idir
351
352 ds = m_zero
353 message(1) = 'Surface element with gygi curvilinear coordinates not implemented'
354 call messages_fatal(1)
355
356 end function curv_gygi_surface_element
357
358 ! ---------------------------------------------------------
359 function curv_gygi_jacobian(this, chi) result(jacobian)
360 class(curv_gygi_t), intent(in) :: this
361 real(real64), intent(in) :: chi(:)
362 real(real64) :: jacobian(1:this%dim, 1:this%dim)
363
364 jacobian = this%jacobian_inverse(chi)
365 call lalg_inverse(this%dim, jacobian, "dir")
366 end function curv_gygi_jacobian
367
368 ! ---------------------------------------------------------
369 function curv_gygi_jacobian_inverse(this, chi) result(jacobian_inverse)
370 class(curv_gygi_t), intent(in) :: this
371 real(real64), intent(in) :: chi(:)
372 real(real64) :: jacobian_inverse(1:this%dim, 1:this%dim)
373
374 real(real64) :: xx(1:this%dim)
375
376 ! no PUSH_SUB, called too often
377
378 xx(:) = this%to_cartesian(chi)
379 jacobian_inverse = curv_gygi_jacobian_inverse_cartesian(this, xx)
380
381 end function curv_gygi_jacobian_inverse
382
383 ! ---------------------------------------------------------
384 function curv_gygi_jacobian_inverse_cartesian(this, xx) result(jacobian_inverse)
385 class(curv_gygi_t), intent(in) :: this
386 real(real64), intent(in) :: xx(:)
387 real(real64) :: jacobian_inverse(1:this%dim, 1:this%dim)
388
389 integer :: i, ix, iy
390 real(real64) :: r, diff(1:this%dim), dfdr_over_r, f
391
392 ! no PUSH_SUB, called too often
393
394 jacobian_inverse(1:this%dim, 1:this%dim) = diagonal_matrix(this%dim, m_one)
395
396 do i = 1, this%npos
397 diff = xx(1:this%dim) - this%pos(1:this%dim, i)
398 r = norm2(diff)
400 f = gygi_f(this, r)
401 dfdr_over_r = gygi_dfdr_over_r(this, r)
402
403 do ix = 1, this%dim
404 jacobian_inverse(ix, ix) = jacobian_inverse(ix, ix) + f
405 do iy = 1, this%dim
406 jacobian_inverse(ix, iy) = jacobian_inverse(ix, iy) + diff(ix)*diff(iy)*dfdr_over_r
407 end do
408 end do
409 end do
410
412
413
414 ! ---------------------------------------------------------
415 pure subroutine curv_gygi_hessian(this, xx, hessian, natoms)
416 class(curv_gygi_t), intent(in) :: this
417 real(real64), intent(in) :: xx(:)
418 real(real64), intent(out) :: hessian(:, :, :)
419 integer, optional, intent(in) :: natoms
420
421 integer :: ia, i, j, k, natoms_
422 real(real64) :: r, diff(this%dim), dfdr_over_r, d2fdr2_combination
423
424 ! no PUSH_SUB, called too often
425
426 hessian(1:this%dim, 1:this%dim, 1:this%dim) = m_zero
427
428 natoms_ = this%npos
429 if (present(natoms)) natoms_ = natoms
430
431 do ia = 1, natoms_
432 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
433 r = norm2(diff)
434
435 dfdr_over_r = gygi_dfdr_over_r(this, r)
436 d2fdr2_combination = gygi_d2fdr2_combination(this, r)
437
438 do i = 1, this%dim
439 do j = 1, this%dim
440 do k = 1, this%dim
441 if (i == j) then
442 hessian(i, j, k) = hessian(i, j, k) + diff(k)*dfdr_over_r
443 end if
444 if (i == k) then
445 hessian(i, j, k) = hessian(i, j, k) + diff(j)*dfdr_over_r
446 end if
447 if (j == k) then
448 hessian(i, j, k) = hessian(i, j, k) + diff(i)*dfdr_over_r
449 end if
450 hessian(i, j, k) = hessian(i, j, k) + diff(i)*diff(j)*diff(k)*d2fdr2_combination
451 end do
452 end do
453 end do
454 end do
455
456 end subroutine curv_gygi_hessian
457
458 ! ---------------------------------------------------------
459 function curv_gygi_trace_hessian(this, chi) result(trace_hessian)
460 class(curv_gygi_t), intent(in) :: this
461 real(real64), intent(in) :: chi(:)
462 real(real64) :: trace_hessian(1:this%dim)
463
464 integer :: ia, i, k
465 real(real64) :: r, xx(this%dim), diff(this%dim), dfdr_over_r, d2fdr2_combination
466
467 ! no PUSH_SUB, called too often
468
469 xx(:) = this%to_cartesian(chi)
470 trace_hessian(:) = m_zero
471
472 do ia = 1, this%npos
473 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
474 r = norm2(diff)
475
476 dfdr_over_r = gygi_dfdr_over_r(this, r)
477 d2fdr2_combination = gygi_d2fdr2_combination(this, r)
478
479 do i = 1, this%dim
480 do k = 1, this%dim
481 if (i == k) then
482 trace_hessian(i) = trace_hessian(i) + m_two*diff(k)*dfdr_over_r
483 end if
484 trace_hessian(i) = trace_hessian(i) + diff(i)*dfdr_over_r
485 trace_hessian(i) = trace_hessian(i) + diff(i)*diff(k)*diff(k)*d2fdr2_combination
486 end do
487 end do
488 end do
489
490 end function curv_gygi_trace_hessian
491
492 ! ---------------------------------------------------------
493 subroutine getf(y, f, jf)
494 real(real64), intent(in) :: y(:)
495 real(real64), intent(out) :: f(:), jf(:, :)
496
497 ! no PUSH_SUB, called too often
498
500 f = gygi_p%from_cartesian(y)
501 f(1:gygi_p%dim) = f(1:gygi_p%dim) - chi_p(1:gygi_p%dim)
502
503 end subroutine getf
504
505end module curv_gygi_oct_m
506
507!! Local Variables:
508!! mode: f90
509!! coding: utf-8
510!! 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:118
type(curv_gygi_t), target gygi_global
Definition: curv_gygi.F90:166
subroutine curv_gygi_write_info(this, iunit, namespace)
Definition: curv_gygi.F90:422
real(real64), dimension(:), allocatable chi_p
Definition: curv_gygi.F90:167
real(real64) pure function gygi_f(this, r)
Definition: curv_gygi.F90:275
real(real64) pure function gygi_d2fdr2_combination(this, r)
Definition: curv_gygi.F90:317
real(real64) pure function gygi_dfdr_over_r(this, r)
Definition: curv_gygi.F90:293
pure subroutine, public curv_gygi_hessian(this, xx, hessian, natoms)
Definition: curv_gygi.F90:509
real(real64) function, dimension(1:this%dim, 1:this%dim) curv_gygi_jacobian_inverse_cartesian(this, xx)
Definition: curv_gygi.F90:478
real(real64) function, dimension(1:this%dim) curv_gygi_trace_hessian(this, chi)
Definition: curv_gygi.F90:553
pure real(real64) function, dimension(1:this%dim) curv_gygi_from_cartesian(this, xx)
Definition: curv_gygi.F90:400
subroutine, public curv_gygi_copy(this_out, this_in)
Definition: curv_gygi.F90:242
real(real64) function, dimension(1:this%dim) curv_gygi_to_cartesian(this, chi)
Definition: curv_gygi.F90:344
subroutine getf(y, f, jf)
Definition: curv_gygi.F90:587
real(real64) function curv_gygi_surface_element(this, idir)
Definition: curv_gygi.F90:442
class(curv_gygi_t) function, pointer curv_gygi_constructor(namespace, dim, npos, pos)
Definition: curv_gygi.F90:173
class(curv_gygi_t), pointer gygi_p
Definition: curv_gygi.F90:165
real(real64) function, dimension(1:this%dim, 1:this%dim) curv_gygi_jacobian(this, chi)
Definition: curv_gygi.F90:453
real(real64) function, dimension(1:this%dim, 1:this%dim) curv_gygi_jacobian_inverse(this, chi)
Definition: curv_gygi.F90:463
subroutine curv_gygi_finalize(this)
Definition: curv_gygi.F90:264
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_min_exp_arg
Definition: global.F90:207
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
type(namespace_t), public global_namespace
Definition: namespace.F90:132
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)