Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_hamiltonian_splitting_base.F90
Go to the documentation of this file.
1 
7 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_working_precision.h"
9 
10  implicit none
11 
12  public :: &
14 
15  private
16 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
17 
20 
21  contains
22  procedure(splitting), deferred :: lie_splitting
23  procedure(splitting), deferred :: lie_splitting_back
24  procedure(splitting), deferred :: strang_splitting
25  procedure :: splitting_fourth
28  procedure(empty), deferred :: free
29  procedure :: reinit_fields
31 
32 
33  abstract interface
34  subroutine splitting(self, dt, number_steps)
37  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
38  sll_real64, intent(in) :: dt
39  sll_int32, intent(in) :: number_steps
40  end subroutine splitting
41  end interface
42 
43  abstract interface
44  subroutine empty(self)
46  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
47  end subroutine empty
48  end interface
49 
50 contains
51 
52  subroutine splitting_fourth( self, dt, number_steps )
53  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
54  sll_real64, intent(in) :: dt
55  sll_int32, intent(in) :: number_steps
56 
57  sll_real64 :: dt1, dt2
58  sll_int32 :: j
59 
60  dt1 = dt / (2.0_f64 - 2.0_f64**(1.0_f64/3.0_f64))
61  dt2 = - 2.0_f64**(1.0_f64/3.0_f64)* dt1
62 
63  do j=1,number_steps
64 
65  call self%strang_splitting( dt1, 1 )
66  call self%strang_splitting( dt2, 1 )
67  call self%strang_splitting( dt1, 1 )
68 
69  end do
70 
71 
72  end subroutine splitting_fourth
73 
74 
75 
76  subroutine splitting_fourth_10steps( self, dt, number_steps )
77  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
78  sll_real64, intent(in) :: dt
79  sll_int32, intent(in) :: number_steps
80 
81  sll_real64 :: a1, a2, a3, a4, a5
82  sll_int32 :: j
83 
84 
85  a1 = (146.0_f64 + 5.0_f64*sqrt(19.0_f64))/540.0_f64
86  a2 = (-2.0_f64+10.0_f64*sqrt(19.0_f64))/135.0_f64
87  a3 = 0.2_f64
88  a4 = (-23.0_f64-20.0_f64*sqrt(19.0_f64))/270.0_f64
89  a5 = (14.0_f64-sqrt(19.0_f64))/108.0_f64
90 
91  do j=1,number_steps
92 
93  call self%lie_splitting_back( a5*dt,1 )
94  call self%lie_splitting( a1*dt,1 )
95  call self%lie_splitting_back( a4*dt,1 )
96  call self%lie_splitting( a2*dt,1 )
97  call self%lie_splitting_back( a3*dt,1 )
98  call self%lie_splitting( a3*dt,1 )
99  call self%lie_splitting_back( a2*dt,1 )
100  call self%lie_splitting( a4*dt,1 )
101  call self%lie_splitting_back( a1*dt,1 )
102  call self%lie_splitting( a5*dt,1 )
103  end do
104 
105 
106  end subroutine splitting_fourth_10steps
107 
108 
109  subroutine splitting_second_4steps( self, dt, number_steps )
110  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
111  sll_real64, intent(in) :: dt
112  sll_int32, intent(in) :: number_steps
113 
114  sll_real64 :: a1, a2
115  sll_int32 :: j
116 
117 
118  a1 = 0.1932_f64!0.25_f64
119  a2 = 0.5_f64 - a1
120 
121  do j=1,number_steps
122 
123  call self%lie_splitting_back( a1*dt,1 )
124  call self%lie_splitting( a2*dt,1 )
125  call self%lie_splitting_back( a2*dt,1 )
126  call self%lie_splitting( a1*dt,1 )
127  end do
128 
129  end subroutine splitting_second_4steps
130 
131  subroutine reinit_fields( self )
132  class(sll_c_hamiltonian_splitting_base), intent(inout) :: self
133 
134  end subroutine reinit_fields
135 
136 
Base class for Hamiltonian splittings.
subroutine splitting_fourth_10steps(self, dt, number_steps)
subroutine splitting_fourth(self, dt, number_steps)
subroutine splitting_second_4steps(self, dt, number_steps)
Module to select the kind parameter.
    Report Typos and Errors