Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
Derived types and interfaces | Functions/Subroutines
sll_m_advection_2d_oblic Module Reference

Description

use oblic interpolation

data are on uniform (fine) grid in x1

Derived types and interfaces

type  sll_t_oblic_2d_advector
 

Functions/Subroutines

type(sll_t_oblic_2d_advector) function, pointer, public sll_f_new_oblic_2d_advector (Nc_x1, adv_x1, Nc_x2, x2_min, x2_max, stencil_r, stencil_s)
 
subroutine initialize_oblic_2d_advector (adv, Nc_x1, adv_x1, Nc_x2, x2_min, x2_max, stencil_r, stencil_s)
 
subroutine, public sll_s_oblic_advect_2d_constant (adv, A1, A2, dt, input, output)
 solves \partial_t f + \nabla A \cdot f = 0, \ A = (A1,A2) for time step dt interpolation in aligned along iota = A1/A2 iota is the number of tours that has been done at x2_max the direction of iota is the line passing through (x1_min,x2_min) to (x1_max,x2_min+iota*(x2_max-x2_min)) here iota is real number for specific cases where x2_min+iota*(x2_max-x2_min) is a mesh point we refer to sll_advection_2d_integer_oblic here A1 and A2 are real numbers we refer to (future) oblic_advect_2d, when (A1,A2) are 2D arrays that is A1 and A2 depend both on x1 and x2 we use here constant advection in the x1 direction with an abstract 1d advector that can do the constant advection in the x2 direction, we use (for the moment) Lagrange interpolation with stencil (r,s) (r,s) = (0,1) : LAG1 (r,s) = (-1,2) : LAG3 (r,s) = (-2,3) : LAG5 periodic conditions are used in both x1 and x2 directions More...
 
subroutine oblic_advect_2d (adv, A1, A2, dt, input, output)
 

Function/Subroutine Documentation

◆ initialize_oblic_2d_advector()

subroutine sll_m_advection_2d_oblic::initialize_oblic_2d_advector ( type(sll_t_oblic_2d_advector), intent(inout)  adv,
integer(kind=i32), intent(in)  Nc_x1,
class(sll_c_advector_1d), pointer  adv_x1,
integer(kind=i32), intent(in)  Nc_x2,
real(kind=f64), intent(in)  x2_min,
real(kind=f64), intent(in)  x2_max,
integer(kind=i32), intent(in)  stencil_r,
integer(kind=i32), intent(in)  stencil_s 
)
private

Definition at line 89 of file sll_m_advection_2d_oblic.F90.

Here is the caller graph for this function:

◆ oblic_advect_2d()

subroutine sll_m_advection_2d_oblic::oblic_advect_2d ( type(sll_t_oblic_2d_advector adv,
real(kind=f64), dimension(:, :), intent(in)  A1,
real(kind=f64), dimension(:, :), intent(in)  A2,
real(kind=f64), intent(in)  dt,
real(kind=f64), dimension(:, :), intent(in)  input,
real(kind=f64), dimension(:, :), intent(out)  output 
)
private

Definition at line 216 of file sll_m_advection_2d_oblic.F90.

◆ sll_f_new_oblic_2d_advector()

type(sll_t_oblic_2d_advector) function, pointer, public sll_m_advection_2d_oblic::sll_f_new_oblic_2d_advector ( integer(kind=i32), intent(in)  Nc_x1,
class(sll_c_advector_1d), pointer  adv_x1,
integer(kind=i32), intent(in)  Nc_x2,
real(kind=f64), intent(in)  x2_min,
real(kind=f64), intent(in)  x2_max,
integer(kind=i32), intent(in)  stencil_r,
integer(kind=i32), intent(in)  stencil_s 
)

Definition at line 56 of file sll_m_advection_2d_oblic.F90.

Here is the call graph for this function:

◆ sll_s_oblic_advect_2d_constant()

subroutine, public sll_m_advection_2d_oblic::sll_s_oblic_advect_2d_constant ( type(sll_t_oblic_2d_advector adv,
real(kind=f64), intent(in)  A1,
real(kind=f64), intent(in)  A2,
real(kind=f64), intent(in)  dt,
real(kind=f64), dimension(:, :), intent(in)  input,
real(kind=f64), dimension(:, :), intent(out)  output 
)

solves \partial_t f + \nabla A \cdot f = 0, \ A = (A1,A2) for time step dt interpolation in aligned along iota = A1/A2 iota is the number of tours that has been done at x2_max the direction of iota is the line passing through (x1_min,x2_min) to (x1_max,x2_min+iota*(x2_max-x2_min)) here iota is real number for specific cases where x2_min+iota*(x2_max-x2_min) is a mesh point we refer to sll_advection_2d_integer_oblic here A1 and A2 are real numbers we refer to (future) oblic_advect_2d, when (A1,A2) are 2D arrays that is A1 and A2 depend both on x1 and x2 we use here constant advection in the x1 direction with an abstract 1d advector that can do the constant advection in the x2 direction, we use (for the moment) Lagrange interpolation with stencil (r,s) (r,s) = (0,1) : LAG1 (r,s) = (-1,2) : LAG3 (r,s) = (-2,3) : LAG5 periodic conditions are used in both x1 and x2 directions

Definition at line 151 of file sll_m_advection_2d_oblic.F90.

Here is the call graph for this function:
    Report Typos and Errors