DDC 0.1.0
Loading...
Searching...
No Matches
ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, Solver, IDimX > Class Template Reference

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

#include <spline_builder_2d.hpp>

Collaboration diagram for ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, 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 builder_type1 = ddc::SplineBuilder< ExecSpace, MemorySpace, BSpline1, IDimI1, BcLower1, BcUpper1, Solver, IDimX... >
 The type of the SplineBuilder used by this class to spline-approximate along first dimension.
 
using builder_type2 = ddc::SplineBuilder< ExecSpace, MemorySpace, BSpline2, IDimI2, BcLower2, BcUpper2, Solver, std::conditional_t< std::is_same_v< IDimX, IDimI1 >, BSpline1, IDimX >... >
 The type of the SplineBuilder used by this class to spline-approximate along second dimension.
 
using builder_deriv_type1 = ddc::SplineBuilder< ExecSpace, MemorySpace, BSpline1, IDimI1, BcLower1, BcUpper1, Solver, std::conditional_t< std::is_same_v< IDimX, IDimI2 >, typename builder_type2::deriv_type, IDimX >... >
 The type of the SplineBuilder used by this class to spline-approximate the second-dimension-derivatives along first dimension.
 
using continuous_dimension_type1 = typename builder_type1::continuous_dimension_type
 The type of the first interpolation continuous dimension.
 
using continuous_dimension_type2 = typename builder_type2::continuous_dimension_type
 The type of the second interpolation continuous dimension.
 
using interpolation_discrete_dimension_type1 = typename builder_type1::interpolation_discrete_dimension_type
 The type of the first interpolation discrete dimension.
 
using interpolation_discrete_dimension_type2 = typename builder_type2::interpolation_discrete_dimension_type
 The type of the second interpolation discrete dimension.
 
using bsplines_type1 = typename builder_type1::bsplines_type
 The type of the B-splines in the first dimension.
 
using bsplines_type2 = typename builder_type2::bsplines_type
 The type of the B-splines in the second dimension.
 
using deriv_type1 = typename builder_type1::deriv_type
 The type of the Deriv domain on boundaries in the first dimension.
 
using deriv_type2 = typename builder_type2::deriv_type
 The type of the Deriv domain on boundaries in the second dimension.
 
using interpolation_domain_type1 = typename builder_type1::interpolation_discrete_dimension_type
 The type of the domain for the interpolation mesh in the first dimension.
 
using interpolation_domain_type2 = typename builder_type2::interpolation_discrete_dimension_type
 The type of the domain for the interpolation mesh in the second dimension.
 
using interpolation_domain_type = ddc::DiscreteDomain< interpolation_discrete_dimension_type1, interpolation_discrete_dimension_type2 >
 The type of the domain for the interpolation mesh in the 2D dimension.
 
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_type1, interpolation_discrete_dimension_type2 >
 The type of the batch domain (obtained by removing the dimensions of interest from the whole domain).
 
using batched_spline_domain_type = ddc::detail::convert_type_seq_to_discrete_domain_t< ddc::type_seq_replace_t< ddc::detail::TypeSeq< IDimX... >, ddc::detail::TypeSeq< interpolation_discrete_dimension_type1, interpolation_discrete_dimension_type2 >, ddc::detail::TypeSeq< bsplines_type1, bsplines_type2 > > >
 The type of the whole spline domain (cartesian product of 2D spline domain and batch domain) preserving the order of dimensions.
 
using batched_derivs_domain_type1 = typename builder_type1::batched_derivs_domain_type
 The type of the whole Derivs domain (cartesian product of the 1D Deriv domain and the associated batch domain) in the first dimension, preserving the order of dimensions.
 
using batched_derivs_domain_type2 = ddc::replace_dim_of_t< batched_interpolation_domain_type, interpolation_discrete_dimension_type2, deriv_type2 >
 The type of the whole Derivs domain (cartesian product of the 1D Deriv domain and the associated batch domain) in the second dimension, preserving the order of dimensions.
 
using batched_derivs_domain_type = ddc::detail::convert_type_seq_to_discrete_domain_t< ddc::type_seq_replace_t< ddc::detail::TypeSeq< IDimX... >, ddc::detail::TypeSeq< interpolation_discrete_dimension_type1, interpolation_discrete_dimension_type2 >, ddc::detail::TypeSeq< deriv_type1, deriv_type2 > > >
 The type of the whole Derivs domain (cartesian product of the 2D Deriv domain and the batch domain) in the second dimension, preserving the order of dimensions.
 

