1 #ifndef DOXYGEN_SHOULD_SKIP_THIS
6 module sll_m_time_splitting
8 #include "sll_memory.h"
9 #include "sll_working_precision.h"
19 type,
abstract :: sll_c_time_splitting
20 sll_real64 :: current_time = 0.0_f64
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
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
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
49 do i = 1, number_steps
50 call this%operator1(dt)
51 call this%operator2(dt)
52 this%current_time = this%current_time + dt
54 end subroutine lie_splitting
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
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)
73 this%current_time = this%current_time + dt
75 call this%operator2(dt)
76 call this%operator1(0.5_f64*dt)
77 end subroutine strang_splitting
79 end module sll_m_time_splitting
81 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
Module to select the kind parameter.