Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
|
Provides capabilities for data and derivative interpolation with cubic B-splines and different boundary conditions.
Supported boundary conditions presently: periodic, hermite. User needs to declare a pointer to the spline object. The pointer needs to be allocated and initialized with the appropriate 'new_' function. Afterwards this initialized object can be used to compute the spline coefficients and interpolate data at a given point or collections of points. Derivatives can also be interpolated. After usage, resources allocated for the spline object can be liberated by a call to 'sll_o_delete()'. More details by following the link sll_m_cubic_splines
Derived types and interfaces | |
type | sll_t_cubic_spline_1d |
basic type for one-dimensional cubic spline data. More... | |
type | sll_t_cubic_spline_2d |
Basic type for one-dimensional cubic spline data. More... | |
Functions/Subroutines | |
subroutine, public | sll_s_cubic_spline_1d_init (self, num_points, xmin, xmax, bc_type, sl, sr, fast_algorithm) |
Returns a pointer to a heap-allocated cubic spline object. More... | |
subroutine, public | sll_s_cubic_spline_1d_compute_interpolant (f, spline) |
Computes the cubic spline coefficients corresponding to the 1D data array 'f'. More... | |
subroutine | compute_spline_1d_periodic_aux (f, num_pts, d, coeffs) |
subroutine | compute_spline_1d_hermite_aux (f, num_pts, d, slope_l, slope_r, delta, coeffs) |
subroutine | compute_spline_1d_periodic (f, spline) |
subroutine | compute_spline_1d_hermite (f, spline) |
real(kind=f64) function | interpolate_value_aux (x, xmin, rh, coeffs) |
real(kind=f64) function, public | sll_f_cubic_spline_1d_eval (x, spline) |
returns the value of the interpolated image of the abscissa x, using the spline coefficients stored in the spline object. More... | |
subroutine, public | sll_s_cubic_spline_1d_eval_array (a_in, a_out, n, spline) |
returns the values of the images of a collection of abscissae, represented by a 1D array in another output array. The spline coefficients used are stored in the spline object pointer. More... | |
subroutine | interpolate_pointer_values (ptr_in, ptr_out, n, spline) |
returns the values of the images of a collection of abscissae, More... | |
real(kind=f64) function | interpolate_derivative_aux (x, xmin, rh, coeffs) |
real(kind=f64) function, public | sll_f_cubic_spline_1d_eval_deriv (x, spline) |
returns the value of the derivative at the image of the abscissa 'x', using the spline coefficients stored in the spline object. More... | |
subroutine, public | sll_s_cubic_spline_1d_eval_array_deriv (array_in, array_out, num_pts, spline) |
returns the values of the derivatives evaluated at a collection of abscissae stored in the input 1D array and stores the results in a different output array. The spline coefficients used are stored in the spline object pointer. More... | |
subroutine | interpolate_pointer_derivatives (ptr_in, ptr_out, num_pts, spline) |
analogous to the sll_s_cubic_spline_1d_eval_array_deriv() subroutine but its input and output arrays are pointers. More... | |
subroutine, public | sll_s_cubic_spline_1d_free (spline) |
deallocate the sll_t_cubic_spline_1d object More... | |
type(sll_t_cubic_spline_2d) function, pointer | sll_f_new_cubic_spline_2d (num_pts_x1, num_pts_x2, x1_min, x1_max, x2_min, x2_max, x1_bc_type, x2_bc_type, const_slope_x1_min, const_slope_x1_max, const_slope_x2_min, const_slope_x2_max, x1_min_slopes, x1_max_slopes, x2_min_slopes, x2_max_slopes) |
Returns a pointer to a heap-allocated 2D cubic spline object. More... | |
subroutine, public | sll_s_cubic_spline_2d_init (this, num_pts_x1, num_pts_x2, x1_min, x1_max, x2_min, x2_max, x1_bc_type, x2_bc_type, const_slope_x1_min, const_slope_x1_max, const_slope_x2_min, const_slope_x2_max, x1_min_slopes, x1_max_slopes, x2_min_slopes, x2_max_slopes) |
subroutine | compute_spline_2d_prdc_prdc (data, spline) |
subroutine | compute_spline_2d_hrmt_prdc (data, spline) |
subroutine | compute_spline_2d_prdc_hrmt (data, spline) |
subroutine | compute_spline_2d_hrmt_hrmt (data, spline) |
subroutine, public | sll_s_cubic_spline_2d_compute_interpolant (data, spline) |
Computes the spline coefficients for the given data. The coeffcients are first computed in the second direction of the array (i.e. the x2 direction) and then in the first (x1) direction. More... | |
subroutine, public | sll_s_cubic_spline_2d_deposit_value (x1, x2, spline, a_out) |
Updated distribution function at time \( t^{n+1} \). More... | |
real(kind=f64) function, public | sll_f_cubic_spline_2d_eval (x1, x2, spline) |
Returns the interpolated value of the image of the point (x1,x2) using the spline decomposition stored in the spline object. More... | |
real(kind=f64) function, public | sll_f_cubic_spline_2d_eval_deriv_x1 (x1, x2, spline) |
Returns the interpolated value of the derivative in the x1 direction at the point (x1,x2) using the spline decomposition stored in the spline object. More... | |
real(kind=f64) function, public | sll_f_cubic_spline_2d_eval_deriv_x2 (x1, x2, spline) |
Returns the interpolated value of the derivative. More... | |
subroutine, public | sll_s_cubic_spline_2d_get_coeff (spline, coeff) |
subroutine, public | sll_s_cubic_spline_2d_free (spline) |
subroutine, public | sll_s_cubic_spline_1d_eval_disp (spline, alpha, output_array) |
Computes the interpolated values at each grid point replaced by alpha for the precomputed spline coefficients. More... | |
subroutine | spline_interpolate_from_interpolant_cell_dx (spline, cell, dx, out) |
Helper function for sll_s_cubic_spline_1d_eval_disp: evaluate spline in given cell and normalized displacement dx. More... | |
Variables | |
real(kind=f64), parameter | inv_6 = 1._f64/6._f64 |
|
private |
Definition at line 692 of file sll_m_cubic_splines.F90.
|
private |
|
private |
Definition at line 654 of file sll_m_cubic_splines.F90.
|
private |
|
private |
Definition at line 1795 of file sll_m_cubic_splines.F90.
|
private |
Definition at line 1582 of file sll_m_cubic_splines.F90.
|
private |
Definition at line 1699 of file sll_m_cubic_splines.F90.
|
private |
Definition at line 1526 of file sll_m_cubic_splines.F90.
|
private |
Definition at line 1025 of file sll_m_cubic_splines.F90.
|
private |
analogous to the sll_s_cubic_spline_1d_eval_array_deriv() subroutine but its input and output arrays are pointers.
[in] | ptr_in | input double-precison element array pointer containing the abscissae at which the derivatives are wanted. |
[out] | ptr_out | output double-precision element array pointer containing the results. |
[in] | num_pts | the number of elements of the input array pointer. |
[in,out] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant() subroutine. |
Definition at line 1128 of file sll_m_cubic_splines.F90.
|
private |
returns the values of the images of a collection of abscissae,
represented by a 1D array pointer, in another output array pointer. The spline coefficients used are stored in the spline object pointer.
[in] | ptr_in | input double-precison element pointer containing the abscissae to be interpolated. |
[out] | ptr_out | output double-precision element pointer containing the results of the interpolation. |
[in] | n | the number of elements of the input pointer which are to be interpolated. |
[in,out] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant() subroutine. |
Definition at line 971 of file sll_m_cubic_splines.F90.
|
private |
real(kind=f64) function, public sll_m_cubic_splines::sll_f_cubic_spline_1d_eval | ( | real(kind=f64), intent(in) | x, |
type(sll_t_cubic_spline_1d) | spline | ||
) |
returns the value of the interpolated image of the abscissa x, using the spline coefficients stored in the spline object.
[in] | x | the value of the abscissa where the data should be interpolated. |
[in] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant subroutine. |
Definition at line 874 of file sll_m_cubic_splines.F90.
real(kind=f64) function, public sll_m_cubic_splines::sll_f_cubic_spline_1d_eval_deriv | ( | real(kind=f64), intent(in) | x, |
type(sll_t_cubic_spline_1d) | spline | ||
) |
returns the value of the derivative at the image of the abscissa 'x', using the spline coefficients stored in the spline object.
[in] | x | the value of the abscissa at which the derivative should be interpolated. |
[in,out] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant() subroutine. |
Definition at line 1065 of file sll_m_cubic_splines.F90.
real(kind=f64) function, public sll_m_cubic_splines::sll_f_cubic_spline_2d_eval | ( | real(kind=f64), intent(in) | x1, |
real(kind=f64), intent(in) | x2, | ||
type(sll_t_cubic_spline_2d), intent(in) | spline | ||
) |
Returns the interpolated value of the image of the point (x1,x2) using the spline decomposition stored in the spline object.
[in] | x1 | first coordinate. |
[in] | x2 | second coordinate. |
[in] | spline | pointer to spline object. |
Definition at line 2328 of file sll_m_cubic_splines.F90.
real(kind=f64) function, public sll_m_cubic_splines::sll_f_cubic_spline_2d_eval_deriv_x1 | ( | real(kind=f64), intent(in) | x1, |
real(kind=f64), intent(in) | x2, | ||
type(sll_t_cubic_spline_2d) | spline | ||
) |
Returns the interpolated value of the derivative in the x1 direction at the point (x1,x2) using the spline decomposition stored in the spline object.
[in] | x1 | first coordinate. |
[in] | x2 | second coordinate. |
[in] | spline | pointer to spline object. |
Definition at line 2412 of file sll_m_cubic_splines.F90.
real(kind=f64) function, public sll_m_cubic_splines::sll_f_cubic_spline_2d_eval_deriv_x2 | ( | real(kind=f64), intent(in) | x1, |
real(kind=f64), intent(in) | x2, | ||
type(sll_t_cubic_spline_2d) | spline | ||
) |
Returns the interpolated value of the derivative.
in the x2 direction at the point(x1,x2) using the spline decomposition stored in the spline object.
[in] | x1 | first coordinate. |
[in] | x2 | second coordinate. |
[in] | spline | pointer to spline object. |
Definition at line 2498 of file sll_m_cubic_splines.F90.
|
private |
Returns a pointer to a heap-allocated 2D cubic spline object.
[in] | num_pts_x1 | Dimension in the x1 direction of the 2D array that stores the data to be interpolated. |
[in] | num_pts_x2 | Dimension in the x2 direction of the 2D array that stores the data to be interpolated. |
[in] | x1_min | In the x1 direction, the minimum value of the domain where the data to be interpolated are represented. to be interpolated. |
[in] | x1_max | In the x1 direction, the maximum value of the domain where the data to be interpolated are represented. to be interpolated. |
[in] | x2_min | In the x2 direction, the minimum value of the domain where the data to be interpolated are represented. to be interpolated. |
[in] | x2_max | In the x2 direction, the maximum value of the domain where the data to be interpolated are represented. to be interpolated. |
[in] | x1_bc_type | Boundary condition specifier in x1 direction. Must be one of the symbols defined in the SLL_BOUNDARY_CONDITION_DESCRIPTORS module. |
[in] | x2_bc_type | Boundary condition specifier in x2 direction. Must be one of the symbols defined in the SLL_BOUNDARY_CONDITION_DESCRIPTORS module. |
[in] | const_slope_x1_min | OPTIONAL: The value of the slope at x1_min, in case that a constant slope value should to be imposed. |
[in] | const_slope_x1_max | OPTIONAL: The value of the slope at x1_max, in case that a constant slope value should to be imposed. |
[in] | const_slope_x2_min | OPTIONAL: The value of the slope at x2_min, in case that a constant slope value should to be imposed. |
[in] | const_slope_x2_max | OPTIONAL: The value of the slope at x2_max, in case that a constant slope value should to be imposed. of hermite boundary conditions. |
[in] | x1_min_slopes | OPTIONAL: The values of the slopes at x1_min, in case that a specific slope value should to be imposed at each border point. Default behavior is to compute the slope consistent with the given data. |
[in] | x1_max_slopes | OPTIONAL: The values of the slopes at x1_max, in case that a specific slope value should to be imposed at each border point. Default behavior is to compute the slope consistent with the given data. |
[in] | x2_min_slopes | OPTIONAL: The values of the slopes at x2_min, in case that a specific slope value should to be imposed at each border point. Default behavior is to compute the slope consistent with the given data. |
[in] | x2_max_slopes | OPTIONAL: The values of the slopes at x2_max, in case that a specific slope value should to be imposed at each border point. Default behavior is to compute the slope consistent with the given data. |
Definition at line 1229 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_compute_interpolant | ( | intent(in) | f, |
type(sll_t_cubic_spline_1d) | spline | ||
) |
Computes the cubic spline coefficients corresponding to the 1D data array 'f'.
[in] | f | a 1D array with the data to interpolate. Its domain of definition and number of points should conform to the values given when the spline object was initialized. |
[in,out] | spline | a 1D cubic spline object pointer, previously initialized. |
Definition at line 498 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_eval_array | ( | real(kind=f64), dimension(1:n), intent(in) | a_in, |
real(kind=f64), dimension(1:n), intent(out) | a_out, | ||
integer(kind=i32), intent(in) | n, | ||
type(sll_t_cubic_spline_1d) | spline | ||
) |
returns the values of the images of a collection of abscissae, represented by a 1D array in another output array. The spline coefficients used are stored in the spline object pointer.
[in] | a_in | input double-precison element array containing the abscissae to be interpolated. |
[out] | a_out | output double-precision element array containing the results of the interpolation. |
[in] | n | the number of elements of the input array which are to be interpolated. |
[in,out] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant() subroutine. |
Definition at line 903 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_eval_array_deriv | ( | real(kind=f64), dimension(:), intent(in) | array_in, |
real(kind=f64), dimension(:), intent(out) | array_out, | ||
integer(kind=i32), intent(in) | num_pts, | ||
type(sll_t_cubic_spline_1d) | spline | ||
) |
returns the values of the derivatives evaluated at a collection of abscissae stored in the input 1D array and stores the results in a different output array. The spline coefficients used are stored in the spline object pointer.
[in] | array_in | input double-precison element array containing the abscissae at which the derivatives are wanted. |
[out] | array_out | output double-precision element array containing the results. |
[in] | num_pts | the number of elements of the input array. |
[in,out] | spline | the spline object pointer, duly initialized and already operated on by the sll_s_cubic_spline_1d_compute_interpolant() subroutine. |
Definition at line 1094 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_eval_disp | ( | type(sll_t_cubic_spline_1d), intent(in) | spline, |
real(kind=f64), intent(in) | alpha, | ||
real(kind=f64), dimension(:), intent(out) | output_array | ||
) |
Computes the interpolated values at each grid point replaced by alpha for the precomputed spline coefficients.
[in] | spline | spline object |
[in] | alpha | displacement |
[out] | output_array | output_array(i) holds the interpolated value at x(i)+alpha on output |
Definition at line 2616 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_free | ( | type(sll_t_cubic_spline_1d) | spline | ) |
deallocate the sll_t_cubic_spline_1d object
Definition at line 1155 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_1d_init | ( | type(sll_t_cubic_spline_1d) | self, |
integer(kind=i32), intent(in) | num_points, | ||
real(kind=f64), intent(in) | xmin, | ||
real(kind=f64), intent(in) | xmax, | ||
integer(kind=i32), intent(in) | bc_type, | ||
real(kind=f64), intent(in), optional | sl, | ||
real(kind=f64), intent(in), optional | sr, | ||
logical, intent(in), optional | fast_algorithm | ||
) |
Returns a pointer to a heap-allocated cubic spline object.
[in] | num_points | Number of points where the data to be interpolated are represented. |
[in] | xmin | Minimum value of the abscissae where the data are meant to be interpolated. |
[in] | xmax | Maximum value of the abscissae where the data are meant to be interpolated. |
[in] | bc_type | A boundary condition specifier. Must be one of the symbols defined in the SLL_BOUNDARY_CONDITION_DESCRIPTORS module. |
[in] | sl | OPTIONAL: The value of the slope at xmin, for use in the case of hermite boundary conditions. |
[in] | sr | OPTIONAL: The value of the slope at xmin, for use in the case of hermite boundary conditions. |
[in] | sr | OPTIONAL: Flag to change between fast algorithm and LU-based algorithm (for num_points>NUM_TERMS). Default: true for num_points>NUM_TERMS. |
Definition at line 249 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_2d_compute_interpolant | ( | real(kind=f64), dimension(:, :), intent(in), target | data, |
type(sll_t_cubic_spline_2d) | spline | ||
) |
Computes the spline coefficients for the given data. The coeffcients are first computed in the second direction of the array (i.e. the x2 direction) and then in the first (x1) direction.
[in] | data | The 2D array with the data for which the cubic spline decomposition is sought. |
[in,out] | spline | a pointer to an initialized spline object. |
Definition at line 1964 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_2d_deposit_value | ( | real(kind=f64), dimension(1:, 1:), intent(in) | x1, |
real(kind=f64), dimension(1:, 1:), intent(in) | x2, | ||
type(sll_t_cubic_spline_2d) | spline, | ||
real(kind=f64), dimension(:, :), intent(out) | a_out | ||
) |
Updated distribution function at time \( t^{n+1} \).
sll_s_cubic_spline_2d_deposit_value(): given a spline that describes the decomposition of the distribution function at time \( t^n\) , and two 2D arrays x1 and x2 where the foot of the forward characteristics are stored, returns a 2D array a_out which is the updated distribution function at time \( t^{n+1} \).
the boundary conditions are taken into account and any type of BC are allowed
Definition at line 2036 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_2d_free | ( | type(sll_t_cubic_spline_2d) | spline | ) |
Definition at line 2582 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_2d_get_coeff | ( | type(sll_t_cubic_spline_2d) | spline, |
real(kind=f64), dimension(:), intent(out) | coeff | ||
) |
Definition at line 2562 of file sll_m_cubic_splines.F90.
subroutine, public sll_m_cubic_splines::sll_s_cubic_spline_2d_init | ( | type(sll_t_cubic_spline_2d) | this, |
integer(kind=i32), intent(in) | num_pts_x1, | ||
integer(kind=i32), intent(in) | num_pts_x2, | ||
real(kind=f64), intent(in) | x1_min, | ||
real(kind=f64), intent(in) | x1_max, | ||
real(kind=f64), intent(in) | x2_min, | ||
real(kind=f64), intent(in) | x2_max, | ||
integer(kind=i32), intent(in) | x1_bc_type, | ||
integer(kind=i32), intent(in) | x2_bc_type, | ||
real(kind=f64), intent(in), optional | const_slope_x1_min, | ||
real(kind=f64), intent(in), optional | const_slope_x1_max, | ||
real(kind=f64), intent(in), optional | const_slope_x2_min, | ||
real(kind=f64), intent(in), optional | const_slope_x2_max, | ||
real(kind=f64), dimension(:), intent(in), optional | x1_min_slopes, | ||
real(kind=f64), dimension(:), intent(in), optional | x1_max_slopes, | ||
real(kind=f64), dimension(:), intent(in), optional | x2_min_slopes, | ||
real(kind=f64), dimension(:), intent(in), optional | x2_max_slopes | ||
) |
Definition at line 1287 of file sll_m_cubic_splines.F90.
|
private |
Helper function for sll_s_cubic_spline_1d_eval_disp: evaluate spline in given cell and normalized displacement dx.
[in] | spline | spline object |
[in] | dx | normalized displacement |
Definition at line 2685 of file sll_m_cubic_splines.F90.
|
private |
Definition at line 82 of file sll_m_cubic_splines.F90.