Public Member Functions

 SplineBuilder2D (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 SplineBuilder2D acting on batched_interpolation_domain.
 
 SplineBuilder2D (SplineBuilder2D const &x)=delete
 Copy-constructor is deleted.
 
 SplineBuilder2D (SplineBuilder2D &&x)=default
 Move-constructs.
 
 ~SplineBuilder2D ()=default
 Destructs.
 
SplineBuilder2Doperator= (SplineBuilder2D const &x)=delete
 Copy-assignment is deleted.
 
SplineBuilder2Doperator= (SplineBuilder2D &&x)=default
 Move-assigns.
 
interpolation_domain_type interpolation_domain () const noexcept
 Get the domain for the 2D 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_type1, bsplines_type2spline_domain () const noexcept
 Get the 2D 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.
 
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_type1, Layout, memory_space > > derivs_min1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1, Layout, memory_space > > derivs_max1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2, Layout, memory_space > > derivs_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2, Layout, memory_space > > derivs_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > mixed_derivs_min1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > mixed_derivs_max1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > mixed_derivs_min1_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > > mixed_derivs_max1_max2=std::nullopt) const
 Compute a 2D spline approximation of a function.
 

Detailed Description

template<class ExecSpace, class MemorySpace, class BSpline1, class BSpline2, class IDimI1, class IDimI2, ddc::BoundCond BcLower1, ddc::BoundCond BcUpper1, ddc::BoundCond BcLower2, ddc::BoundCond BcUpper2, ddc::SplineSolver Solver, class... IDimX>
class ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, Solver, IDimX >

A class for creating a 2D spline approximation of a function.

A class which contains an operator () which can be used to build a 2D spline approximation of a function. A 2D spline approximation uses a cross-product between two 1D SplineBuilder.

See also
SplineBuilder

Definition at line 38 of file spline_builder_2d.hpp.

Member Typedef Documentation

◆ exec_space

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

Definition at line 42 of file spline_builder_2d.hpp.

◆ memory_space

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

Definition at line 45 of file spline_builder_2d.hpp.

◆ builder_type1

The type of the SplineBuilder used by this class to spline-approximate along first dimension.

Definition at line 48 of file spline_builder_2d.hpp.

◆ builder_type2

The type of the SplineBuilder used by this class to spline-approximate along second dimension.

Definition at line 59 of file spline_builder_2d.hpp.

◆ builder_deriv_type1

The type of the SplineBuilder used by this class to spline-approximate the second-dimension-derivatives along first dimension.

Definition at line 70 of file spline_builder_2d.hpp.

◆ continuous_dimension_type1

The type of the first interpolation continuous dimension.

Definition at line 84 of file spline_builder_2d.hpp.

◆ continuous_dimension_type2

The type of the second interpolation continuous dimension.

Definition at line 87 of file spline_builder_2d.hpp.

◆ interpolation_discrete_dimension_type1

The type of the first interpolation discrete dimension.

Definition at line 90 of file spline_builder_2d.hpp.

◆ interpolation_discrete_dimension_type2

The type of the second interpolation discrete dimension.

Definition at line 94 of file spline_builder_2d.hpp.

◆ bsplines_type1

The type of the B-splines in the first dimension.

Definition at line 98 of file spline_builder_2d.hpp.

◆ bsplines_type2

The type of the B-splines in the second dimension.

Definition at line 101 of file spline_builder_2d.hpp.

◆ deriv_type1

The type of the Deriv domain on boundaries in the first dimension.

Definition at line 104 of file spline_builder_2d.hpp.

◆ deriv_type2

The type of the Deriv domain on boundaries in the second dimension.

Definition at line 107 of file spline_builder_2d.hpp.

◆ interpolation_domain_type1

The type of the domain for the interpolation mesh in the first dimension.

Definition at line 110 of file spline_builder_2d.hpp.

◆ interpolation_domain_type2

The type of the domain for the interpolation mesh in the second dimension.

Definition at line 114 of file spline_builder_2d.hpp.

◆ interpolation_domain_type

The type of the domain for the interpolation mesh in the 2D dimension.

Definition at line 118 of file spline_builder_2d.hpp.

◆ batched_interpolation_domain_type

The type of the whole domain representing interpolation points.

Definition at line 123 of file spline_builder_2d.hpp.

◆ batch_domain_type

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

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

Definition at line 132 of file spline_builder_2d.hpp.

◆ batched_spline_domain_type

template<class ExecSpace , class MemorySpace , class BSpline1 , class BSpline2 , class IDimI1 , class IDimI2 , ddc::BoundCond BcLower1, ddc::BoundCond BcUpper1, ddc::BoundCond BcLower2, ddc::BoundCond BcUpper2, ddc::SplineSolver Solver, class... IDimX>
using ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, Solver, IDimX >::batched_spline_domain_type = ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t< ddc::detail::TypeSeq<IDimX...>, ddc::detail::TypeSeq< interpolation_discrete_dimension_type1, interpolation_discrete_dimension_type2>, ddc::detail::TypeSeq<bsplines_type1, bsplines_type2> >>

The type of the whole spline domain (cartesian product of 2D spline domain and batch domain) preserving the order of dimensions.

Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y (associated to B-splines tags BSplinesX and BSplinesY), this is DiscreteDomain<BSplinesX, BSplinesY, Z>

Definition at line 144 of file spline_builder_2d.hpp.

◆ batched_derivs_domain_type1

The type of the whole Derivs domain (cartesian product of the 1D Deriv domain and the associated batch domain) in the first dimension, preserving the order of dimensions.

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

