DDC 0.10.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 "deriv.hpp"
17#include "integrals.hpp"
18#include "periodic_extrapolation_rule.hpp"
19
20namespace ddc {
21
22/**
23 * @brief A class to evaluate, differentiate or integrate a spline function.
24 *
25 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a spline function.
26 *
27 * @tparam ExecSpace The Kokkos execution space on which the spline evaluation is performed.
28 * @tparam MemorySpace The Kokkos memory space on which the data (spline coefficients and evaluation) is stored.
29 * @tparam BSplines The discrete dimension representing the B-splines.
30 * @tparam EvaluationDDim The discrete dimension on which evaluation points are defined.
31 * @tparam LowerExtrapolationRule The lower extrapolation rule type.
32 * @tparam UpperExtrapolationRule The upper extrapolation rule type.
33 */
34template <
35 class ExecSpace,
36 class MemorySpace,
37 class BSplines,
38 class EvaluationDDim,
39 class LowerExtrapolationRule,
40 class UpperExtrapolationRule>
42{
43public:
44 /// @brief The type of the Kokkos execution space used by this class.
45 using exec_space = ExecSpace;
46
47 /// @brief The type of the Kokkos memory space used by this class.
48 using memory_space = MemorySpace;
49
50 /// @brief The type of the evaluation continuous dimension (continuous dimension of interest) used by this class.
51 using continuous_dimension_type = typename BSplines::continuous_dimension_type;
52
53 /// @brief The type of the evaluation discrete dimension (discrete dimension of interest) used by this class.
54 using evaluation_discrete_dimension_type = EvaluationDDim;
55
56 /// @brief The discrete dimension representing the B-splines.
57 using bsplines_type = BSplines;
58
59 /// @brief The type of the domain for the 1D evaluation mesh used by this class.
60 using evaluation_domain_type = ddc::DiscreteDomain<evaluation_discrete_dimension_type>;
61
62 /**
63 * @brief The type of the whole domain representing evaluation points.
64 *
65 * @tparam The batched discrete domain on which the interpolation points are defined.
66 */
67 template <
68 class BatchedInterpolationDDom,
69 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
70 using batched_evaluation_domain_type = BatchedInterpolationDDom;
71
72 /// @brief The type of the 1D spline domain corresponding to the dimension of interest.
73 using spline_domain_type = ddc::DiscreteDomain<bsplines_type>;
74
75 /**
76 * @brief The type of the batch domain (obtained by removing the dimension of interest
77 * from the whole domain).
78 *
79 * @tparam The batched discrete domain on which the interpolation points are defined.
80 */
81 template <
82 class BatchedInterpolationDDom,
83 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
84 using batch_domain_type
85 = ddc::remove_dims_of_t<BatchedInterpolationDDom, evaluation_discrete_dimension_type>;
86
87 /**
88 * @brief The type of the whole spline domain (cartesian product of 1D spline domain
89 * and batch domain) preserving the order of dimensions.
90 *
91 * @tparam The batched discrete domain on which the interpolation points are defined.
92 */
93 template <
94 class BatchedInterpolationDDom,
95 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
96 using batched_spline_domain_type = ddc::replace_dim_of_t<
97 BatchedInterpolationDDom,
98 evaluation_discrete_dimension_type,
99 bsplines_type>;
100
101 /// @brief The type of the extrapolation rule at the lower boundary.
102 using lower_extrapolation_rule_type = LowerExtrapolationRule;
103
104 /// @brief The type of the extrapolation rule at the upper boundary.
105 using upper_extrapolation_rule_type = UpperExtrapolationRule;
106
107
108private:
109 LowerExtrapolationRule m_lower_extrap_rule;
110
111 UpperExtrapolationRule m_upper_extrap_rule;
112
113public:
114 static_assert(
115 std::is_same_v<LowerExtrapolationRule,
116 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
117 == bsplines_type::is_periodic()
118 && std::is_same_v<
119 UpperExtrapolationRule,
120 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
121 == bsplines_type::is_periodic(),
122 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
123 static_assert(
124 std::is_invocable_r_v<
125 double,
126 LowerExtrapolationRule,
127 ddc::Coordinate<continuous_dimension_type>,
128 ddc::ChunkSpan<
129 double const,
130 spline_domain_type,
131 Kokkos::layout_right,
132 memory_space>>,
133 "LowerExtrapolationRule::operator() has to be callable with usual arguments.");
134 static_assert(
135 std::is_invocable_r_v<
136 double,
137 UpperExtrapolationRule,
138 ddc::Coordinate<continuous_dimension_type>,
139 ddc::ChunkSpan<
140 double const,
141 spline_domain_type,
142 Kokkos::layout_right,
143 memory_space>>,
144 "UpperExtrapolationRule::operator() has to be callable with usual arguments.");
145
146 /**
147 * @brief Build a SplineEvaluator acting on batched_spline_domain.
148 *
149 * @param lower_extrap_rule The extrapolation rule at the lower boundary.
150 * @param upper_extrap_rule The extrapolation rule at the upper boundary.
151 *
152 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
153 */
154 explicit SplineEvaluator(
155 LowerExtrapolationRule const& lower_extrap_rule,
156 UpperExtrapolationRule const& upper_extrap_rule)
157 : m_lower_extrap_rule(lower_extrap_rule)
158 , m_upper_extrap_rule(upper_extrap_rule)
159 {
160 }
161
162 /**
163 * @brief Copy-constructs.
164 *
165 * @param x A reference to another SplineEvaluator.
166 */
167 SplineEvaluator(SplineEvaluator const& x) = default;
168
169 /**
170 * @brief Move-constructs.
171 *
172 * @param x An rvalue to another SplineEvaluator.
173 */
174 SplineEvaluator(SplineEvaluator&& x) = default;
175
176 /// @brief Destructs
177 ~SplineEvaluator() = default;
178
179 /**
180 * @brief Copy-assigns.
181 *
182 * @param x A reference to another SplineEvaluator.
183 * @return A reference to this object.
184 */
185 SplineEvaluator& operator=(SplineEvaluator const& x) = default;
186
187 /**
188 * @brief Move-assigns.
189 *
190 * @param x An rvalue to another SplineEvaluator.
191 * @return A reference to this object.
192 */
194
195 /**
196 * @brief Get the lower extrapolation rule.
197 *
198 * 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.
199 *
200 * @return The lower extrapolation rule.
201 *
202 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
203 */
204 lower_extrapolation_rule_type lower_extrapolation_rule() const
205 {
206 return m_lower_extrap_rule;
207 }
208
209 /**
210 * @brief Get the upper extrapolation rule.
211 *
212 * 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.
213 *
214 * @return The upper extrapolation rule.
215 *
216 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
217 */
218 upper_extrapolation_rule_type upper_extrapolation_rule() const
219 {
220 return m_upper_extrap_rule;
221 }
222
223 /**
224 * @brief Evaluate 1D spline function (described by its spline coefficients) at a given coordinate.
225 *
226 * 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.
227 *
228 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
229 *
230 * @param coord_eval The coordinate where the spline is evaluated. Note that only the component along the dimension of interest is used.
231 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
232 *
233 * @return The value of the spline function at the desired coordinate.
234 */
235 template <class Layout, class... CoordsDims>
236 KOKKOS_FUNCTION double operator()(
237 ddc::Coordinate<CoordsDims...> const& coord_eval,
238 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
239 spline_coef) const
240 {
241 return eval(coord_eval, spline_coef);
242 }
243
244 /**
245 * @brief Evaluate spline function (described by its spline coefficients) on a mesh.
246 *
247 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
248 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
249 *
250 * This is not a multidimensional evaluation. This is a batched 1D evaluation. This means that for each slice of coordinates
251 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 1D set of
252 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
253 *
254 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
255 *
256 * @param[out] spline_eval The values of the spline function at the desired coordinates. For practical reasons those are
257 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
258 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
259 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
260 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
261 * the set of 1D spline coefficients retained to perform the evaluation).
262 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
263 */
264 template <
265 class Layout1,
266 class Layout2,
267 class Layout3,
268 class BatchedInterpolationDDom,
269 class... CoordsDims>
270 void operator()(
271 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
272 spline_eval,
273 ddc::ChunkSpan<
274 ddc::Coordinate<CoordsDims...> const,
275 BatchedInterpolationDDom,
276 Layout2,
277 memory_space> const coords_eval,
278 ddc::ChunkSpan<
279 double const,
280 batched_spline_domain_type<BatchedInterpolationDDom>,
281 Layout3,
282 memory_space> const spline_coef) const
283 {
284 evaluation_domain_type const evaluation_domain(spline_eval.domain());
285 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
286
287 ddc::parallel_for_each(
288 "ddc_splines_evaluate",
289 exec_space(),
290 batch_domain,
291 KOKKOS_CLASS_LAMBDA(
292 typename batch_domain_type<
293 BatchedInterpolationDDom>::discrete_element_type const j) {
294 auto const spline_eval_1D = spline_eval[j];
295 auto const coords_eval_1D = coords_eval[j];
296 auto const spline_coef_1D = spline_coef[j];
297 for (auto const i : evaluation_domain) {
298 spline_eval_1D(i) = eval(coords_eval_1D(i), spline_coef_1D);
299 }
300 });
301 }
302
303 /**
304 * @brief Evaluate a spline function (described by its spline coefficients) on a mesh.
305 *
306 * The spline coefficients represent a spline function defined on a cartesian
307 * product of batch_domain and B-splines (basis splines). They can be obtained
308 * via various methods, such as using a SplineBuilder.
309 *
310 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
311 * This means that for each slice of spline_eval the evaluation is performed with
312 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
313 *
314 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline
315 * interpolation.
316 *
317 * @param[out] spline_eval The values of the spline function at the coordinates
318 * of the mesh.
319 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
320 */
321 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
322 void operator()(
323 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
324 spline_eval,
325 ddc::ChunkSpan<
326 double const,
327 batched_spline_domain_type<BatchedInterpolationDDom>,
328 Layout2,
329 memory_space> const spline_coef) const
330 {
331 evaluation_domain_type const evaluation_domain(spline_eval.domain());
332 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
333
334 ddc::parallel_for_each(
335 "ddc_splines_evaluate",
336 exec_space(),
337 batch_domain,
338 KOKKOS_CLASS_LAMBDA(
339 typename batch_domain_type<
340 BatchedInterpolationDDom>::discrete_element_type const j) {
341 auto const spline_eval_1D = spline_eval[j];
342 auto const spline_coef_1D = spline_coef[j];
343 for (auto const i : evaluation_domain) {
344 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
345 = ddc::coordinate(i);
346 spline_eval_1D(i) = eval(coord_eval_1D, spline_coef_1D);
347 }
348 });
349 }
350
351#if defined(DDC_BUILD_DEPRECATED_CODE)
352 /**
353 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
354 *
355 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
356 * obtained via various methods, such as using a SplineBuilder.
357 *
358 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
359 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
360 *
361 * @return The derivative of the spline function at the desired coordinate.
362 */
363 template <class Layout, class... CoordsDims>
364 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] KOKKOS_FUNCTION double
365 deriv(ddc::Coordinate<CoordsDims...> const& coord_eval,
366 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const spline_coef)
367 const
368 {
369 return deriv(
370 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
371 coord_eval,
372 spline_coef);
373 }
374
375 /**
376 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
377 *
378 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
379 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
380 *
381 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
382 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
383 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
384 *
385 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
386 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
387 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
388 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
389 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
390 * the set of 1D spline coefficients retained to perform the evaluation).
391 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
392 */
393 template <
394 class Layout1,
395 class Layout2,
396 class Layout3,
397 class BatchedInterpolationDDom,
398 class... CoordsDims>
399 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
400 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
401 spline_eval,
402 ddc::ChunkSpan<
403 ddc::Coordinate<CoordsDims...> const,
404 BatchedInterpolationDDom,
405 Layout2,
406 memory_space> const coords_eval,
407 ddc::ChunkSpan<
408 double const,
409 batched_spline_domain_type<BatchedInterpolationDDom>,
410 Layout3,
411 memory_space> const spline_coef) const
412 {
413 return deriv(
414 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
415 spline_eval,
416 coords_eval,
417 spline_coef);
418 }
419
420 /**
421 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
422 *
423 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
424 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
425 *
426 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
427 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
428 * This means that for each slice of spline_eval the evaluation is performed with
429 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
430 *
431 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
432 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
433 */
434 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
435 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
436 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
437 spline_eval,
438 ddc::ChunkSpan<
439 double const,
440 batched_spline_domain_type<BatchedInterpolationDDom>,
441 Layout2,
442 memory_space> const spline_coef) const
443 {
444 return deriv(
445 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
446 spline_eval,
447 spline_coef);
448 }
449#endif
450
451 /**
452 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
453 *
454 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
455 * obtained via various methods, such as using a SplineBuilder.
456 *
457 * @param deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
458 * If the dimension is not present, the order of derivation is considered to be 0.
459 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
460 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
461 *
462 * @return The derivative of the spline function at the desired coordinate.
463 */
464 template <class DElem, class Layout, class... CoordsDims>
465 KOKKOS_FUNCTION double deriv(
466 DElem const& deriv_order,
467 ddc::Coordinate<CoordsDims...> const& coord_eval,
468 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
469 spline_coef) const
470 {
471 static_assert(is_discrete_element_v<DElem>);
472 return eval_no_bc(deriv_order, coord_eval, spline_coef);
473 }
474
475 /**
476 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
477 *
478 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
479 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
480 *
481 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
482 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
483 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
484 *
485 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
486 * If the dimension is not present, the order of derivation is considered to be 0.
487 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
488 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
489 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
490 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
491 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
492 * the set of 1D spline coefficients retained to perform the evaluation).
493 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
494 */
495 template <
496 class DElem,
497 class Layout1,
498 class Layout2,
499 class Layout3,
500 class BatchedInterpolationDDom,
501 class... CoordsDims>
502 void deriv(
503 DElem const& deriv_order,
504 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
505 spline_eval,
506 ddc::ChunkSpan<
507 ddc::Coordinate<CoordsDims...> const,
508 BatchedInterpolationDDom,
509 Layout2,
510 memory_space> const coords_eval,
511 ddc::ChunkSpan<
512 double const,
513 batched_spline_domain_type<BatchedInterpolationDDom>,
514 Layout3,
515 memory_space> const spline_coef) const
516 {
517 static_assert(is_discrete_element_v<DElem>);
518
519 evaluation_domain_type const evaluation_domain(spline_eval.domain());
520 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
521
522 ddc::parallel_for_each(
523 "ddc_splines_differentiate",
524 exec_space(),
525 batch_domain,
526 KOKKOS_CLASS_LAMBDA(
527 typename batch_domain_type<
528 BatchedInterpolationDDom>::discrete_element_type const j) {
529 auto const spline_eval_1D = spline_eval[j];
530 auto const coords_eval_1D = coords_eval[j];
531 auto const spline_coef_1D = spline_coef[j];
532 for (auto const i : evaluation_domain) {
533 spline_eval_1D(i)
534 = eval_no_bc(deriv_order, coords_eval_1D(i), spline_coef_1D);
535 }
536 });
537 }
538
539 /**
540 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
541 *
542 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
543 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
544 *
545 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
546 * This means that for each slice of spline_eval the evaluation is performed with
547 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
548 *
549 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
550 * If the dimension is not present, the order of derivation is considered to be 0.
551 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
552 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
553 */
554 template <class DElem, class Layout1, class Layout2, class BatchedInterpolationDDom>
555 void deriv(
556 DElem const& deriv_order,
557 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
558 spline_eval,
559 ddc::ChunkSpan<
560 double const,
561 batched_spline_domain_type<BatchedInterpolationDDom>,
562 Layout2,
563 memory_space> const spline_coef) const
564 {
565 static_assert(is_discrete_element_v<DElem>);
566
567 evaluation_domain_type const evaluation_domain(spline_eval.domain());
568 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
569
570 ddc::parallel_for_each(
571 "ddc_splines_differentiate",
572 exec_space(),
573 batch_domain,
574 KOKKOS_CLASS_LAMBDA(
575 typename batch_domain_type<
576 BatchedInterpolationDDom>::discrete_element_type const j) {
577 auto const spline_eval_1D = spline_eval[j];
578 auto const spline_coef_1D = spline_coef[j];
579 for (auto const i : evaluation_domain) {
580 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
581 = ddc::coordinate(i);
582 spline_eval_1D(i) = eval_no_bc(deriv_order, coord_eval_1D, spline_coef_1D);
583 }
584 });
585 }
586
587 /** @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.
588 *
589 * The spline coefficients represent a spline function defined on a B-splines (basis splines). They can be obtained via the SplineBuilder.
590 *
591 * The integration is not performed in a multidimensional way (in any sense). This is a batched 1D integration.
592 * This means that for each element of integrals, the integration is performed with the 1D set of
593 * spline coefficients identified by the same DiscreteElement.
594 *
595 * @param[out] integrals The integrals of the spline function on the subdomain of batch_domain. For practical reasons those are
596 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
597 * points represented by this domain are unused and irrelevant.
598 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
599 */
600 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
601 void integrate(
602 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
603 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
604 spline_coef) const
605 {
606 static_assert(
607 ddc::in_tags_v<bsplines_type, to_type_seq_t<BatchedSplineDDom>>,
608 "The spline coefficients domain must contain the bsplines dimension");
609 using batch_domain_type = ddc::remove_dims_of_t<BatchedSplineDDom, bsplines_type>;
610 static_assert(
611 std::is_same_v<batch_domain_type, BatchedDDom>,
612 "The integrals domain must only contain the batch dimensions");
613
614 BatchedDDom const batch_domain(integrals.domain());
615 ddc::Chunk values_alloc(
616 ddc::DiscreteDomain<bsplines_type>(spline_coef.domain()),
617 ddc::KokkosAllocator<double, memory_space>());
618 ddc::ChunkSpan const values = values_alloc.span_view();
619 ddc::integrals(exec_space(), values);
620
621 ddc::parallel_for_each(
622 "ddc_splines_integrate",
623 exec_space(),
624 batch_domain,
625 KOKKOS_LAMBDA(typename BatchedDDom::discrete_element_type const j) {
626 integrals(j) = 0;
627 for (typename spline_domain_type::discrete_element_type const i :
628 values.domain()) {
629 integrals(j) += spline_coef(i, j) * values(i);
630 }
631 });
632 }
633
634private:
635 template <class Layout, class... CoordsDims>
636 KOKKOS_INLINE_FUNCTION double eval(
637 ddc::Coordinate<CoordsDims...> const& coord_eval,
638 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
639 spline_coef) const
640 {
641 ddc::Coordinate<continuous_dimension_type> coord_eval_interest(coord_eval);
642 if constexpr (bsplines_type::is_periodic()) {
643 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()
644 || coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
645 coord_eval_interest -= Kokkos::floor(
646 (coord_eval_interest
647 - ddc::discrete_space<bsplines_type>().rmin())
648 / ddc::discrete_space<bsplines_type>().length())
649 * ddc::discrete_space<bsplines_type>().length();
650 }
651 } else {
652 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()) {
653 return m_lower_extrap_rule(coord_eval_interest, spline_coef);
654 }
655 if (coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
656 return m_upper_extrap_rule(coord_eval_interest, spline_coef);
657 }
658 }
659 return eval_no_bc(ddc::DiscreteElement<>(), coord_eval_interest, spline_coef);
660 }
661
662 template <class... DerivDims, class Layout, class... CoordsDims>
663 KOKKOS_INLINE_FUNCTION double eval_no_bc(
664 ddc::DiscreteElement<DerivDims...> const& deriv_order,
665 ddc::Coordinate<CoordsDims...> const& coord_eval,
666 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
667 spline_coef) const
668 {
669 using deriv_dim = ddc::Deriv<continuous_dimension_type>;
670 static_assert(
671 sizeof...(DerivDims) == 0
672 || type_seq_same_v<
673 detail::TypeSeq<DerivDims...>,
674 detail::TypeSeq<deriv_dim>>,
675 "The only valid dimension for deriv_order is Deriv<Dim>");
676
677 ddc::DiscreteElement<bsplines_type> jmin;
678 std::array<double, bsplines_type::degree() + 1> vals_ptr;
679 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
680 vals(vals_ptr.data());
681 ddc::Coordinate<continuous_dimension_type> const coord_eval_interest(coord_eval);
682
683 if constexpr (sizeof...(DerivDims) == 0) {
684 jmin = ddc::discrete_space<bsplines_type>().eval_basis(vals, coord_eval_interest);
685 } else {
686 auto const order = deriv_order.uid();
687 KOKKOS_ASSERT(order > 0 && order <= bsplines_type::degree())
688
689 std::array<double, (bsplines_type::degree() + 1) * (bsplines_type::degree() + 1)>
690 derivs_ptr;
691 Kokkos::mdspan<
692 double,
693 Kokkos::extents<
694 std::size_t,
695 bsplines_type::degree() + 1,
696 Kokkos::dynamic_extent>> const derivs(derivs_ptr.data(), order + 1);
697
698 jmin = ddc::discrete_space<bsplines_type>()
699 .eval_basis_and_n_derivs(derivs, coord_eval_interest, order);
700
701 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
702 vals[i] = DDC_MDSPAN_ACCESS_OP(derivs, i, order);
703 }
704 }
705
706 double y = 0.0;
707 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
708 y += spline_coef(ddc::DiscreteElement<bsplines_type>(jmin + i)) * vals[i];
709 }
710 return y;
711 }
712};
713
714} // namespace ddc
friend class ChunkSpan
friend class Chunk
Definition chunk.hpp:81
friend class DiscreteDomain
KOKKOS_DEFAULTED_FUNCTION constexpr DiscreteElement()=default
A class to evaluate, differentiate or integrate a spline function.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Evaluate spline function (described by its spline coefficients) on a mesh.
upper_extrapolation_rule_type upper_extrapolation_rule() const
Get the upper extrapolation rule.
SplineEvaluator & operator=(SplineEvaluator const &x)=default
Copy-assigns.
SplineEvaluator & operator=(SplineEvaluator &&x)=default
Move-assigns.
SplineEvaluator(LowerExtrapolationRule const &lower_extrap_rule, UpperExtrapolationRule const &upper_extrap_rule)
Build a SplineEvaluator acting on batched_spline_domain.
lower_extrapolation_rule_type lower_extrapolation_rule() const
Get the lower extrapolation rule.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) on a mesh.
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.
KOKKOS_FUNCTION double deriv(DElem const &deriv_order, 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.
SplineEvaluator(SplineEvaluator const &x)=default
Copy-constructs.
SplineEvaluator(SplineEvaluator &&x)=default
Move-constructs.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) on a mesh.
void integrate(ddc::ChunkSpan< double, BatchedDDom, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, BatchedSplineDDom, 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, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Evaluate a spline function (described by its spline coefficients) on a mesh.
~SplineEvaluator()=default
Destructs.
The top-level namespace of DDC.
A templated struct representing a discrete dimension storing the derivatives of a function along a co...
Definition deriv.hpp:15