DDC 0.10.0
Loading...
Searching...
No Matches
spline_evaluator_3d.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 3D spline function.
23 *
24 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a 3D 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 BSplines1 The discrete dimension representing the B-splines along the first dimension of interest.
29 * @tparam BSplines2 The discrete dimension representing the B-splines along the second dimension of interest.
30 * @tparam BSplines3 The discrete dimension representing the B-splines along the third dimension of interest.
31 * @tparam EvaluationDDim1 The first discrete dimension on which evaluation points are defined.
32 * @tparam EvaluationDDim2 The second discrete dimension on which evaluation points are defined.
33 * @tparam EvaluationDDim3 The third discrete dimension on which evaluation points are defined.
34 * @tparam LowerExtrapolationRule1 The lower extrapolation rule type along first dimension of interest.
35 * @tparam UpperExtrapolationRule1 The upper extrapolation rule type along first dimension of interest.
36 * @tparam LowerExtrapolationRule2 The lower extrapolation rule type along second dimension of interest.
37 * @tparam UpperExtrapolationRule2 The upper extrapolation rule type along second dimension of interest.
38 * @tparam LowerExtrapolationRule3 The lower extrapolation rule type along third dimension of interest.
39 * @tparam UpperExtrapolationRule3 The upper extrapolation rule type along third dimension of interest.
40 */
41template <
42 class ExecSpace,
43 class MemorySpace,
44 class BSplines1,
45 class BSplines2,
46 class BSplines3,
47 class EvaluationDDim1,
48 class EvaluationDDim2,
49 class EvaluationDDim3,
50 class LowerExtrapolationRule1,
51 class UpperExtrapolationRule1,
52 class LowerExtrapolationRule2,
53 class UpperExtrapolationRule2,
54 class LowerExtrapolationRule3,
55 class UpperExtrapolationRule3>
57{
58public:
59 /// @brief The type of the first evaluation continuous dimension used by this class.
60 using continuous_dimension_type1 = typename BSplines1::continuous_dimension_type;
61
62 /// @brief The type of the second evaluation continuous dimension used by this class.
63 using continuous_dimension_type2 = typename BSplines2::continuous_dimension_type;
64
65 /// @brief The type of the third evaluation continuous dimension used by this class.
66 using continuous_dimension_type3 = typename BSplines3::continuous_dimension_type;
67
68 /// @brief The type of the Kokkos execution space used by this class.
69 using exec_space = ExecSpace;
70
71 /// @brief The type of the Kokkos memory space used by this class.
72 using memory_space = MemorySpace;
73
74 /// @brief The type of the first discrete dimension of interest used by this class.
75 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
76
77 /// @brief The type of the second discrete dimension of interest used by this class.
78 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
79
80 /// @brief The type of the third discrete dimension of interest used by this class.
81 using evaluation_discrete_dimension_type3 = EvaluationDDim3;
82
83 /// @brief The discrete dimension representing the B-splines along first dimension.
84 using bsplines_type1 = BSplines1;
85
86 /// @brief The discrete dimension representing the B-splines along second dimension.
87 using bsplines_type2 = BSplines2;
88
89 /// @brief The discrete dimension representing the B-splines along third dimension.
90 using bsplines_type3 = BSplines3;
91
92 /// @brief The type of the domain for the 1D evaluation mesh along first dimension used by this class.
93 using evaluation_domain_type1 = ddc::DiscreteDomain<evaluation_discrete_dimension_type1>;
94
95 /// @brief The type of the domain for the 1D evaluation mesh along second dimension used by this class.
96 using evaluation_domain_type2 = ddc::DiscreteDomain<evaluation_discrete_dimension_type2>;
97
98 /// @brief The type of the domain for the 1D evaluation mesh along third dimension used by this class.
99 using evaluation_domain_type3 = ddc::DiscreteDomain<evaluation_discrete_dimension_type3>;
100
101 /// @brief The type of the domain for the 3D evaluation mesh used by this class.
102 using evaluation_domain_type = ddc::DiscreteDomain<
103 evaluation_discrete_dimension_type1,
104 evaluation_discrete_dimension_type2,
105 evaluation_discrete_dimension_type3>;
106
107 /**
108 * @brief The type of the whole domain representing evaluation points.
109 *
110 * @tparam The batched discrete domain on which the interpolation points are defined.
111 */
112 template <
113 class BatchedInterpolationDDom,
114 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
115 using batched_evaluation_domain_type = BatchedInterpolationDDom;
116
117 /// @brief The type of the 1D spline domain corresponding to the first dimension of interest.
118 using spline_domain_type1 = ddc::DiscreteDomain<bsplines_type1>;
119
120 /// @brief The type of the 1D spline domain corresponding to the second dimension of interest.
121 using spline_domain_type2 = ddc::DiscreteDomain<bsplines_type2>;
122
123 /// @brief The type of the 1D spline domain corresponding to the third dimension of interest.
124 using spline_domain_type3 = ddc::DiscreteDomain<bsplines_type3>;
125
126 /// @brief The type of the 3D spline domain corresponding to the dimensions of interest.
127 using spline_domain_type = ddc::DiscreteDomain<bsplines_type1, bsplines_type2, bsplines_type3>;
128
129 /**
130 * @brief The type of the batch domain (obtained by removing the dimensions of interest
131 * from the whole domain).
132 *
133 * @tparam The batched discrete domain on which the interpolation points are defined.
134 */
135 template <
136 class BatchedInterpolationDDom,
137 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
138 using batch_domain_type = typename ddc::remove_dims_of_t<
139 BatchedInterpolationDDom,
140 evaluation_discrete_dimension_type1,
141 evaluation_discrete_dimension_type2,
142 evaluation_discrete_dimension_type3>;
143
144 /**
145 * @brief The type of the whole spline domain (cartesian product of 3D spline domain
146 * and batch domain) preserving the underlying memory layout (order of dimensions).
147 *
148 * @tparam The batched discrete domain on which the interpolation points are defined.
149 */
150 template <
151 class BatchedInterpolationDDom,
152 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
153 using batched_spline_domain_type =
154 typename ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
155 ddc::to_type_seq_t<BatchedInterpolationDDom>,
156 ddc::detail::TypeSeq<
157 evaluation_discrete_dimension_type1,
158 evaluation_discrete_dimension_type2,
159 evaluation_discrete_dimension_type3>,
160 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>>>;
161
162 /// @brief The type of the extrapolation rule at the lower boundary along the first dimension.
163 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
164
165 /// @brief The type of the extrapolation rule at the upper boundary along the first dimension.
166 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
167
168 /// @brief The type of the extrapolation rule at the lower boundary along the second dimension.
169 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
170
171 /// @brief The type of the extrapolation rule at the upper boundary along the second dimension.
172 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
173
174 /// @brief The type of the extrapolation rule at the lower boundary along the third dimension.
175 using lower_extrapolation_rule_3_type = LowerExtrapolationRule3;
176
177 /// @brief The type of the extrapolation rule at the upper boundary along the third dimension.
178 using upper_extrapolation_rule_3_type = UpperExtrapolationRule3;
179
180private:
181 LowerExtrapolationRule1 m_lower_extrap_rule_1;
182
183 UpperExtrapolationRule1 m_upper_extrap_rule_1;
184
185 LowerExtrapolationRule2 m_lower_extrap_rule_2;
186
187 UpperExtrapolationRule2 m_upper_extrap_rule_2;
188
189 LowerExtrapolationRule3 m_lower_extrap_rule_3;
190
191 UpperExtrapolationRule3 m_upper_extrap_rule_3;
192
193public:
194 static_assert(
195 std::is_same_v<LowerExtrapolationRule1,
196 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
197 == bsplines_type1::is_periodic()
198 && std::is_same_v<
199 UpperExtrapolationRule1,
200 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
201 == bsplines_type1::is_periodic()
202 && std::is_same_v<
203 LowerExtrapolationRule2,
204 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
205 == bsplines_type2::is_periodic()
206 && std::is_same_v<
207 UpperExtrapolationRule2,
208 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
209 == bsplines_type2::is_periodic()
210 && std::is_same_v<
211 LowerExtrapolationRule3,
212 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
213 == bsplines_type3::is_periodic()
214 && std::is_same_v<
215 UpperExtrapolationRule3,
216 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
217 == bsplines_type3::is_periodic(),
218 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
219 static_assert(
220 std::is_invocable_r_v<
221 double,
222 LowerExtrapolationRule1,
223 ddc::Coordinate<continuous_dimension_type1>,
224 ddc::ChunkSpan<
225 double const,
226 spline_domain_type,
227 Kokkos::layout_right,
228 memory_space>>,
229 "LowerExtrapolationRule1::operator() has to be callable "
230 "with usual arguments.");
231 static_assert(
232 std::is_invocable_r_v<
233 double,
234 UpperExtrapolationRule1,
235 ddc::Coordinate<continuous_dimension_type1>,
236 ddc::ChunkSpan<
237 double const,
238 spline_domain_type,
239 Kokkos::layout_right,
240 memory_space>>,
241 "UpperExtrapolationRule1::operator() has to be callable "
242 "with usual arguments.");
243 static_assert(
244 std::is_invocable_r_v<
245 double,
246 LowerExtrapolationRule2,
247 ddc::Coordinate<continuous_dimension_type2>,
248 ddc::ChunkSpan<
249 double const,
250 spline_domain_type,
251 Kokkos::layout_right,
252 memory_space>>,
253 "LowerExtrapolationRule2::operator() has to be callable "
254 "with usual arguments.");
255 static_assert(
256 std::is_invocable_r_v<
257 double,
258 UpperExtrapolationRule2,
259 ddc::Coordinate<continuous_dimension_type2>,
260 ddc::ChunkSpan<
261 double const,
262 spline_domain_type,
263 Kokkos::layout_right,
264 memory_space>>,
265 "UpperExtrapolationRule2::operator() has to be callable "
266 "with usual arguments.");
267 static_assert(
268 std::is_invocable_r_v<
269 double,
270 LowerExtrapolationRule3,
271 ddc::Coordinate<continuous_dimension_type3>,
272 ddc::ChunkSpan<
273 double const,
274 spline_domain_type,
275 Kokkos::layout_right,
276 memory_space>>,
277 "LowerExtrapolationRule3::operator() has to be callable "
278 "with usual arguments.");
279 static_assert(
280 std::is_invocable_r_v<
281 double,
282 UpperExtrapolationRule3,
283 ddc::Coordinate<continuous_dimension_type3>,
284 ddc::ChunkSpan<
285 double const,
286 spline_domain_type,
287 Kokkos::layout_right,
288 memory_space>>,
289 "UpperExtrapolationRule3::operator() has to be callable "
290 "with usual arguments.");
291
292 /**
293 * @brief Build a SplineEvaluator3D acting on batched_spline_domain.
294 *
295 * @param lower_extrap_rule1 The extrapolation rule at the lower boundary along the first dimension.
296 * @param upper_extrap_rule1 The extrapolation rule at the upper boundary along the first dimension.
297 * @param lower_extrap_rule2 The extrapolation rule at the lower boundary along the second dimension.
298 * @param upper_extrap_rule2 The extrapolation rule at the upper boundary along the second dimension.
299 * @param lower_extrap_rule3 The extrapolation rule at the lower boundary along the third dimension.
300 * @param upper_extrap_rule3 The extrapolation rule at the upper boundary along the third dimension.
301 *
302 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
303 */
304 explicit SplineEvaluator3D(
305 LowerExtrapolationRule1 const& lower_extrap_rule1,
306 UpperExtrapolationRule1 const& upper_extrap_rule1,
307 LowerExtrapolationRule2 const& lower_extrap_rule2,
308 UpperExtrapolationRule2 const& upper_extrap_rule2,
309 LowerExtrapolationRule3 const& lower_extrap_rule3,
310 UpperExtrapolationRule3 const& upper_extrap_rule3)
311 : m_lower_extrap_rule_1(lower_extrap_rule1)
312 , m_upper_extrap_rule_1(upper_extrap_rule1)
313 , m_lower_extrap_rule_2(lower_extrap_rule2)
314 , m_upper_extrap_rule_2(upper_extrap_rule2)
315 , m_lower_extrap_rule_3(lower_extrap_rule3)
316 , m_upper_extrap_rule_3(upper_extrap_rule3)
317 {
318 }
319
320 /**
321 * @brief Copy-constructs.
322 *
323 * @param x A reference to another SplineEvaluator.
324 */
325 SplineEvaluator3D(SplineEvaluator3D const& x) = default;
326
327 /**
328 * @brief Move-constructs.
329 *
330 * @param x An rvalue to another SplineEvaluator.
331 */
333
334 /// @brief Destructs.
335 ~SplineEvaluator3D() = default;
336
337 /**
338 * @brief Copy-assigns.
339 *
340 * @param x A reference to another SplineEvaluator.
341 * @return A reference to this object.
342 */
343 SplineEvaluator3D& operator=(SplineEvaluator3D const& x) = default;
344
345 /**
346 * @brief Move-assigns.
347 *
348 * @param x An rvalue to another SplineEvaluator.
349 * @return A reference to this object.
350 */
352
353 /**
354 * @brief Get the lower extrapolation rule along the first dimension.
355 *
356 * 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.
357 *
358 * @return The lower extrapolation rule along the first dimension.
359 *
360 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
361 */
362 lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
363 {
364 return m_lower_extrap_rule_1;
365 }
366
367 /**
368 * @brief Get the upper extrapolation rule along the first dimension.
369 *
370 * 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.
371 *
372 * @return The upper extrapolation rule along the first dimension.
373 *
374 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
375 */
376 upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
377 {
378 return m_upper_extrap_rule_1;
379 }
380
381 /**
382 * @brief Get the lower extrapolation rule along the second dimension.
383 *
384 * 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.
385 *
386 * @return The lower extrapolation rule along the second dimension.
387 *
388 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
389 */
390 lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
391 {
392 return m_lower_extrap_rule_2;
393 }
394
395 /**
396 * @brief Get the upper extrapolation rule along the second dimension.
397 *
398 * 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.
399 *
400 * @return The upper extrapolation rule along the second dimension.
401 *
402 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
403 */
404 upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
405 {
406 return m_upper_extrap_rule_2;
407 }
408
409 /**
410 * @brief Get the lower extrapolation rule along the third dimension.
411 *
412 * 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.
413 *
414 * @return The lower extrapolation rule along the third dimension.
415 *
416 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
417 */
418 lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
419 {
420 return m_lower_extrap_rule_3;
421 }
422
423 /**
424 * @brief Get the upper extrapolation rule along the third dimension.
425 *
426 * 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.
427 *
428 * @return The upper extrapolation rule along the third dimension.
429 *
430 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
431 */
432 upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
433 {
434 return m_upper_extrap_rule_3;
435 }
436
437 /**
438 * @brief Evaluate 3D spline function (described by its spline coefficients) at a given coordinate.
439 *
440 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
441 *
442 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
443 *
444 * @param coord_eval The coordinate where the spline is evaluated. Note that only the components along the dimensions of interest are used.
445 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
446 *
447 * @return The value of the spline function at the desired coordinate.
448 */
449 template <class Layout, class... CoordsDims>
450 KOKKOS_FUNCTION double operator()(
451 ddc::Coordinate<CoordsDims...> const& coord_eval,
452 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
453 spline_coef) const
454 {
455 return eval(coord_eval, spline_coef);
456 }
457
458 /**
459 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
460 *
461 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
462 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
463 *
464 * This is not a nD evaluation. This is a batched 3D evaluation. This means that for each slice of coordinates
465 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 3D set of
466 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
467 *
468 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
469 *
470 * @param[out] spline_eval The values of the 3D spline function at the desired coordinates. For practical reasons those are
471 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
472 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
473 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
474 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
475 * the set of 3D spline coefficients retained to perform the evaluation).
476 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
477 */
478 template <
479 class Layout1,
480 class Layout2,
481 class Layout3,
482 class BatchedInterpolationDDom,
483 class... CoordsDims>
484 void operator()(
485 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
486 spline_eval,
487 ddc::ChunkSpan<
488 ddc::Coordinate<CoordsDims...> const,
489 BatchedInterpolationDDom,
490 Layout2,
491 memory_space> const coords_eval,
492 ddc::ChunkSpan<
493 double const,
494 batched_spline_domain_type<BatchedInterpolationDDom>,
495 Layout3,
496 memory_space> const spline_coef) const
497 {
498 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
499 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
500 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
501 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
502 ddc::parallel_for_each(
503 "ddc_splines_evaluate_3d",
504 exec_space(),
505 batch_domain,
506 KOKKOS_CLASS_LAMBDA(
507 typename batch_domain_type<
508 BatchedInterpolationDDom>::discrete_element_type const j) {
509 auto const spline_eval_3D = spline_eval[j];
510 auto const coords_eval_3D = coords_eval[j];
511 auto const spline_coef_3D = spline_coef[j];
512 for (auto const i1 : evaluation_domain1) {
513 for (auto const i2 : evaluation_domain2) {
514 for (auto const i3 : evaluation_domain3) {
515 spline_eval_3D(i1, i2, i3)
516 = eval(coords_eval_3D(i1, i2, i3), spline_coef_3D);
517 }
518 }
519 }
520 });
521 }
522
523 /**
524 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
525 *
526 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
527 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
528 *
529 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
530 * This means that for each slice of spline_eval the evaluation is performed with
531 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
532 *
533 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
534 *
535 * @param[out] spline_eval The values of the 3D spline function at their coordinates.
536 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
537 */
538 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
539 void operator()(
540 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
541 spline_eval,
542 ddc::ChunkSpan<
543 double const,
544 batched_spline_domain_type<BatchedInterpolationDDom>,
545 Layout2,
546 memory_space> const spline_coef) const
547 {
548 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
549 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
550 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
551 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
552 ddc::parallel_for_each(
553 "ddc_splines_evaluate_3d",
554 exec_space(),
555 batch_domain,
556 KOKKOS_CLASS_LAMBDA(
557 typename batch_domain_type<
558 BatchedInterpolationDDom>::discrete_element_type const j) {
559 auto const spline_eval_3D = spline_eval[j];
560 auto const spline_coef_3D = spline_coef[j];
561 for (auto const i1 : evaluation_domain1) {
562 for (auto const i2 : evaluation_domain2) {
563 for (auto const i3 : evaluation_domain3) {
564 ddc::Coordinate<
565 continuous_dimension_type1,
566 continuous_dimension_type2,
567 continuous_dimension_type3>
568 coord_eval_3D(
569 ddc::coordinate(i1),
570 ddc::coordinate(i2),
571 ddc::coordinate(i3));
572 spline_eval_3D(i1, i2, i3)
573 = eval(coord_eval_3D(i1, i2, i3), spline_coef_3D);
574 }
575 }
576 }
577 });
578 }
579
580 /**
581 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the dimensions of interest.
582 *
583 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
584 * obtained via various methods, such as using a SplineBuilder3D.
585 *
586 * @param deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
587 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
588 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
589 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
590 *
591 * @return The derivative of the spline function at the desired coordinate.
592 */
593 template <class DElem, class Layout, class... CoordsDims>
594 KOKKOS_FUNCTION double deriv(
595 DElem const& deriv_order,
596 ddc::Coordinate<CoordsDims...> const& coord_eval,
597 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
598 spline_coef) const
599 {
600 static_assert(is_discrete_element_v<DElem>);
601 return eval_no_bc(deriv_order, coord_eval, spline_coef);
602 }
603
604 /**
605 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
606 *
607 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
608 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
609 *
610 * This is not a nD differentiation. This is a batched 3D differentiation.
611 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
612 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
613 *
614 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
615 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
616 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
617 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
618 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
619 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
620 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
621 * the set of 3D spline coefficients retained to perform the evaluation).
622 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
623 */
624 template <
625 class DElem,
626 class Layout1,
627 class Layout2,
628 class Layout3,
629 class BatchedInterpolationDDom,
630 class... CoordsDims>
631 void deriv(
632 DElem const& deriv_order,
633 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
634 spline_eval,
635 ddc::ChunkSpan<
636 ddc::Coordinate<CoordsDims...> const,
637 BatchedInterpolationDDom,
638 Layout2,
639 memory_space> const coords_eval,
640 ddc::ChunkSpan<
641 double const,
642 batched_spline_domain_type<BatchedInterpolationDDom>,
643 Layout3,
644 memory_space> const spline_coef) const
645 {
646 static_assert(is_discrete_element_v<DElem>);
647
648 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
649 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
650 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
651 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
652 ddc::parallel_for_each(
653 "ddc_splines_cross_differentiate_3d",
654 exec_space(),
655 batch_domain,
656 KOKKOS_CLASS_LAMBDA(
657 typename batch_domain_type<
658 BatchedInterpolationDDom>::discrete_element_type const j) {
659 auto const spline_eval_3D = spline_eval[j];
660 auto const coords_eval_3D = coords_eval[j];
661 auto const spline_coef_3D = spline_coef[j];
662 for (auto const i1 : evaluation_domain1) {
663 for (auto const i2 : evaluation_domain2) {
664 for (auto const i3 : evaluation_domain3) {
665 spline_eval_3D(i1, i2, i3) = eval_no_bc(
666 deriv_order,
667 coords_eval_3D(i1, i2, i3),
668 spline_coef_3D);
669 }
670 }
671 }
672 });
673 }
674
675 /**
676 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
677 *
678 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
679 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
680 *
681 * This is not a multidimensional differentiation. This is a batched 3D differentiation.
682 * This means that for each slice of spline_eval the differentiation is performed with
683 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
684 *
685 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
686 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
687 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
688 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
689 */
690 template <class DElem, class Layout1, class Layout2, class BatchedInterpolationDDom>
691 void deriv(
692 DElem const& deriv_order,
693 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
694 spline_eval,
695 ddc::ChunkSpan<
696 double const,
697 batched_spline_domain_type<BatchedInterpolationDDom>,
698 Layout2,
699 memory_space> const spline_coef) const
700 {
701 static_assert(is_discrete_element_v<DElem>);
702
703 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
704 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
705 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
706 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
707 ddc::parallel_for_each(
708 "ddc_splines_cross_differentiate_3d",
709 exec_space(),
710 batch_domain,
711 KOKKOS_CLASS_LAMBDA(
712 typename batch_domain_type<
713 BatchedInterpolationDDom>::discrete_element_type const j) {
714 auto const spline_eval_3D = spline_eval[j];
715 auto const spline_coef_3D = spline_coef[j];
716 for (auto const i1 : evaluation_domain1) {
717 for (auto const i2 : evaluation_domain2) {
718 for (auto const i3 : evaluation_domain3) {
719 ddc::Coordinate<
720 continuous_dimension_type1,
721 continuous_dimension_type2,
722 continuous_dimension_type3>
723 coord_eval_3D(
724 ddc::coordinate(i1),
725 ddc::coordinate(i2),
726 ddc::coordinate(i3));
727 spline_eval_3D(i1, i2, i3)
728 = eval_no_bc(deriv_order, coord_eval_3D, spline_coef_3D);
729 }
730 }
731 }
732 });
733 }
734
735 /** @brief Perform batched 3D integrations of a spline function (described by its spline coefficients) along the dimensions of interest and store results on a subdomain of batch_domain.
736 *
737 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
738 *
739 * This is not a nD integration. This is a batched 3D integration.
740 * This means that for each element of integrals, the integration is performed with the 3D set of
741 * spline coefficients identified by the same DiscreteElement.
742 *
743 * @param[out] integrals The integrals of the 3D spline function on the subdomain of batch_domain. For practical reasons those are
744 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
745 * points represented by this domain are unused and irrelevant.
746 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
747 */
748 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
749 void integrate(
750 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
751 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
752 spline_coef) const
753 {
754 static_assert(
755 ddc::type_seq_contains_v<
756 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>,
757 to_type_seq_t<BatchedSplineDDom>>,
758 "The spline coefficients domain must contain the bsplines dimensions");
759 using batch_domain_type = ddc::
760 remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2, bsplines_type3>;
761 static_assert(
762 std::is_same_v<batch_domain_type, BatchedDDom>,
763 "The integrals domain must only contain the batch dimensions");
764
765 batch_domain_type batch_domain(integrals.domain());
766 ddc::Chunk values1_alloc(
767 ddc::DiscreteDomain<bsplines_type1>(spline_coef.domain()),
768 ddc::KokkosAllocator<double, memory_space>());
769 ddc::ChunkSpan values1 = values1_alloc.span_view();
770 ddc::integrals(exec_space(), values1);
771 ddc::Chunk values2_alloc(
772 ddc::DiscreteDomain<bsplines_type2>(spline_coef.domain()),
773 ddc::KokkosAllocator<double, memory_space>());
774 ddc::ChunkSpan values2 = values2_alloc.span_view();
775 ddc::integrals(exec_space(), values2);
776 ddc::Chunk values3_alloc(
777 ddc::DiscreteDomain<bsplines_type3>(spline_coef.domain()),
778 ddc::KokkosAllocator<double, memory_space>());
779 ddc::ChunkSpan values3 = values3_alloc.span_view();
780 ddc::integrals(exec_space(), values3);
781
782 ddc::parallel_for_each(
783 "ddc_splines_integrate_bsplines",
784 exec_space(),
785 batch_domain,
786 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
787 integrals(j) = 0;
788 for (typename spline_domain_type1::discrete_element_type const i1 :
789 values1.domain()) {
790 for (typename spline_domain_type2::discrete_element_type const i2 :
791 values2.domain()) {
792 for (typename spline_domain_type3::discrete_element_type const i3 :
793 values3.domain()) {
794 integrals(j) += spline_coef(i1, i2, i3, j) * values1(i1)
795 * values2(i2) * values3(i3);
796 }
797 }
798 }
799 });
800 }
801
802private:
803 /**
804 * @brief Evaluate the function on B-splines at the coordinate given.
805 *
806 * This function firstly deals with the boundary conditions and calls the SplineEvaluator3D::eval_no_bc function
807 * to evaluate.
808 *
809 * @param[in] coord_eval The 3D coordinate where we want to evaluate.
810 * @param[in] spline_coef The B-splines coefficients of the function we want to evaluate.
811 * @param[out] vals1 A ChunkSpan with the not-null values of each function of the spline in the first dimension.
812 * @param[out] vals2 A ChunkSpan with the not-null values of each function of the spline in the second dimension.
813 *
814 * @return A double with the value of the function at the coordinate given.
815 *
816 * @see SplineBoundaryValue
817 */
818 template <class Layout, class... CoordsDims>
819 KOKKOS_INLINE_FUNCTION double eval(
820 ddc::Coordinate<CoordsDims...> coord_eval,
821 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
822 spline_coef) const
823 {
824 using Dim1 = continuous_dimension_type1;
825 using Dim2 = continuous_dimension_type2;
826 using Dim3 = continuous_dimension_type3;
827 if constexpr (bsplines_type1::is_periodic()) {
828 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()
829 || ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
830 ddc::get<Dim1>(coord_eval)
831 -= Kokkos::floor(
832 (ddc::get<Dim1>(coord_eval)
833 - ddc::discrete_space<bsplines_type1>().rmin())
834 / ddc::discrete_space<bsplines_type1>().length())
835 * ddc::discrete_space<bsplines_type1>().length();
836 }
837 }
838 if constexpr (bsplines_type2::is_periodic()) {
839 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()
840 || ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
841 ddc::get<Dim2>(coord_eval)
842 -= Kokkos::floor(
843 (ddc::get<Dim2>(coord_eval)
844 - ddc::discrete_space<bsplines_type2>().rmin())
845 / ddc::discrete_space<bsplines_type2>().length())
846 * ddc::discrete_space<bsplines_type2>().length();
847 }
848 }
849 if constexpr (bsplines_type3::is_periodic()) {
850 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()
851 || ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
852 ddc::get<Dim3>(coord_eval)
853 -= Kokkos::floor(
854 (ddc::get<Dim3>(coord_eval)
855 - ddc::discrete_space<bsplines_type3>().rmin())
856 / ddc::discrete_space<bsplines_type3>().length())
857 * ddc::discrete_space<bsplines_type3>().length();
858 }
859 }
860 if constexpr (!bsplines_type1::is_periodic()) {
861 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()) {
862 return m_lower_extrap_rule_1(coord_eval, spline_coef);
863 }
864 if (ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
865 return m_upper_extrap_rule_1(coord_eval, spline_coef);
866 }
867 }
868 if constexpr (!bsplines_type2::is_periodic()) {
869 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()) {
870 return m_lower_extrap_rule_2(coord_eval, spline_coef);
871 }
872 if (ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
873 return m_upper_extrap_rule_2(coord_eval, spline_coef);
874 }
875 }
876 if constexpr (!bsplines_type3::is_periodic()) {
877 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()) {
878 return m_lower_extrap_rule_3(coord_eval, spline_coef);
879 }
880 if (ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
881 return m_upper_extrap_rule_3(coord_eval, spline_coef);
882 }
883 }
884 return eval_no_bc(
885 ddc::DiscreteElement<>(),
886 ddc::Coordinate<
887 continuous_dimension_type1,
888 continuous_dimension_type2,
889 continuous_dimension_type3>(
890 ddc::get<Dim1>(coord_eval),
891 ddc::get<Dim2>(coord_eval),
892 ddc::get<Dim3>(coord_eval)),
893 spline_coef);
894 }
895
896 /**
897 * @brief Evaluate the function or its derivative at the coordinate given.
898 *
899 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
900 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
901 * @param[in] coord_eval The coordinate where we want to evaluate.
902 * @param[in] splne_coef The B-splines coefficients of the function we want to evaluate.
903 */
904 template <class... DerivDims, class Layout, class... CoordsDims>
905 KOKKOS_INLINE_FUNCTION double eval_no_bc(
906 ddc::DiscreteElement<DerivDims...> const& deriv_order,
907 ddc::Coordinate<CoordsDims...> const& coord_eval,
908 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
909 spline_coef) const
910 {
911 using deriv_dim1 = Deriv<continuous_dimension_type1>;
912 using deriv_dim2 = Deriv<continuous_dimension_type2>;
913 using deriv_dim3 = Deriv<continuous_dimension_type3>;
914 using deriv_dims = detail::TypeSeq<DerivDims...>;
915
916 // Check that the tags are valid
917 static_assert(
918 (in_tags_v<DerivDims, ddc::detail::TypeSeq<deriv_dim1, deriv_dim2, deriv_dim3>>
919 && ...),
920 "The only valid dimensions for deriv_order are Deriv<Dim1>, Deriv<Dim2> and "
921 "Deriv<Dim3>");
922
923 ddc::DiscreteElement<bsplines_type1> jmin1;
924 ddc::DiscreteElement<bsplines_type2> jmin2;
925 ddc::DiscreteElement<bsplines_type3> jmin3;
926
927 std::array<double, bsplines_type1::degree() + 1> vals1_ptr;
928 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>> const
929 vals1(vals1_ptr.data());
930 std::array<double, bsplines_type2::degree() + 1> vals2_ptr;
931 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>> const
932 vals2(vals2_ptr.data());
933 std::array<double, bsplines_type3::degree() + 1> vals3_ptr;
934 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type3::degree() + 1>> const
935 vals3(vals3_ptr.data());
936 ddc::Coordinate<continuous_dimension_type1> const coord_eval_interest1(coord_eval);
937 ddc::Coordinate<continuous_dimension_type2> const coord_eval_interest2(coord_eval);
938 ddc::Coordinate<continuous_dimension_type3> const coord_eval_interest3(coord_eval);
939
940 if constexpr (!in_tags_v<deriv_dim1, deriv_dims>) {
941 jmin1 = ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
942 } else {
943 auto const order1 = deriv_order.template uid<deriv_dim1>();
944 KOKKOS_ASSERT(order1 > 0 && order1 <= bsplines_type1::degree())
945
946 std::array<double, (bsplines_type1::degree() + 1) * (bsplines_type1::degree() + 1)>
947 derivs1_ptr;
948 Kokkos::mdspan<
949 double,
950 Kokkos::extents<
951 std::size_t,
952 bsplines_type1::degree() + 1,
953 Kokkos::dynamic_extent>> const derivs1(derivs1_ptr.data(), order1 + 1);
954
955 jmin1 = ddc::discrete_space<bsplines_type1>()
956 .eval_basis_and_n_derivs(derivs1, coord_eval_interest1, order1);
957
958 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
959 vals1[i] = DDC_MDSPAN_ACCESS_OP(derivs1, i, order1);
960 }
961 }
962
963 if constexpr (!in_tags_v<deriv_dim2, deriv_dims>) {
964 jmin2 = ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
965 } else {
966 auto const order2 = deriv_order.template uid<deriv_dim2>();
967 KOKKOS_ASSERT(order2 > 0 && order2 <= bsplines_type2::degree())
968
969 std::array<double, (bsplines_type2::degree() + 1) * (bsplines_type2::degree() + 1)>
970 derivs2_ptr;
971 Kokkos::mdspan<
972 double,
973 Kokkos::extents<
974 std::size_t,
975 bsplines_type2::degree() + 1,
976 Kokkos::dynamic_extent>> const derivs2(derivs2_ptr.data(), order2 + 1);
977
978 jmin2 = ddc::discrete_space<bsplines_type2>()
979 .eval_basis_and_n_derivs(derivs2, coord_eval_interest2, order2);
980
981 for (std::size_t i = 0; i < bsplines_type2::degree() + 1; ++i) {
982 vals2[i] = DDC_MDSPAN_ACCESS_OP(derivs2, i, order2);
983 }
984 }
985
986 if constexpr (!in_tags_v<deriv_dim3, deriv_dims>) {
987 jmin3 = ddc::discrete_space<bsplines_type3>().eval_basis(vals3, coord_eval_interest3);
988 } else {
989 auto const order3 = deriv_order.template uid<deriv_dim3>();
990 KOKKOS_ASSERT(order3 > 0 && order3 <= bsplines_type3::degree())
991
992 std::array<double, (bsplines_type3::degree() + 1) * (bsplines_type3::degree() + 1)>
993 derivs3_ptr;
994 Kokkos::mdspan<
995 double,
996 Kokkos::extents<
997 std::size_t,
998 bsplines_type3::degree() + 1,
999 Kokkos::dynamic_extent>> const derivs3(derivs3_ptr.data(), order3 + 1);
1000
1001 jmin3 = ddc::discrete_space<bsplines_type3>()
1002 .eval_basis_and_n_derivs(derivs3, coord_eval_interest3, order3);
1003
1004 for (std::size_t i = 0; i < bsplines_type3::degree() + 1; ++i) {
1005 vals3[i] = DDC_MDSPAN_ACCESS_OP(derivs3, i, order3);
1006 }
1007 }
1008
1009 double y = 0.0;
1010 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
1011 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
1012 for (std::size_t k = 0; k < bsplines_type3::degree() + 1; ++k) {
1013 y += spline_coef(
1014 ddc::DiscreteElement<
1015 bsplines_type1,
1016 bsplines_type2,
1017 bsplines_type3>(jmin1 + i, jmin2 + j, jmin3 + k))
1018 * vals1[i] * vals2[j] * vals3[k];
1019 }
1020 }
1021 }
1022 return y;
1023 }
1024};
1025
1026} // 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 3D spline function.
upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
Get the upper extrapolation rule along the first dimension.
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 3D spline function (described by its spline coefficients) on a mesh along the dimension...
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 3D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator3D & operator=(SplineEvaluator3D &&x)=default
Move-assigns.
upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
Get the upper extrapolation rule along the third dimension.
SplineEvaluator3D(LowerExtrapolationRule1 const &lower_extrap_rule1, UpperExtrapolationRule1 const &upper_extrap_rule1, LowerExtrapolationRule2 const &lower_extrap_rule2, UpperExtrapolationRule2 const &upper_extrap_rule2, LowerExtrapolationRule3 const &lower_extrap_rule3, UpperExtrapolationRule3 const &upper_extrap_rule3)
Build a SplineEvaluator3D acting on batched_spline_domain.
~SplineEvaluator3D()=default
Destructs.
lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
Get the lower extrapolation rule along the first dimension.
upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
Get the upper extrapolation rule along the second dimension.
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 3D spline function (described by its spline coefficients) on a mesh along the dimension...
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 3D spline function (described by its spline coefficients) at a given coordinate along t...
SplineEvaluator3D(SplineEvaluator3D const &x)=default
Copy-constructs.
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 3D spline function (described by its spline coefficients) at a given coordinate.
lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
Get the lower extrapolation rule along the third dimension.
SplineEvaluator3D(SplineEvaluator3D &&x)=default
Move-constructs.
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 3D integrations of a spline function (described by its spline coefficients) along the...
lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
Get the lower extrapolation rule along the second dimension.
SplineEvaluator3D & operator=(SplineEvaluator3D const &x)=default
Copy-assigns.
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 3D spline function (described by its spline coefficients) on a mesh.
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