28 use,
intrinsic :: iso_fortran_env
51 real(real64),
public :: A
52 real(real64),
public :: alpha
53 real(real64),
public :: beta
54 real(real64),
allocatable :: pos(:, :)
56 type(root_solver_t) :: rs
71 procedure curv_gygi_constructor
75 class(curv_gygi_t),
pointer :: gygi_p
77 real(real64),
allocatable :: chi_p(:)
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
94 gygi%local_basis = .
true.
95 gygi%orthogonal = .
true.
98 safe_allocate(gygi%pos(1:dim, 1:gygi%npos))
99 gygi%pos(:, :) = pos(:, :)
142 gygi%min_mesh_scaling_product = (
m_one / (
m_one + gygi%A))**gygi%dim
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
172 safe_deallocate_a(this%pos)
178 real(real64)
pure function
gygi_f(this, r)
180 real(real64),
intent(in) :: r
184 if (r < 1.0e-4_real64)
then
188 gygi_f = this%A * this%alpha/r *
tanh(r/this%alpha) *
exp(-(r/this%beta)**2)
196 real(real64),
intent(in) :: r
200 if (r < 3.0e-3_real64)
then
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)
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)
218 real(real64),
intent(in) :: r
222 if (r < 3.0e-3_real64)
then
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 &
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
244 real(real64),
intent(in) :: chi(:)
245 real(real64) :: xx(1:this%dim)
254 safe_allocate(
chi_p(1:this%dim))
257 call droot_solver_run(this%rs,
getf, xx, conv, startval = chi)
263 call droot_solver_run(this%rs,
getf, xx, conv, startval = xx(1:this%dim))
268 safe_deallocate_a(
chi_p)
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)
284 real(real64),
intent(in) :: xx(:)
285 real(real64) :: chi(1:this%dim)
288 real(real64) :: diff(this%dim),
f
292 chi(1:this%dim) = xx(1:this%dim)
294 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
297 chi(i) = chi(i) + diff(i) *
f
306 real(real64),
intent(in) :: xx(:)
307 real(real64),
intent(in) :: chi(:)
309 real(real64) :: dummy(this%dim)
310 real(real64) :: jac(1:this%dim, 1:this%dim)
315 jdet = m_one/lalg_determinant(this%dim, jac, preserve_mat = .false.)
322 integer,
optional,
intent(in) :: iunit
323 type(namespace_t),
optional,
intent(in) :: namespace
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)
342 integer,
intent(in) :: idir
345 message(1) =
'Surface element with gygi curvilinear coordinates not implemented'
346 call messages_fatal(1)
353 real(real64),
intent(in) :: xx(:)
354 real(real64),
intent(out) :: chi(:)
355 real(real64),
intent(out) :: jac(:, :)
356 integer,
optional,
intent(in) :: natoms
358 integer :: i, ix, iy, natoms_
359 real(real64) :: r, diff(this%dim), dfdr_over_r,
f
363 jac(1:this%dim, 1:this%dim) = diagonal_matrix(this%dim, m_one)
364 chi(1:this%dim) = xx(1:this%dim)
367 if (
present(natoms)) natoms_ = natoms
370 diff = xx(1:this%dim) - this%pos(1:this%dim, i)
377 chi(ix) = chi(ix) + diff(ix) *
f
379 jac(ix, ix) = jac(ix, ix) +
f
381 jac(ix, iy) = jac(ix, iy) + diff(ix)*diff(iy)*dfdr_over_r
391 real(real64),
intent(in) :: xx(:)
392 real(real64),
intent(out) :: hessian(:, :, :)
393 integer,
optional,
intent(in) :: natoms
395 integer :: ia, i, j, k, natoms_
396 real(real64) :: r, diff(this%dim), dfdr_over_r, d2fdr2_combination
400 hessian(1:this%dim, 1:this%dim, 1:this%dim) = m_zero
403 if (
present(natoms)) natoms_ = natoms
406 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
416 hessian(i, j, k) = hessian(i, j, k) + diff(k)*dfdr_over_r
419 hessian(i, j, k) = hessian(i, j, k) + diff(j)*dfdr_over_r
422 hessian(i, j, k) = hessian(i, j, k) + diff(i)*dfdr_over_r
424 hessian(i, j, k) = hessian(i, j, k) + diff(i)*diff(j)*diff(k)*d2fdr2_combination
435 real(real64),
intent(in) :: xx(:)
436 real(real64),
intent(out) :: hessian_trace(:)
437 integer,
optional,
intent(in) :: natoms
439 integer :: ia, i, k, natoms_
440 real(real64) :: r, diff(this%dim), dfdr_over_r, d2fdr2_combination
444 hessian_trace(1:this%dim) = m_zero
447 if (
present(natoms)) natoms_ = natoms
450 diff = xx(1:this%dim) - this%pos(1:this%dim, ia)
459 hessian_trace(i) = hessian_trace(i) + m_two*diff(k)*dfdr_over_r
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
470 pure subroutine getf(y, f, jf)
471 real(real64),
intent(in) :: y(:)
472 real(real64),
intent(out) ::
f(:), jf(:, :)
483#include "curv_gygi_inc.F90"
486#include "complex.F90"
487#include "curv_gygi_inc.F90"
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)...
pure subroutine, public curv_gygi_jacobian(this, xx, chi, jac, natoms)
pure subroutine zcurv_gygi_vector_from_cartesian(this, xx, vv, src)
subroutine curv_gygi_write_info(this, iunit, namespace)
real(real64), dimension(:), allocatable chi_p
real(real64) pure function gygi_f(this, r)
real(real64) pure function gygi_d2fdr2_combination(this, r)
pure subroutine zcurv_gygi_covector_to_cartesian(this, xx, cv, src)
real(real64) pure function gygi_dfdr_over_r(this, r)
pure subroutine getf(y, f, jf)
real(real64) function curv_gygi_det_jac(this, xx, chi)
pure subroutine, public curv_gygi_hessian(this, xx, hessian, natoms)
pure subroutine dcurv_gygi_vector_from_cartesian(this, xx, vv, src)
pure real(real64) function, dimension(1:this%dim) curv_gygi_from_cartesian(this, xx)
subroutine, public curv_gygi_copy(this_out, this_in)
real(real64) function, dimension(1:this%dim) curv_gygi_to_cartesian(this, chi)
real(real64) function curv_gygi_surface_element(this, idir)
class(curv_gygi_t) function, pointer curv_gygi_constructor(namespace, dim, npos, pos)
class(curv_gygi_t), pointer gygi_p
pure subroutine dcurv_gygi_covector_to_cartesian(this, xx, cv, src)
subroutine curv_gygi_finalize(this)
pure subroutine, public curv_gygi_hessian_trace(this, xx, hessian_trace, natoms)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public messages_input_error(namespace, var, details, row, column)
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.
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