9 #include "sll_assert.h"
31 real(
wp) :: x0(2) = [0.0_wp, 0.0_wp]
32 real(
wp) :: delta = 0.0_wp
33 real(
wp) :: kappa = 0.0_wp
38 procedure :: eval => f_singular_mapping_analytic_target__eval
39 procedure :: jmat => f_singular_mapping_analytic_target__jmat
40 procedure :: jmat_comp => f_singular_mapping_analytic_target__jmat_comp
50 real(wp),
optional,
intent(in) :: x0(2)
51 real(wp),
optional,
intent(in) :: Delta
52 real(wp),
optional,
intent(in) :: kappa
55 if (
present(x0)) self%x0 = x0
56 if (
present(delta)) self%Delta = delta
57 if (
present(kappa)) self%kappa = kappa
62 sll_pure
function f_singular_mapping_analytic_target__eval(self, eta)
result(x)
64 real(wp),
intent(in) :: eta(2)
67 associate(s => eta(1), t => eta(2), x0 => self%x0, delta => self%Delta, kappa => self%kappa)
69 x(1) = x0(1) + (1.0_wp - kappa)*s*cos(t) - delta*s**2
70 x(2) = x0(2) + (1.0_wp + kappa)*s*sin(t)
74 end function f_singular_mapping_analytic_target__eval
77 sll_pure
function f_singular_mapping_analytic_target__jmat(self, eta)
result(jmat)
79 real(wp),
intent(in) :: eta(2)
80 real(wp) :: jmat(2, 2)
82 associate(s => eta(1), t => eta(2), delta => self%Delta, kappa => self%kappa)
88 jmat(1, 1) = -2.0_wp*delta*s + (1.0_wp - kappa)*cos(t)
89 jmat(1, 2) = -(1.0_wp - kappa)*s*sin(t)
90 jmat(2, 1) = (1.0_wp + kappa)*sin(t)
91 jmat(2, 2) = (1.0_wp + kappa)*s*cos(t)
95 end function f_singular_mapping_analytic_target__jmat
98 sll_pure
function f_singular_mapping_analytic_target__jmat_comp(self, eta)
result(jmat_comp)
100 real(wp),
intent(in) :: eta(2)
101 real(wp) :: jmat_comp(2, 2)
103 associate(s => eta(1), t => eta(2), delta => self%Delta, kappa => self%kappa)
105 jmat_comp(1, 1) = (1.0_wp + kappa)
106 jmat_comp(1, 2) = 2.0_wp*delta*s*sin(t)
107 jmat_comp(2, 1) = 0.0_wp
108 jmat_comp(2, 2) = (1.0_wp - kappa) - 2.0_wp*delta*s*cos(t)
110 jmat_comp = jmat_comp/((1.0_wp + kappa)*((1.0_wp - kappa) - 2.0_wp*delta*s*cos(t)))
114 end function f_singular_mapping_analytic_target__jmat_comp
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_twopi
integer, parameter wp
Working precision.
subroutine s_singular_mapping_analytic_target__init(self, x0, Delta, kappa)
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
Abstract type, analytical singular mapping.
Concrete type, analytical singular mapping.