A class for creating a spline approximation of a function. More...
#include <spline_builder.hpp>
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). | |
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. | |
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.
ExecSpace | The Kokkos execution space on which the spline approximation is performed. |
MemorySpace | The Kokkos memory space on which the data (interpolation function and splines coefficients) is stored. |
BSplines | The discrete dimension representing the B-splines. |
InterpolationDDim | The discrete dimension on which interpolation points are defined. |
BcLower | The lower boundary condition. |
BcUpper | The upper boundary condition. |
Solver | The SplineSolver giving the backend used to perform the spline approximation. |
IDimX | A variadic template of all the discrete dimensions forming the full space (InterpolationDDim + batched dimensions). |
Definition at line 63 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::exec_space = ExecSpace |
The type of the Kokkos execution space used by this class.
Definition at line 73 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::memory_space = MemorySpace |
The type of the Kokkos memory space used by this class.
Definition at line 76 of file spline_builder.hpp.
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.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::interpolation_discrete_dimension_type = InterpolationDDim |
The type of the interpolation discrete dimension (discrete dimension of interest) used by this class.
Definition at line 82 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::bsplines_type = BSplines |
The discrete dimension representing the B-splines.
Definition at line 85 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::deriv_type = ddc::Deriv<continuous_dimension_type> |
The type of the Deriv dimension at the boundaries.
Definition at line 88 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::interpolation_domain_type = ddc::DiscreteDomain<interpolation_discrete_dimension_type> |
The type of the domain for the 1D interpolation mesh used by this class.
Definition at line 91 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::batched_interpolation_domain_type = ddc::DiscreteDomain<IDimX...> |
The type of the whole domain representing interpolation points.
Definition at line 94 of file spline_builder.hpp.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::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).
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.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::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).
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.
using ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::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).
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.
|
inlineexplicit |
Build a SplineBuilder acting on batched_interpolation_domain.
batched_interpolation_domain | The domain on which the interpolation points are defined. |
cols_per_chunk | A 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_size | A 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. |
Definition at line 195 of file spline_builder.hpp.
|
delete |
Copy-constructor is deleted.
|
default |
Move-constructs.
x | An rvalue to another SplineBuilder. |
|
default |
Destructs.
delete |
Copy-assignment is deleted.
|
default |
|
inlinenoexcept |
Get the domain for the 1D interpolation mesh used by this class.
This is 1D because it is defined along the dimension of interest.
Definition at line 254 of file spline_builder.hpp.
|
inlinenoexcept |
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).
Definition at line 267 of file spline_builder.hpp.
|
inlinenoexcept |
Get the batch domain.
Obtained by removing the dimension of interest from the whole interpolation domain.
Definition at line 279 of file spline_builder.hpp.
|
inlinenoexcept |
Get the 1D domain on which spline coefficients are defined.
The 1D spline domain corresponding to the dimension of interest.
Definition at line 291 of file spline_builder.hpp.
|
inlinenoexcept |
Get the whole domain on which spline coefficients are defined.
Spline approximations (spline-transformed functions) are computed on this domain.
Definition at line 303 of file spline_builder.hpp.
|
inlinenoexcept |
Get the whole domain on which derivatives on lower boundary are defined.
This is only used with BoundCond::HERMITE boundary conditions.
Definition at line 336 of file spline_builder.hpp.
|
inlinenoexcept |
Get the whole domain on which derivatives on upper boundary are defined.
This is only used with BoundCond::HERMITE boundary conditions.
Definition at line 352 of file spline_builder.hpp.
|
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.
Definition at line 377 of file spline_builder.hpp.
void ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::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.
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.
[out] | spline | The coefficients of the spline computed by this SplineBuilder. |
[in] | vals | The values of the function on the interpolation mesh. |
[in] | derivs_xmin | The values of the derivatives at the lower boundary (used only with BoundCond::HERMITE lower boundary condition). |
[in] | derivs_xmax | The 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.
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 > > > ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver, IDimX >::quadrature_coefficients | ( | ) | const |
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.
OutMemorySpace | The Kokkos::MemorySpace on which the quadrature coefficients are be returned (but they are computed on ExecSpace then copied). |
Definition at line 957 of file spline_builder.hpp.
|
staticconstexpr |
Indicates if the degree of the splines is odd or even.
Definition at line 148 of file spline_builder.hpp.
|
staticconstexpr |
The number of equations defining the boundary condition at the lower bound.
Definition at line 151 of file spline_builder.hpp.
|
staticconstexpr |
The number of equations defining the boundary condition at the upper bound.
Definition at line 154 of file spline_builder.hpp.
|
staticconstexpr |
The boundary condition implemented at the lower bound.
Definition at line 157 of file spline_builder.hpp.
|
staticconstexpr |
The boundary condition implemented at the upper bound.
Definition at line 160 of file spline_builder.hpp.
|
staticconstexpr |
The SplineSolver giving the backend used to perform the spline approximation.
Definition at line 163 of file spline_builder.hpp.