DDC 0.10.0
Loading...
Searching...
No Matches
knots_as_interpolation_points.hpp
1// Copyright (C) The DDC development team, see COPYRIGHT.md file
2//
3// SPDX-License-Identifier: MIT
4
5#pragma once
6
7#include <tuple>
8#include <type_traits>
9#include <vector>
10
11#include <ddc/ddc.hpp>
12
13#include <Kokkos_Core.hpp>
14
18
19namespace ddc {
20
21/**
22 * @brief Helper class for the initialisation of the mesh of interpolation points.
23 *
24 * A helper class for the initialisation of the mesh of interpolation points. This
25 * class should be used when the interpolation points should be located at the
26 * knots of the spline. This is possible with any kind of boundary condition except
27 * Greville boundary conditions (as there will not be enough interpolation points).
28 * In the case of strongly non-uniform splines this choice may result in a less
29 * well conditioned problem, however most mathematical stability results are proven
30 * with this choice of interpolation points.
31 *
32 * @tparam BSplines The type of the uniform or non-uniform spline basis whose knots are used as interpolation points.
33 * @tparam BcLower The lower boundary condition.
34 * @tparam BcLower The upper boundary condition.
35 */
36template <class BSplines, ddc::BoundCond BcLower, ddc::BoundCond BcUpper>
38{
39 static_assert(BcLower != ddc::BoundCond::GREVILLE);
40 static_assert(BcUpper != ddc::BoundCond::GREVILLE);
41
42 using continuous_dimension_type = typename BSplines::continuous_dimension_type;
43
44public:
45 /**
46 * Get the sampling of interpolation points.
47 *
48 * @return sampling The DDC point sampling of the interpolation points.
49 */
50 template <typename Sampling, typename U = BSplines>
51 static auto get_sampling()
52 {
53 if constexpr (U::is_uniform()) {
54 return std::get<0>(Sampling::
55 init(ddc::discrete_space<BSplines>().rmin(),
56 ddc::discrete_space<BSplines>().rmax(),
57 ddc::DiscreteVector<Sampling>(
58 ddc::discrete_space<BSplines>().ncells() + 1)));
59 } else {
60 using SamplingImpl = typename Sampling::template Impl<Sampling, Kokkos::HostSpace>;
61 std::vector<double> knots(ddc::discrete_space<BSplines>().npoints());
62 ddc::DiscreteDomain<knot_discrete_dimension_t<BSplines>> break_point_domain(
63 ddc::discrete_space<BSplines>().break_point_domain());
64 ddc::host_for_each(
65 break_point_domain,
66 [&](ddc::DiscreteElement<knot_discrete_dimension_t<BSplines>> ik) {
67 knots[ik - break_point_domain.front()] = ddc::coordinate(ik);
68 });
69 return SamplingImpl(knots);
70 }
71 }
72
73 /// The DDC type of the sampling for the interpolation points.
74 using interpolation_discrete_dimension_type = std::conditional_t<
75 is_uniform_bsplines_v<BSplines>,
76 ddc::UniformPointSampling<continuous_dimension_type>,
77 ddc::NonUniformPointSampling<continuous_dimension_type>>;
78
79 /**
80 * Get the domain which can be used to access the interpolation points in the sampling.
81 *
82 * @return domain The discrete domain which maps to the sampling of interpolation points.
83 */
84 template <typename Sampling>
85 static ddc::DiscreteDomain<Sampling> get_domain()
86 {
87 int const npoints = ddc::discrete_space<BSplines>().ncells() + !BSplines::is_periodic();
88 return ddc::DiscreteDomain<Sampling>(
89 ddc::DiscreteElement<Sampling>(0),
90 ddc::DiscreteVector<Sampling>(npoints));
91 }
92};
93
94} // namespace ddc
friend class ChunkSpan
friend class DiscreteDomain
KOKKOS_FUNCTION constexpr bool operator!=(DiscreteVector< OTags... > const &rhs) const noexcept
A class which provides helper functions to initialise the Greville points from a B-Spline definition.
static ddc::DiscreteDomain< Sampling > get_domain()
Get the domain which gives us access to all of the Greville points.
static auto get_sampling()
Get the UniformPointSampling defining the Greville points.
Helper class for the initialisation of the mesh of interpolation points.
static auto get_sampling()
Get the sampling of interpolation points.
static ddc::DiscreteDomain< Sampling > get_domain()
Get the domain which can be used to access the interpolation points in the sampling.
Storage class of the static attributes of the discrete dimension.
Impl & operator=(Impl &&x)=default
Move-assigns.
Impl(RandomIt breaks_begin, RandomIt breaks_end)
Constructs an Impl by iterating over a range of break points from begin to end.
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmin() const noexcept
Returns the coordinate of the first break point of the domain on which the B-splines are defined.
Impl(std::vector< ddc::Coordinate< CDim > > const &breaks)
Constructs an Impl using a std::vector.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis(DSpan1D values, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-splines at a given coordinate.
KOKKOS_INLINE_FUNCTION std::size_t size() const noexcept
Returns the number of elements necessary to construct a spline representation of a function.
Impl(Impl< DDim, OriginMemorySpace > const &impl)
Copy-constructs from another Impl with a different Kokkos memory space.
~Impl()=default
Destructs.
KOKKOS_INLINE_FUNCTION ddc::DiscreteDomain< knot_discrete_dimension_type > break_point_domain() const
Returns the discrete domain which describes the break points.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_last_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the last support knot associated to a DiscreteElement identifying a B-splin...
Impl(Impl &&x)=default
Move-constructs.
Impl(std::initializer_list< ddc::Coordinate< CDim > > breaks)
Constructs an Impl using a brace-list, i.e.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis_and_n_derivs(ddc::DSpan2D derivs, ddc::Coordinate< CDim > const &x, std::size_t n) const
Evaluates non-zero B-spline values and derivatives at a given coordinate.
KOKKOS_INLINE_FUNCTION std::size_t ncells() const noexcept
Returns the number of cells over which the B-splines are defined.
KOKKOS_INLINE_FUNCTION discrete_domain_type full_domain() const
Returns the discrete domain including eventual additional B-splines in the periodic case.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_first_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the first support knot associated to a DiscreteElement identifying a B-spli...
KOKKOS_INLINE_FUNCTION std::size_t npoints() const noexcept
The number of break points.
KOKKOS_INLINE_FUNCTION std::size_t nbasis() const noexcept
Returns the number of basis functions.
Impl(Impl const &x)=default
Copy-constructs.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_deriv(DSpan1D derivs, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-spline derivatives at a given coordinate.
KOKKOS_INLINE_FUNCTION double length() const noexcept
Returns the length of the domain.
Impl & operator=(Impl const &x)=default
Copy-assigns.
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmax() const noexcept
Returns the coordinate of the last break point of the domain on which the B-splines are defined.
The type of a non-uniform 1D spline basis (B-spline).
static constexpr std::size_t degree() noexcept
The degree of B-splines.
static constexpr bool is_periodic() noexcept
Indicates if the B-splines are periodic or not.
static constexpr bool is_uniform() noexcept
Indicates if the B-splines are uniform or not (this is not the case here).
NonUniformPointSampling models a non-uniform discretization of the CDim segment .
Storage class of the static attributes of the discrete dimension.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_last_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the last support knot associated to a DiscreteElement identifying a B-splin...
Impl(ddc::Coordinate< CDim > rmin, ddc::Coordinate< CDim > rmax, std::size_t ncells)
Constructs a spline basis (B-splines) with n equidistant knots over .
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmax() const noexcept
Returns the coordinate of the upper bound of the domain on which the B-splines are defined.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis(DSpan1D values, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-splines at a given coordinate.
KOKKOS_INLINE_FUNCTION ddc::DiscreteDomain< knot_discrete_dimension_type > break_point_domain() const
Returns the discrete domain which describes the break points.
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmin() const noexcept
Returns the coordinate of the lower bound of the domain on which the B-splines are defined.
~Impl()=default
Destructs.
KOKKOS_INLINE_FUNCTION std::size_t nbasis() const noexcept
Returns the number of basis functions.
Impl(Impl const &x)=default
Copy-constructs.
KOKKOS_INLINE_FUNCTION std::size_t size() const noexcept
Returns the number of elements necessary to construct a spline representation of a function.
Impl(Impl &&x)=default
Move-constructs.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis_and_n_derivs(ddc::DSpan2D derivs, ddc::Coordinate< CDim > const &x, std::size_t n) const
Evaluates non-zero B-spline values and derivatives at a given coordinate.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_first_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the first support knot associated to a DiscreteElement identifying a B-spli...
KOKKOS_INLINE_FUNCTION double length() const noexcept
Returns the length of the domain.
Impl(Impl< DDim, OriginMemorySpace > const &impl)
Copy-constructs from another Impl with a different Kokkos memory space.
KOKKOS_INLINE_FUNCTION std::size_t ncells() const noexcept
Returns the number of cells over which the B-splines are defined.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_deriv(DSpan1D derivs, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-spline derivatives at a given coordinate.
Impl & operator=(Impl &&x)=default
Move-assigns.
KOKKOS_INLINE_FUNCTION discrete_domain_type full_domain() const
Returns the discrete domain including eventual additional B-splines in the periodic case.
Impl & operator=(Impl const &x)=default
Copy-assigns.
The type of a uniform 1D spline basis (B-spline).
static constexpr bool is_uniform() noexcept
Indicates if the B-splines are uniform or not (this is the case here).
static constexpr bool is_periodic() noexcept
Indicates if the B-splines are periodic or not.
static constexpr std::size_t degree() noexcept
The degree of B-splines.
UniformPointSampling models a uniform discretization of the provided continuous dimension.
The top-level namespace of DDC.
constexpr int n_boundary_equations(ddc::BoundCond const bc, std::size_t const degree)
Return the number of equations needed to describe a given boundary condition.
constexpr bool is_uniform_bsplines_v
Indicates if a tag corresponds to uniform B-splines or not.
BoundCond
An enum representing a spline boundary condition.
@ GREVILLE
Use Greville points instead of conditions on derivative for B-Spline interpolation.
@ HERMITE
Hermite boundary condition.
@ PERIODIC
Periodic boundary condition u(1)=u(n)
ddc::ChunkSpan< double, ddc::DiscreteDomain< DDim >, Layout, MemorySpace > integrals(ExecSpace const &execution_space, ddc::ChunkSpan< double, ddc::DiscreteDomain< DDim >, Layout, MemorySpace > int_vals)
Compute the integrals of the B-splines.
constexpr bool is_non_uniform_bsplines_v
Indicates if a tag corresponds to non-uniform B-splines or not.
A templated struct representing a discrete dimension storing the derivatives of a function along a co...
Definition deriv.hpp:15
If the type DDim is a B-spline, defines type to the discrete dimension of the associated knots.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos, ddc::Coordinate< DimNI > eval_pos_not_interest_min, ddc::Coordinate< DimNI > eval_pos_not_interest_max)
Instantiate a ConstantExtrapolationRule.
KOKKOS_FUNCTION double operator()(CoordType coord_extrap, ddc::ChunkSpan< double const, ddc::DiscreteDomain< BSplines1, BSplines2 >, Layout, MemorySpace > const spline_coef) const
Get the value of the function on B-splines at a coordinate outside the domain.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos)
Instantiate a ConstantExtrapolationRule.
KOKKOS_FUNCTION double operator()(CoordType pos, ddc::ChunkSpan< double const, ddc::DiscreteDomain< BSplines >, Layout, MemorySpace > const spline_coef) const
Get the value of the function on B-splines at a coordinate outside the domain.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos)
Instantiate a ConstantExtrapolationRule.