DDC 0.1.0
Loading...
Searching...
No Matches
ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX > Class Template Reference

A class for creating a spline approximation of a function. More...

#include <spline_builder.hpp>

Collaboration diagram for ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >:
Collaboration graph

Public Types

using exec_space = ExecSpace
 The type of the Kokkos execution space used by this class.
 
using memory_space = MemorySpace
 The type of the Kokkos memory space used by this class.
 
using continuous_dimension_type = typename InterpolationDDim::continuous_dimension_type
 The type of the interpolation continuous dimension (continuous dimension of interest) used by this class.
 
using interpolation_discrete_dimension_type = InterpolationDDim
 The type of the interpolation discrete dimension (discrete dimension of interest) used by this class.
 
using bsplines_type = BSplines
 The discrete dimension representing the B-splines.
 
using deriv_type = ddc::Deriv< continuous_dimension_type >
 The type of the Deriv dimension at the boundaries.
 
using interpolation_domain_type = ddc::DiscreteDomain< interpolation_discrete_dimension_type >
 The type of the domain for the 1D interpolation mesh used by this class.
 
using batched_interpolation_domain_type = ddc::DiscreteDomain< IDimX... >
 The type of the whole domain representing interpolation points.
 
using batch_domain_type = ddc::remove_dims_of_t< batched_interpolation_domain_type, interpolation_discrete_dimension_type >
 The type of the batch domain (obtained by removing the dimension of interest from the whole domain).
 
using batched_spline_domain_type = ddc::replace_dim_of_t< batched_interpolation_domain_type, interpolation_discrete_dimension_type, bsplines_type >
 The type of the whole spline domain (cartesian product of 1D spline domain and batch domain) preserving the underlying memory layout (order of dimensions).
 
using batched_derivs_domain_type = ddc::replace_dim_of_t< batched_interpolation_domain_type, interpolation_discrete_dimension_type, deriv_type >
 The type of the whole Deriv domain (cartesian product of 1D Deriv domain and batch domain) preserving the underlying memory layout (order of dimensions).
 

Public Member Functions

 SplineBuilder (batched_interpolation_domain_type const &batched_interpolation_domain, std::optional< std::size_t > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditioner_max_block_size=std::nullopt)
 Build a SplineBuilder acting on batched_interpolation_domain.
 
 SplineBuilder (SplineBuilder const &x)=delete
 Copy-constructor is deleted.
 
 SplineBuilder (SplineBuilder &&x)=default
 Move-constructs.
 
 ~SplineBuilder ()=default
 Destructs.
 
SplineBuilderoperator= (SplineBuilder const &x)=delete
 Copy-assignment is deleted.
 
SplineBuilderoperator= (SplineBuilder &&x)=default
 Move-assigns.
 
interpolation_domain_type interpolation_domain () const noexcept
 Get the domain for the 1D interpolation mesh used by this class.
 
batched_interpolation_domain_type batched_interpolation_domain () const noexcept
 Get the whole domain representing interpolation points.
 
batch_domain_type batch_domain () const noexcept
 Get the batch domain.
 
ddc::DiscreteDomain< bsplines_typespline_domain () const noexcept
 Get the 1D domain on which spline coefficients are defined.
 
batched_spline_domain_type batched_spline_domain () const noexcept
 Get the whole domain on which spline coefficients are defined.
 
batched_derivs_domain_type batched_derivs_xmin_domain () const noexcept
 Get the whole domain on which derivatives on lower boundary are defined.
 
batched_derivs_domain_type batched_derivs_xmax_domain () const noexcept
 Get the whole domain on which derivatives on upper boundary are defined.
 
const ddc::detail::SplinesLinearProblem< exec_space > & get_interpolation_matrix () const noexcept
 Get the interpolation matrix.
 
template<class Layout >
void operator() (ddc::ChunkSpan< double, batched_spline_domain_type, Layout, memory_space > spline, ddc::ChunkSpan< double const, batched_interpolation_domain_type, Layout, memory_space > vals, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > derivs_xmin=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > derivs_xmax=std::nullopt) const
 Compute a spline approximation of a function.
 
