5 #include "sll_memory.h"
6 #include "sll_working_precision.h"
35 sll_real64,
dimension(:, :),
pointer ::
data
45 procedure, pass(this) :: operatort =>
adv1
46 procedure, pass(this) :: operatorv =>
adv2
51 split_case, split_step, nb_split_step, split_begin_T, dt) &
54 sll_real64,
dimension(:, :),
pointer,
intent(in) ::
data
55 sll_int32,
intent(in) :: n1
56 sll_int32,
intent(in) :: n2
57 sll_real64,
intent(in) :: a1
58 sll_real64,
intent(in) :: a2
61 sll_int32,
intent(in) :: split_case
62 sll_real64,
dimension(:),
intent(in),
optional :: split_step
63 sll_int32,
intent(in),
optional :: nb_split_step
64 logical,
intent(in),
optional :: split_begin_t
65 sll_real64,
intent(in),
optional :: dt
69 sll_allocate(this, ierr)
71 split_case, split_step, nb_split_step, split_begin_t, dt)
77 split_case, split_step, nb_split_step, split_begin_T, dt)
79 sll_real64,
dimension(:, :),
pointer,
intent(in) ::
data
80 sll_int32,
intent(in) :: n1
81 sll_int32,
intent(in) :: n2
82 sll_real64,
intent(in) :: a1
83 sll_real64,
intent(in) :: a2
86 sll_int32,
intent(in) :: split_case
87 sll_real64,
dimension(:),
intent(in),
optional :: split_step
88 sll_int32,
intent(in),
optional :: nb_split_step
89 logical,
intent(in),
optional :: split_begin_T
90 sll_real64,
intent(in),
optional :: dt
97 this%interp1 => interp1
98 this%interp2 => interp2
112 sll_real64,
intent(in) :: dt
114 sll_real64,
dimension(:),
pointer :: f1d
115 sll_real64 :: displacement
119 displacement = -this%a1*dt
120 f1d => this%data(:, j)
121 call this%interp1%interpolate_array_disp_inplace(this%n1, f1d, displacement)
128 sll_real64,
intent(in) :: dt
130 sll_real64,
dimension(:),
pointer :: f1d
131 sll_real64 :: displacement
135 displacement = -this%a2*dt
136 f1d => this%data(i, :)
137 call this%interp2%interpolate_array_disp_inplace(this%n2, f1d, displacement)
Implements split operators for constant coefficient advection.
subroutine initialize_const_coef_advection_2d(this, data, n1, n2, a1, a2, interp1, interp2, split_case, split_step, nb_split_step, split_begin_T, dt)
Initialise sll_t_const_coef_advection_2d object.
class(sll_t_const_coef_advection_2d) function, pointer, public sll_f_new_const_coef_advection_2d(data, n1, n2, a1, a2, interp1, interp2, split_case, split_step, nb_split_step, split_begin_T, dt)
subroutine adv2(this, dt)
Constant coefficient advection operator in second direction.
subroutine adv1(this, dt)
Constant coefficient advection operator in first direction.
Module for 1D interpolation and reconstruction.
Base class of operator splitting library. It is only used with particle-in-cell method.
subroutine, public sll_s_operator_splitting_init(split, split_case, split_step, nb_split_step, split_begin_T, dt)
Initialises data for operator splitting.
Simple operator splitting type for 2D constant coefficient advection Extends operator splitting.
Abstract class for 1D interpolation and reconstruction.
operator splitting object