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_hermite_interpolator_2d Module Reference

Description

Hermite interpolation in 2d.

derivatives are given with finite stencil formulae of order p which can be arbitrary in each direction

do not hesitate to take large odd p, like the favourite p=17

Todo:
for the moment only in implementation for the case DIRICHLET x PERIODIC

Derived types and interfaces

type  sll_hermite_interpolator_2d
 The hermite-based interpolator is only a wrapper around the capabilities of the hermite interpolation. More...
 

Functions/Subroutines

type(sll_hermite_interpolator_2d) function, pointer, public sll_f_new_hermite_interpolator_2d (npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, degree1, degree2, eta1_hermite_continuity, eta2_hermite_continuity, eta1_bc_type, eta2_bc_type, const_eta1_min_slope, const_eta1_max_slope, const_eta2_min_slope, const_eta2_max_slope, eta1_min_slopes, eta1_max_slopes, eta2_min_slopes, eta2_max_slopes)
 
subroutine initialize_hermite_interpolator_2d (interpolator, npts1, npts2, eta1_min, eta1_max, eta2_min, eta2_max, degree1, degree2, eta1_hermite_continuity, eta2_hermite_continuity, eta1_bc_type, eta2_bc_type, const_eta1_min_slope, const_eta1_max_slope, const_eta2_min_slope, const_eta2_max_slope, eta1_min_slopes, eta1_max_slopes, eta2_min_slopes, eta2_max_slopes)
 
subroutine wrap_compute_interpolants_hermite_2d (interpolator, data_array, eta1_coords, size_eta1_coords, eta2_coords, size_eta2_coords)
 
real(kind=f64) function wrap_interpolate_value_hermite_2d (interpolator, eta1, eta2)
 
real(kind=f64) function wrap_interpolate_deriv1_hermite_2d (interpolator, eta1, eta2)
 
real(kind=f64) function wrap_interpolate_deriv2_hermite_2d (interpolator, eta1, eta2)
 
subroutine wrap_interpolate_array_hermite_2d (this, num_points1, num_points2, data_in, eta1, eta2, data_out)
 
subroutine wrap_interpolate2d_disp_hermite_2d (this, num_points1, num_points2, data_in, alpha1, alpha2, data_out)
 
subroutine wrap_set_coefficients_hermite_2d (interpolator, coeffs_1d, coeffs_2d, coeff2d_size1, coeff2d_size2, knots1, size_knots1, knots2, size_knots2)
 
real(kind=f64) function, dimension(:, :), pointer wrap_get_coefficients_hermite_2d (interpolator)
 
logical function wrap_coefficients_are_set_hermite_2d (interpolator)
 
subroutine delete_sll_hermite_interpolator_2d (interpolator)
 

Function/Subroutine Documentation

◆ delete_sll_hermite_interpolator_2d()

subroutine sll_m_hermite_interpolator_2d::delete_sll_hermite_interpolator_2d ( class(sll_hermite_interpolator_2d), intent(inout)  interpolator)
private

Definition at line 423 of file sll_m_hermite_interpolator_2d.F90.

◆ initialize_hermite_interpolator_2d()

subroutine sll_m_hermite_interpolator_2d::initialize_hermite_interpolator_2d ( class(sll_hermite_interpolator_2d), intent(inout)  interpolator,
integer(kind=i32), intent(in)  npts1,
integer(kind=i32), intent(in)  npts2,
real(kind=f64), intent(in)  eta1_min,
real(kind=f64), intent(in)  eta1_max,
real(kind=f64), intent(in)  eta2_min,
real(kind=f64), intent(in)  eta2_max,
integer(kind=i32), intent(in)  degree1,
integer(kind=i32), intent(in)  degree2,
integer(kind=i32), intent(in)  eta1_hermite_continuity,
integer(kind=i32), intent(in)  eta2_hermite_continuity,
integer(kind=i32), intent(in)  eta1_bc_type,
integer(kind=i32), intent(in)  eta2_bc_type,
real(kind=f64), intent(in), optional  const_eta1_min_slope,
real(kind=f64), intent(in), optional  const_eta1_max_slope,
real(kind=f64), intent(in), optional  const_eta2_min_slope,
real(kind=f64), intent(in), optional  const_eta2_max_slope,
real(kind=f64), dimension(:), intent(in), optional  eta1_min_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta1_max_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta2_min_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta2_max_slopes 
)
private

Definition at line 174 of file sll_m_hermite_interpolator_2d.F90.

Here is the call graph for this function:

◆ sll_f_new_hermite_interpolator_2d()

type(sll_hermite_interpolator_2d) function, pointer, public sll_m_hermite_interpolator_2d::sll_f_new_hermite_interpolator_2d ( integer(kind=i32), intent(in)  npts1,
integer(kind=i32), intent(in)  npts2,
real(kind=f64), intent(in)  eta1_min,
real(kind=f64), intent(in)  eta1_max,
real(kind=f64), intent(in)  eta2_min,
real(kind=f64), intent(in)  eta2_max,
integer(kind=i32), intent(in)  degree1,
integer(kind=i32), intent(in)  degree2,
integer(kind=i32), intent(in)  eta1_hermite_continuity,
integer(kind=i32), intent(in)  eta2_hermite_continuity,
integer(kind=i32), intent(in)  eta1_bc_type,
integer(kind=i32), intent(in)  eta2_bc_type,
real(kind=f64), intent(in), optional  const_eta1_min_slope,
real(kind=f64), intent(in), optional  const_eta1_max_slope,
real(kind=f64), intent(in), optional  const_eta2_min_slope,
real(kind=f64), intent(in), optional  const_eta2_max_slope,
real(kind=f64), dimension(:), intent(in), optional  eta1_min_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta1_max_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta2_min_slopes,
real(kind=f64), dimension(:), intent(in), optional  eta2_max_slopes 
)