template<class OutMemorySpace = MemorySpace>
std::tuple< ddc::Chunk< double, ddc::DiscreteDomain< ddc::Deriv< typename InterpolationDDim::continuous_dimension_type > >, ddc::KokkosAllocator< double, OutMemorySpace > >, ddc::Chunk< double, ddc::DiscreteDomain< InterpolationDDim >, ddc::KokkosAllocator< double, OutMemorySpace > >, ddc::Chunk< double, ddc::DiscreteDomain< ddc::Deriv< typename InterpolationDDim::continuous_dimension_type > >, ddc::KokkosAllocator< double, OutMemorySpace > > > quadrature_coefficients () const
 Compute the quadrature coefficients associated to the b-splines used by this SplineBuilder.
 

Static Public Attributes

static constexpr bool s_odd = BSplines::degree() % 2
 Indicates if the degree of the splines is odd or even.
 
static constexpr int s_nbc_xmin = n_boundary_equations(BcLower, BSplines::degree())
 The number of equations defining the boundary condition at the lower bound.
 
static constexpr int s_nbc_xmax = n_boundary_equations(BcUpper, BSplines::degree())
 The number of equations defining the boundary condition at the upper bound.
 
static constexpr ddc::BoundCond s_bc_xmin = BcLower
 The boundary condition implemented at the lower bound.
 
static constexpr ddc::BoundCond s_bc_xmax = BcUpper
 The boundary condition implemented at the upper bound.
 
static constexpr SplineSolver s_spline_solver = Solver
 The SplineSolver giving the backend used to perform the spline approximation.
 

Detailed Description

template<class ExecSpace, class MemorySpace, class BSplines, class InterpolationDDim, ddc::BoundCond BcLower, ddc::BoundCond BcUpper, SplineSolver Solver, class... IDimX>
class ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >

A class for creating a spline approximation of a function.

A class which contains an operator () which can be used to build a spline approximation of a function. A spline approximation is represented by coefficients stored in a Chunk of B-splines. The spline is constructed such that it respects the boundary conditions BcLower and BcUpper, and it interpolates the function at the points on the interpolation_discrete_dimension associated with interpolation_discrete_dimension_type.

Template Parameters
ExecSpaceThe Kokkos execution space on which the spline approximation is performed.
MemorySpaceThe Kokkos memory space on which the data (interpolation function and splines coefficients) is stored.
BSplinesThe discrete dimension representing the B-splines.
InterpolationDDimThe discrete dimension on which interpolation points are defined.
BcLowerThe lower boundary condition.
BcUpperThe upper boundary condition.
SolverThe SplineSolver giving the backend used to perform the spline approximation.
IDimXA variadic template of all the discrete dimensions forming the full space (InterpolationDDim + batched dimensions).

Definition at line 63 of file spline_builder.hpp.

Member Typedef Documentation

◆ exec_space

The type of the Kokkos execution space used by this class.

Definition at line 73 of file spline_builder.hpp.

◆ memory_space

The type of the Kokkos memory space used by this class.

Definition at line 76 of file spline_builder.hpp.

◆ continuous_dimension_type

template<class ExecSpace , class MemorySpace , class BSplines , class InterpolationDDim , ddc::BoundCond BcLower, ddc::BoundCond BcUpper, SplineSolver Solver, class... IDimX>
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::continuous_dimension_type = typename InterpolationDDim::continuous_dimension_type

The type of the interpolation continuous dimension (continuous dimension of interest) used by this class.

Definition at line 79 of file spline_builder.hpp.

◆ interpolation_discrete_dimension_type

The type of the interpolation discrete dimension (discrete dimension of interest) used by this class.

Definition at line 82 of file spline_builder.hpp.

◆ bsplines_type

The discrete dimension representing the B-splines.

Definition at line 85 of file spline_builder.hpp.

◆ deriv_type

The type of the Deriv dimension at the boundaries.

Definition at line 88 of file spline_builder.hpp.

◆ interpolation_domain_type

The type of the domain for the 1D interpolation mesh used by this class.

Definition at line 91 of file spline_builder.hpp.

◆ batched_interpolation_domain_type

The type of the whole domain representing interpolation points.

Definition at line 94 of file spline_builder.hpp.

