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_Core.hpp>
14
15#include "deriv.hpp"
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 */
33template <
34 class ExecSpace,
35 class MemorySpace,
36 class BSplines,
37 class EvaluationDDim,
38 class LowerExtrapolationRule,
39 class UpperExtrapolationRule>
41{
42public:
43 /// @brief The type of the Kokkos execution space used by this class.
44 using exec_space = ExecSpace;
45
46 /// @brief The type of the Kokkos memory space used by this class.
47 using memory_space = MemorySpace;
48
49 /// @brief The type of the evaluation continuous dimension (continuous dimension of interest) used by this class.
50 using continuous_dimension_type = typename BSplines::continuous_dimension_type;
51
52 /// @brief The type of the evaluation discrete dimension (discrete dimension of interest) used by this class.
53 using evaluation_discrete_dimension_type = EvaluationDDim;
54
55 /// @brief The discrete dimension representing the B-splines.
56 using bsplines_type = BSplines;
57
58 /// @brief The type of the domain for the 1D evaluation mesh used by this class.
59 using evaluation_domain_type = ddc::DiscreteDomain<evaluation_discrete_dimension_type>;
60
61 /**
62 * @brief The type of the whole domain representing evaluation points.
63 *
64 * @tparam The batched discrete domain on which the interpolation points are defined.
65 */
66 template <
67 class BatchedInterpolationDDom,
68 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
69 using batched_evaluation_domain_type = BatchedInterpolationDDom;
70
71 /// @brief The type of the 1D spline domain corresponding to the dimension of interest.
72 using spline_domain_type = ddc::DiscreteDomain<bsplines_type>;
73
74 /**
75 * @brief The type of the batch domain (obtained by removing the dimension of interest
76 * from the whole domain).
77 *
78 * @tparam The batched discrete domain on which the interpolation points are defined.
79 */
80 template <
81 class BatchedInterpolationDDom,
82 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
83 using batch_domain_type
84 = ddc::remove_dims_of_t<BatchedInterpolationDDom, evaluation_discrete_dimension_type>;
85
86 /**
87 * @brief The type of the whole spline domain (cartesian product of 1D spline domain
88 * and batch domain) preserving the order of dimensions.
89 *
90 * @tparam The batched discrete domain on which the interpolation points are defined.
91 */
92 template <
93 class BatchedInterpolationDDom,
94 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
95 using batched_spline_domain_type = ddc::replace_dim_of_t<
96 BatchedInterpolationDDom,
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 <
264 class Layout1,
265 class Layout2,
266 class Layout3,
267 class BatchedInterpolationDDom,
268 class... CoordsDims>
269 void operator()(
270 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
271 spline_eval,
272 ddc::ChunkSpan<
273 ddc::Coordinate<CoordsDims...> const,
274 BatchedInterpolationDDom,
275 Layout2,
276 memory_space> const coords_eval,
277 ddc::ChunkSpan<
278 double const,
279 batched_spline_domain_type<BatchedInterpolationDDom>,
280 Layout3,
281 memory_space> const spline_coef) const
282 {
283 evaluation_domain_type const evaluation_domain(spline_eval.domain());
284 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
285
286 ddc::parallel_for_each(
287 "ddc_splines_evaluate",
288 exec_space(),
289 batch_domain,
290 KOKKOS_CLASS_LAMBDA(
291 typename batch_domain_type<
292 BatchedInterpolationDDom>::discrete_element_type const j) {
293 auto const spline_eval_1D = spline_eval[j];
294 auto const coords_eval_1D = coords_eval[j];
295 auto const spline_coef_1D = spline_coef[j];
296 for (auto const i : evaluation_domain) {
297 spline_eval_1D(i) = eval(coords_eval_1D(i), spline_coef_1D);
298 }
299 });
300 }
301
302 /**
303 * @brief Evaluate a spline function (described by its spline coefficients) on a mesh.
304 *
305 * The spline coefficients represent a spline function defined on a cartesian
306 * product of batch_domain and B-splines (basis splines). They can be obtained
307 * via various methods, such as using a SplineBuilder.
308 *
309 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
310 * This means that for each slice of spline_eval the evaluation is performed with
311 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
312 *
313 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline
314 * interpolation.
315 *
316 * @param[out] spline_eval The values of the spline function at the coordinates
317 * of the mesh.
318 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
319 */
320 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
321 void operator()(
322 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
323 spline_eval,
324 ddc::ChunkSpan<
325 double const,
326 batched_spline_domain_type<BatchedInterpolationDDom>,
327 Layout2,
328 memory_space> const spline_coef) const
329 {
330 evaluation_domain_type const evaluation_domain(spline_eval.domain());
331 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
332
333 ddc::parallel_for_each(
334 "ddc_splines_evaluate",
335 exec_space(),
336 batch_domain,
337 KOKKOS_CLASS_LAMBDA(
338 typename batch_domain_type<
339 BatchedInterpolationDDom>::discrete_element_type const j) {
340 auto const spline_eval_1D = spline_eval[j];
341 auto const spline_coef_1D = spline_coef[j];
342 for (auto const i : evaluation_domain) {
343 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
344 = ddc::coordinate(i);
345 spline_eval_1D(i) = eval(coord_eval_1D, spline_coef_1D);
346 }
347 });
348 }
349
350#if defined(DDC_BUILD_DEPRECATED_CODE)
351 /**
352 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
353 *
354 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
355 * obtained via various methods, such as using a SplineBuilder.
356 *
357 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
358 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
359 *
360 * @return The derivative of the spline function at the desired coordinate.
361 */
362 template <class Layout, class... CoordsDims>
363 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] KOKKOS_FUNCTION double
364 deriv(ddc::Coordinate<CoordsDims...> const& coord_eval,
365 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const spline_coef)
366 const
367 {
368 return deriv(
369 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
370 coord_eval,
371 spline_coef);
372 }
373
374 /**
375 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
376 *
377 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
378 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
379 *
380 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
381 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
382 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
383 *
384 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
385 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
386 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
387 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
388 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
389 * the set of 1D spline coefficients retained to perform the evaluation).
390 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
391 */
392 template <
393 class Layout1,
394 class Layout2,
395 class Layout3,
396 class BatchedInterpolationDDom,
397 class... CoordsDims>
398 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
399 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
400 spline_eval,
401 ddc::ChunkSpan<
402 ddc::Coordinate<CoordsDims...> const,
403 BatchedInterpolationDDom,
404 Layout2,
405 memory_space> const coords_eval,
406 ddc::ChunkSpan<
407 double const,
408 batched_spline_domain_type<BatchedInterpolationDDom>,
409 Layout3,
410 memory_space> const spline_coef) const
411 {
412 return deriv(
413 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
414 spline_eval,
415 coords_eval,
416 spline_coef);
417 }
418
419 /**
420 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
421 *
422 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
423 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
424 *
425 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
426 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
427 * This means that for each slice of spline_eval the evaluation is performed with
428 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
429 *
430 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
431 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
432 */
433 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
434 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
435 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
436 spline_eval,
437 ddc::ChunkSpan<
438 double const,
439 batched_spline_domain_type<BatchedInterpolationDDom>,
440 Layout2,
441 memory_space> const spline_coef) const
442 {
443 return deriv(
444 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
445 spline_eval,
446 spline_coef);
447 }
448#endif
449
450 /**
451 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
452 *
453 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
454 * obtained via various methods, such as using a SplineBuilder.
455 *
456 * @param deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
457 * If the dimension is not present, the order of derivation is considered to be 0.
458 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
459 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
460 *
461 * @return The derivative of the spline function at the desired coordinate.
462 */
463 template <class DElem, class Layout, class... CoordsDims>
464 KOKKOS_FUNCTION double deriv(
465 DElem const& deriv_order,
466 ddc::Coordinate<CoordsDims...> const& coord_eval,
467 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
468 spline_coef) const
469 {
470 static_assert(is_discrete_element_v<DElem>);
471 return eval_no_bc(deriv_order, coord_eval, spline_coef);
472 }
473
474 /**
475 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
476 *
477 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
478 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
479 *
480 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
481 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
482 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
483 *
484 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
485 * If the dimension is not present, the order of derivation is considered to be 0.
486 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
487 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
488 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
489 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
490 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
491 * the set of 1D spline coefficients retained to perform the evaluation).
492 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
493 */
494 template <
495 class DElem,
496 class Layout1,
497 class Layout2,
498 class Layout3,
499 class BatchedInterpolationDDom,
500 class... CoordsDims>
501 void deriv(
502 DElem const& deriv_order,
503 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
504 spline_eval,
505 ddc::ChunkSpan<
506 ddc::Coordinate<CoordsDims...> const,
507 BatchedInterpolationDDom,
508 Layout2,
509 memory_space> const coords_eval,
510 ddc::ChunkSpan<
511 double const,
512 batched_spline_domain_type<BatchedInterpolationDDom>,
513 Layout3,
514 memory_space> const spline_coef) const
515 {
516 static_assert(is_discrete_element_v<DElem>);
517
518 evaluation_domain_type const evaluation_domain(spline_eval.domain());
519 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
520
521 ddc::parallel_for_each(
522 "ddc_splines_differentiate",
523 exec_space(),
524 batch_domain,
525 KOKKOS_CLASS_LAMBDA(
526 typename batch_domain_type<
527 BatchedInterpolationDDom>::discrete_element_type const j) {
528 auto const spline_eval_1D = spline_eval[j];
529 auto const coords_eval_1D = coords_eval[j];
530 auto const spline_coef_1D = spline_coef[j];
531 for (auto const i : evaluation_domain) {
532 spline_eval_1D(i)
533 = eval_no_bc(deriv_order, coords_eval_1D(i), spline_coef_1D);
534 }
535 });
536 }
537
538 /**
539 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
540 *
541 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
542 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
543 *
544 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
545 * This means that for each slice of spline_eval the evaluation is performed with
546 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
547 *
548 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
549 * If the dimension is not present, the order of derivation is considered to be 0.
550 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
551 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
552 */
553 template <class DElem, class Layout1, class Layout2, class BatchedInterpolationDDom>
554 void deriv(
555 DElem const& deriv_order,
556 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
557 spline_eval,
558 ddc::ChunkSpan<
559 double const,
560 batched_spline_domain_type<BatchedInterpolationDDom>,
561 Layout2,
562 memory_space> const spline_coef) const
563 {
564 static_assert(is_discrete_element_v<DElem>);
565
566 evaluation_domain_type const evaluation_domain(spline_eval.domain());
567 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
568
569 ddc::parallel_for_each(
570 "ddc_splines_differentiate",
571 exec_space(),
572 batch_domain,
573 KOKKOS_CLASS_LAMBDA(
574 typename batch_domain_type<
575 BatchedInterpolationDDom>::discrete_element_type const j) {
576 auto const spline_eval_1D = spline_eval[j];
577 auto const spline_coef_1D = spline_coef[j];
578 for (auto const i : evaluation_domain) {
579 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
580 = ddc::coordinate(i);
581 spline_eval_1D(i) = eval_no_bc(deriv_order, coord_eval_1D, spline_coef_1D);
582 }
583 });
584 }
585
586 /** @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.
587 *
588 * The spline coefficients represent a spline function defined on a B-splines (basis splines). They can be obtained via the SplineBuilder.
589 *
590 * The integration is not performed in a multidimensional way (in any sense). This is a batched 1D integration.
591 * This means that for each element of integrals, the integration is performed with the 1D set of
592 * spline coefficients identified by the same DiscreteElement.
593 *
594 * @param[out] integrals The integrals of the spline function on the subdomain of batch_domain. For practical reasons those are
595 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
596 * points represented by this domain are unused and irrelevant.
597 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
598 */
599 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
600 void integrate(
601 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
602 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
603 spline_coef) const
604 {
605 static_assert(
606 ddc::in_tags_v<bsplines_type, to_type_seq_t<BatchedSplineDDom>>,
607 "The spline coefficients domain must contain the bsplines dimension");
608 using batch_domain_type = ddc::remove_dims_of_t<BatchedSplineDDom, bsplines_type>;
609 static_assert(
610 std::is_same_v<batch_domain_type, BatchedDDom>,
611 "The integrals domain must only contain the batch dimensions");
612
613 BatchedDDom const batch_domain(integrals.domain());
614 ddc::Chunk values_alloc(
615 ddc::DiscreteDomain<bsplines_type>(spline_coef.domain()),
616 ddc::KokkosAllocator<double, memory_space>());
617 ddc::ChunkSpan const values = values_alloc.span_view();
618 ddc::integrals(exec_space(), values);
619
620 ddc::parallel_for_each(
621 "ddc_splines_integrate",
622 exec_space(),
623 batch_domain,
624 KOKKOS_LAMBDA(typename BatchedDDom::discrete_element_type const j) {
625 integrals(j) = 0;
626 for (typename spline_domain_type::discrete_element_type const i :
627 values.domain()) {
628 integrals(j) += spline_coef(i, j) * values(i);
629 }
630 });
631 }
632
633private:
634 template <class Layout, class... CoordsDims>
635 KOKKOS_INLINE_FUNCTION double eval(
636 ddc::Coordinate<CoordsDims...> const& coord_eval,
637 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
638 spline_coef) const
639 {
640 ddc::Coordinate<continuous_dimension_type> coord_eval_interest(coord_eval);
641 if constexpr (bsplines_type::is_periodic()) {
642 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()
643 || coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
644 coord_eval_interest -= Kokkos::floor(
645 (coord_eval_interest
646 - ddc::discrete_space<bsplines_type>().rmin())
647 / ddc::discrete_space<bsplines_type>().length())
648 * ddc::discrete_space<bsplines_type>().length();
649 }
650 } else {
651 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()) {
652 return m_lower_extrap_rule(coord_eval_interest, spline_coef);
653 }
654 if (coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
655 return m_upper_extrap_rule(coord_eval_interest, spline_coef);
656 }
657 }
658 return eval_no_bc(ddc::DiscreteElement<>(), coord_eval_interest, spline_coef);
659 }
660
661 template <class... DerivDims, class Layout, class... CoordsDims>
662 KOKKOS_INLINE_FUNCTION double eval_no_bc(
663 ddc::DiscreteElement<DerivDims...> const& deriv_order,
664 ddc::Coordinate<CoordsDims...> const& coord_eval,
665 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
666 spline_coef) const
667 {
668 using deriv_dim = ddc::Deriv<continuous_dimension_type>;
669 static_assert(
670 sizeof...(DerivDims) == 0
671 || type_seq_same_v<
672 detail::TypeSeq<DerivDims...>,
673 detail::TypeSeq<deriv_dim>>,
674 "The only valid dimension for deriv_order is Deriv<Dim>");
675
676 ddc::DiscreteElement<bsplines_type> jmin;
677 std::array<double, bsplines_type::degree() + 1> vals_ptr;
678 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
679 vals(vals_ptr.data());
680 ddc::Coordinate<continuous_dimension_type> const coord_eval_interest(coord_eval);
681
682 if constexpr (sizeof...(DerivDims) == 0) {
683 jmin = ddc::discrete_space<bsplines_type>().eval_basis(vals, coord_eval_interest);
684 } else {
685 auto const order = deriv_order.uid();
686 KOKKOS_ASSERT(order > 0 && order <= bsplines_type::degree())
687
688 std::array<double, (bsplines_type::degree() + 1) * (bsplines_type::degree() + 1)>
689 derivs_ptr;
690 Kokkos::mdspan<
691 double,
692 Kokkos::extents<
693 std::size_t,
694 bsplines_type::degree() + 1,
695 Kokkos::dynamic_extent>> const derivs(derivs_ptr.data(), order + 1);
696
697 jmin = ddc::discrete_space<bsplines_type>()
698 .eval_basis_and_n_derivs(derivs, coord_eval_interest, order);
699
700 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
701 vals[i] = DDC_MDSPAN_ACCESS_OP(derivs, i, order);
702 }
703 }
704
705 double y = 0.0;
706 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
707 y += spline_coef(ddc::DiscreteElement<bsplines_type>(jmin + i)) * vals[i];
708 }
709 return y;
710 }
711};
712
713} // namespace ddc
friend class ChunkSpan
friend class Chunk
Definition chunk.hpp:83
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