32 use,
intrinsic :: iso_fortran_env
51 real(real64),
allocatable :: lsize(:)
55 real(real64),
allocatable :: Jlocal(:)
56 real(real64),
allocatable :: Jrange(:)
58 real(real64),
allocatable :: chi_atoms(:,:)
59 real(real64),
allocatable :: csi(:,:)
62 type(root_solver_t) :: rs
73 procedure curv_modine_constructor
77 class(curv_modine_t),
pointer :: modine_p
78 real(real64),
allocatable :: x_p(:)
84 type(namespace_t),
intent(in) :: namespace
85 integer,
intent(in) :: dim
86 integer,
intent(in) :: npos
87 real(real64),
intent(in) :: pos(1:dim,1:npos)
88 real(real64),
intent(in) :: lsize(1:dim)
89 real(real64),
intent(in) :: spacing(1:dim)
90 class(curv_modine_t),
pointer :: modine
97 modine%local_basis = .
true.
98 modine%orthogonal = .
true.
122 safe_allocate(modine%lsize(1:dim))
123 modine%lsize = lsize / modine%Jbar
125 if (modine%xbar <
m_zero .or. modine%xbar >
m_one)
then
126 message(1) =
'The parameter "CurvModineXBar" must lie between 0 and 1.'
130 safe_allocate(modine%Jlocal(1:npos))
131 safe_allocate(modine%Jrange(1:npos))
143 call parse_variable(namespace,
'CurvModineJlocal', 0.25_real64, modine%Jlocal(1))
156 message(1) =
'The parameter "CurvModineJlocal" must lie between 0 and 1.'
160 modine%Jlocal(:) = modine%Jlocal(1)
161 modine%Jrange(:) = modine%Jrange(1)
165 solver_type =
root_newton, maxiter = 500, abs_tolerance = 1.0e-10_real64)
170 modine%min_mesh_scaling_product = modine%Jbar**modine%dim
181 safe_allocate(modine%csi(1:modine%dim, 1:modine%npos))
185 safe_allocate(modine%chi_atoms(1:modine%dim, 1:modine%npos))
187 do iatom = 1, modine%npos
188 modine%chi_atoms(:, iatom) = modine%from_cartesian(pos(:, iatom))
190 modine%csi(:,:) = modine%chi_atoms(:,:)
193 do iatom = 1, modine%npos
195 modine%chi_atoms(:, iatom) = nint(modine%chi_atoms(:, iatom) / spacing(:)) * spacing(:)
203 integer :: iatom, idim, index
204 real(real64),
allocatable :: my_csi(:), start_csi(:)
210 safe_allocate(x_p(1:modine%dim*modine%npos))
211 safe_allocate(my_csi(1:modine%dim*modine%npos))
212 safe_allocate(start_csi(1:modine%dim*modine%npos))
214 do iatom = 1, modine%npos
215 do idim = 1, modine%dim
216 index = (iatom-1)*dim + idim
217 x_p(index) = pos(idim, iatom)
218 start_csi(index) = modine%chi_atoms(idim, iatom)
225 message(1) =
"During the construction of the adaptive grid, the Newton-Raphson"
226 message(2) =
"method did not converge."
231 do iatom = 1, modine%npos
232 do idim = 1, modine%dim
233 index = (iatom-1)*modine%dim + idim
234 modine_p%csi(idim, iatom) = my_csi(index)
238 safe_deallocate_a(x_p)
239 safe_deallocate_a(my_csi)
240 safe_deallocate_a(start_csi)
256 this_out%dim = this_in%dim
257 this_out%local_basis = this_in%local_basis
258 safe_allocate_source_a(this_out%lsize, this_in%lsize)
259 this_out%xbar = this_in%xbar
260 this_out%Jbar = this_in%Jbar
261 safe_allocate_source_a(this_out%Jlocal, this_in%Jlocal)
262 safe_allocate_source_a(this_out%Jrange, this_in%Jrange)
263 safe_allocate_source_a(this_out%chi_atoms, this_in%chi_atoms)
264 safe_allocate_source_a(this_out%csi, this_in%csi)
265 this_out%npos = this_in%npos
276 safe_deallocate_a(this%lsize)
277 safe_deallocate_a(this%Jlocal)
278 safe_deallocate_a(this%Jrange)
279 safe_deallocate_a(this%chi_atoms)
280 safe_deallocate_a(this%csi)
288 real(real64),
intent(in) :: chi(:)
289 real(real64) :: xx(1:this%dim)
291 real(real64) :: chi2(this%dim), rr2, dd
299 do iatom = 1, this%npos
300 rr2 = sum((chi2(:) - this%csi(:, iatom))**2)
301 dd =
exp(-rr2/(
m_two*this%Jrange(iatom)**2))
303 xx(:) = xx(:) - this%Jlocal(iatom)*(chi2(:) - this%csi(:, iatom)) * dd
311 real(real64),
intent(in) :: xx(:)
312 real(real64) :: chi(1:this%dim)
319 safe_allocate(x_p(1:this%dim))
324 safe_deallocate_a(x_p)
328 message(1) =
"During the construction of the adaptive grid, the Newton-Raphson"
329 message(2) =
"method did not converge for point:"
330 write(
message(3),
'(3f14.6)') xx(1:this%dim)
337 real(real64) function curv_modine_det_jac(this, xx, chi) result(jdet)
339 real(real64),
intent(in) :: xx(:)
340 real(real64),
intent(in) :: chi(:)
342 real(real64) :: dummy(this%dim), jac(1:this%dim, 1:this%dim)
354 integer,
optional,
intent(in) :: iunit
355 type(namespace_t),
optional,
intent(in) :: namespace
359 write(message(1),
'(a)')
' Curvilinear Method = modine'
360 call messages_info(1, iunit=iunit, namespace=namespace)
366 real(real64) function curv_modine_surface_element(this, idir) result(ds)
368 integer,
intent(in) :: idir
371 message(1) =
'Surface element with modine curvilinear coordinates not implemented'
372 call messages_fatal(1)
379 real(real64),
intent(in) :: chi_(:)
380 real(real64),
intent(out) :: chi2(:)
381 real(real64),
optional,
intent(out) :: jac(:)
383 integer,
parameter :: qq = 3
384 real(real64) :: chibar(this%dim), rr, chi
390 chibar = this%xbar*this%lsize
396 chi2(i) = this%Jbar * chi
397 if (
present(jac)) jac(i) = this%Jbar
399 if (chi > chibar(i))
then
400 rr = (chi - chibar(i))/(this%lsize(i) - chibar(i))
402 chi2(i) = chi2(i) + this%lsize(i)/m_two*(1 - this%Jbar) * rr**qq * &
403 (qq + m_one - (qq - m_one)*rr)
405 if (
present(jac))
then
406 jac(i) = jac(i) + this%lsize(i)/m_two*(m_one - this%Jbar) * rr**(qq - 1)/(this%lsize(i) - chibar(i)) * &
407 (qq*(qq + m_one) - (qq**2 - m_one)*rr)
411 if (neg) chi2(i) = -chi2(i)
420 real(real64),
intent(in) :: chi(:)
421 real(real64),
intent(out) :: xx(:)
422 real(real64),
intent(out) :: jac(:, :)
424 real(real64) :: chi2(this%dim), rr2, dd, j2(this%dim)
425 integer :: iatom, idim, idim2
434 do idim = 1, this%dim
435 jac(idim, idim) = m_one
438 do iatom = 1, this%npos
439 rr2 = sum((chi2(:) - this%csi(:, iatom))**2)
440 dd =
exp(-rr2/(m_two*this%Jrange(iatom)**2))
442 xx(:) = xx(:) - this%Jlocal(iatom)*(chi2(:) - this%csi(:, iatom)) * dd
444 do idim = 1, this%dim
445 jac(idim, idim) = jac(idim, idim) - this%Jlocal(iatom) * dd
446 do idim2 = 1, this%dim
447 jac(idim, idim2) = jac(idim, idim2) + &
448 this%Jlocal(iatom)*(chi2(idim) - this%csi(idim, iatom))*(chi2(idim2) - this%csi(idim2, iatom)) * &
449 m_two/(m_two*this%Jrange(iatom)**2) * dd
454 do idim = 1, this%dim
455 jac(idim, :) = jac(idim, :) * j2(:)
461 pure subroutine getf(yy, ff, jf)
462 real(real64),
intent(in) :: yy(:)
463 real(real64),
intent(out) :: ff(:), jf(:, :)
468 ff(:) = ff(:) -
x_p(:)
473 subroutine getf2(csi, ff, jf)
474 real(real64),
intent(in) :: csi(:)
475 real(real64),
intent(out) :: ff(:), jf(:, :)
477 integer :: i1, j1, i2, j2, index1, index2
478 real(real64) :: rr2, dd, dd2
479 real(real64),
allocatable :: xx(:), chi2(:)
502 rr2 = sum((chi2 -
modine_p%csi(:,i2))**2)
510 ff(index1) = xx(j1) -
x_p(index1)
513 rr2 = sum((chi2 -
modine_p%csi(:,i2))**2)
515 dd2 = -m_two/(m_two*
modine_p%Jrange(i2)**2)*dd
518 jf(index1, index2) =
modine_p%Jlocal(i2) * dd
523 jf(index1, index2) = jf(index1, index2) +
modine_p%Jlocal(i2) * dd2 * &
530 safe_deallocate_a(xx)
531 safe_deallocate_a(chi2)
subroutine find_atom_points()
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
double exp(double __x) __attribute__((__nothrow__
This module implements the curvilinear coordinates given in N. A. Modine, G. Zumbach,...
real(real64) function curv_modine_det_jac(this, xx, chi)
class(curv_modine_t), pointer modine_p
real(real64) function, dimension(1:this%dim) curv_modine_from_cartesian(this, xx)
real(real64) function curv_modine_surface_element(this, idir)
pure subroutine curv_modine_jacobian_inv(this, chi, xx, Jac)
class(curv_modine_t) function, pointer curv_modine_constructor(namespace, dim, npos, pos, lsize, spacing)
subroutine getf2(csi, ff, jf)
pure subroutine getf(yy, ff, jf)
subroutine curv_modine_finalize(this)
real(real64), dimension(:), allocatable x_p
pure subroutine curv_modine_chi2chi2(this, chi_, chi2, Jac)
subroutine, public curv_modine_copy(this_out, this_in)
subroutine curv_modine_write_info(this, iunit, namespace)
real(real64) function, dimension(1:this%dim) curv_modine_to_cartesian(this, chi)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_half
real(real64), parameter, public m_one
real(real64), parameter, public m_three
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
integer, parameter, public root_newton
subroutine, public root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, rel_tolerance, abs_tolerance)
subroutine, public droot_solver_run(rs, func, root, success, startval)
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
abstract class to describe coordinate systems