Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_singular_mapping_analytic_target.F90
Go to the documentation of this file.
1 ! This module implements the analytical mapping given by
2 ! equation (2) of https://doi.org/10.1016/j.jcp.2019.108889:
3 !
4 ! x(s,theta)=x0+(1-kappa)*s*cos(theta)-Delta*s**2
5 ! y(s,theta)=y0+(1+kappa)*s*sin(theta)
6 
8 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9 #include "sll_assert.h"
10 
11  use sll_m_working_precision, only: f64
12 
13  use sll_m_constants, only: sll_p_twopi
14 
16 
17  implicit none
18 
20 
21  private
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 
25  integer, parameter :: wp = f64
26 
29 
30  ! Default values (circle), can be overwritten from 'init' method
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
34 
35  contains
36 
38  procedure :: eval => f_singular_mapping_analytic_target__eval
39  procedure :: jmat => f_singular_mapping_analytic_target__jmat ! Jacobian matrix
40  procedure :: jmat_comp => f_singular_mapping_analytic_target__jmat_comp
41 
43 
44 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
45 contains
46 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 
48  subroutine s_singular_mapping_analytic_target__init(self, x0, Delta, kappa)
49  class(sll_t_singular_mapping_analytic_target), intent(inout) :: self
50  real(wp), optional, intent(in) :: x0(2)
51  real(wp), optional, intent(in) :: Delta
52  real(wp), optional, intent(in) :: kappa
53 
54  ! Overwrite parameters
55  if (present(x0)) self%x0 = x0
56  if (present(delta)) self%Delta = delta
57  if (present(kappa)) self%kappa = kappa
58 
60 
61  !-----------------------------------------------------------------------------
62  sll_pure function f_singular_mapping_analytic_target__eval(self, eta) result(x)
63  class(sll_t_singular_mapping_analytic_target), intent(in) :: self
64  real(wp), intent(in) :: eta(2)
65  real(wp) :: x(2)
66 
67  associate(s => eta(1), t => eta(2), x0 => self%x0, delta => self%Delta, kappa => self%kappa)
68 
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)
71 
72  end associate
73 
74  end function f_singular_mapping_analytic_target__eval
75 
76  !-----------------------------------------------------------------------------
77  sll_pure function f_singular_mapping_analytic_target__jmat(self, eta) result(jmat)
78  class(sll_t_singular_mapping_analytic_target), intent(in) :: self
79  real(wp), intent(in) :: eta(2)
80  real(wp) :: jmat(2, 2)
81 
82  associate(s => eta(1), t => eta(2), delta => self%Delta, kappa => self%kappa)
83 
84  ! J_11 = d(x1)/d(eta1)
85  ! J_12 = d(x1)/d(eta2)
86  ! J_21 = d(x2)/d(eta1)
87  ! J_22 = d(x2)/d(eta2)
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)
92 
93  end associate
94 
95  end function f_singular_mapping_analytic_target__jmat
96 
97  !-----------------------------------------------------------------------------
98  sll_pure function f_singular_mapping_analytic_target__jmat_comp(self, eta) result(jmat_comp)
99  class(sll_t_singular_mapping_analytic_target), intent(in) :: self
100  real(wp), intent(in) :: eta(2)
101  real(wp) :: jmat_comp(2, 2)
102 
103  associate(s => eta(1), t => eta(2), delta => self%Delta, kappa => self%kappa)
104 
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)
109 
110  jmat_comp = jmat_comp/((1.0_wp + kappa)*((1.0_wp - kappa) - 2.0_wp*delta*s*cos(t)))
111 
112  end associate
113 
114  end function f_singular_mapping_analytic_target__jmat_comp
115 
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_twopi
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)
    Report Typos and Errors