DDC 0.1.0
Loading...
Searching...
No Matches
spline_evaluator.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 <array>
8#include <cstddef>
9#include <type_traits>
10
11#include <ddc/ddc.hpp>
12
13#include <Kokkos_Macros.hpp>
14#include <Kokkos_MathematicalFunctions.hpp>
15
16#include "integrals.hpp"
17#include "periodic_extrapolation_rule.hpp"
18
19namespace ddc {
20
21/**
22 * @brief A class to evaluate, differentiate or integrate a spline function.
23 *
24 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a spline function.
25 *
26 * @tparam ExecSpace The Kokkos execution space on which the spline evaluation is performed.
27 * @tparam MemorySpace The Kokkos memory space on which the data (spline coefficients and evaluation) is stored.
28 * @tparam BSplines The discrete dimension representing the B-splines.
29 * @tparam EvaluationDDim The discrete dimension on which evaluation points are defined.
30 * @tparam LowerExtrapolationRule The lower extrapolation rule type.
31 * @tparam UpperExtrapolationRule The upper extrapolation rule type.
32 * @tparam IDimX A variadic template of all the discrete dimensions forming the full space (EvaluationDDim + batched dimensions).
33 */
34template <
35 class ExecSpace,
36 class MemorySpace,
37 class BSplines,
38 class EvaluationDDim,
39 class LowerExtrapolationRule,
40 class UpperExtrapolationRule,
41 class... IDimX>
43{
44private:
45 /**
46 * @brief Tag to indicate that the value of the spline should be evaluated.
47 */
48 struct eval_type
49 {
50 };
51
52 /**
53 * @brief Tag to indicate that derivative of the spline should be evaluated.
54 */
55 struct eval_deriv_type
56 {
57 };
58
59public:
60 /// @brief The type of the Kokkos execution space used by this class.
61 using exec_space = ExecSpace;
62
63 /// @brief The type of the Kokkos memory space used by this class.
64 using memory_space = MemorySpace;
65
66 /// @brief The type of the evaluation continuous dimension (continuous dimension of interest) used by this class.
67 using continuous_dimension_type = typename BSplines::continuous_dimension_type;
68
69 /// @brief The type of the evaluation discrete dimension (discrete dimension of interest) used by this class.
70 using evaluation_discrete_dimension_type = EvaluationDDim;
71
72 /// @brief The discrete dimension representing the B-splines.
73 using bsplines_type = BSplines;
74
75 /// @brief The type of the domain for the 1D evaluation mesh used by this class.
76 using evaluation_domain_type = ddc::DiscreteDomain<evaluation_discrete_dimension_type>;
77
78 /// @brief The type of the whole domain representing evaluation points.
79 using batched_evaluation_domain_type = ddc::DiscreteDomain<IDimX...>;
80
81 /// @brief The type of the 1D spline domain corresponding to the dimension of interest.
82 using spline_domain_type = ddc::DiscreteDomain<bsplines_type>;
83
84 /**
85 * @brief The type of the batch domain (obtained by removing the dimension of interest
86 * from the whole domain).
87 */
88 using batch_domain_type = ddc::
89 remove_dims_of_t<batched_evaluation_domain_type, evaluation_discrete_dimension_type>;
90
91 /**
92 * @brief The type of the whole spline domain (cartesian product of 1D spline domain
93 * and batch domain) preserving the order of dimensions.
94 */
95 using batched_spline_domain_type = ddc::replace_dim_of_t<
96 batched_evaluation_domain_type,
97 evaluation_discrete_dimension_type,
98 bsplines_type>;
99
100 /// @brief The type of the extrapolation rule at the lower boundary.
101 using lower_extrapolation_rule_type = LowerExtrapolationRule;
102
103 /// @brief The type of the extrapolation rule at the upper boundary.
104 using upper_extrapolation_rule_type = UpperExtrapolationRule;
105
106
107private:
108 LowerExtrapolationRule m_lower_extrap_rule;
109
110 UpperExtrapolationRule m_upper_extrap_rule;
111
112public:
113 static_assert(
114 std::is_same_v<LowerExtrapolationRule,
115 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
116 == bsplines_type::is_periodic()
117 && std::is_same_v<
118 UpperExtrapolationRule,
119 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
120 == bsplines_type::is_periodic(),
121 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
122 static_assert(
123 std::is_invocable_r_v<
124 double,
125 LowerExtrapolationRule,
126 ddc::Coordinate<continuous_dimension_type>,
127 ddc::ChunkSpan<
128 double const,
129 spline_domain_type,
130 Kokkos::layout_right,
131 memory_space>>,
132 "LowerExtrapolationRule::operator() has to be callable with usual arguments.");
133 static_assert(
134 std::is_invocable_r_v<
135 double,
136 UpperExtrapolationRule,
137 ddc::Coordinate<continuous_dimension_type>,
138 ddc::ChunkSpan<
139 double const,
140 spline_domain_type,
141 Kokkos::layout_right,
142 memory_space>>,
143 "UpperExtrapolationRule::operator() has to be callable with usual arguments.");
144
145 /**
146 * @brief Build a SplineEvaluator acting on batched_spline_domain.
147 *
148 * @param lower_extrap_rule The extrapolation rule at the lower boundary.
149 * @param upper_extrap_rule The extrapolation rule at the upper boundary.
150 *
151 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
152 */
153 explicit SplineEvaluator(
154 LowerExtrapolationRule const& lower_extrap_rule,
155 UpperExtrapolationRule const& upper_extrap_rule)
156 : m_lower_extrap_rule(lower_extrap_rule)
157 , m_upper_extrap_rule(upper_extrap_rule)
158 {
159 }
160
161 /**
162 * @brief Copy-constructs.
163 *
164 * @param x A reference to another SplineEvaluator.
165 */
166 SplineEvaluator(SplineEvaluator const& x) = default;
167
168 /**
169 * @brief Move-constructs.
170 *
171 * @param x An rvalue to another SplineEvaluator.
172 */
173 SplineEvaluator(SplineEvaluator&& x) = default;
174
175 /// @brief Destructs
176 ~SplineEvaluator() = default;
177
178 /**
179 * @brief Copy-assigns.
180 *
181 * @param x A reference to another SplineEvaluator.
182 * @return A reference to this object.
183 */
184 SplineEvaluator& operator=(SplineEvaluator const& x) = default;
185
186 /**
187 * @brief Move-assigns.
188 *
189 * @param x An rvalue to another SplineEvaluator.
190 * @return A reference to this object.
191 */
193
194 /**
195 * @brief Get the lower extrapolation rule.
196 *
197 * Extrapolation rules are functors used to define the behavior of the SplineEvaluator out of the domain where the break points of the B-splines are defined.
198 *
199 * @return The lower extrapolation rule.
200 *
201 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
202 */
203 lower_extrapolation_rule_type lower_extrapolation_rule() const
204 {
205 return m_lower_extrap_rule;
206 }
207
208 /**
209 * @brief Get the upper extrapolation rule.
210 *
211 * Extrapolation rules are functors used to define the behavior of the SplineEvaluator out of the domain where the break points of the B-splines are defined.
212 *
213 * @return The upper extrapolation rule.
214 *
215 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
216 */
217 upper_extrapolation_rule_type upper_extrapolation_rule() const
218 {
219 return m_upper_extrap_rule;
220 }
221
222 /**
223 * @brief Evaluate 1D spline function (described by its spline coefficients) at a given coordinate.
224 *
225 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
226 *
227 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
228 *
229 * @param coord_eval The coordinate where the spline is evaluated. Note that only the component along the dimension of interest is used.
230 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
231 *
232 * @return The value of the spline function at the desired coordinate.
233 */
234 template <class Layout, class... CoordsDims>
235 KOKKOS_FUNCTION double operator()(
236 ddc::Coordinate<CoordsDims...> const& coord_eval,
237 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
238 spline_coef) const
239 {
240 return eval(coord_eval, spline_coef);
241 }
242
243 /**
244 * @brief Evaluate spline function (described by its spline coefficients) on a mesh.
245 *
246 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
247 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
248 *
249 * This is not a multidimensional evaluation. This is a batched 1D evaluation. This means that for each slice of coordinates
250 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 1D set of
251 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
252 *
253 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
254 *
255 * @param[out] spline_eval The values of the spline function at the desired coordinates. For practical reasons those are
256 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
257 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
258 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
259 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
260 * the set of 1D spline coefficients retained to perform the evaluation).
261 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
262 */
263 template <class Layout1, class Layout2, class Layout3, class... CoordsDims>
264 void operator()(
265 ddc::ChunkSpan<double, batched_evaluation_domain_type, Layout1, memory_space> const
266 spline_eval,
267 ddc::ChunkSpan<
268 ddc::Coordinate<CoordsDims...> const,
269 batched_evaluation_domain_type,
270 Layout2,
271 memory_space> const coords_eval,
272 ddc::ChunkSpan<double const, batched_spline_domain_type, Layout3, memory_space> const
273 spline_coef) const
274 {
275 evaluation_domain_type const evaluation_domain(spline_eval.domain());
276 batch_domain_type const batch_domain(spline_eval.domain());
277
278 ddc::parallel_for_each(
279 "ddc_splines_evaluate",
280 exec_space(),
281 batch_domain,
282 KOKKOS_CLASS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
283 const auto spline_eval_1D = spline_eval[j];
284 const auto coords_eval_1D = coords_eval[j];
285 const auto spline_coef_1D = spline_coef[j];
286 for (auto const i : evaluation_domain) {
287 spline_eval_1D(i) = eval(coords_eval_1D(i), spline_coef_1D);
288 }
289 });
290 }
291
292 /**
293 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
294 *
295 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
296 * obtained via various methods, such as using a SplineBuilder.
297 *
298 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
299 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
300 *
301 * @return The derivative of the spline function at the desired coordinate.
302 */
303 template <class Layout, class... CoordsDims>
304 KOKKOS_FUNCTION double deriv(
305 ddc::Coordinate<CoordsDims...> const& coord_eval,
306 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
307 spline_coef) const
308 {
309 return eval_no_bc<eval_deriv_type>(coord_eval, spline_coef);
310 }
311
312 /**
313 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
314 *
315 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
316 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
317 *
318 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
319 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
320 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
321 *
322 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
323 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
324 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
325 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
326 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
327 * the set of 1D spline coefficients retained to perform the evaluation).
328 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
329 */
330 template <class Layout1, class Layout2, class Layout3, class... CoordsDims>
331 void deriv(
332 ddc::ChunkSpan<double, batched_evaluation_domain_type, Layout1, memory_space> const
333 spline_eval,
334 ddc::ChunkSpan<
335 ddc::Coordinate<CoordsDims...> const,
336 batched_evaluation_domain_type,
337 Layout2,
338 memory_space> const coords_eval,
339 ddc::ChunkSpan<double const, batched_spline_domain_type, Layout3, memory_space> const
340 spline_coef) const
341 {
342 evaluation_domain_type const evaluation_domain(spline_eval.domain());
343 batch_domain_type const batch_domain(spline_eval.domain());
344
345 ddc::parallel_for_each(
346 "ddc_splines_differentiate",
347 exec_space(),
348 batch_domain,
349 KOKKOS_CLASS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
350 const auto spline_eval_1D = spline_eval[j];
351 const auto coords_eval_1D = coords_eval[j];
352 const auto spline_coef_1D = spline_coef[j];
353 for (auto const i : evaluation_domain) {
354 spline_eval_1D(i)
355 = eval_no_bc<eval_deriv_type>(coords_eval_1D(i), spline_coef_1D);
356 }
357 });
358 }
359
360 /** @brief Perform batched 1D integrations of a spline function (described by its spline coefficients) along the dimension of interest and store results on a subdomain of batch_domain.
361 *
362 * The spline coefficients represent a spline function defined on a B-splines (basis splines). They can be obtained via the SplineBuilder.
363 *
364 * The integration is not performed in a multidimensional way (in any sense). This is a batched 1D integration.
365 * This means that for each element of integrals, the integration is performed with the 1D set of
366 * spline coefficients identified by the same DiscreteElement.
367 *
368 * @param[out] integrals The integrals of the spline function on the subdomain of batch_domain. For practical reasons those are
369 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
370 * points represented by this domain are unused and irrelevant.
371 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
372 */
373 template <class Layout1, class Layout2>
374 void integrate(
375 ddc::ChunkSpan<double, batch_domain_type, Layout1, memory_space> const integrals,
376 ddc::ChunkSpan<double const, batched_spline_domain_type, Layout2, memory_space> const
377 spline_coef) const
378 {
379 batch_domain_type const batch_domain(integrals.domain());
380 ddc::Chunk values_alloc(
381 ddc::DiscreteDomain<bsplines_type>(spline_coef.domain()),
382 ddc::KokkosAllocator<double, memory_space>());
383 ddc::ChunkSpan const values = values_alloc.span_view();
384 ddc::integrals(exec_space(), values);
385
386 ddc::parallel_for_each(
387 "ddc_splines_integrate",
388 exec_space(),
389 batch_domain,
390 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
391 integrals(j) = 0;
392 for (typename spline_domain_type::discrete_element_type const i :
393 values.domain()) {
394 integrals(j) += spline_coef(i, j) * values(i);
395 }
396 });
397 }
398
399private:
400 template <class Layout, class... CoordsDims>
401 KOKKOS_INLINE_FUNCTION double eval(
402 ddc::Coordinate<CoordsDims...> const& coord_eval,
403 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
404 spline_coef) const
405 {
406 ddc::Coordinate<continuous_dimension_type> coord_eval_interest
407 = ddc::select<continuous_dimension_type>(coord_eval);
408 if constexpr (bsplines_type::is_periodic()) {
409 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()
410 || coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
411 coord_eval_interest -= Kokkos::floor(
412 (coord_eval_interest
413 - ddc::discrete_space<bsplines_type>().rmin())
414 / ddc::discrete_space<bsplines_type>().length())
415 * ddc::discrete_space<bsplines_type>().length();
416 }
417 } else {
418 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()) {
419 return m_lower_extrap_rule(coord_eval_interest, spline_coef);
420 }
421 if (coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
422 return m_upper_extrap_rule(coord_eval_interest, spline_coef);
423 }
424 }
425 return eval_no_bc<eval_type>(coord_eval_interest, spline_coef);
426 }
427
428 template <class EvalType, class Layout, class... CoordsDims>
429 KOKKOS_INLINE_FUNCTION double eval_no_bc(
430 ddc::Coordinate<CoordsDims...> const& coord_eval,
431 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
432 spline_coef) const
433 {
434 static_assert(
435 std::is_same_v<EvalType, eval_type> || std::is_same_v<EvalType, eval_deriv_type>);
436 ddc::DiscreteElement<bsplines_type> jmin;
437 std::array<double, bsplines_type::degree() + 1> vals_ptr;
438 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
439 vals(vals_ptr.data());
440 ddc::Coordinate<continuous_dimension_type> const coord_eval_interest
441 = ddc::select<continuous_dimension_type>(coord_eval);
442 if constexpr (std::is_same_v<EvalType, eval_type>) {
443 jmin = ddc::discrete_space<bsplines_type>().eval_basis(vals, coord_eval_interest);
444 } else if constexpr (std::is_same_v<EvalType, eval_deriv_type>) {
445 jmin = ddc::discrete_space<bsplines_type>().eval_deriv(vals, coord_eval_interest);
446 }
447 double y = 0.0;
448 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
449 y += spline_coef(ddc::DiscreteElement<bsplines_type>(jmin + i)) * vals[i];
450 }
451 return y;
452 }
453};
454
455} // namespace ddc
friend class DiscreteDomain
A class to evaluate, differentiate or integrate a spline function.
void deriv(ddc::ChunkSpan< double, batched_evaluation_domain_type, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, batched_evaluation_domain_type, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type, Layout3, memory_space > const spline_coef) const
Differentiate spline function (described by its spline coefficients) on a mesh.
KOKKOS_FUNCTION double deriv(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
lower_extrapolation_rule_type lower_extrapolation_rule() const
Get the lower extrapolation rule.
SplineEvaluator(LowerExtrapolationRule const &lower_extrap_rule, UpperExtrapolationRule const &upper_extrap_rule)
Build a SplineEvaluator acting on batched_spline_domain.
SplineEvaluator(SplineEvaluator &&x)=default
Move-constructs.
SplineEvaluator & operator=(SplineEvaluator &&x)=default
Move-assigns.
SplineEvaluator(SplineEvaluator const &x)=default
Copy-constructs.
~SplineEvaluator()=default
Destructs.
void integrate(ddc::ChunkSpan< double, batch_domain_type, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, batched_spline_domain_type, Layout2, memory_space > const spline_coef) const
Perform batched 1D integrations of a spline function (described by its spline coefficients) along the...
void operator()(ddc::ChunkSpan< double, batched_evaluation_domain_type, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, batched_evaluation_domain_type, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type, Layout3, memory_space > const spline_coef) const
Evaluate spline function (described by its spline coefficients) on a mesh.
SplineEvaluator & operator=(SplineEvaluator const &x)=default
Copy-assigns.
upper_extrapolation_rule_type upper_extrapolation_rule() const
Get the upper extrapolation rule.
KOKKOS_FUNCTION double operator()(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Evaluate 1D spline function (described by its spline coefficients) at a given coordinate.
The top-level namespace of DDC.