Definition at line 99 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_coefficients_are_set_hermite_2d()

logical function sll_m_hermite_interpolator_2d::wrap_coefficients_are_set_hermite_2d ( class(sll_hermite_interpolator_2d), intent(in)  interpolator)
private

Definition at line 413 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_compute_interpolants_hermite_2d()

subroutine sll_m_hermite_interpolator_2d::wrap_compute_interpolants_hermite_2d ( class(sll_hermite_interpolator_2d), intent(inout)  interpolator,
real(kind=f64), dimension(:, :), intent(in)  data_array,
real(kind=f64), dimension(:), intent(in), optional  eta1_coords,
integer(kind=i32), intent(in), optional  size_eta1_coords,
real(kind=f64), dimension(:), intent(in), optional  eta2_coords,
integer(kind=i32), intent(in), optional  size_eta2_coords 
)
private

Definition at line 243 of file sll_m_hermite_interpolator_2d.F90.

Here is the call graph for this function:

◆ wrap_get_coefficients_hermite_2d()

real(kind=f64) function, dimension(:, :), pointer sll_m_hermite_interpolator_2d::wrap_get_coefficients_hermite_2d ( class(sll_hermite_interpolator_2d), intent(in)  interpolator)
private

Definition at line 402 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_interpolate2d_disp_hermite_2d()

subroutine sll_m_hermite_interpolator_2d::wrap_interpolate2d_disp_hermite_2d ( class(sll_hermite_interpolator_2d), intent(inout)  this,
integer(kind=i32), intent(in)  num_points1,
integer(kind=i32), intent(in)  num_points2,
real(kind=f64), dimension(:, :), intent(in)  data_in,
real(kind=f64), dimension(:, :), intent(in)  alpha1,
real(kind=f64), dimension(:, :), intent(in)  alpha2,
real(kind=f64), dimension(num_points1, num_points2), intent(out)  data_out 
)
private

Definition at line 329 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_interpolate_array_hermite_2d()

subroutine sll_m_hermite_interpolator_2d::wrap_interpolate_array_hermite_2d ( class(sll_hermite_interpolator_2d), intent(inout)  this,
integer(kind=i32), intent(in)  num_points1,
integer(kind=i32), intent(in)  num_points2,
real(kind=f64), dimension(:, :), intent(in)  data_in,
real(kind=f64), dimension(:, :), intent(in)  eta1,
real(kind=f64), dimension(:, :), intent(in)  eta2,
real(kind=f64), dimension(num_points1, num_points2), intent(out)  data_out 
)
private

Definition at line 301 of file sll_m_hermite_interpolator_2d.F90.

Here is the call graph for this function:

◆ wrap_interpolate_deriv1_hermite_2d()

real(kind=f64) function sll_m_hermite_interpolator_2d::wrap_interpolate_deriv1_hermite_2d ( class(sll_hermite_interpolator_2d), intent(in)  interpolator,
real(kind=f64), intent(in)  eta1,
real(kind=f64), intent(in)  eta2 
)
private

Definition at line 275 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_interpolate_deriv2_hermite_2d()

real(kind=f64) function sll_m_hermite_interpolator_2d::wrap_interpolate_deriv2_hermite_2d ( class(sll_hermite_interpolator_2d), intent(in)  interpolator,
real(kind=f64), intent(in)  eta1,
real(kind=f64), intent(in)  eta2 
)
private

Definition at line 288 of file sll_m_hermite_interpolator_2d.F90.

◆ wrap_interpolate_value_hermite_2d()

real(kind=f64) function sll_m_hermite_interpolator_2d::wrap_interpolate_value_hermite_2d ( class(sll_hermite_interpolator_2d), intent(in)  interpolator,
real(kind=f64), intent(in)  eta1,
real(kind=f64), intent(in)  eta2 
)
private

Definition at line 266 of file sll_m_hermite_interpolator_2d.F90.

Here is the call graph for this function:

◆ wrap_set_coefficients_hermite_2d()

subroutine sll_m_hermite_interpolator_2d::wrap_set_coefficients_hermite_2d ( class(sll_hermite_interpolator_2d), intent(inout)  interpolator,
real(kind=f64), dimension(:), intent(in), optional  coeffs_1d,
real(kind=f64), dimension(:, :), intent(in), optional  coeffs_2d,
integer(kind=i32), intent(in), optional  coeff2d_size1,
integer(kind=i32), intent(in), optional  coeff2d_size2,
real(kind=f64), dimension(:), intent(in), optional  knots1,
integer(kind=i32), intent(in), optional  size_knots1,
real(kind=f64), dimension(:), intent(in), optional  knots2,
integer(kind=i32), intent(in), optional  size_knots2 
)
private

Definition at line 353 of file sll_m_hermite_interpolator_2d.F90.

    Report Typos and Errors