DDC 0.0.0

a discrete domain computation library

ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationMesh, BcXmin, BcXmax, Solver, IDimX > Class Template Reference

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

Public Types

using exec_space = ExecSpace
 
using memory_space = MemorySpace
 
using interpolation_mesh_type = InterpolationMesh
 The type of the interpolation mesh used by this class.
 
using bsplines_type = BSplines
 The type of the BSplines which are compatible with this class.
 
using deriv_type = ddc::Deriv<tag_type>
 
using interpolation_domain_type = ddc::DiscreteDomain<interpolation_mesh_type>
 The type of the domain for the interpolation mesh used by this class.
 
using batched_interpolation_domain_type = ddc::DiscreteDomain<IDimX...>
 
using batch_domain_type
 
using batched_spline_domain_type
 
using batched_spline_tr_domain_type
 
using batched_derivs_domain_type
 

Public Member Functions

 SplineBuilder (batched_interpolation_domain_type const &batched_interpolation_domain, std::optional< int > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditionner_max_block_size=std::nullopt)
 
 SplineBuilder (SplineBuilder const &x)=delete
 
 SplineBuilder (SplineBuilder &&x)=default
 Create a new SplineBuilder by copy.
 
 ~SplineBuilder ()=default
 
SplineBuilderoperator= (SplineBuilder const &x)=delete
 
SplineBuilderoperator= (SplineBuilder &&x)=default
 Copy a SplineBuilder.
 
batched_interpolation_domain_type batched_interpolation_domain () const noexcept
 
interpolation_domain_type interpolation_domain () const noexcept
 Get the domain from which the approximation is defined.
 
batch_domain_type batch_domain () const noexcept
 
ddc::DiscreteDomain< bsplines_typespline_domain () const noexcept
 
batched_spline_domain_type batched_spline_domain () const noexcept
 Get the domain on which the approximation is defined.
 
batched_spline_tr_domain_type batched_spline_tr_domain () const noexcept
 
batched_derivs_domain_type batched_derivs_xmin_domain () const noexcept
 
batched_derivs_domain_type batched_derivs_xmax_domain () const noexcept
 
const ddc::detail::Matrix & 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 > > const derivs_xmin=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > const derivs_xmax=std::nullopt) const
 Build a spline approximation of a function.
 

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(BcXmin, BSplines::degree())
 The number of equations which define the boundary conditions at the lower bound.
 
static constexpr int s_nbc_xmax = n_boundary_equations(BcXmax, BSplines::degree())
 The number of equations which define the boundary conditions at the upper bound.
 
static constexpr ddc::BoundCond s_bc_xmin = BcXmin
 The boundary condition implemented at the lower bound.
 
static constexpr ddc::BoundCond s_bc_xmax = BcXmax
 The boundary condition implemented at the upper bound.
 

Detailed Description

template<class ExecSpace, class MemorySpace, class BSplines, class InterpolationMesh, ddc::BoundCond BcXmin, ddc::BoundCond BcXmax, SplineSolver Solver, class... IDimX>
class ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationMesh, BcXmin, BcXmax, 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 BSplines. The spline is constructed such that it respects the boundary conditions BcXmin and BcXmax, and it interpolates the function at the points on the interpolation_mesh associated with interpolation_mesh_type.

Member Typedef Documentation

◆ exec_space

◆ memory_space

◆ interpolation_mesh_type

The type of the interpolation mesh used by this class.

◆ bsplines_type

The type of the BSplines which are compatible with this class.

◆ deriv_type

◆ interpolation_domain_type

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

◆ batched_interpolation_domain_type

◆ batch_domain_type

Initial value:
typename ddc::detail::convert_type_seq_to_discrete_domain<ddc::type_seq_remove_t<
ddc::detail::TypeSeq<IDimX...>,
ddc::detail::TypeSeq<interpolation_mesh_type>>>
constexpr bool enable_chunk
Definition chunk_traits.hpp:16

◆ batched_spline_domain_type

Initial value:
typename ddc::detail::convert_type_seq_to_discrete_domain<ddc::type_seq_replace_t<
ddc::detail::TypeSeq<IDimX...>,
ddc::detail::TypeSeq<interpolation_mesh_type>,
ddc::detail::TypeSeq<bsplines_type>>>

◆ batched_spline_tr_domain_type

Initial value:
typename ddc::detail::convert_type_seq_to_discrete_domain<ddc::type_seq_merge_t<
ddc::detail::TypeSeq<bsplines_type>,
ddc::detail::TypeSeq<IDimX...>,
ddc::detail::TypeSeq<interpolation_mesh_type>>>>

◆ batched_derivs_domain_type

Initial value:
typename ddc::detail::convert_type_seq_to_discrete_domain<ddc::type_seq_replace_t<
ddc::detail::TypeSeq<IDimX...>,
ddc::detail::TypeSeq<interpolation_mesh_type>,
ddc::detail::TypeSeq<deriv_type>>>

Constructor & Destructor Documentation

◆ SplineBuilder() [1/3]

ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationMesh, BcXmin, BcXmax, Solver, IDimX >::SplineBuilder ( batched_interpolation_domain_type const & batched_interpolation_domain,
std::optional< int > cols_per_chunk = std::nullopt,
std::optional< unsigned int > preconditionner_max_block_size = std::nullopt )
inlineexplicit

◆ SplineBuilder() [2/3]

◆ SplineBuilder() [3/3]

◆ ~SplineBuilder()

Member Function Documentation

◆ operator=() [1/2]

◆ operator=() [2/2]

◆ batched_interpolation_domain()

◆ interpolation_domain()

Get the domain from which the approximation is defined.

Get the domain on which values of the function must be provided in order to build a spline approximation of the function.

Returns
The domain for the grid points.

◆ batch_domain()

◆ spline_domain()

◆ batched_spline_domain()

Get the domain on which the approximation is defined.

Get the domain of the basis-splines for which the coefficients of the spline approximation must be calculated.

Returns
The domain for the splines.

◆ batched_spline_tr_domain()

◆ batched_derivs_xmin_domain()

◆ batched_derivs_xmax_domain()

◆ get_interpolation_matrix()

const ddc::detail::Matrix & ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationMesh, BcXmin, BcXmax, Solver, IDimX >::get_interpolation_matrix ( ) const
inlinenoexcept

Get the interpolation matrix.

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

Returns
A reference to the interpolation matrix.

◆ operator()()

Build a spline approximation of a function.

Use the values of a function at known grid points (as specified by SplineBuilder::interpolation_domain) and the derivatives of the function at the boundaries (if necessary for the chosen boundary conditions) to calculate a spline approximation of a function.

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

Parameters
[out]splineThe coefficients of the spline calculated by the function.
[in]valsThe values of the function at the grid points.
[in]derivs_xminThe values of the derivatives at the lower boundary.
[in]derivs_xmaxThe values of the derivatives at the upper boundary.

Member Data Documentation

◆ s_odd

constexpr bool ddc::SplineBuilder< ExecSpace, MemorySpace, BSplines, InterpolationMesh, BcXmin, BcXmax, Solver, IDimX >::s_odd = BSplines::degree() % 2
staticconstexpr

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

◆ s_nbc_xmin

The number of equations which define the boundary conditions at the lower bound.

◆ s_nbc_xmax

The number of equations which define the boundary conditions at the upper bound.

◆ s_bc_xmin

The boundary condition implemented at the lower bound.

◆ s_bc_xmax

The boundary condition implemented at the upper bound.


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