◆ batch_domain_type

The type of the batch domain (obtained by removing the dimension of interest from the whole domain).

Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and a dimension of interest Y, this is DiscreteDomain<X,Z>

Definition at line 103 of file spline_builder.hpp.

◆ batched_spline_domain_type

The type of the whole spline domain (cartesian product of 1D spline domain and batch domain) preserving the underlying memory layout (order of dimensions).

Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and a dimension of interest Y (associated to a B-splines tag BSplinesY), this is DiscreteDomain<X,BSplinesY,Z>.

Definition at line 114 of file spline_builder.hpp.

◆ batched_derivs_domain_type

The type of the whole Deriv domain (cartesian product of 1D Deriv domain and batch domain) preserving the underlying memory layout (order of dimensions).

Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and a dimension of interest Y, this is DiscreteDomain<X,Deriv<Y>,Z>

Definition at line 142 of file spline_builder.hpp.

Constructor & Destructor Documentation

◆ SplineBuilder() [1/3]

template<class ExecSpace , class MemorySpace , class BSplines , class InterpolationDDim , ddc::BoundCond BcLower, ddc::BoundCond BcUpper, SplineSolver Solver, class... IDimX>
ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::SplineBuilder ( batched_interpolation_domain_type const batched_interpolation_domain,
std::optional< std::size_t >  cols_per_chunk = std::nullopt,
std::optional< unsigned int preconditioner_max_block_size = std::nullopt 
)
inlineexplicit

Build a SplineBuilder acting on batched_interpolation_domain.

Parameters
batched_interpolation_domainThe domain on which the interpolation points are defined.
cols_per_chunkA parameter used by the slicer (internal to the solver) to define the size of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated by the linear solver one-after-the-other). This value is optional. If no value is provided then the default value is chosen by the requested solver.
preconditioner_max_block_sizeA parameter used by the slicer (internal to the solver) to define the size of a block used by the Block-Jacobi preconditioner. This value is optional. If no value is provided then the default value is chosen by the requested solver.
See also
MatrixSparse

Definition at line 195 of file spline_builder.hpp.

◆ SplineBuilder() [2/3]

◆ SplineBuilder() [3/3]

Move-constructs.

Parameters
xAn rvalue to another SplineBuilder.

◆ ~SplineBuilder()

Member Function Documentation

◆ operator=() [1/2]

◆ operator=() [2/2]

Move-assigns.

Parameters
xAn rvalue to another SplineBuilder.
Returns
A reference to this object.

◆ interpolation_domain()

Get the domain for the 1D interpolation mesh used by this class.

This is 1D because it is defined along the dimension of interest.

Returns
The 1D domain for the interpolation mesh.

Definition at line 254 of file spline_builder.hpp.

◆ batched_interpolation_domain()

Get the whole domain representing interpolation points.

Values of the function must be provided on this domain in order to build a spline representation of the function (cartesian product of 1D interpolation_domain and batch_domain).

Returns
The domain for the interpolation mesh.

Definition at line 267 of file spline_builder.hpp.

◆ batch_domain()

Get the batch domain.

Obtained by removing the dimension of interest from the whole interpolation domain.

Returns
The batch domain.

Definition at line 279 of file spline_builder.hpp.

◆ spline_domain()

Get the 1D domain on which spline coefficients are defined.

The 1D spline domain corresponding to the dimension of interest.

Returns
The 1D domain for the spline coefficients.

Definition at line 291 of file spline_builder.hpp.

◆ batched_spline_domain()

Get the whole domain on which spline coefficients are defined.

Spline approximations (spline-transformed functions) are computed on this domain.

Returns
The domain for the spline coefficients.

Definition at line 303 of file spline_builder.hpp.

◆ batched_derivs_xmin_domain()

Get the whole domain on which derivatives on lower boundary are defined.

This is only used with BoundCond::HERMITE boundary conditions.

Returns
The domain for the Derivs values.

Definition at line 336 of file spline_builder.hpp.

◆ batched_derivs_xmax_domain()

Get the whole domain on which derivatives on upper boundary are defined.

This is only used with BoundCond::HERMITE boundary conditions.

Returns
The domain for the Derivs values.

