Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_time_splitting.F90
Go to the documentation of this file.
1 #ifndef DOXYGEN_SHOULD_SKIP_THIS
2 
6 module sll_m_time_splitting
7 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8 #include "sll_memory.h"
9 #include "sll_working_precision.h"
10 
11  implicit none
12 
13  public :: &
14  sll_c_time_splitting
15 
16  private
17 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
18 
19  type, abstract :: sll_c_time_splitting
20  sll_real64 :: current_time = 0.0_f64
21  contains
22  procedure(time_splitting_operator), pass(this), deferred :: operator1
23  procedure(time_splitting_operator), pass(this), deferred :: operator2
24  procedure, pass(this) :: lie_splitting
25  procedure, pass(this) :: strang_splitting
26  end type sll_c_time_splitting
27 
28  abstract interface
29  subroutine time_splitting_operator(this, dt)
31  import sll_c_time_splitting
32  class(sll_c_time_splitting) :: this
33  sll_real64, intent(in) :: dt
34  end subroutine time_splitting_operator
35  end interface
36 
37 contains
42  subroutine lie_splitting(this, dt, number_steps)
43  class(sll_c_time_splitting) :: this
44  sll_real64, intent(in) :: dt
45  sll_int32, intent(in) :: number_steps
46  ! local variables
47  sll_int32 :: i
48 
49  do i = 1, number_steps
50  call this%operator1(dt)
51  call this%operator2(dt)
52  this%current_time = this%current_time + dt
53  end do
54  end subroutine lie_splitting
55 
60  subroutine strang_splitting(this, dt, number_steps)
61  class(sll_c_time_splitting) :: this
62  sll_real64, intent(in) :: dt
63  sll_int32, intent(in) :: number_steps
64  ! local variables
65  sll_int32 :: i
66 
67  call this%operator1(0.5_f64*dt)
68  do i = 1, number_steps - 2
69  call this%operator2(dt)
70  call this%operator1(dt)
71  ! warning this implies that operator 1 does not depend on value of current_time
72  ! if this is not true the two steps implying operator one should not be combined
73  this%current_time = this%current_time + dt
74  end do
75  call this%operator2(dt)
76  call this%operator1(0.5_f64*dt)
77  end subroutine strang_splitting
78 
79 end module sll_m_time_splitting
80 
81 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
Module to select the kind parameter.
    Report Typos and Errors