Definition at line 159 of file spline_builder_2d.hpp.

◆ batched_derivs_domain_type2

The type of the whole Derivs domain (cartesian product of the 1D Deriv domain and the associated batch domain) in the second dimension, preserving the order of dimensions.

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

Definition at line 168 of file spline_builder_2d.hpp.

◆ batched_derivs_domain_type

template<class ExecSpace , class MemorySpace , class BSpline1 , class BSpline2 , class IDimI1 , class IDimI2 , ddc::BoundCond BcLower1, ddc::BoundCond BcUpper1, ddc::BoundCond BcLower2, ddc::BoundCond BcUpper2, ddc::SplineSolver Solver, class... IDimX>
using ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, Solver, IDimX >::batched_derivs_domain_type = ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t< ddc::detail::TypeSeq<IDimX...>, ddc::detail::TypeSeq< interpolation_discrete_dimension_type1, interpolation_discrete_dimension_type2>, ddc::detail::TypeSeq<deriv_type1, deriv_type2> >>

The type of the whole Derivs domain (cartesian product of the 2D Deriv domain and the batch domain) in the second dimension, preserving the order of dimensions.

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

Definition at line 180 of file spline_builder_2d.hpp.

Constructor & Destructor Documentation

◆ SplineBuilder2D() [1/3]

template<class ExecSpace , class MemorySpace , class BSpline1 , class BSpline2 , class IDimI1 , class IDimI2 , ddc::BoundCond BcLower1, ddc::BoundCond BcUpper1, ddc::BoundCond BcLower2, ddc::BoundCond BcUpper2, ddc::SplineSolver Solver, class... IDimX>
ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, Solver, IDimX >::SplineBuilder2D ( 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 SplineBuilder2D 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
SplinesLinearProblemSparse

Definition at line 210 of file spline_builder_2d.hpp.

◆ SplineBuilder2D() [2/3]

◆ SplineBuilder2D() [3/3]

◆ ~SplineBuilder2D()

Member Function Documentation

◆ operator=() [1/2]

◆ operator=() [2/2]

◆ interpolation_domain()

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

This is 2D because it is defined along the dimensions of interest.

Returns
The 2D domain for the interpolation mesh.

Definition at line 261 of file spline_builder_2d.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 2D interpolation_domain and batch_domain).

Returns
The domain for the interpolation mesh.

Definition at line 276 of file spline_builder_2d.hpp.

◆ batch_domain()

Get the batch domain.

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

Returns
The batch domain.

Definition at line 288 of file spline_builder_2d.hpp.

◆ spline_domain()

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

The 2D spline domain corresponding to the dimensions of interest.

Returns
The 2D domain for the spline coefficients.

Definition at line 300 of file spline_builder_2d.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 314 of file spline_builder_2d.hpp.

◆ operator()()

template<class ExecSpace , class MemorySpace , class BSpline1 , class BSpline2 , class IDimI1 , class IDimI2 , ddc::BoundCond BcLower1, ddc::BoundCond BcUpper1, ddc::BoundCond BcLower2, ddc::BoundCond BcUpper2, ddc::SplineSolver Solver, class... IDimX>
template<class Layout >
void ddc::SplineBuilder2D< ExecSpace, MemorySpace, BSpline1, BSpline2, IDimI1, IDimI2, BcLower1, BcUpper1, BcLower2, BcUpper2, 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_type1, Layout, memory_space > >  derivs_min1 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1, Layout, memory_space > >  derivs_max1 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2, Layout, memory_space > >  derivs_min2 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2, Layout, memory_space > >  derivs_max2 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > >  mixed_derivs_min1_min2 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > >  mixed_derivs_max1_min2 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > >  mixed_derivs_min1_max2 = std::nullopt,
std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type, Layout, memory_space > >  mixed_derivs_max1_max2 = std::nullopt 
) const

Compute a 2D spline approximation of a function.

Use the values of a function (defined on SplineBuilder2D::batched_interpolation_domain) and the derivatives of the function at the boundaries (in the case of BoundCond::HERMITE only) to calculate a 2D 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 at the interpolation mesh.
[in]derivs_min1The values of the derivatives at the lower boundary in the first dimension.
[in]derivs_max1The values of the derivatives at the upper boundary in the first dimension.
[in]derivs_min2The values of the derivatives at the lower boundary in the second dimension.
[in]derivs_max2The values of the derivatives at the upper boundary in the second dimension.
[in]mixed_derivs_min1_min2The values of the the cross-derivatives at the lower boundary in the first dimension and the lower boundary in the second dimension.
[in]mixed_derivs_max1_min2The values of the the cross-derivatives at the upper boundary in the first dimension and the lower boundary in the second dimension.
[in]mixed_derivs_min1_max2The values of the the cross-derivatives at the lower boundary in the first dimension and the upper boundary in the second dimension.
[in]mixed_derivs_max1_max2The values of the the cross-derivatives at the upper boundary in the first dimension and the upper boundary in the second dimension.

Definition at line 425 of file spline_builder_2d.hpp.


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