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

Description

Interpolator 1d using cubic splines on regular mesh with halo cells.

Author
Katharina Kormann, IPP

The module provides an optimized implementation of the interpolator function interpolate_array_disp for a domain with halo cells, i.e. we do not have explicit boundary conditions but compute the Lagrange interpolation for a number of cells provided that enough cells around this domain are present such that no boundary conditions need to be imposed. The module also provides a version for periodic boundary conditions. The implementation is based on the algorithms described in Section 5.4.4 (Fast local spline interpolation) of Sonnendrücker, Numerical Methods for the Vlasov-Maxwell equations, to appear.

Functions/Subroutines

subroutine, public sll_s_cubic_spline_halo_1d_prepare_exchange (fdata, si, num_points, d_0, c_np2)
 Compute the part of d(0) and c(num_points+2) that need to be send to neighboring processors. More...
 
subroutine, public sll_s_cubic_spline_halo_1d_finish_boundary_conditions (fdata, si, num_points, d_0, c_np2)
 Complete d(0) and c(num_points+2) with local data after their values have been received from the neighboring processors. More...
 
subroutine, public sll_s_cubic_spline_halo_1d_compute_interpolant (f, num_points, d, coeffs)
 Compute the coefficients of the local interpolating spline (after d(0) and c(num_points+2) have been computed. More...
 
subroutine, public sll_s_cubic_spline_halo_1d_eval_disp (coeffs, alpha, num_points, fout)
 This function corresponds to the interpolate_array_disp function of the interpolators but the displacement is normalized and between [0,1]. More...
 
subroutine, public sll_s_cubic_spline_halo_1d_periodic (fin, alpha, num_cells, fout)
 
subroutine sll_s_cubic_spline_halo_1d_compute_periodic (num_points, d, f_coeffs)
 Compute the coefficients of the local interpolating spline (after d(0) and c(num_points+2) have been computed. More...
 

Variables

real(kind=f64), parameter p_a = sqrt((2.0_f64 + sqrt(3.0_f64))/6.0_f64)
 
real(kind=f64), parameter p_r_a = 1.0_f64/p_a
 
real(kind=f64), parameter p_b = sqrt((2.0_f64 - sqrt(3.0_f64))/6.0_f64)
 
real(kind=f64), parameter p_b_a = p_b/p_a
 
real(kind=f64), parameter p_sqrt3 = sqrt(3.0_f64)
 
real(kind=f64), parameter p_inv_6 = 1._f64/6._f64
 
real(kind=f64), dimension(0:27), parameter pba_pow = [ (-p_b_a)**0, (-p_b_a)**1, (-p_b_a)**2, (-p_b_a)**3, (-p_b_a)**4, (-p_b_a)**5, (-p_b_a)**6, (-p_b_a)**7, (-p_b_a)**8, (-p_b_a)**9, (-p_b_a)**10, (-p_b_a)**11, (-p_b_a)**12, (-p_b_a)**13, (-p_b_a)**14, (-p_b_a)**15, (-p_b_a)**16, (-p_b_a)**17, (-p_b_a)**18, (-p_b_a)**19, (-p_b_a)**20, (-p_b_a)**21, (-p_b_a)**22, (-p_b_a)**23, (-p_b_a)**24, (-p_b_a)**25, (-p_b_a)**26, (-p_b_a)**27 ]
 

Function/Subroutine Documentation

◆ sll_s_cubic_spline_halo_1d_compute_interpolant()

subroutine, public sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_compute_interpolant ( real(kind=f64), dimension(0:), intent(in)  f,
integer(kind=i32), intent(in)  num_points,
real(kind=f64), dimension(0:), intent(out)  d,
real(kind=f64), dimension(0:), intent(out)  coeffs 
)

Compute the coefficients of the local interpolating spline (after d(0) and c(num_points+2) have been computed.

Parameters
[in]fdata values including all points to be interpolated at
[in]num_pointsnumber of local data points
[out]dhelper variable for spline coefficients (vales from forward recursion)
[out]coeffsspline coefficients

Definition at line 143 of file sll_m_cubic_spline_halo_1d.F90.

Here is the caller graph for this function:

◆ sll_s_cubic_spline_halo_1d_compute_periodic()

subroutine sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_compute_periodic ( integer(kind=i32), intent(in)  num_points,
real(kind=f64), dimension(0:), intent(out)  d,
real(kind=f64), dimension(0:), intent(inout)  f_coeffs 
)
private

Compute the coefficients of the local interpolating spline (after d(0) and c(num_points+2) have been computed.

Parameters
[in]num_pointsnumber of local data points
[out]dhelper variable for spline coefficients (vales from forward recursion)
[in,out]f_coeffson input: data values including all points to be interpolated at; on output: spline coefficients

Definition at line 216 of file sll_m_cubic_spline_halo_1d.F90.

Here is the caller graph for this function:

◆ sll_s_cubic_spline_halo_1d_eval_disp()

subroutine, public sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_eval_disp ( real(kind=f64), dimension(0:), intent(in)  coeffs,
real(kind=f64), intent(in)  alpha,
integer(kind=i32), intent(in)  num_points,
real(kind=f64), dimension(:), intent(out)  fout 
)

This function corresponds to the interpolate_array_disp function of the interpolators but the displacement is normalized and between [0,1].

Parameters
[in]coeffsSpline coefficients (centered around the cell into which we displace, i.e. integer part of displacement is already build in here)
[in]alphaDisplacement normalized by dx and only remainder of modulo 1 (in [0,1])
[in]num_pointsNumber of local data points
[out]foutInterpolated values

Definition at line 179 of file sll_m_cubic_spline_halo_1d.F90.

Here is the caller graph for this function:

◆ sll_s_cubic_spline_halo_1d_finish_boundary_conditions()

subroutine, public sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_finish_boundary_conditions ( real(kind=f64), dimension(:), intent(in)  fdata,
integer(kind=i32), intent(in)  si,
integer(kind=i32), intent(in)  num_points,
intent(out)  d_0,
intent(inout)  c_np2 
)

Complete d(0) and c(num_points+2) with local data after their values have been received from the neighboring processors.

Parameters
[in]fdataLocal data values
[in]siInteger part of the shift
[in]num_pointsnumber of local data values

Definition at line 110 of file sll_m_cubic_spline_halo_1d.F90.

Here is the caller graph for this function:

◆ sll_s_cubic_spline_halo_1d_periodic()

subroutine, public sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_periodic ( real(kind=f64), dimension(0:), intent(inout)  fin,
real(kind=f64), intent(in)  alpha,
integer(kind=i32), intent(in)  num_cells,
real(kind=f64), dimension(:), intent(out)  fout 
)

Definition at line 204 of file sll_m_cubic_spline_halo_1d.F90.

Here is the call graph for this function:

◆ sll_s_cubic_spline_halo_1d_prepare_exchange()

subroutine, public sll_m_cubic_spline_halo_1d::sll_s_cubic_spline_halo_1d_prepare_exchange ( real(kind=f64), dimension(:), intent(in)  fdata,
integer(kind=i32), intent(in)  si,
integer(kind=i32), intent(in)  num_points,
intent(out)  d_0,
intent(out)  c_np2 
)

Compute the part of d(0) and c(num_points+2) that need to be send to neighboring processors.

Parameters
[in]fdataLocal data values
[in]siInteger part of the shift
[in]num_pointsnumber of local data values

Definition at line 71 of file sll_m_cubic_spline_halo_1d.F90.

Here is the caller graph for this function:

Variable Documentation

◆ p_a

real(kind=f64), parameter p_a = sqrt((2.0_f64 + sqrt(3.0_f64))/6.0_f64)
private

Definition at line 51 of file sll_m_cubic_spline_halo_1d.F90.

◆ p_b

real(kind=f64), parameter p_b = sqrt((2.0_f64 - sqrt(3.0_f64))/6.0_f64)
private

Definition at line 53 of file sll_m_cubic_spline_halo_1d.F90.

◆ p_b_a

real(kind=f64), parameter p_b_a = p_b/p_a
private

Definition at line 54 of file sll_m_cubic_spline_halo_1d.F90.

◆ p_inv_6

real(kind=f64), parameter p_inv_6 = 1._f64/6._f64
private

Definition at line 56 of file sll_m_cubic_spline_halo_1d.F90.

◆ p_r_a

real(kind=f64), parameter p_r_a = 1.0_f64/p_a
private

Definition at line 52 of file sll_m_cubic_spline_halo_1d.F90.

◆ p_sqrt3

real(kind=f64), parameter p_sqrt3 = sqrt(3.0_f64)
private

Definition at line 55 of file sll_m_cubic_spline_halo_1d.F90.

◆ pba_pow

real(kind=f64), dimension(0:27), parameter pba_pow = [ (-p_b_a)**0, (-p_b_a)**1, (-p_b_a)**2, (-p_b_a)**3, (-p_b_a)**4, (-p_b_a)**5, (-p_b_a)**6, (-p_b_a)**7, (-p_b_a)**8, (-p_b_a)**9, (-p_b_a)**10, (-p_b_a)**11, (-p_b_a)**12, (-p_b_a)**13, (-p_b_a)**14, (-p_b_a)**15, (-p_b_a)**16, (-p_b_a)**17, (-p_b_a)**18, (-p_b_a)**19, (-p_b_a)**20, (-p_b_a)**21, (-p_b_a)**22, (-p_b_a)**23, (-p_b_a)**24, (-p_b_a)**25, (-p_b_a)**26, (-p_b_a)**27 ]
private

Definition at line 60 of file sll_m_cubic_spline_halo_1d.F90.

    Report Typos and Errors