DDC 0.9.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 "integrals.hpp"
16#include "periodic_extrapolation_rule.hpp"
17
18namespace ddc {
19
20/**
21 * @brief A class to evaluate, differentiate or integrate a 3D spline function.
22 *
23 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a 3D spline function.
24 *
25 * @tparam ExecSpace The Kokkos execution space on which the spline evaluation is performed.
26 * @tparam MemorySpace The Kokkos memory space on which the data (spline coefficients and evaluation) is stored.
27 * @tparam BSplines1 The discrete dimension representing the B-splines along the first dimension of interest.
28 * @tparam BSplines2 The discrete dimension representing the B-splines along the second dimension of interest.
29 * @tparam BSplines3 The discrete dimension representing the B-splines along the third dimension of interest.
30 * @tparam EvaluationDDim1 The first discrete dimension on which evaluation points are defined.
31 * @tparam EvaluationDDim2 The second discrete dimension on which evaluation points are defined.
32 * @tparam EvaluationDDim3 The third discrete dimension on which evaluation points are defined.
33 * @tparam LowerExtrapolationRule1 The lower extrapolation rule type along first dimension of interest.
34 * @tparam UpperExtrapolationRule1 The upper extrapolation rule type along first dimension of interest.
35 * @tparam LowerExtrapolationRule2 The lower extrapolation rule type along second dimension of interest.
36 * @tparam UpperExtrapolationRule2 The upper extrapolation rule type along second dimension of interest.
37 * @tparam LowerExtrapolationRule3 The lower extrapolation rule type along third dimension of interest.
38 * @tparam UpperExtrapolationRule3 The upper extrapolation rule type along third dimension of interest.
39 */
40template <
41 class ExecSpace,
42 class MemorySpace,
43 class BSplines1,
44 class BSplines2,
45 class BSplines3,
46 class EvaluationDDim1,
47 class EvaluationDDim2,
48 class EvaluationDDim3,
49 class LowerExtrapolationRule1,
50 class UpperExtrapolationRule1,
51 class LowerExtrapolationRule2,
52 class UpperExtrapolationRule2,
53 class LowerExtrapolationRule3,
54 class UpperExtrapolationRule3>
56{
57private:
58 /**
59 * @brief Tag to indicate that the value of the spline should be evaluated.
60 */
61 struct eval_type
62 {
63 };
64
65 /**
66 * @brief Tag to indicate that derivative of the spline should be evaluated.
67 */
68 struct eval_deriv_type
69 {
70 };
71
72public:
73 /// @brief The type of the first evaluation continuous dimension used by this class.
74 using continuous_dimension_type1 = typename BSplines1::continuous_dimension_type;
75
76 /// @brief The type of the second evaluation continuous dimension used by this class.
77 using continuous_dimension_type2 = typename BSplines2::continuous_dimension_type;
78
79 /// @brief The type of the third evaluation continuous dimension used by this class.
80 using continuous_dimension_type3 = typename BSplines3::continuous_dimension_type;
81
82 /// @brief The type of the Kokkos execution space used by this class.
83 using exec_space = ExecSpace;
84
85 /// @brief The type of the Kokkos memory space used by this class.
86 using memory_space = MemorySpace;
87
88 /// @brief The type of the first discrete dimension of interest used by this class.
89 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
90
91 /// @brief The type of the second discrete dimension of interest used by this class.
92 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
93
94 /// @brief The type of the third discrete dimension of interest used by this class.
95 using evaluation_discrete_dimension_type3 = EvaluationDDim3;
96
97 /// @brief The discrete dimension representing the B-splines along first dimension.
98 using bsplines_type1 = BSplines1;
99
100 /// @brief The discrete dimension representing the B-splines along second dimension.
101 using bsplines_type2 = BSplines2;
102
103 /// @brief The discrete dimension representing the B-splines along third dimension.
104 using bsplines_type3 = BSplines3;
105
106 /// @brief The type of the domain for the 1D evaluation mesh along first dimension used by this class.
107 using evaluation_domain_type1 = ddc::DiscreteDomain<evaluation_discrete_dimension_type1>;
108
109 /// @brief The type of the domain for the 1D evaluation mesh along second dimension used by this class.
110 using evaluation_domain_type2 = ddc::DiscreteDomain<evaluation_discrete_dimension_type2>;
111
112 /// @brief The type of the domain for the 1D evaluation mesh along third dimension used by this class.
113 using evaluation_domain_type3 = ddc::DiscreteDomain<evaluation_discrete_dimension_type3>;
114
115 /// @brief The type of the domain for the 3D evaluation mesh used by this class.
116 using evaluation_domain_type = ddc::DiscreteDomain<
117 evaluation_discrete_dimension_type1,
118 evaluation_discrete_dimension_type2,
119 evaluation_discrete_dimension_type3>;
120
121 /**
122 * @brief The type of the whole domain representing evaluation points.
123 *
124 * @tparam The batched discrete domain on which the interpolation points are defined.
125 */
126 template <
127 class BatchedInterpolationDDom,
128 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
129 using batched_evaluation_domain_type = BatchedInterpolationDDom;
130
131 /// @brief The type of the 1D spline domain corresponding to the first dimension of interest.
132 using spline_domain_type1 = ddc::DiscreteDomain<bsplines_type1>;
133
134 /// @brief The type of the 1D spline domain corresponding to the second dimension of interest.
135 using spline_domain_type2 = ddc::DiscreteDomain<bsplines_type2>;
136
137 /// @brief The type of the 1D spline domain corresponding to the third dimension of interest.
138 using spline_domain_type3 = ddc::DiscreteDomain<bsplines_type3>;
139
140 /// @brief The type of the 3D spline domain corresponding to the dimensions of interest.
141 using spline_domain_type = ddc::DiscreteDomain<bsplines_type1, bsplines_type2, bsplines_type3>;
142
143 /**
144 * @brief The type of the batch domain (obtained by removing the dimensions of interest
145 * from the whole domain).
146 *
147 * @tparam The batched discrete domain on which the interpolation points are defined.
148 */
149 template <
150 class BatchedInterpolationDDom,
151 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
152 using batch_domain_type = typename ddc::remove_dims_of_t<
153 BatchedInterpolationDDom,
154 evaluation_discrete_dimension_type1,
155 evaluation_discrete_dimension_type2,
156 evaluation_discrete_dimension_type3>;
157
158 /**
159 * @brief The type of the whole spline domain (cartesian product of 3D spline domain
160 * and batch domain) preserving the underlying memory layout (order of dimensions).
161 *
162 * @tparam The batched discrete domain on which the interpolation points are defined.
163 */
164 template <
165 class BatchedInterpolationDDom,
166 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
167 using batched_spline_domain_type =
168 typename ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
169 ddc::to_type_seq_t<BatchedInterpolationDDom>,
170 ddc::detail::TypeSeq<
171 evaluation_discrete_dimension_type1,
172 evaluation_discrete_dimension_type2,
173 evaluation_discrete_dimension_type3>,
174 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>>>;
175
176 /// @brief The type of the extrapolation rule at the lower boundary along the first dimension.
177 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
178
179 /// @brief The type of the extrapolation rule at the upper boundary along the first dimension.
180 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
181
182 /// @brief The type of the extrapolation rule at the lower boundary along the second dimension.
183 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
184
185 /// @brief The type of the extrapolation rule at the upper boundary along the second dimension.
186 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
187
188 /// @brief The type of the extrapolation rule at the lower boundary along the third dimension.
189 using lower_extrapolation_rule_3_type = LowerExtrapolationRule3;
190
191 /// @brief The type of the extrapolation rule at the upper boundary along the third dimension.
192 using upper_extrapolation_rule_3_type = UpperExtrapolationRule3;
193
194private:
195 LowerExtrapolationRule1 m_lower_extrap_rule_1;
196
197 UpperExtrapolationRule1 m_upper_extrap_rule_1;
198
199 LowerExtrapolationRule2 m_lower_extrap_rule_2;
200
201 UpperExtrapolationRule2 m_upper_extrap_rule_2;
202
203 LowerExtrapolationRule3 m_lower_extrap_rule_3;
204
205 UpperExtrapolationRule3 m_upper_extrap_rule_3;
206
207public:
208 static_assert(
209 std::is_same_v<LowerExtrapolationRule1,
210 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
211 == bsplines_type1::is_periodic()
212 && std::is_same_v<
213 UpperExtrapolationRule1,
214 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
215 == bsplines_type1::is_periodic()
216 && std::is_same_v<
217 LowerExtrapolationRule2,
218 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
219 == bsplines_type2::is_periodic()
220 && std::is_same_v<
221 UpperExtrapolationRule2,
222 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
223 == bsplines_type2::is_periodic()
224 && std::is_same_v<
225 LowerExtrapolationRule3,
226 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
227 == bsplines_type3::is_periodic()
228 && std::is_same_v<
229 UpperExtrapolationRule3,
230 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
231 == bsplines_type3::is_periodic(),
232 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
233 static_assert(
234 std::is_invocable_r_v<
235 double,
236 LowerExtrapolationRule1,
237 ddc::Coordinate<continuous_dimension_type1>,
238 ddc::ChunkSpan<
239 double const,
240 spline_domain_type,
241 Kokkos::layout_right,
242 memory_space>>,
243 "LowerExtrapolationRule1::operator() has to be callable "
244 "with usual arguments.");
245 static_assert(
246 std::is_invocable_r_v<
247 double,
248 UpperExtrapolationRule1,
249 ddc::Coordinate<continuous_dimension_type1>,
250 ddc::ChunkSpan<
251 double const,
252 spline_domain_type,
253 Kokkos::layout_right,
254 memory_space>>,
255 "UpperExtrapolationRule1::operator() has to be callable "
256 "with usual arguments.");
257 static_assert(
258 std::is_invocable_r_v<
259 double,
260 LowerExtrapolationRule2,
261 ddc::Coordinate<continuous_dimension_type2>,
262 ddc::ChunkSpan<
263 double const,
264 spline_domain_type,
265 Kokkos::layout_right,
266 memory_space>>,
267 "LowerExtrapolationRule2::operator() has to be callable "
268 "with usual arguments.");
269 static_assert(
270 std::is_invocable_r_v<
271 double,
272 UpperExtrapolationRule2,
273 ddc::Coordinate<continuous_dimension_type2>,
274 ddc::ChunkSpan<
275 double const,
276 spline_domain_type,
277 Kokkos::layout_right,
278 memory_space>>,
279 "UpperExtrapolationRule2::operator() has to be callable "
280 "with usual arguments.");
281 static_assert(
282 std::is_invocable_r_v<
283 double,
284 LowerExtrapolationRule3,
285 ddc::Coordinate<continuous_dimension_type3>,
286 ddc::ChunkSpan<
287 double const,
288 spline_domain_type,
289 Kokkos::layout_right,
290 memory_space>>,
291 "LowerExtrapolationRule3::operator() has to be callable "
292 "with usual arguments.");
293 static_assert(
294 std::is_invocable_r_v<
295 double,
296 UpperExtrapolationRule3,
297 ddc::Coordinate<continuous_dimension_type3>,
298 ddc::ChunkSpan<
299 double const,
300 spline_domain_type,
301 Kokkos::layout_right,
302 memory_space>>,
303 "UpperExtrapolationRule3::operator() has to be callable "
304 "with usual arguments.");
305
306 /**
307 * @brief Build a SplineEvaluator3D acting on batched_spline_domain.
308 *
309 * @param lower_extrap_rule1 The extrapolation rule at the lower boundary along the first dimension.
310 * @param upper_extrap_rule1 The extrapolation rule at the upper boundary along the first dimension.
311 * @param lower_extrap_rule2 The extrapolation rule at the lower boundary along the second dimension.
312 * @param upper_extrap_rule2 The extrapolation rule at the upper boundary along the second dimension.
313 * @param lower_extrap_rule3 The extrapolation rule at the lower boundary along the third dimension.
314 * @param upper_extrap_rule3 The extrapolation rule at the upper boundary along the third dimension.
315 *
316 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
317 */
318 explicit SplineEvaluator3D(
319 LowerExtrapolationRule1 const& lower_extrap_rule1,
320 UpperExtrapolationRule1 const& upper_extrap_rule1,
321 LowerExtrapolationRule2 const& lower_extrap_rule2,
322 UpperExtrapolationRule2 const& upper_extrap_rule2,
323 LowerExtrapolationRule3 const& lower_extrap_rule3,
324 UpperExtrapolationRule3 const& upper_extrap_rule3)
325 : m_lower_extrap_rule_1(lower_extrap_rule1)
326 , m_upper_extrap_rule_1(upper_extrap_rule1)
327 , m_lower_extrap_rule_2(lower_extrap_rule2)
328 , m_upper_extrap_rule_2(upper_extrap_rule2)
329 , m_lower_extrap_rule_3(lower_extrap_rule3)
330 , m_upper_extrap_rule_3(upper_extrap_rule3)
331 {
332 }
333
334 /**
335 * @brief Copy-constructs.
336 *
337 * @param x A reference to another SplineEvaluator.
338 */
339 SplineEvaluator3D(SplineEvaluator3D const& x) = default;
340
341 /**
342 * @brief Move-constructs.
343 *
344 * @param x An rvalue to another SplineEvaluator.
345 */
347
348 /// @brief Destructs.
349 ~SplineEvaluator3D() = default;
350
351 /**
352 * @brief Copy-assigns.
353 *
354 * @param x A reference to another SplineEvaluator.
355 * @return A reference to this object.
356 */
357 SplineEvaluator3D& operator=(SplineEvaluator3D const& x) = default;
358
359 /**
360 * @brief Move-assigns.
361 *
362 * @param x An rvalue to another SplineEvaluator.
363 * @return A reference to this object.
364 */
366
367 /**
368 * @brief Get the lower 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 lower extrapolation rule along the first dimension.
373 *
374 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
375 */
376 lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
377 {
378 return m_lower_extrap_rule_1;
379 }
380
381 /**
382 * @brief Get the upper extrapolation rule along the first 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 upper extrapolation rule along the first dimension.
387 *
388 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
389 */
390 upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
391 {
392 return m_upper_extrap_rule_1;
393 }
394
395 /**
396 * @brief Get the lower 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 lower extrapolation rule along the second dimension.
401 *
402 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
403 */
404 lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
405 {
406 return m_lower_extrap_rule_2;
407 }
408
409 /**
410 * @brief Get the upper extrapolation rule along the second 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 upper extrapolation rule along the second dimension.
415 *
416 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
417 */
418 upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
419 {
420 return m_upper_extrap_rule_2;
421 }
422
423 /**
424 * @brief Get the lower 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 lower extrapolation rule along the third dimension.
429 *
430 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
431 */
432 lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
433 {
434 return m_lower_extrap_rule_3;
435 }
436
437 /**
438 * @brief Get the upper extrapolation rule along the third dimension.
439 *
440 * 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.
441 *
442 * @return The upper extrapolation rule along the third dimension.
443 *
444 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
445 */
446 upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
447 {
448 return m_upper_extrap_rule_3;
449 }
450
451 /**
452 * @brief Evaluate 3D spline function (described by its spline coefficients) at a given coordinate.
453 *
454 * 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.
455 *
456 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
457 *
458 * @param coord_eval The coordinate where the spline is evaluated. Note that only the components along the dimensions of interest are used.
459 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
460 *
461 * @return The value of the spline function at the desired coordinate.
462 */
463 template <class Layout, class... CoordsDims>
464 KOKKOS_FUNCTION double operator()(
465 ddc::Coordinate<CoordsDims...> const& coord_eval,
466 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
467 spline_coef) const
468 {
469 return eval(coord_eval, spline_coef);
470 }
471
472 /**
473 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
474 *
475 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
476 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
477 *
478 * This is not a nD evaluation. This is a batched 3D evaluation. This means that for each slice of coordinates
479 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 3D set of
480 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
481 *
482 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
483 *
484 * @param[out] spline_eval The values of the 3D spline function at the desired coordinates. For practical reasons those are
485 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
486 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
487 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
488 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
489 * the set of 3D spline coefficients retained to perform the evaluation).
490 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
491 */
492 template <
493 class Layout1,
494 class Layout2,
495 class Layout3,
496 class BatchedInterpolationDDom,
497 class... CoordsDims>
498 void operator()(
499 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
500 spline_eval,
501 ddc::ChunkSpan<
502 ddc::Coordinate<CoordsDims...> const,
503 BatchedInterpolationDDom,
504 Layout2,
505 memory_space> const coords_eval,
506 ddc::ChunkSpan<
507 double const,
508 batched_spline_domain_type<BatchedInterpolationDDom>,
509 Layout3,
510 memory_space> const spline_coef) const
511 {
512 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
513 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
514 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
515 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
516 ddc::parallel_for_each(
517 "ddc_splines_evaluate_3d",
518 exec_space(),
519 batch_domain,
520 KOKKOS_CLASS_LAMBDA(
521 typename batch_domain_type<
522 BatchedInterpolationDDom>::discrete_element_type const j) {
523 auto const spline_eval_3D = spline_eval[j];
524 auto const coords_eval_3D = coords_eval[j];
525 auto const spline_coef_3D = spline_coef[j];
526 for (auto const i1 : evaluation_domain1) {
527 for (auto const i2 : evaluation_domain2) {
528 for (auto const i3 : evaluation_domain3) {
529 spline_eval_3D(i1, i2, i3)
530 = eval(coords_eval_3D(i1, i2, i3), spline_coef_3D);
531 }
532 }
533 }
534 });
535 }
536
537 /**
538 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
539 *
540 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
541 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
542 *
543 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
544 * This means that for each slice of spline_eval the evaluation is performed with
545 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
546 *
547 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
548 *
549 * @param[out] spline_eval The values of the 3D spline function at their coordinates.
550 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
551 */
552 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
553 void operator()(
554 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
555 spline_eval,
556 ddc::ChunkSpan<
557 double const,
558 batched_spline_domain_type<BatchedInterpolationDDom>,
559 Layout2,
560 memory_space> const spline_coef) const
561 {
562 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
563 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
564 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
565 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
566 ddc::parallel_for_each(
567 "ddc_splines_evaluate_3d",
568 exec_space(),
569 batch_domain,
570 KOKKOS_CLASS_LAMBDA(
571 typename batch_domain_type<
572 BatchedInterpolationDDom>::discrete_element_type const j) {
573 auto const spline_eval_3D = spline_eval[j];
574 auto const spline_coef_3D = spline_coef[j];
575 for (auto const i1 : evaluation_domain1) {
576 for (auto const i2 : evaluation_domain2) {
577 for (auto const i3 : evaluation_domain3) {
578 ddc::Coordinate<
579 continuous_dimension_type1,
580 continuous_dimension_type2,
581 continuous_dimension_type3>
582 coord_eval_3D(
583 ddc::coordinate(i1),
584 ddc::coordinate(i2),
585 ddc::coordinate(i3));
586 spline_eval_3D(i1, i2, i3)
587 = eval(coord_eval_3D(i1, i2, i3), spline_coef_3D);
588 }
589 }
590 }
591 });
592 }
593
594 /**
595 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along first dimension of interest.
596 *
597 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
598 * obtained via various methods, such as using a SplineBuilder3D.
599 *
600 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
601 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
602 *
603 * @return The derivative of the spline function at the desired coordinate.
604 */
605 template <class Layout, class... CoordsDims>
607 ddc::Coordinate<CoordsDims...> const& coord_eval,
608 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
609 spline_coef) const
610 {
611 return eval_no_bc<eval_deriv_type, eval_type, eval_type>(coord_eval, spline_coef);
612 }
613
614 /**
615 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along second dimension of interest.
616 *
617 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
618 * obtained via various methods, such as using a SplineBuilder3D.
619 *
620 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
621 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
622 *
623 * @return The derivative of the spline function at the desired coordinate.
624 */
625 template <class Layout, class... CoordsDims>
627 ddc::Coordinate<CoordsDims...> const& coord_eval,
628 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
629 spline_coef) const
630 {
631 return eval_no_bc<eval_type, eval_deriv_type, eval_type>(coord_eval, spline_coef);
632 }
633
634 /**
635 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along third dimension of interest.
636 *
637 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
638 * obtained via various methods, such as using a SplineBuilder3D.
639 *
640 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
641 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
642 *
643 * @return The derivative of the spline function at the desired coordinate.
644 */
645 template <class Layout, class... CoordsDims>
647 ddc::Coordinate<CoordsDims...> const& coord_eval,
648 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
649 spline_coef) const
650 {
651 return eval_no_bc<eval_type, eval_type, eval_deriv_type>(coord_eval, spline_coef);
652 }
653
654 /**
655 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the first and second dimensions.
656 *
657 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
658 * obtained via various methods, such as using a SplineBuilder3D.
659 *
660 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
661 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
662 *
663 * @return The derivative of the spline function at the desired coordinate.
664 */
665 template <class Layout, class... CoordsDims>
667 ddc::Coordinate<CoordsDims...> const& coord_eval,
668 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
669 spline_coef) const
670 {
671 return eval_no_bc<eval_deriv_type, eval_deriv_type, eval_type>(coord_eval, spline_coef);
672 }
673
674 /**
675 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the second and third dimensions.
676 *
677 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
678 * obtained via various methods, such as using a SplineBuilder3D.
679 *
680 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
681 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
682 *
683 * @return The derivative of the spline function at the desired coordinate.
684 */
685 template <class Layout, class... CoordsDims>
687 ddc::Coordinate<CoordsDims...> const& coord_eval,
688 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
689 spline_coef) const
690 {
691 return eval_no_bc<eval_type, eval_deriv_type, eval_deriv_type>(coord_eval, spline_coef);
692 }
693
694 /**
695 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the first and third dimensions.
696 *
697 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
698 * obtained via various methods, such as using a SplineBuilder3D.
699 *
700 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
701 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
702 *
703 * @return The derivative of the spline function at the desired coordinate.
704 */
705 template <class Layout, class... CoordsDims>
707 ddc::Coordinate<CoordsDims...> const& coord_eval,
708 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
709 spline_coef) const
710 {
711 return eval_no_bc<eval_deriv_type, eval_type, eval_deriv_type>(coord_eval, spline_coef);
712 }
713
714 /**
715 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the dimensions of interest.
716 *
717 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
718 * obtained via various methods, such as using a SplineBuilder3D.
719 *
720 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
721 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
722 *
723 * @return The derivative of the spline function at the desired coordinate.
724 */
725 template <class Layout, class... CoordsDims>
727 ddc::Coordinate<CoordsDims...> const& coord_eval,
728 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
729 spline_coef) const
730 {
731 return eval_no_bc<
732 eval_deriv_type,
733 eval_deriv_type,
734 eval_deriv_type>(coord_eval, spline_coef);
735 }
736
737 /**
738 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along a specified dimension of interest.
739 *
740 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
741 * obtained via various methods, such as using a SplineBuilder3D.
742 *
743 * @tparam InterestDim Dimension along which differentiation is performed.
744 *
745 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
746 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
747 *
748 * @return The derivative of the spline function at the desired coordinate.
749 */
750 template <class InterestDim, class Layout, class... CoordsDims>
751 KOKKOS_FUNCTION double deriv(
752 ddc::Coordinate<CoordsDims...> const& coord_eval,
753 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
754 spline_coef) const
755 {
756 static_assert(
757 std::is_same_v<InterestDim, continuous_dimension_type1>
758 || std::is_same_v<InterestDim, continuous_dimension_type2>
759 || std::is_same_v<InterestDim, continuous_dimension_type3>);
760 if constexpr (std::is_same_v<InterestDim, continuous_dimension_type1>) {
761 return deriv_dim_1(coord_eval, spline_coef);
762 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type2>) {
763 return deriv_dim_2(coord_eval, spline_coef);
764 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type3>) {
765 return deriv_dim_3(coord_eval, spline_coef);
766 }
767 }
768
769 /**
770 * @brief Double-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along specified dimensions of interest.
771 *
772 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
773 * obtained via various methods, such as using a SplineBuilder3D.
774 *
775 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
776 *
777 * @tparam InterestDim1 First dimension along which differentiation is performed.
778 * @tparam InterestDim2 Second dimension along which differentiation is performed.
779 *
780 * @param coord_eval The coordinate where the spline is double-differentiated. Note that only the components along the dimensions of interest are used.
781 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
782 *
783 * @return The derivative of the spline function at the desired coordinate.
784 */
785 template <class InterestDim1, class InterestDim2, class Layout, class... CoordsDims>
786 KOKKOS_FUNCTION double deriv2(
787 ddc::Coordinate<CoordsDims...> const& coord_eval,
788 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
789 spline_coef) const
790 {
791 static_assert(
792 (std::is_same_v<InterestDim1, continuous_dimension_type1>
793 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
794 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
795 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
796 || (std::is_same_v<InterestDim1, continuous_dimension_type2>
797 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
798 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
799 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
800 || (std::is_same_v<InterestDim1, continuous_dimension_type1>
801 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
802 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
803 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
804
805 if constexpr (
806 (std::is_same_v<InterestDim1, continuous_dimension_type1>
807 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
808 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
809 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
810 return deriv_1_and_2(coord_eval, spline_coef);
811 } else if constexpr (
812 (std::is_same_v<InterestDim1, continuous_dimension_type2>
813 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
814 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
815 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
816 return deriv_2_and_3(coord_eval, spline_coef);
817 } else if constexpr (
818 (std::is_same_v<InterestDim1, continuous_dimension_type1>
819 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
820 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
821 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
822 return deriv_1_and_3(coord_eval, spline_coef);
823 }
824 }
825
826 /**
827 * @brief Triple-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along specified dimensions of interest.
828 *
829 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
830 * obtained via various methods, such as using a SplineBuilder3D.
831 *
832 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
833 *
834 * @tparam InterestDim1 First dimension along which differentiation is performed.
835 * @tparam InterestDim2 Second dimension along which differentiation is performed.
836 * @tparam InterestDim3 Third dimension along which differentiation is performed.
837 *
838 * @param coord_eval The coordinate where the spline is triple-differentiated. Note that only the components along the dimensions of interest are used.
839 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
840 *
841 * @return The derivative of the spline function at the desired coordinate.
842 */
843 template <
844 class InterestDim1,
845 class InterestDim2,
846 class InterestDim3,
847 class Layout,
848 class... CoordsDims>
849 KOKKOS_FUNCTION double deriv3(
850 ddc::Coordinate<CoordsDims...> const& coord_eval,
851 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
852 spline_coef) const
853 {
854 static_assert(
855 (std::is_same_v<InterestDim1, continuous_dimension_type1>
856 && std::is_same_v<InterestDim2, continuous_dimension_type2>
857 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
858 || (std::is_same_v<InterestDim3, continuous_dimension_type1>
859 && std::is_same_v<InterestDim1, continuous_dimension_type2>
860 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
861 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
862 && std::is_same_v<InterestDim3, continuous_dimension_type2>
863 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
864
865 return deriv_1_2_3(coord_eval, spline_coef);
866 }
867
868 /**
869 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
870 *
871 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
872 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
873 *
874 * This is not a nD evaluation. This is a batched 3D differentiation.
875 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
876 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
877 *
878 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
879 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
880 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
881 * the set of 3D spline coefficients retained to perform the evaluation).
882 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
883 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
884 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
885 * the set of 3D spline coefficients retained to perform the evaluation).
886 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
887 */
888 template <
889 class Layout1,
890 class Layout2,
891 class Layout3,
892 class BatchedInterpolationDDom,
893 class... CoordsDims>
894 void deriv_dim_1(
895 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
896 spline_eval,
897 ddc::ChunkSpan<
898 ddc::Coordinate<CoordsDims...> const,
899 BatchedInterpolationDDom,
900 Layout2,
901 memory_space> const coords_eval,
902 ddc::ChunkSpan<
903 double const,
904 batched_spline_domain_type<BatchedInterpolationDDom>,
905 Layout3,
906 memory_space> const spline_coef) const
907 {
908 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
909 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
910 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
911 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
912 ddc::parallel_for_each(
913 "ddc_splines_differentiate_3d_dim_1",
914 exec_space(),
915 batch_domain,
916 KOKKOS_CLASS_LAMBDA(
917 typename batch_domain_type<
918 BatchedInterpolationDDom>::discrete_element_type const j) {
919 auto const spline_eval_3D = spline_eval[j];
920 auto const coords_eval_3D = coords_eval[j];
921 auto const spline_coef_3D = spline_coef[j];
922 for (auto const i1 : evaluation_domain1) {
923 for (auto const i2 : evaluation_domain2) {
924 for (auto const i3 : evaluation_domain3) {
925 spline_eval_3D(i1, i2, i3) = eval_no_bc<
926 eval_deriv_type,
927 eval_type,
928 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
929 }
930 }
931 }
932 });
933 }
934
935 /**
936 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
937 *
938 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
939 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
940 *
941 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
942 * This means that for each slice of spline_eval the evaluation is performed with
943 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
944 *
945 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
946 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
947 */
948 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
949 void deriv_dim_1(
950 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
951 spline_eval,
952 ddc::ChunkSpan<
953 double const,
954 batched_spline_domain_type<BatchedInterpolationDDom>,
955 Layout2,
956 memory_space> const spline_coef) const
957 {
958 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
959 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
960 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
961 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
962 ddc::parallel_for_each(
963 "ddc_splines_differentiate_3d_dim_1",
964 exec_space(),
965 batch_domain,
966 KOKKOS_CLASS_LAMBDA(
967 typename batch_domain_type<
968 BatchedInterpolationDDom>::discrete_element_type const j) {
969 auto const spline_eval_3D = spline_eval[j];
970 auto const spline_coef_3D = spline_coef[j];
971 for (auto const i1 : evaluation_domain1) {
972 for (auto const i2 : evaluation_domain2) {
973 for (auto const i3 : evaluation_domain3) {
974 ddc::Coordinate<
975 continuous_dimension_type1,
976 continuous_dimension_type2,
977 continuous_dimension_type3>
978 coord_eval_3D(
979 ddc::coordinate(i1),
980 ddc::coordinate(i2),
981 ddc::coordinate(i3));
982 spline_eval_3D(i1, i2, i3) = eval_no_bc<
983 eval_deriv_type,
984 eval_type,
985 eval_type>(coord_eval_3D, spline_coef_3D);
986 }
987 }
988 }
989 });
990 }
991
992 /**
993 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
994 *
995 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
996 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
997 *
998 * This is not a nD differentiation. This is a batched 3D differentiation.
999 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1000 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1001 *
1002 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1003 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1004 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1005 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1006 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1007 * the set of 3D spline coefficients retained to perform the evaluation).
1008 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1009 */
1010 template <
1011 class Layout1,
1012 class Layout2,
1013 class Layout3,
1014 class BatchedInterpolationDDom,
1015 class... CoordsDims>
1016 void deriv_dim_2(
1017 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1018 spline_eval,
1019 ddc::ChunkSpan<
1020 ddc::Coordinate<CoordsDims...> const,
1021 BatchedInterpolationDDom,
1022 Layout2,
1023 memory_space> const coords_eval,
1024 ddc::ChunkSpan<
1025 double const,
1026 batched_spline_domain_type<BatchedInterpolationDDom>,
1027 Layout3,
1028 memory_space> const spline_coef) const
1029 {
1030 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1031 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1032 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1033 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1034 ddc::parallel_for_each(
1035 "ddc_splines_differentiate_3d_dim_2",
1036 exec_space(),
1037 batch_domain,
1038 KOKKOS_CLASS_LAMBDA(
1039 typename batch_domain_type<
1040 BatchedInterpolationDDom>::discrete_element_type const j) {
1041 auto const spline_eval_3D = spline_eval[j];
1042 auto const coords_eval_3D = coords_eval[j];
1043 auto const spline_coef_3D = spline_coef[j];
1044 for (auto const i1 : evaluation_domain1) {
1045 for (auto const i2 : evaluation_domain2) {
1046 for (auto const i3 : evaluation_domain3) {
1047 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1048 eval_type,
1049 eval_deriv_type,
1050 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
1051 }
1052 }
1053 }
1054 });
1055 }
1056
1057 /**
1058 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
1059 *
1060 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1061 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1062 *
1063 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1064 * This means that for each slice of spline_eval the evaluation is performed with
1065 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1066 *
1067 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1068 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1069 */
1070 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1071 void deriv_dim_2(
1072 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1073 spline_eval,
1074 ddc::ChunkSpan<
1075 double const,
1076 batched_spline_domain_type<BatchedInterpolationDDom>,
1077 Layout2,
1078 memory_space> const spline_coef) const
1079 {
1080 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1081 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1082 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1083 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1084 ddc::parallel_for_each(
1085 "ddc_splines_differentiate_3d_dim_2",
1086 exec_space(),
1087 batch_domain,
1088 KOKKOS_CLASS_LAMBDA(
1089 typename batch_domain_type<
1090 BatchedInterpolationDDom>::discrete_element_type const j) {
1091 auto const spline_eval_3D = spline_eval[j];
1092 auto const spline_coef_3D = spline_coef[j];
1093 for (auto const i1 : evaluation_domain1) {
1094 for (auto const i2 : evaluation_domain2) {
1095 for (auto const i3 : evaluation_domain3) {
1096 ddc::Coordinate<
1097 continuous_dimension_type1,
1098 continuous_dimension_type2,
1099 continuous_dimension_type3>
1100 coord_eval_3D(
1101 ddc::coordinate(i1),
1102 ddc::coordinate(i2),
1103 ddc::coordinate(i3));
1104 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1105 eval_type,
1106 eval_deriv_type,
1107 eval_type>(coord_eval_3D, spline_coef_3D);
1108 }
1109 }
1110 }
1111 });
1112 }
1113
1114 /**
1115 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along third dimension of interest.
1116 *
1117 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1118 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1119 *
1120 * This is not a nD differentiation. This is a batched 3D differentiation.
1121 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1122 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1123 *
1124 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1125 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1126 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1127 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1128 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1129 * the set of 3D spline coefficients retained to perform the evaluation).
1130 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1131 */
1132 template <
1133 class Layout1,
1134 class Layout2,
1135 class Layout3,
1136 class BatchedInterpolationDDom,
1137 class... CoordsDims>
1138 void deriv_dim_3(
1139 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1140 spline_eval,
1141 ddc::ChunkSpan<
1142 ddc::Coordinate<CoordsDims...> const,
1143 BatchedInterpolationDDom,
1144 Layout2,
1145 memory_space> const coords_eval,
1146 ddc::ChunkSpan<
1147 double const,
1148 batched_spline_domain_type<BatchedInterpolationDDom>,
1149 Layout3,
1150 memory_space> const spline_coef) const
1151 {
1152 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1153 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1154 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1155 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1156 ddc::parallel_for_each(
1157 "ddc_splines_differentiate_3d_dim_3",
1158 exec_space(),
1159 batch_domain,
1160 KOKKOS_CLASS_LAMBDA(
1161 typename batch_domain_type<
1162 BatchedInterpolationDDom>::discrete_element_type const j) {
1163 auto const spline_eval_3D = spline_eval[j];
1164 auto const coords_eval_3D = coords_eval[j];
1165 auto const spline_coef_3D = spline_coef[j];
1166 for (auto const i1 : evaluation_domain1) {
1167 for (auto const i2 : evaluation_domain2) {
1168 for (auto const i3 : evaluation_domain3) {
1169 spline_eval_3D(i1, i2, i3)
1170 = eval_no_bc<eval_type, eval_type, eval_deriv_type>(
1171 coords_eval_3D(i1, i2, i3),
1172 spline_coef_3D);
1173 }
1174 }
1175 }
1176 });
1177 }
1178
1179 /**
1180 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along third dimension of interest.
1181 *
1182 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1183 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1184 *
1185 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1186 * This means that for each slice of spline_eval the evaluation is performed with
1187 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1188 *
1189 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1190 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1191 */
1192 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1193 void deriv_dim_3(
1194 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1195 spline_eval,
1196 ddc::ChunkSpan<
1197 double const,
1198 batched_spline_domain_type<BatchedInterpolationDDom>,
1199 Layout2,
1200 memory_space> const spline_coef) const
1201 {
1202 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1203 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1204 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1205 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1206 ddc::parallel_for_each(
1207 "ddc_splines_differentiate_3d_dim_3",
1208 exec_space(),
1209 batch_domain,
1210 KOKKOS_CLASS_LAMBDA(
1211 typename batch_domain_type<
1212 BatchedInterpolationDDom>::discrete_element_type const j) {
1213 auto const spline_eval_3D = spline_eval[j];
1214 auto const spline_coef_3D = spline_coef[j];
1215 for (auto const i1 : evaluation_domain1) {
1216 for (auto const i2 : evaluation_domain2) {
1217 for (auto const i3 : evaluation_domain3) {
1218 ddc::Coordinate<
1219 continuous_dimension_type1,
1220 continuous_dimension_type2,
1221 continuous_dimension_type3>
1222 coord_eval_3D(
1223 ddc::coordinate(i1),
1224 ddc::coordinate(i2),
1225 ddc::coordinate(i3));
1226 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1227 eval_type,
1228 eval_type,
1229 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1230 }
1231 }
1232 }
1233 });
1234 }
1235
1236 /**
1237 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and second dimensions of interest.
1238 *
1239 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1240 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1241 *
1242 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1243 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1244 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1245 *
1246 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1247 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1248 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1249 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1250 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1251 * the set of 3D spline coefficients retained to perform the evaluation).
1252 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1253 */
1254 template <
1255 class Layout1,
1256 class Layout2,
1257 class Layout3,
1258 class BatchedInterpolationDDom,
1259 class... CoordsDims>
1260 void deriv_1_and_2(
1261 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1262 spline_eval,
1263 ddc::ChunkSpan<
1264 ddc::Coordinate<CoordsDims...> const,
1265 BatchedInterpolationDDom,
1266 Layout2,
1267 memory_space> const coords_eval,
1268 ddc::ChunkSpan<
1269 double const,
1270 batched_spline_domain_type<BatchedInterpolationDDom>,
1271 Layout3,
1272 memory_space> const spline_coef) const
1273 {
1274 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1275 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1276 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1277 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1278 ddc::parallel_for_each(
1279 "ddc_splines_cross_differentiate_3d_dim1_2",
1280 exec_space(),
1281 batch_domain,
1282 KOKKOS_CLASS_LAMBDA(
1283 typename batch_domain_type<
1284 BatchedInterpolationDDom>::discrete_element_type const j) {
1285 auto const spline_eval_3D = spline_eval[j];
1286 auto const coords_eval_3D = coords_eval[j];
1287 auto const spline_coef_3D = spline_coef[j];
1288 for (auto const i1 : evaluation_domain1) {
1289 for (auto const i2 : evaluation_domain2) {
1290 for (auto const i3 : evaluation_domain3) {
1291 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1292 eval_deriv_type,
1293 eval_deriv_type,
1294 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
1295 }
1296 }
1297 }
1298 });
1299 }
1300
1301 /**
1302 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and second dimensions of interest.
1303 *
1304 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1305 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1306 *
1307 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1308 * This means that for each slice of spline_eval the evaluation is performed with
1309 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1310 *
1311 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1312 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1313 */
1314 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1315 void deriv_1_and_2(
1316 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1317 spline_eval,
1318 ddc::ChunkSpan<
1319 double const,
1320 batched_spline_domain_type<BatchedInterpolationDDom>,
1321 Layout2,
1322 memory_space> const spline_coef) const
1323 {
1324 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1325 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1326 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1327 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1328 ddc::parallel_for_each(
1329 "ddc_splines_cross_differentiate_3d_dim1_2",
1330 exec_space(),
1331 batch_domain,
1332 KOKKOS_CLASS_LAMBDA(
1333 typename batch_domain_type<
1334 BatchedInterpolationDDom>::discrete_element_type const j) {
1335 auto const spline_eval_3D = spline_eval[j];
1336 auto const spline_coef_3D = spline_coef[j];
1337 for (auto const i1 : evaluation_domain1) {
1338 for (auto const i2 : evaluation_domain2) {
1339 for (auto const i3 : evaluation_domain3) {
1340 ddc::Coordinate<
1341 continuous_dimension_type1,
1342 continuous_dimension_type2,
1343 continuous_dimension_type3>
1344 coord_eval_3D(
1345 ddc::coordinate(i1),
1346 ddc::coordinate(i2),
1347 ddc::coordinate(i3));
1348 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1349 eval_deriv_type,
1350 eval_deriv_type,
1351 eval_type>(coord_eval_3D, spline_coef_3D);
1352 }
1353 }
1354 }
1355 });
1356 }
1357
1358 /**
1359 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the second and third dimensions of interest.
1360 *
1361 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1362 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1363 *
1364 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1365 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1366 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1367 *
1368 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1369 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1370 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1371 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1372 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1373 * the set of 3D spline coefficients retained to perform the evaluation).
1374 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1375 */
1376 template <
1377 class Layout1,
1378 class Layout2,
1379 class Layout3,
1380 class BatchedInterpolationDDom,
1381 class... CoordsDims>
1382 void deriv_2_and_3(
1383 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1384 spline_eval,
1385 ddc::ChunkSpan<
1386 ddc::Coordinate<CoordsDims...> const,
1387 BatchedInterpolationDDom,
1388 Layout2,
1389 memory_space> const coords_eval,
1390 ddc::ChunkSpan<
1391 double const,
1392 batched_spline_domain_type<BatchedInterpolationDDom>,
1393 Layout3,
1394 memory_space> const spline_coef) const
1395 {
1396 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1397 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1398 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1399 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1400 ddc::parallel_for_each(
1401 "ddc_splines_cross_differentiate_3d_dim2_3",
1402 exec_space(),
1403 batch_domain,
1404 KOKKOS_CLASS_LAMBDA(
1405 typename batch_domain_type<
1406 BatchedInterpolationDDom>::discrete_element_type const j) {
1407 auto const spline_eval_3D = spline_eval[j];
1408 auto const coords_eval_3D = coords_eval[j];
1409 auto const spline_coef_3D = spline_coef[j];
1410 for (auto const i1 : evaluation_domain1) {
1411 for (auto const i2 : evaluation_domain2) {
1412 for (auto const i3 : evaluation_domain3) {
1413 spline_eval_3D(i1, i2, i3)
1414 = eval_no_bc<eval_type, eval_deriv_type, eval_deriv_type>(
1415 coords_eval_3D(i1, i2, i3),
1416 spline_coef_3D);
1417 }
1418 }
1419 }
1420 });
1421 }
1422
1423 /**
1424 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the second and third dimensions of interest.
1425 *
1426 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1427 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1428 *
1429 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1430 * This means that for each slice of spline_eval the evaluation is performed with
1431 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1432 *
1433 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1434 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1435 */
1436 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1437 void deriv_2_and_3(
1438 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1439 spline_eval,
1440 ddc::ChunkSpan<
1441 double const,
1442 batched_spline_domain_type<BatchedInterpolationDDom>,
1443 Layout2,
1444 memory_space> const spline_coef) const
1445 {
1446 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1447 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1448 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1449 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1450 ddc::parallel_for_each(
1451 "ddc_splines_cross_differentiate_3d_dim2_3",
1452 exec_space(),
1453 batch_domain,
1454 KOKKOS_CLASS_LAMBDA(
1455 typename batch_domain_type<
1456 BatchedInterpolationDDom>::discrete_element_type const j) {
1457 auto const spline_eval_3D = spline_eval[j];
1458 auto const spline_coef_3D = spline_coef[j];
1459 for (auto const i1 : evaluation_domain1) {
1460 for (auto const i2 : evaluation_domain2) {
1461 for (auto const i3 : evaluation_domain3) {
1462 ddc::Coordinate<
1463 continuous_dimension_type1,
1464 continuous_dimension_type2,
1465 continuous_dimension_type3>
1466 coord_eval_3D(
1467 ddc::coordinate(i1),
1468 ddc::coordinate(i2),
1469 ddc::coordinate(i3));
1470 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1471 eval_type,
1472 eval_deriv_type,
1473 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1474 }
1475 }
1476 }
1477 });
1478 }
1479
1480 /**
1481 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and third dimensions of interest.
1482 *
1483 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1484 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1485 *
1486 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1487 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1488 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1489 *
1490 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1491 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1492 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1493 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1494 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1495 * the set of 3D spline coefficients retained to perform the evaluation).
1496 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1497 */
1498 template <
1499 class Layout1,
1500 class Layout2,
1501 class Layout3,
1502 class BatchedInterpolationDDom,
1503 class... CoordsDims>
1504 void deriv_1_and_3(
1505 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1506 spline_eval,
1507 ddc::ChunkSpan<
1508 ddc::Coordinate<CoordsDims...> const,
1509 BatchedInterpolationDDom,
1510 Layout2,
1511 memory_space> const coords_eval,
1512 ddc::ChunkSpan<
1513 double const,
1514 batched_spline_domain_type<BatchedInterpolationDDom>,
1515 Layout3,
1516 memory_space> const spline_coef) const
1517 {
1518 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1519 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1520 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1521 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1522 ddc::parallel_for_each(
1523 "ddc_splines_cross_differentiate_3d_dim1_3",
1524 exec_space(),
1525 batch_domain,
1526 KOKKOS_CLASS_LAMBDA(
1527 typename batch_domain_type<
1528 BatchedInterpolationDDom>::discrete_element_type const j) {
1529 auto const spline_eval_3D = spline_eval[j];
1530 auto const coords_eval_3D = coords_eval[j];
1531 auto const spline_coef_3D = spline_coef[j];
1532 for (auto const i1 : evaluation_domain1) {
1533 for (auto const i2 : evaluation_domain2) {
1534 for (auto const i3 : evaluation_domain3) {
1535 spline_eval_3D(i1, i2, i3)
1536 = eval_no_bc<eval_deriv_type, eval_type, eval_deriv_type>(
1537 coords_eval_3D(i1, i2, i3),
1538 spline_coef_3D);
1539 }
1540 }
1541 }
1542 });
1543 }
1544
1545 /**
1546 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and third dimensions of interest.
1547 *
1548 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1549 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1550 *
1551 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1552 * This means that for each slice of spline_eval the evaluation is performed with
1553 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1554 *
1555 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1556 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1557 */
1558 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1559 void deriv_1_and_3(
1560 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1561 spline_eval,
1562 ddc::ChunkSpan<
1563 double const,
1564 batched_spline_domain_type<BatchedInterpolationDDom>,
1565 Layout2,
1566 memory_space> const spline_coef) const
1567 {
1568 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1569 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1570 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1571 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1572 ddc::parallel_for_each(
1573 "ddc_splines_cross_differentiate_3d_dim1_3",
1574 exec_space(),
1575 batch_domain,
1576 KOKKOS_CLASS_LAMBDA(
1577 typename batch_domain_type<
1578 BatchedInterpolationDDom>::discrete_element_type const j) {
1579 auto const spline_eval_3D = spline_eval[j];
1580 auto const spline_coef_3D = spline_coef[j];
1581 for (auto const i1 : evaluation_domain1) {
1582 for (auto const i2 : evaluation_domain2) {
1583 for (auto const i3 : evaluation_domain3) {
1584 ddc::Coordinate<
1585 continuous_dimension_type1,
1586 continuous_dimension_type2,
1587 continuous_dimension_type3>
1588 coord_eval_3D(
1589 ddc::coordinate(i1),
1590 ddc::coordinate(i2),
1591 ddc::coordinate(i3));
1592 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1593 eval_deriv_type,
1594 eval_type,
1595 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1596 }
1597 }
1598 }
1599 });
1600 }
1601
1602 /**
1603 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
1604 *
1605 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1606 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1607 *
1608 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1609 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1610 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1611 *
1612 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1613 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1614 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1615 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1616 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1617 * the set of 3D spline coefficients retained to perform the evaluation).
1618 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1619 */
1620 template <
1621 class Layout1,
1622 class Layout2,
1623 class Layout3,
1624 class BatchedInterpolationDDom,
1625 class... CoordsDims>
1626 void deriv_1_2_3(
1627 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1628 spline_eval,
1629 ddc::ChunkSpan<
1630 ddc::Coordinate<CoordsDims...> const,
1631 BatchedInterpolationDDom,
1632 Layout2,
1633 memory_space> const coords_eval,
1634 ddc::ChunkSpan<
1635 double const,
1636 batched_spline_domain_type<BatchedInterpolationDDom>,
1637 Layout3,
1638 memory_space> const spline_coef) const
1639 {
1640 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1641 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1642 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1643 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1644 ddc::parallel_for_each(
1645 "ddc_splines_cross_differentiate_3d_dim1_2_3",
1646 exec_space(),
1647 batch_domain,
1648 KOKKOS_CLASS_LAMBDA(
1649 typename batch_domain_type<
1650 BatchedInterpolationDDom>::discrete_element_type const j) {
1651 auto const spline_eval_3D = spline_eval[j];
1652 auto const coords_eval_3D = coords_eval[j];
1653 auto const spline_coef_3D = spline_coef[j];
1654 for (auto const i1 : evaluation_domain1) {
1655 for (auto const i2 : evaluation_domain2) {
1656 for (auto const i3 : evaluation_domain3) {
1657 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1658 eval_deriv_type,
1659 eval_deriv_type,
1660 eval_deriv_type>(
1661 coords_eval_3D(i1, i2, i3),
1662 spline_coef_3D);
1663 }
1664 }
1665 }
1666 });
1667 }
1668
1669 /**
1670 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
1671 *
1672 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1673 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1674 *
1675 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1676 * This means that for each slice of spline_eval the evaluation is performed with
1677 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1678 *
1679 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1680 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1681 */
1682 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1683 void deriv_1_2_3(
1684 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1685 spline_eval,
1686 ddc::ChunkSpan<
1687 double const,
1688 batched_spline_domain_type<BatchedInterpolationDDom>,
1689 Layout2,
1690 memory_space> const spline_coef) const
1691 {
1692 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1693 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1694 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1695 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1696 ddc::parallel_for_each(
1697 "ddc_splines_cross_differentiate_3d_dim1_2_3",
1698 exec_space(),
1699 batch_domain,
1700 KOKKOS_CLASS_LAMBDA(
1701 typename batch_domain_type<
1702 BatchedInterpolationDDom>::discrete_element_type const j) {
1703 auto const spline_eval_3D = spline_eval[j];
1704 auto const spline_coef_3D = spline_coef[j];
1705 for (auto const i1 : evaluation_domain1) {
1706 for (auto const i2 : evaluation_domain2) {
1707 for (auto const i3 : evaluation_domain3) {
1708 ddc::Coordinate<
1709 continuous_dimension_type1,
1710 continuous_dimension_type2,
1711 continuous_dimension_type3>
1712 coord_eval_3D(
1713 ddc::coordinate(i1),
1714 ddc::coordinate(i2),
1715 ddc::coordinate(i3));
1716 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1717 eval_deriv_type,
1718 eval_deriv_type,
1719 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1720 }
1721 }
1722 }
1723 });
1724 }
1725
1726 /**
1727 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
1728 *
1729 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1730 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1731 *
1732 * This is not a nD evaluation. This is a batched 3D differentiation.
1733 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1734 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1735 *
1736 * @tparam InterestDim Dimension along which differentiation is performed.
1737 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1738 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1739 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1740 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1741 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1742 * the set of 3D spline coefficients retained to perform the evaluation).
1743 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1744 */
1745 template <
1746 class InterestDim,
1747 class Layout1,
1748 class Layout2,
1749 class Layout3,
1750 class BatchedInterpolationDDom,
1751 class... CoordsDims>
1752 void deriv(
1753 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1754 spline_eval,
1755 ddc::ChunkSpan<
1756 ddc::Coordinate<CoordsDims...> const,
1757 BatchedInterpolationDDom,
1758 Layout2,
1759 memory_space> const coords_eval,
1760 ddc::ChunkSpan<
1761 double const,
1762 batched_spline_domain_type<BatchedInterpolationDDom>,
1763 Layout3,
1764 memory_space> const spline_coef) const
1765 {
1766 static_assert(
1767 std::is_same_v<InterestDim, continuous_dimension_type1>
1768 || std::is_same_v<InterestDim, continuous_dimension_type2>
1769 || std::is_same_v<InterestDim, continuous_dimension_type3>);
1770 if constexpr (std::is_same_v<InterestDim, continuous_dimension_type1>) {
1771 return deriv_dim_1(spline_eval, coords_eval, spline_coef);
1772 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type2>) {
1773 return deriv_dim_2(spline_eval, coords_eval, spline_coef);
1774 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type3>) {
1775 return deriv_dim_3(spline_eval, coords_eval, spline_coef);
1776 }
1777 }
1778
1779 /**
1780 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
1781 *
1782 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1783 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1784 *
1785 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1786 * This means that for each slice of spline_eval the evaluation is performed with
1787 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1788 *
1789 * @tparam InterestDim Dimension along which differentiation is performed.
1790 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1791 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1792 */
1793 template <class InterestDim, class Layout1, class Layout2, class BatchedInterpolationDDom>
1794 void deriv(
1795 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1796 spline_eval,
1797 ddc::ChunkSpan<
1798 double const,
1799 batched_spline_domain_type<BatchedInterpolationDDom>,
1800 Layout2,
1801 memory_space> const spline_coef) const
1802 {
1803 static_assert(
1804 std::is_same_v<InterestDim, continuous_dimension_type1>
1805 || std::is_same_v<InterestDim, continuous_dimension_type2>
1806 || std::is_same_v<InterestDim, continuous_dimension_type3>);
1807 if constexpr (std::is_same_v<InterestDim, continuous_dimension_type1>) {
1808 return deriv_dim_1(spline_eval, spline_coef);
1809 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type2>) {
1810 return deriv_dim_2(spline_eval, spline_coef);
1811 } else if constexpr (std::is_same_v<InterestDim, continuous_dimension_type3>) {
1812 return deriv_dim_3(spline_eval, spline_coef);
1813 }
1814 }
1815
1816 /**
1817 * @brief Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1818 *
1819 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1820 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1821 *
1822 * This is not a nD evaluation. This is a batched 3D differentiation.
1823 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1824 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1825 *
1826 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1827 *
1828 * @tparam InterestDim1 First dimension along which differentiation is performed.
1829 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1830 *
1831 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1832 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1833 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1834 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1835 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1836 * the set of 3D spline coefficients retained to perform the evaluation).
1837 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1838 */
1839 template <
1840 class InterestDim1,
1841 class InterestDim2,
1842 class Layout1,
1843 class Layout2,
1844 class Layout3,
1845 class BatchedInterpolationDDom,
1846 class... CoordsDims>
1847 void deriv2(
1848 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1849 spline_eval,
1850 ddc::ChunkSpan<
1851 ddc::Coordinate<CoordsDims...> const,
1852 BatchedInterpolationDDom,
1853 Layout2,
1854 memory_space> const coords_eval,
1855 ddc::ChunkSpan<
1856 double const,
1857 batched_spline_domain_type<BatchedInterpolationDDom>,
1858 Layout3,
1859 memory_space> const spline_coef) const
1860 {
1861 static_assert(
1862 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1863 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1864 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1865 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
1866 || (std::is_same_v<InterestDim1, continuous_dimension_type2>
1867 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1868 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
1869 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
1870 || (std::is_same_v<InterestDim1, continuous_dimension_type1>
1871 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1872 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1873 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
1874
1875 if constexpr (
1876 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1877 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1878 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1879 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
1880 return deriv_1_and_2(spline_eval, coords_eval, spline_coef);
1881 } else if constexpr (
1882 (std::is_same_v<InterestDim1, continuous_dimension_type2>
1883 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1884 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
1885 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1886 return deriv_2_and_3(spline_eval, coords_eval, spline_coef);
1887 } else if constexpr (
1888 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1889 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1890 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1891 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1892 return deriv_1_and_3(spline_eval, coords_eval, spline_coef);
1893 }
1894 }
1895
1896 /**
1897 * @brief Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1898 *
1899 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1900 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1901 *
1902 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1903 * This means that for each slice of spline_eval the evaluation is performed with
1904 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1905 *
1906 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1907 *
1908 * @tparam InterestDim1 First dimension along which differentiation is performed.
1909 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1910 *
1911 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1912 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1913 */
1914 template <
1915 class InterestDim1,
1916 class InterestDim2,
1917 class Layout1,
1918 class Layout2,
1919 class BatchedInterpolationDDom>
1920 void deriv2(
1921 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1922 spline_eval,
1923 ddc::ChunkSpan<
1924 double const,
1925 batched_spline_domain_type<BatchedInterpolationDDom>,
1926 Layout2,
1927 memory_space> const spline_coef) const
1928 {
1929 static_assert(
1930 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1931 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1932 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1933 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
1934 || (std::is_same_v<InterestDim1, continuous_dimension_type2>
1935 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1936 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
1937 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
1938 || (std::is_same_v<InterestDim1, continuous_dimension_type1>
1939 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1940 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1941 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
1942
1943 if constexpr (
1944 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1945 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1946 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1947 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
1948 return deriv_1_and_2(spline_eval, spline_coef);
1949 } else if constexpr (
1950 (std::is_same_v<InterestDim1, continuous_dimension_type2>
1951 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1952 || (std::is_same_v<InterestDim2, continuous_dimension_type2>
1953 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1954 return deriv_2_and_3(spline_eval, spline_coef);
1955 } else if constexpr (
1956 (std::is_same_v<InterestDim1, continuous_dimension_type1>
1957 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1958 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
1959 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1960 return deriv_1_and_3(spline_eval, spline_coef);
1961 }
1962 }
1963
1964 /**
1965 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimensions of interest.
1966 *
1967 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1968 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1969 *
1970 * This is not a nD evaluation. This is a batched 3D differentiation.
1971 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1972 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1973 *
1974 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
1975 *
1976 * @tparam InterestDim1 First dimension along which differentiation is performed.
1977 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1978 * @tparam InterestDim3 Third dimension along which differentiation is performed.
1979 *
1980 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1981 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1982 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1983 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1984 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1985 * the set of 3D spline coefficients retained to perform the evaluation).
1986 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1987 */
1988 template <
1989 class InterestDim1,
1990 class InterestDim2,
1991 class InterestDim3,
1992 class Layout1,
1993 class Layout2,
1994 class Layout3,
1995 class BatchedInterpolationDDom,
1996 class... CoordsDims>
1997 void deriv3(
1998 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1999 spline_eval,
2000 ddc::ChunkSpan<
2001 ddc::Coordinate<CoordsDims...> const,
2002 BatchedInterpolationDDom,
2003 Layout2,
2004 memory_space> const coords_eval,
2005 ddc::ChunkSpan<
2006 double const,
2007 batched_spline_domain_type<BatchedInterpolationDDom>,
2008 Layout3,
2009 memory_space> const spline_coef) const
2010 {
2011 static_assert(
2012 (std::is_same_v<InterestDim1, continuous_dimension_type1>
2013 && std::is_same_v<InterestDim2, continuous_dimension_type2>
2014 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
2015 || (std::is_same_v<InterestDim3, continuous_dimension_type1>
2016 && std::is_same_v<InterestDim1, continuous_dimension_type2>
2017 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2018 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
2019 && std::is_same_v<InterestDim3, continuous_dimension_type2>
2020 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
2021
2022 return deriv_1_2_3(spline_eval, coords_eval, spline_coef);
2023 }
2024
2025 /**
2026 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
2027 *
2028 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
2029 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
2030 *
2031 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
2032 * This means that for each slice of spline_eval the evaluation is performed with
2033 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
2034 *
2035 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
2036 *
2037 * @tparam InterestDim1 First dimension along which differentiation is performed.
2038 * @tparam InterestDim2 Second dimension along which differentiation is performed.
2039 * @tparam InterestDim3 Third dimension along which differentiation is performed.
2040 *
2041 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
2042 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
2043 */
2044 template <
2045 class InterestDim1,
2046 class InterestDim2,
2047 class InterestDim3,
2048 class Layout1,
2049 class Layout2,
2050 class BatchedInterpolationDDom>
2051 void deriv3(
2052 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
2053 spline_eval,
2054 ddc::ChunkSpan<
2055 double const,
2056 batched_spline_domain_type<BatchedInterpolationDDom>,
2057 Layout2,
2058 memory_space> const spline_coef) const
2059 {
2060 static_assert(
2061 (std::is_same_v<InterestDim1, continuous_dimension_type1>
2062 && std::is_same_v<InterestDim2, continuous_dimension_type2>
2063 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
2064 || (std::is_same_v<InterestDim3, continuous_dimension_type1>
2065 && std::is_same_v<InterestDim1, continuous_dimension_type2>
2066 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2067 || (std::is_same_v<InterestDim2, continuous_dimension_type1>
2068 && std::is_same_v<InterestDim3, continuous_dimension_type2>
2069 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
2070
2071 return deriv_1_2_3(spline_eval, spline_coef);
2072 }
2073
2074 /** @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.
2075 *
2076 * 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.
2077 *
2078 * This is not a nD integration. This is a batched 3D integration.
2079 * This means that for each element of integrals, the integration is performed with the 3D set of
2080 * spline coefficients identified by the same DiscreteElement.
2081 *
2082 * @param[out] integrals The integrals of the 3D spline function on the subdomain of batch_domain. For practical reasons those are
2083 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
2084 * points represented by this domain are unused and irrelevant.
2085 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
2086 */
2087 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
2088 void integrate(
2089 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
2090 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
2091 spline_coef) const
2092 {
2093 static_assert(
2094 ddc::type_seq_contains_v<
2095 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>,
2096 to_type_seq_t<BatchedSplineDDom>>,
2097 "The spline coefficients domain must contain the bsplines dimensions");
2098 using batch_domain_type = ddc::
2099 remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2, bsplines_type3>;
2100 static_assert(
2101 std::is_same_v<batch_domain_type, BatchedDDom>,
2102 "The integrals domain must only contain the batch dimensions");
2103
2104 batch_domain_type batch_domain(integrals.domain());
2105 ddc::Chunk values1_alloc(
2106 ddc::DiscreteDomain<bsplines_type1>(spline_coef.domain()),
2107 ddc::KokkosAllocator<double, memory_space>());
2108 ddc::ChunkSpan values1 = values1_alloc.span_view();
2109 ddc::integrals(exec_space(), values1);
2110 ddc::Chunk values2_alloc(
2111 ddc::DiscreteDomain<bsplines_type2>(spline_coef.domain()),
2112 ddc::KokkosAllocator<double, memory_space>());
2113 ddc::ChunkSpan values2 = values2_alloc.span_view();
2114 ddc::integrals(exec_space(), values2);
2115 ddc::Chunk values3_alloc(
2116 ddc::DiscreteDomain<bsplines_type3>(spline_coef.domain()),
2117 ddc::KokkosAllocator<double, memory_space>());
2118 ddc::ChunkSpan values3 = values3_alloc.span_view();
2119 ddc::integrals(exec_space(), values3);
2120
2121 ddc::parallel_for_each(
2122 "ddc_splines_integrate_bsplines",
2123 exec_space(),
2124 batch_domain,
2125 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
2126 integrals(j) = 0;
2127 for (typename spline_domain_type1::discrete_element_type const i1 :
2128 values1.domain()) {
2129 for (typename spline_domain_type2::discrete_element_type const i2 :
2130 values2.domain()) {
2131 for (typename spline_domain_type3::discrete_element_type const i3 :
2132 values3.domain()) {
2133 integrals(j) += spline_coef(i1, i2, i3, j) * values1(i1)
2134 * values2(i2) * values3(i3);
2135 }
2136 }
2137 }
2138 });
2139 }
2140
2141private:
2142 /**
2143 * @brief Evaluate the function on B-splines at the coordinate given.
2144 *
2145 * This function firstly deals with the boundary conditions and calls the SplineEvaluator3D::eval_no_bc function
2146 * to evaluate.
2147 *
2148 * @param[in] coord_eval The 3D coordinate where we want to evaluate.
2149 * @param[in] spline_coef The B-splines coefficients of the function we want to evaluate.
2150 * @param[out] vals1 A ChunkSpan with the not-null values of each function of the spline in the first dimension.
2151 * @param[out] vals2 A ChunkSpan with the not-null values of each function of the spline in the second dimension.
2152 *
2153 * @return A double with the value of the function at the coordinate given.
2154 *
2155 * @see SplineBoundaryValue
2156 */
2157 template <class Layout, class... CoordsDims>
2158 KOKKOS_INLINE_FUNCTION double eval(
2159 ddc::Coordinate<CoordsDims...> coord_eval,
2160 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
2161 spline_coef) const
2162 {
2163 using Dim1 = continuous_dimension_type1;
2164 using Dim2 = continuous_dimension_type2;
2165 using Dim3 = continuous_dimension_type3;
2166 if constexpr (bsplines_type1::is_periodic()) {
2167 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()
2168 || ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
2169 ddc::get<Dim1>(coord_eval)
2170 -= Kokkos::floor(
2171 (ddc::get<Dim1>(coord_eval)
2172 - ddc::discrete_space<bsplines_type1>().rmin())
2173 / ddc::discrete_space<bsplines_type1>().length())
2174 * ddc::discrete_space<bsplines_type1>().length();
2175 }
2176 }
2177 if constexpr (bsplines_type2::is_periodic()) {
2178 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()
2179 || ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
2180 ddc::get<Dim2>(coord_eval)
2181 -= Kokkos::floor(
2182 (ddc::get<Dim2>(coord_eval)
2183 - ddc::discrete_space<bsplines_type2>().rmin())
2184 / ddc::discrete_space<bsplines_type2>().length())
2185 * ddc::discrete_space<bsplines_type2>().length();
2186 }
2187 }
2188 if constexpr (bsplines_type3::is_periodic()) {
2189 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()
2190 || ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
2191 ddc::get<Dim3>(coord_eval)
2192 -= Kokkos::floor(
2193 (ddc::get<Dim3>(coord_eval)
2194 - ddc::discrete_space<bsplines_type3>().rmin())
2195 / ddc::discrete_space<bsplines_type3>().length())
2196 * ddc::discrete_space<bsplines_type3>().length();
2197 }
2198 }
2199 if constexpr (!bsplines_type1::is_periodic()) {
2200 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()) {
2201 return m_lower_extrap_rule_1(coord_eval, spline_coef);
2202 }
2203 if (ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
2204 return m_upper_extrap_rule_1(coord_eval, spline_coef);
2205 }
2206 }
2207 if constexpr (!bsplines_type2::is_periodic()) {
2208 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()) {
2209 return m_lower_extrap_rule_2(coord_eval, spline_coef);
2210 }
2211 if (ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
2212 return m_upper_extrap_rule_2(coord_eval, spline_coef);
2213 }
2214 }
2215 if constexpr (!bsplines_type3::is_periodic()) {
2216 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()) {
2217 return m_lower_extrap_rule_3(coord_eval, spline_coef);
2218 }
2219 if (ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
2220 return m_upper_extrap_rule_3(coord_eval, spline_coef);
2221 }
2222 }
2223 return eval_no_bc<eval_type, eval_type, eval_type>(
2224 ddc::Coordinate<
2225 continuous_dimension_type1,
2226 continuous_dimension_type2,
2227 continuous_dimension_type3>(
2228 ddc::get<Dim1>(coord_eval),
2229 ddc::get<Dim2>(coord_eval),
2230 ddc::get<Dim3>(coord_eval)),
2231 spline_coef);
2232 }
2233
2234 /**
2235 * @brief Evaluate the function or its derivative at the coordinate given.
2236 *
2237 * @param[in] coord_eval The coordinate where we want to evaluate.
2238 * @param[in] splne_coef The B-splines coefficients of the function we want to evaluate.
2239 * @tparam EvalType1 A flag indicating if we evaluate the function or its derivative in the first dimension. The type of this object is either `eval_type` or `eval_deriv_type`.
2240 * @tparam EvalType2 A flag indicating if we evaluate the function or its derivative in the second dimension. The type of this object is either `eval_type` or `eval_deriv_type`.
2241 */
2242 template <class EvalType1, class EvalType2, class EvalType3, class Layout, class... CoordsDims>
2243 KOKKOS_INLINE_FUNCTION double eval_no_bc(
2244 ddc::Coordinate<CoordsDims...> const& coord_eval,
2245 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
2246 spline_coef) const
2247 {
2248 static_assert(
2249 std::is_same_v<EvalType1, eval_type> || std::is_same_v<EvalType1, eval_deriv_type>);
2250 static_assert(
2251 std::is_same_v<EvalType2, eval_type> || std::is_same_v<EvalType2, eval_deriv_type>);
2252 static_assert(
2253 std::is_same_v<EvalType3, eval_type> || std::is_same_v<EvalType3, eval_deriv_type>);
2254 ddc::DiscreteElement<bsplines_type1> jmin1;
2255 ddc::DiscreteElement<bsplines_type2> jmin2;
2256 ddc::DiscreteElement<bsplines_type3> jmin3;
2257
2258 std::array<double, bsplines_type1::degree() + 1> vals1_ptr;
2259 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>> const
2260 vals1(vals1_ptr.data());
2261 std::array<double, bsplines_type2::degree() + 1> vals2_ptr;
2262 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>> const
2263 vals2(vals2_ptr.data());
2264 std::array<double, bsplines_type3::degree() + 1> vals3_ptr;
2265 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type3::degree() + 1>> const
2266 vals3(vals3_ptr.data());
2267 ddc::Coordinate<continuous_dimension_type1> const coord_eval_interest1(coord_eval);
2268 ddc::Coordinate<continuous_dimension_type2> const coord_eval_interest2(coord_eval);
2269 ddc::Coordinate<continuous_dimension_type3> const coord_eval_interest3(coord_eval);
2270
2271 if constexpr (std::is_same_v<EvalType1, eval_type>) {
2272 jmin1 = ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
2273 } else if constexpr (std::is_same_v<EvalType1, eval_deriv_type>) {
2274 jmin1 = ddc::discrete_space<bsplines_type1>().eval_deriv(vals1, coord_eval_interest1);
2275 }
2276 if constexpr (std::is_same_v<EvalType2, eval_type>) {
2277 jmin2 = ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
2278 } else if constexpr (std::is_same_v<EvalType2, eval_deriv_type>) {
2279 jmin2 = ddc::discrete_space<bsplines_type2>().eval_deriv(vals2, coord_eval_interest2);
2280 }
2281 if constexpr (std::is_same_v<EvalType3, eval_type>) {
2282 jmin3 = ddc::discrete_space<bsplines_type3>().eval_basis(vals3, coord_eval_interest3);
2283 } else if constexpr (std::is_same_v<EvalType3, eval_deriv_type>) {
2284 jmin3 = ddc::discrete_space<bsplines_type3>().eval_deriv(vals3, coord_eval_interest3);
2285 }
2286
2287 double y = 0.0;
2288 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
2289 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
2290 for (std::size_t k = 0; k < bsplines_type3::degree() + 1; ++k) {
2291 y += spline_coef(
2292 ddc::DiscreteElement<
2293 bsplines_type1,
2294 bsplines_type2,
2295 bsplines_type3>(jmin1 + i, jmin2 + j, jmin3 + k))
2296 * vals1[i] * vals2[j] * vals3[k];
2297 }
2298 }
2299 }
2300 return y;
2301 }
2302};
2303
2304} // namespace ddc
friend class ChunkSpan
friend class Chunk
Definition chunk.hpp:81
friend class DiscreteDomain
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.
KOKKOS_FUNCTION double deriv_dim_2(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 s...
void deriv(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 spline function (described by its spline coefficients) on a mesh along a specified dime...
void deriv_dim_2(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 second dimens...
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.
void deriv_2_and_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the sec...
void deriv_dim_1(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 first dimensi...
upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
Get the upper extrapolation rule along the third dimension.
void deriv_1_2_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dim...
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.
void deriv_dim_2(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 second dimens...
KOKKOS_FUNCTION double deriv_dim_1(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 f...
void deriv(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 spline function (described by its spline coefficients) on a mesh along a specified dime...
KOKKOS_FUNCTION double deriv_1_and_3(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate a...
~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.
KOKKOS_FUNCTION double deriv(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along a...
KOKKOS_FUNCTION double deriv_dim_3(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...
KOKKOS_FUNCTION double deriv2(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Double-differentiate 3D spline function (described by its spline coefficients) at a given coordinate ...
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.
void deriv2(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
Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specif...
KOKKOS_FUNCTION double deriv_1_2_3(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate a...
void deriv_dim_3(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 third dimensi...
void deriv_dim_1(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 first dimensi...
KOKKOS_FUNCTION double deriv_1_and_2(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate a...
SplineEvaluator3D(SplineEvaluator3D &&x)=default
Move-constructs.
void deriv3(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 spline function (described by its spline coefficients) on a mesh along a specified dime...
void deriv_dim_3(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 third dimensi...
void deriv_2_and_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the sec...
void deriv_1_2_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dim...
void deriv2(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
Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specif...
void deriv3(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 spline function (described by its spline coefficients) on a mesh along specified dimens...
void deriv_1_and_2(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the fir...
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...
void deriv_1_and_2(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the fir...
void deriv_1_and_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the fir...
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 deriv_1_and_3(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
Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the fir...
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.
KOKKOS_FUNCTION double deriv3(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Triple-differentiate 3D spline function (described by its spline coefficients) at a given coordinate ...
KOKKOS_FUNCTION double deriv_2_and_3(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Cross-differentiate 3D spline function (described by its spline coefficients) at a given coordinate a...
The top-level namespace of DDC.