![]() |
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Implements split operators for constant coefficient advection.
Derived types and interfaces | |
| type | sll_t_const_coef_advection_2d |
| Simple operator splitting type for 2D constant coefficient advection Extends operator splitting. More... | |
Functions/Subroutines | |
| 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 | 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. More... | |
| subroutine | adv1 (this, dt) |
| Constant coefficient advection operator in first direction. More... | |
| subroutine | adv2 (this, dt) |
| Constant coefficient advection operator in second direction. More... | |
|
private |
Constant coefficient advection operator in first direction.
| [in,out] | this | object |
| [in] | dt | time step |
Definition at line 110 of file sll_m_const_coef_advection_2d.F90.
|
private |
Constant coefficient advection operator in second direction.
| [in,out] | this | object |
| [in] | dt | time step |
Definition at line 126 of file sll_m_const_coef_advection_2d.F90.
|
private |
Initialise sll_t_const_coef_advection_2d object.
| [in,out] | this | object |
| [in] | data | initial value of function |
| [in] | n1 | dimension in first direction |
| [in] | n2 | dimension in second direction |
| [in] | a1 | advection coefficient in first direction |
| [in] | a2 | advection coefficient in second direction |
| interp1 | interpolator for first direction | |
| interp2 | interpolator for second direction | |
| [in] | split_case | defines splitting method |
| [in] | split_step | coefficients of split step |
| [in] | nb_split_step | number of split steps |
| [in] | split_begin_t | begin with operator T if .true. |
| [in] | dt | time step |
Definition at line 76 of file sll_m_const_coef_advection_2d.F90.
| class(sll_t_const_coef_advection_2d) function, pointer, public sll_m_const_coef_advection_2d::sll_f_new_const_coef_advection_2d | ( | real(kind=f64), dimension(:, :), intent(in), pointer | data, |
| integer(kind=i32), intent(in) | n1, | ||
| integer(kind=i32), intent(in) | n2, | ||
| real(kind=f64), intent(in) | a1, | ||
| real(kind=f64), intent(in) | a2, | ||
| class(sll_c_interpolator_1d), pointer | interp1, | ||
| class(sll_c_interpolator_1d), pointer | interp2, | ||
| integer(kind=i32), intent(in) | split_case, | ||
| real(kind=f64), dimension(:), intent(in), optional | split_step, | ||
| integer(kind=i32), intent(in), optional | nb_split_step, | ||
| logical, intent(in), optional | split_begin_T, | ||
| real(kind=f64), intent(in), optional | dt | ||
| ) |
| [in] | data | initial value of function |
| [in] | n1 | dimension in first direction |
| [in] | n2 | dimension in second direction |
| [in] | a1 | advection coeeficient in first direction |
| [in] | a2 | advection coeeficient in first direction |
| interp1 | interpolator for first direction | |
| interp2 | interpolator for second direction | |
| [in] | split_case | defines splitting method |
| [in] | split_step | coefficients of split step |
| [in] | nb_split_step | number of split steps |
| [in] | split_begin_t | begin with operator T if .true. |
| [in] | dt | time step |
Definition at line 50 of file sll_m_const_coef_advection_2d.F90.
1.9.1