Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_singular_mapping_advector_rotating.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_assert.h"
4 
6 
8 
9  implicit none
10 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
11 
12  ! Working precision
13  integer, parameter :: wp = f64
14 
16 
17  real(wp) :: xc(2) ! center of rotation
18  real(wp) :: omega ! angular velocity
19 
20  contains
21 
23  procedure :: velocity_field => s_singular_mapping_advector_rotating__velocity_field
24  procedure :: flow_field => f_singular_mapping_advector_rotating__flow_field
26 
28 
29 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
30 contains
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 
33  subroutine s_singular_mapping_advector_rotating__init(self, xc, omega)
34  class(sll_t_singular_mapping_advector_rotating), intent(inout) :: self
35  real(wp), intent(in) :: xc(2) ! center of rotation
36  real(wp), intent(in) :: omega ! angular velocity
37 
38  self%xc(:) = xc(:)
39  self%omega = omega
40 
42 
43  !-----------------------------------------------------------------------------
44  sll_pure subroutine s_singular_mapping_advector_rotating__velocity_field(self, x, a)
45  class(sll_t_singular_mapping_advector_rotating), intent(in) :: self
46  real(wp), intent(in) :: x(2)
47  real(wp), intent(out) :: a(2)
48 
49  a(1) = -self%omega*(x(2) - self%xc(2))
50  a(2) = self%omega*(x(1) - self%xc(1))
51 
52  end subroutine s_singular_mapping_advector_rotating__velocity_field
53 
54  !-----------------------------------------------------------------------------
55  sll_pure function f_singular_mapping_advector_rotating__flow_field(self, x, h) result(y)
56  class(sll_t_singular_mapping_advector_rotating), intent(in) :: self
57  real(wp), intent(in) :: x(2)
58  real(wp), intent(in) :: h
59  real(wp) :: y(2)
60 
61  y(1) = self%xc(1) &
62  + (x(1) - self%xc(1))*cos(self%omega*h) &
63  - (x(2) - self%xc(2))*sin(self%omega*h)
64 
65  y(2) = self%xc(2) &
66  + (x(1) - self%xc(1))*sin(self%omega*h) &
67  + (x(2) - self%xc(2))*cos(self%omega*h)
68 
69  end function f_singular_mapping_advector_rotating__flow_field
70 
71  !-----------------------------------------------------------------------------
73  class(sll_t_singular_mapping_advector_rotating), intent(inout) :: self
75 
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