Definition at line 352 of file spline_builder.hpp.

◆ get_interpolation_matrix()

template<class ExecSpace , class MemorySpace , class BSplines , class InterpolationDDim , ddc::BoundCond BcLower, ddc::BoundCond BcUpper, SplineSolver Solver, class... IDimX>
const ddc::detail:: SplinesLinearProblem< exec_space > & ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::get_interpolation_matrix ( ) const
inlinenoexcept

Get the interpolation matrix.

This can be useful for debugging (as it allows one to print the matrix) or for more complex quadrature schemes.

Deprecated:
Use quadrature_coefficients instead.
Warning
the returned detail::Matrix class is not supposed to be exposed to user, which means its usage is not supported out of the scope of current class. Use at your own risk.
Returns
A reference to the interpolation matrix.

Definition at line 377 of file spline_builder.hpp.

◆ operator()()

Compute a spline approximation of a function.

Use the values of a function (defined on SplineBuilder::batched_interpolation_domain) and the derivatives of the function at the boundaries (in the case of BoundCond::HERMITE only, defined on SplineBuilder::batched_derivs_xmin_domain and SplineBuilder::batched_derivs_xmax_domain) to calculate a spline approximation of this function.

The spline approximation is stored as a ChunkSpan of coefficients associated with B-splines.

Parameters
[out]splineThe coefficients of the spline computed by this SplineBuilder.
[in]valsThe values of the function on the interpolation mesh.
[in]derivs_xminThe values of the derivatives at the lower boundary (used only with BoundCond::HERMITE lower boundary condition).
[in]derivs_xmaxThe values of the derivatives at the upper boundary (used only with BoundCond::HERMITE upper boundary condition).

Definition at line 768 of file spline_builder.hpp.

◆ quadrature_coefficients()

Compute the quadrature coefficients associated to the b-splines used by this SplineBuilder.

Those coefficients can be used to perform integration way faster than SplineEvaluator::integrate().

This function solves matrix equation A^t*Q=integral_bsplines. In case of HERMITE boundary conditions, integral_bsplines contains the integral coefficients at the boundaries, and Q thus has to be splitted in three parts (quadrature coefficients for the derivatives at lower boundary, for the values inside the domain and for the derivatives at upper boundary).

A discrete function f can then be integrated using sum_j Q_j*f_j for j in interpolation_domain. If boundary condition is HERMITE, sum_j Qderiv_j*(d^j f/dx^j) for j in derivs_domain must be added at the boundary.

Please refer to section 2.8.1 of Emily's Bourne phd (https://theses.fr/2022AIXM0412) for more information and to the (Non)PeriodicSplineBuilderTest for example usage to compute integrals.

Template Parameters
OutMemorySpaceThe Kokkos::MemorySpace on which the quadrature coefficients are be returned (but they are computed on ExecSpace then copied).
Returns
A tuple containing the three Chunks containing the quadrature coefficients (if HERMITE is not used, first and third are empty).

Definition at line 957 of file spline_builder.hpp.

Member Data Documentation

◆ s_odd

template<class ExecSpace , class MemorySpace , class BSplines , class InterpolationDDim , ddc::BoundCond BcLower, ddc::BoundCond BcUpper, SplineSolver Solver, class... IDimX>
constexpr bool ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::s_odd = BSplines::degree() % 2
staticconstexpr

Indicates if the degree of the splines is odd or even.

Definition at line 148 of file spline_builder.hpp.

◆ s_nbc_xmin

The number of equations defining the boundary condition at the lower bound.

Definition at line 151 of file spline_builder.hpp.

◆ s_nbc_xmax

The number of equations defining the boundary condition at the upper bound.

Definition at line 154 of file spline_builder.hpp.

◆ s_bc_xmin

The boundary condition implemented at the lower bound.

Definition at line 157 of file spline_builder.hpp.

◆ s_bc_xmax

The boundary condition implemented at the upper bound.

Definition at line 160 of file spline_builder.hpp.

◆ s_spline_solver

The SplineSolver giving the backend used to perform the spline approximation.

Definition at line 163 of file spline_builder.hpp.


The documentation for this class was generated from the following file: