Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_const_coef_advection_2d.F90
Go to the documentation of this file.
1 
4 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
5 #include "sll_memory.h"
6 #include "sll_working_precision.h"
7 
8  use sll_m_interpolators_1d_base, only: &
10 
11  use sll_m_operator_splitting, only: &
14 
15  implicit none
16 
17  public :: &
20 
21  private
22 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 
31  class(sll_c_interpolator_1d), pointer :: interp1
33  class(sll_c_interpolator_1d), pointer :: interp2
35  sll_real64, dimension(:, :), pointer :: data
37  sll_int32 :: n1
39  sll_int32 :: n2
41  sll_real64 :: a1
43  sll_real64 :: a2
44  contains
45  procedure, pass(this) :: operatort => adv1
46  procedure, pass(this) :: operatorv => adv2
48 
49 contains
50  function sll_f_new_const_coef_advection_2d(data, n1, n2, a1, a2, interp1, interp2, &
51  split_case, split_step, nb_split_step, split_begin_T, dt) &
52  result(this)
53  class(sll_t_const_coef_advection_2d), pointer :: this
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
59  class(sll_c_interpolator_1d), pointer :: interp1
60  class(sll_c_interpolator_1d), pointer :: interp2
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
66  ! local variable
67  sll_int32 :: ierr
68 
69  sll_allocate(this, ierr)
70  call initialize_const_coef_advection_2d(this, data, n1, n2, a1, a2, interp1, interp2, &
71  split_case, split_step, nb_split_step, split_begin_t, dt)
72 
73  end function
74 
76  subroutine initialize_const_coef_advection_2d(this, data, n1, n2, a1, a2, interp1, interp2, &
77  split_case, split_step, nb_split_step, split_begin_T, dt)
78  class(sll_t_const_coef_advection_2d), intent(inout) :: this
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
84  class(sll_c_interpolator_1d), pointer :: interp1
85  class(sll_c_interpolator_1d), pointer :: interp2
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
91 
92  this%data => data
93  this%n1 = n1
94  this%n2 = n2
95  this%a1 = a1
96  this%a2 = a2
97  this%interp1 => interp1
98  this%interp2 => interp2
99 
101  this, &
102  split_case, &
103  split_step, &
104  nb_split_step, &
105  split_begin_t, &
106  dt)
107  end subroutine
108 
110  subroutine adv1(this, dt)
111  class(sll_t_const_coef_advection_2d), intent(inout) :: this
112  sll_real64, intent(in) :: dt
113  ! local variables
114  sll_real64, dimension(:), pointer :: f1d
115  sll_real64 :: displacement
116  sll_int32 :: j
117 
118  do j = 1, this%n2
119  displacement = -this%a1*dt
120  f1d => this%data(:, j)
121  call this%interp1%interpolate_array_disp_inplace(this%n1, f1d, displacement)
122  end do
123  end subroutine
124 
126  subroutine adv2(this, dt)
127  class(sll_t_const_coef_advection_2d), intent(inout) :: this
128  sll_real64, intent(in) :: dt
129  ! local variables
130  sll_real64, dimension(:), pointer :: f1d
131  sll_real64 :: displacement
132  sll_int32 :: i
133 
134  do i = 1, this%n1
135  displacement = -this%a2*dt
136  f1d => this%data(i, :)
137  call this%interp2%interpolate_array_disp_inplace(this%n2, f1d, displacement)
138  end do
139  end subroutine
140 
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.
    Report Typos and Errors