DDC 0.8.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<
761 InterestDim,
762 typename evaluation_discrete_dimension_type1::
763 continuous_dimension_type>) {
764 return deriv_dim_1(coord_eval, spline_coef);
765 } else if constexpr (std::is_same_v<
766 InterestDim,
767 typename evaluation_discrete_dimension_type2::
768 continuous_dimension_type>) {
769 return deriv_dim_2(coord_eval, spline_coef);
770 } else if constexpr (std::is_same_v<
771 InterestDim,
772 typename evaluation_discrete_dimension_type3::
773 continuous_dimension_type>) {
774 return deriv_dim_3(coord_eval, spline_coef);
775 }
776 }
777
778 /**
779 * @brief Double-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along specified dimensions of interest.
780 *
781 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
782 * obtained via various methods, such as using a SplineBuilder3D.
783 *
784 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
785 *
786 * @tparam InterestDim1 First dimension along which differentiation is performed.
787 * @tparam InterestDim2 Second dimension along which differentiation is performed.
788 *
789 * @param coord_eval The coordinate where the spline is double-differentiated. Note that only the components along the dimensions of interest are used.
790 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
791 *
792 * @return The derivative of the spline function at the desired coordinate.
793 */
794 template <class InterestDim1, class InterestDim2, class Layout, class... CoordsDims>
795 KOKKOS_FUNCTION double deriv2(
796 ddc::Coordinate<CoordsDims...> const& coord_eval,
797 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
798 spline_coef) const
799 {
800 static_assert(
801 (std::is_same_v<
802 InterestDim1,
803 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
804 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
805 || (std::is_same_v<
806 InterestDim2,
807 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
808 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
809 || (std::is_same_v<
810 InterestDim1,
811 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
812 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
813 || (std::is_same_v<
814 InterestDim2,
815 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
816 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
817 || (std::is_same_v<
818 InterestDim1,
819 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
820 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
821 || (std::is_same_v<
822 InterestDim2,
823 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
824 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
825
826 if constexpr (
827 (std::is_same_v<
828 InterestDim1,
829 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
830 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
831 || (std::is_same_v<
832 InterestDim2,
833 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
834 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
835 return deriv_1_and_2(coord_eval, spline_coef);
836 } else if constexpr (
837 (std::is_same_v<
838 InterestDim1,
839 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
840 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
841 || (std::is_same_v<
842 InterestDim2,
843 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
844 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
845 return deriv_2_and_3(coord_eval, spline_coef);
846 } else if constexpr (
847 (std::is_same_v<
848 InterestDim1,
849 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
850 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
851 || (std::is_same_v<
852 InterestDim2,
853 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
854 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
855 return deriv_1_and_3(coord_eval, spline_coef);
856 }
857 }
858
859 /**
860 * @brief Triple-differentiate 3D spline function (described by its spline coefficients) at a given coordinate along specified dimensions of interest.
861 *
862 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
863 * obtained via various methods, such as using a SplineBuilder3D.
864 *
865 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
866 *
867 * @tparam InterestDim1 First dimension along which differentiation is performed.
868 * @tparam InterestDim2 Second dimension along which differentiation is performed.
869 * @tparam InterestDim3 Third dimension along which differentiation is performed.
870 *
871 * @param coord_eval The coordinate where the spline is triple-differentiated. Note that only the components along the dimensions of interest are used.
872 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
873 *
874 * @return The derivative of the spline function at the desired coordinate.
875 */
876 template <
877 class InterestDim1,
878 class InterestDim2,
879 class InterestDim3,
880 class Layout,
881 class... CoordsDims>
882 KOKKOS_FUNCTION double deriv3(
883 ddc::Coordinate<CoordsDims...> const& coord_eval,
884 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
885 spline_coef) const
886 {
887 static_assert(
888 (std::is_same_v<
889 InterestDim1,
890 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
891 && std::is_same_v<InterestDim2, continuous_dimension_type2>
892 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
893 || (std::is_same_v<
894 InterestDim3,
895 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
896 && std::is_same_v<InterestDim1, continuous_dimension_type2>
897 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
898 || (std::is_same_v<
899 InterestDim2,
900 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
901 && std::is_same_v<InterestDim3, continuous_dimension_type2>
902 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
903
904 return deriv_1_2_3(coord_eval, spline_coef);
905 }
906
907 /**
908 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
909 *
910 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
911 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
912 *
913 * This is not a nD evaluation. This is a batched 3D differentiation.
914 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
915 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
916 *
917 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
918 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
919 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
920 * the set of 3D spline coefficients retained to perform the evaluation).
921 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
922 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
923 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
924 * the set of 3D spline coefficients retained to perform the evaluation).
925 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
926 */
927 template <
928 class Layout1,
929 class Layout2,
930 class Layout3,
931 class BatchedInterpolationDDom,
932 class... CoordsDims>
933 void deriv_dim_1(
934 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
935 spline_eval,
936 ddc::ChunkSpan<
937 ddc::Coordinate<CoordsDims...> const,
938 BatchedInterpolationDDom,
939 Layout2,
940 memory_space> const coords_eval,
941 ddc::ChunkSpan<
942 double const,
943 batched_spline_domain_type<BatchedInterpolationDDom>,
944 Layout3,
945 memory_space> const spline_coef) const
946 {
947 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
948 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
949 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
950 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
951 ddc::parallel_for_each(
952 "ddc_splines_differentiate_3d_dim_1",
953 exec_space(),
954 batch_domain,
955 KOKKOS_CLASS_LAMBDA(
956 typename batch_domain_type<
957 BatchedInterpolationDDom>::discrete_element_type const j) {
958 auto const spline_eval_3D = spline_eval[j];
959 auto const coords_eval_3D = coords_eval[j];
960 auto const spline_coef_3D = spline_coef[j];
961 for (auto const i1 : evaluation_domain1) {
962 for (auto const i2 : evaluation_domain2) {
963 for (auto const i3 : evaluation_domain3) {
964 spline_eval_3D(i1, i2, i3) = eval_no_bc<
965 eval_deriv_type,
966 eval_type,
967 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
968 }
969 }
970 }
971 });
972 }
973
974 /**
975 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
976 *
977 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
978 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
979 *
980 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
981 * This means that for each slice of spline_eval the evaluation is performed with
982 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
983 *
984 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
985 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
986 */
987 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
988 void deriv_dim_1(
989 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
990 spline_eval,
991 ddc::ChunkSpan<
992 double const,
993 batched_spline_domain_type<BatchedInterpolationDDom>,
994 Layout2,
995 memory_space> const spline_coef) const
996 {
997 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
998 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
999 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1000 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1001 ddc::parallel_for_each(
1002 "ddc_splines_differentiate_3d_dim_1",
1003 exec_space(),
1004 batch_domain,
1005 KOKKOS_CLASS_LAMBDA(
1006 typename batch_domain_type<
1007 BatchedInterpolationDDom>::discrete_element_type const j) {
1008 auto const spline_eval_3D = spline_eval[j];
1009 auto const spline_coef_3D = spline_coef[j];
1010 for (auto const i1 : evaluation_domain1) {
1011 for (auto const i2 : evaluation_domain2) {
1012 for (auto const i3 : evaluation_domain3) {
1013 ddc::Coordinate<
1014 continuous_dimension_type1,
1015 continuous_dimension_type2,
1016 continuous_dimension_type3>
1017 coord_eval_3D(
1018 ddc::coordinate(i1),
1019 ddc::coordinate(i2),
1020 ddc::coordinate(i3));
1021 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1022 eval_deriv_type,
1023 eval_type,
1024 eval_type>(coord_eval_3D, spline_coef_3D);
1025 }
1026 }
1027 }
1028 });
1029 }
1030
1031 /**
1032 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
1033 *
1034 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1035 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1036 *
1037 * This is not a nD differentiation. This is a batched 3D differentiation.
1038 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1039 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1040 *
1041 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1042 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1043 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1044 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1045 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1046 * the set of 3D spline coefficients retained to perform the evaluation).
1047 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1048 */
1049 template <
1050 class Layout1,
1051 class Layout2,
1052 class Layout3,
1053 class BatchedInterpolationDDom,
1054 class... CoordsDims>
1055 void deriv_dim_2(
1056 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1057 spline_eval,
1058 ddc::ChunkSpan<
1059 ddc::Coordinate<CoordsDims...> const,
1060 BatchedInterpolationDDom,
1061 Layout2,
1062 memory_space> const coords_eval,
1063 ddc::ChunkSpan<
1064 double const,
1065 batched_spline_domain_type<BatchedInterpolationDDom>,
1066 Layout3,
1067 memory_space> const spline_coef) const
1068 {
1069 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1070 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1071 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1072 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1073 ddc::parallel_for_each(
1074 "ddc_splines_differentiate_3d_dim_2",
1075 exec_space(),
1076 batch_domain,
1077 KOKKOS_CLASS_LAMBDA(
1078 typename batch_domain_type<
1079 BatchedInterpolationDDom>::discrete_element_type const j) {
1080 auto const spline_eval_3D = spline_eval[j];
1081 auto const coords_eval_3D = coords_eval[j];
1082 auto const spline_coef_3D = spline_coef[j];
1083 for (auto const i1 : evaluation_domain1) {
1084 for (auto const i2 : evaluation_domain2) {
1085 for (auto const i3 : evaluation_domain3) {
1086 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1087 eval_type,
1088 eval_deriv_type,
1089 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
1090 }
1091 }
1092 }
1093 });
1094 }
1095
1096 /**
1097 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
1098 *
1099 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1100 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1101 *
1102 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1103 * This means that for each slice of spline_eval the evaluation is performed with
1104 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1105 *
1106 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1107 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1108 */
1109 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1110 void deriv_dim_2(
1111 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1112 spline_eval,
1113 ddc::ChunkSpan<
1114 double const,
1115 batched_spline_domain_type<BatchedInterpolationDDom>,
1116 Layout2,
1117 memory_space> const spline_coef) const
1118 {
1119 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1120 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1121 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1122 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1123 ddc::parallel_for_each(
1124 "ddc_splines_differentiate_3d_dim_2",
1125 exec_space(),
1126 batch_domain,
1127 KOKKOS_CLASS_LAMBDA(
1128 typename batch_domain_type<
1129 BatchedInterpolationDDom>::discrete_element_type const j) {
1130 auto const spline_eval_3D = spline_eval[j];
1131 auto const spline_coef_3D = spline_coef[j];
1132 for (auto const i1 : evaluation_domain1) {
1133 for (auto const i2 : evaluation_domain2) {
1134 for (auto const i3 : evaluation_domain3) {
1135 ddc::Coordinate<
1136 continuous_dimension_type1,
1137 continuous_dimension_type2,
1138 continuous_dimension_type3>
1139 coord_eval_3D(
1140 ddc::coordinate(i1),
1141 ddc::coordinate(i2),
1142 ddc::coordinate(i3));
1143 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1144 eval_type,
1145 eval_deriv_type,
1146 eval_type>(coord_eval_3D, spline_coef_3D);
1147 }
1148 }
1149 }
1150 });
1151 }
1152
1153 /**
1154 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along third dimension of interest.
1155 *
1156 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1157 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1158 *
1159 * This is not a nD differentiation. This is a batched 3D differentiation.
1160 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1161 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1162 *
1163 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1164 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1165 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1166 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1167 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1168 * the set of 3D spline coefficients retained to perform the evaluation).
1169 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1170 */
1171 template <
1172 class Layout1,
1173 class Layout2,
1174 class Layout3,
1175 class BatchedInterpolationDDom,
1176 class... CoordsDims>
1177 void deriv_dim_3(
1178 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1179 spline_eval,
1180 ddc::ChunkSpan<
1181 ddc::Coordinate<CoordsDims...> const,
1182 BatchedInterpolationDDom,
1183 Layout2,
1184 memory_space> const coords_eval,
1185 ddc::ChunkSpan<
1186 double const,
1187 batched_spline_domain_type<BatchedInterpolationDDom>,
1188 Layout3,
1189 memory_space> const spline_coef) const
1190 {
1191 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1192 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1193 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1194 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1195 ddc::parallel_for_each(
1196 "ddc_splines_differentiate_3d_dim_3",
1197 exec_space(),
1198 batch_domain,
1199 KOKKOS_CLASS_LAMBDA(
1200 typename batch_domain_type<
1201 BatchedInterpolationDDom>::discrete_element_type const j) {
1202 auto const spline_eval_3D = spline_eval[j];
1203 auto const coords_eval_3D = coords_eval[j];
1204 auto const spline_coef_3D = spline_coef[j];
1205 for (auto const i1 : evaluation_domain1) {
1206 for (auto const i2 : evaluation_domain2) {
1207 for (auto const i3 : evaluation_domain3) {
1208 spline_eval_3D(i1, i2, i3)
1209 = eval_no_bc<eval_type, eval_type, eval_deriv_type>(
1210 coords_eval_3D(i1, i2, i3),
1211 spline_coef_3D);
1212 }
1213 }
1214 }
1215 });
1216 }
1217
1218 /**
1219 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along third dimension of interest.
1220 *
1221 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1222 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1223 *
1224 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1225 * This means that for each slice of spline_eval the evaluation is performed with
1226 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1227 *
1228 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1229 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1230 */
1231 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1232 void deriv_dim_3(
1233 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1234 spline_eval,
1235 ddc::ChunkSpan<
1236 double const,
1237 batched_spline_domain_type<BatchedInterpolationDDom>,
1238 Layout2,
1239 memory_space> const spline_coef) const
1240 {
1241 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1242 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1243 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1244 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1245 ddc::parallel_for_each(
1246 "ddc_splines_differentiate_3d_dim_3",
1247 exec_space(),
1248 batch_domain,
1249 KOKKOS_CLASS_LAMBDA(
1250 typename batch_domain_type<
1251 BatchedInterpolationDDom>::discrete_element_type const j) {
1252 auto const spline_eval_3D = spline_eval[j];
1253 auto const spline_coef_3D = spline_coef[j];
1254 for (auto const i1 : evaluation_domain1) {
1255 for (auto const i2 : evaluation_domain2) {
1256 for (auto const i3 : evaluation_domain3) {
1257 ddc::Coordinate<
1258 continuous_dimension_type1,
1259 continuous_dimension_type2,
1260 continuous_dimension_type3>
1261 coord_eval_3D(
1262 ddc::coordinate(i1),
1263 ddc::coordinate(i2),
1264 ddc::coordinate(i3));
1265 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1266 eval_type,
1267 eval_type,
1268 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1269 }
1270 }
1271 }
1272 });
1273 }
1274
1275 /**
1276 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and second dimensions of interest.
1277 *
1278 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1279 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1280 *
1281 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1282 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1283 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1284 *
1285 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1286 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1287 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1288 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1289 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1290 * the set of 3D spline coefficients retained to perform the evaluation).
1291 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1292 */
1293 template <
1294 class Layout1,
1295 class Layout2,
1296 class Layout3,
1297 class BatchedInterpolationDDom,
1298 class... CoordsDims>
1299 void deriv_1_and_2(
1300 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1301 spline_eval,
1302 ddc::ChunkSpan<
1303 ddc::Coordinate<CoordsDims...> const,
1304 BatchedInterpolationDDom,
1305 Layout2,
1306 memory_space> const coords_eval,
1307 ddc::ChunkSpan<
1308 double const,
1309 batched_spline_domain_type<BatchedInterpolationDDom>,
1310 Layout3,
1311 memory_space> const spline_coef) const
1312 {
1313 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1314 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1315 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1316 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1317 ddc::parallel_for_each(
1318 "ddc_splines_cross_differentiate_3d_dim1_2",
1319 exec_space(),
1320 batch_domain,
1321 KOKKOS_CLASS_LAMBDA(
1322 typename batch_domain_type<
1323 BatchedInterpolationDDom>::discrete_element_type const j) {
1324 auto const spline_eval_3D = spline_eval[j];
1325 auto const coords_eval_3D = coords_eval[j];
1326 auto const spline_coef_3D = spline_coef[j];
1327 for (auto const i1 : evaluation_domain1) {
1328 for (auto const i2 : evaluation_domain2) {
1329 for (auto const i3 : evaluation_domain3) {
1330 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1331 eval_deriv_type,
1332 eval_deriv_type,
1333 eval_type>(coords_eval_3D(i1, i2, i3), spline_coef_3D);
1334 }
1335 }
1336 }
1337 });
1338 }
1339
1340 /**
1341 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and second dimensions of interest.
1342 *
1343 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1344 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1345 *
1346 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1347 * This means that for each slice of spline_eval the evaluation is performed with
1348 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1349 *
1350 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1351 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1352 */
1353 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1354 void deriv_1_and_2(
1355 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1356 spline_eval,
1357 ddc::ChunkSpan<
1358 double const,
1359 batched_spline_domain_type<BatchedInterpolationDDom>,
1360 Layout2,
1361 memory_space> const spline_coef) const
1362 {
1363 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1364 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1365 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1366 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1367 ddc::parallel_for_each(
1368 "ddc_splines_cross_differentiate_3d_dim1_2",
1369 exec_space(),
1370 batch_domain,
1371 KOKKOS_CLASS_LAMBDA(
1372 typename batch_domain_type<
1373 BatchedInterpolationDDom>::discrete_element_type const j) {
1374 auto const spline_eval_3D = spline_eval[j];
1375 auto const spline_coef_3D = spline_coef[j];
1376 for (auto const i1 : evaluation_domain1) {
1377 for (auto const i2 : evaluation_domain2) {
1378 for (auto const i3 : evaluation_domain3) {
1379 ddc::Coordinate<
1380 continuous_dimension_type1,
1381 continuous_dimension_type2,
1382 continuous_dimension_type3>
1383 coord_eval_3D(
1384 ddc::coordinate(i1),
1385 ddc::coordinate(i2),
1386 ddc::coordinate(i3));
1387 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1388 eval_deriv_type,
1389 eval_deriv_type,
1390 eval_type>(coord_eval_3D, spline_coef_3D);
1391 }
1392 }
1393 }
1394 });
1395 }
1396
1397 /**
1398 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the second and third dimensions of interest.
1399 *
1400 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1401 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1402 *
1403 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1404 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1405 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1406 *
1407 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1408 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1409 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1410 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1411 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1412 * the set of 3D spline coefficients retained to perform the evaluation).
1413 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1414 */
1415 template <
1416 class Layout1,
1417 class Layout2,
1418 class Layout3,
1419 class BatchedInterpolationDDom,
1420 class... CoordsDims>
1421 void deriv_2_and_3(
1422 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1423 spline_eval,
1424 ddc::ChunkSpan<
1425 ddc::Coordinate<CoordsDims...> const,
1426 BatchedInterpolationDDom,
1427 Layout2,
1428 memory_space> const coords_eval,
1429 ddc::ChunkSpan<
1430 double const,
1431 batched_spline_domain_type<BatchedInterpolationDDom>,
1432 Layout3,
1433 memory_space> const spline_coef) const
1434 {
1435 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1436 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1437 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1438 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1439 ddc::parallel_for_each(
1440 "ddc_splines_cross_differentiate_3d_dim2_3",
1441 exec_space(),
1442 batch_domain,
1443 KOKKOS_CLASS_LAMBDA(
1444 typename batch_domain_type<
1445 BatchedInterpolationDDom>::discrete_element_type const j) {
1446 auto const spline_eval_3D = spline_eval[j];
1447 auto const coords_eval_3D = coords_eval[j];
1448 auto const spline_coef_3D = spline_coef[j];
1449 for (auto const i1 : evaluation_domain1) {
1450 for (auto const i2 : evaluation_domain2) {
1451 for (auto const i3 : evaluation_domain3) {
1452 spline_eval_3D(i1, i2, i3)
1453 = eval_no_bc<eval_type, eval_deriv_type, eval_deriv_type>(
1454 coords_eval_3D(i1, i2, i3),
1455 spline_coef_3D);
1456 }
1457 }
1458 }
1459 });
1460 }
1461
1462 /**
1463 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the second and third dimensions of interest.
1464 *
1465 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1466 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1467 *
1468 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1469 * This means that for each slice of spline_eval the evaluation is performed with
1470 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1471 *
1472 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1473 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1474 */
1475 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1476 void deriv_2_and_3(
1477 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1478 spline_eval,
1479 ddc::ChunkSpan<
1480 double const,
1481 batched_spline_domain_type<BatchedInterpolationDDom>,
1482 Layout2,
1483 memory_space> const spline_coef) const
1484 {
1485 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1486 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1487 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1488 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1489 ddc::parallel_for_each(
1490 "ddc_splines_cross_differentiate_3d_dim2_3",
1491 exec_space(),
1492 batch_domain,
1493 KOKKOS_CLASS_LAMBDA(
1494 typename batch_domain_type<
1495 BatchedInterpolationDDom>::discrete_element_type const j) {
1496 auto const spline_eval_3D = spline_eval[j];
1497 auto const spline_coef_3D = spline_coef[j];
1498 for (auto const i1 : evaluation_domain1) {
1499 for (auto const i2 : evaluation_domain2) {
1500 for (auto const i3 : evaluation_domain3) {
1501 ddc::Coordinate<
1502 continuous_dimension_type1,
1503 continuous_dimension_type2,
1504 continuous_dimension_type3>
1505 coord_eval_3D(
1506 ddc::coordinate(i1),
1507 ddc::coordinate(i2),
1508 ddc::coordinate(i3));
1509 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1510 eval_type,
1511 eval_deriv_type,
1512 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1513 }
1514 }
1515 }
1516 });
1517 }
1518
1519 /**
1520 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and third dimensions of interest.
1521 *
1522 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1523 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1524 *
1525 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1526 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1527 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1528 *
1529 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1530 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1531 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1532 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1533 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1534 * the set of 3D spline coefficients retained to perform the evaluation).
1535 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1536 */
1537 template <
1538 class Layout1,
1539 class Layout2,
1540 class Layout3,
1541 class BatchedInterpolationDDom,
1542 class... CoordsDims>
1543 void deriv_1_and_3(
1544 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1545 spline_eval,
1546 ddc::ChunkSpan<
1547 ddc::Coordinate<CoordsDims...> const,
1548 BatchedInterpolationDDom,
1549 Layout2,
1550 memory_space> const coords_eval,
1551 ddc::ChunkSpan<
1552 double const,
1553 batched_spline_domain_type<BatchedInterpolationDDom>,
1554 Layout3,
1555 memory_space> const spline_coef) const
1556 {
1557 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1558 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1559 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1560 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1561 ddc::parallel_for_each(
1562 "ddc_splines_cross_differentiate_3d_dim1_3",
1563 exec_space(),
1564 batch_domain,
1565 KOKKOS_CLASS_LAMBDA(
1566 typename batch_domain_type<
1567 BatchedInterpolationDDom>::discrete_element_type const j) {
1568 auto const spline_eval_3D = spline_eval[j];
1569 auto const coords_eval_3D = coords_eval[j];
1570 auto const spline_coef_3D = spline_coef[j];
1571 for (auto const i1 : evaluation_domain1) {
1572 for (auto const i2 : evaluation_domain2) {
1573 for (auto const i3 : evaluation_domain3) {
1574 spline_eval_3D(i1, i2, i3)
1575 = eval_no_bc<eval_deriv_type, eval_type, eval_deriv_type>(
1576 coords_eval_3D(i1, i2, i3),
1577 spline_coef_3D);
1578 }
1579 }
1580 }
1581 });
1582 }
1583
1584 /**
1585 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the first and third dimensions of interest.
1586 *
1587 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1588 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1589 *
1590 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1591 * This means that for each slice of spline_eval the evaluation is performed with
1592 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1593 *
1594 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1595 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1596 */
1597 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1598 void deriv_1_and_3(
1599 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1600 spline_eval,
1601 ddc::ChunkSpan<
1602 double const,
1603 batched_spline_domain_type<BatchedInterpolationDDom>,
1604 Layout2,
1605 memory_space> const spline_coef) const
1606 {
1607 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1608 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1609 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1610 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1611 ddc::parallel_for_each(
1612 "ddc_splines_cross_differentiate_3d_dim1_3",
1613 exec_space(),
1614 batch_domain,
1615 KOKKOS_CLASS_LAMBDA(
1616 typename batch_domain_type<
1617 BatchedInterpolationDDom>::discrete_element_type const j) {
1618 auto const spline_eval_3D = spline_eval[j];
1619 auto const spline_coef_3D = spline_coef[j];
1620 for (auto const i1 : evaluation_domain1) {
1621 for (auto const i2 : evaluation_domain2) {
1622 for (auto const i3 : evaluation_domain3) {
1623 ddc::Coordinate<
1624 continuous_dimension_type1,
1625 continuous_dimension_type2,
1626 continuous_dimension_type3>
1627 coord_eval_3D(
1628 ddc::coordinate(i1),
1629 ddc::coordinate(i2),
1630 ddc::coordinate(i3));
1631 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1632 eval_deriv_type,
1633 eval_type,
1634 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1635 }
1636 }
1637 }
1638 });
1639 }
1640
1641 /**
1642 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
1643 *
1644 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1645 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1646 *
1647 * This is not a nD cross-differentiation. This is a batched 3D cross-differentiation.
1648 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1649 * the cross-differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1650 *
1651 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1652 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1653 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1654 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1655 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1656 * the set of 3D spline coefficients retained to perform the evaluation).
1657 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1658 */
1659 template <
1660 class Layout1,
1661 class Layout2,
1662 class Layout3,
1663 class BatchedInterpolationDDom,
1664 class... CoordsDims>
1665 void deriv_1_2_3(
1666 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1667 spline_eval,
1668 ddc::ChunkSpan<
1669 ddc::Coordinate<CoordsDims...> const,
1670 BatchedInterpolationDDom,
1671 Layout2,
1672 memory_space> const coords_eval,
1673 ddc::ChunkSpan<
1674 double const,
1675 batched_spline_domain_type<BatchedInterpolationDDom>,
1676 Layout3,
1677 memory_space> const spline_coef) const
1678 {
1679 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
1680 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1681 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1682 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1683 ddc::parallel_for_each(
1684 "ddc_splines_cross_differentiate_3d_dim1_2_3",
1685 exec_space(),
1686 batch_domain,
1687 KOKKOS_CLASS_LAMBDA(
1688 typename batch_domain_type<
1689 BatchedInterpolationDDom>::discrete_element_type const j) {
1690 auto const spline_eval_3D = spline_eval[j];
1691 auto const coords_eval_3D = coords_eval[j];
1692 auto const spline_coef_3D = spline_coef[j];
1693 for (auto const i1 : evaluation_domain1) {
1694 for (auto const i2 : evaluation_domain2) {
1695 for (auto const i3 : evaluation_domain3) {
1696 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1697 eval_deriv_type,
1698 eval_deriv_type,
1699 eval_deriv_type>(
1700 coords_eval_3D(i1, i2, i3),
1701 spline_coef_3D);
1702 }
1703 }
1704 }
1705 });
1706 }
1707
1708 /**
1709 * @brief Cross-differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
1710 *
1711 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1712 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1713 *
1714 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1715 * This means that for each slice of spline_eval the evaluation is performed with
1716 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1717 *
1718 * @param[out] spline_eval The cross-derivatives of the 3D spline function at the desired coordinates.
1719 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1720 */
1721 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
1722 void deriv_1_2_3(
1723 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1724 spline_eval,
1725 ddc::ChunkSpan<
1726 double const,
1727 batched_spline_domain_type<BatchedInterpolationDDom>,
1728 Layout2,
1729 memory_space> const spline_coef) const
1730 {
1731 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
1732 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
1733 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
1734 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
1735 ddc::parallel_for_each(
1736 "ddc_splines_cross_differentiate_3d_dim1_2_3",
1737 exec_space(),
1738 batch_domain,
1739 KOKKOS_CLASS_LAMBDA(
1740 typename batch_domain_type<
1741 BatchedInterpolationDDom>::discrete_element_type const j) {
1742 auto const spline_eval_3D = spline_eval[j];
1743 auto const spline_coef_3D = spline_coef[j];
1744 for (auto const i1 : evaluation_domain1) {
1745 for (auto const i2 : evaluation_domain2) {
1746 for (auto const i3 : evaluation_domain3) {
1747 ddc::Coordinate<
1748 continuous_dimension_type1,
1749 continuous_dimension_type2,
1750 continuous_dimension_type3>
1751 coord_eval_3D(
1752 ddc::coordinate(i1),
1753 ddc::coordinate(i2),
1754 ddc::coordinate(i3));
1755 spline_eval_3D(i1, i2, i3) = eval_no_bc<
1756 eval_deriv_type,
1757 eval_deriv_type,
1758 eval_deriv_type>(coord_eval_3D, spline_coef_3D);
1759 }
1760 }
1761 }
1762 });
1763 }
1764
1765 /**
1766 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
1767 *
1768 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1769 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1770 *
1771 * This is not a nD evaluation. This is a batched 3D differentiation.
1772 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1773 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1774 *
1775 * @tparam InterestDim Dimension along which differentiation is performed.
1776 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1777 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1778 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1779 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1780 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1781 * the set of 3D spline coefficients retained to perform the evaluation).
1782 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1783 */
1784 template <
1785 class InterestDim,
1786 class Layout1,
1787 class Layout2,
1788 class Layout3,
1789 class BatchedInterpolationDDom,
1790 class... CoordsDims>
1791 void deriv(
1792 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1793 spline_eval,
1794 ddc::ChunkSpan<
1795 ddc::Coordinate<CoordsDims...> const,
1796 BatchedInterpolationDDom,
1797 Layout2,
1798 memory_space> const coords_eval,
1799 ddc::ChunkSpan<
1800 double const,
1801 batched_spline_domain_type<BatchedInterpolationDDom>,
1802 Layout3,
1803 memory_space> const spline_coef) const
1804 {
1805 static_assert(
1806 std::is_same_v<
1807 InterestDim,
1808 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1809 || std::is_same_v<InterestDim, continuous_dimension_type2>
1810 || std::is_same_v<InterestDim, continuous_dimension_type3>);
1811 if constexpr (std::is_same_v<
1812 InterestDim,
1813 typename evaluation_discrete_dimension_type1::
1814 continuous_dimension_type>) {
1815 return deriv_dim_1(spline_eval, coords_eval, spline_coef);
1816 } else if constexpr (std::is_same_v<
1817 InterestDim,
1818 typename evaluation_discrete_dimension_type2::
1819 continuous_dimension_type>) {
1820 return deriv_dim_2(spline_eval, coords_eval, spline_coef);
1821 } else if constexpr (std::is_same_v<
1822 InterestDim,
1823 typename evaluation_discrete_dimension_type3::
1824 continuous_dimension_type>) {
1825 return deriv_dim_3(spline_eval, coords_eval, spline_coef);
1826 }
1827 }
1828
1829 /**
1830 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
1831 *
1832 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1833 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1834 *
1835 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1836 * This means that for each slice of spline_eval the evaluation is performed with
1837 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1838 *
1839 * @tparam InterestDim Dimension along which differentiation is performed.
1840 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1841 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1842 */
1843 template <class InterestDim, class Layout1, class Layout2, class BatchedInterpolationDDom>
1844 void deriv(
1845 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1846 spline_eval,
1847 ddc::ChunkSpan<
1848 double const,
1849 batched_spline_domain_type<BatchedInterpolationDDom>,
1850 Layout2,
1851 memory_space> const spline_coef) const
1852 {
1853 static_assert(
1854 std::is_same_v<
1855 InterestDim,
1856 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1857 || std::is_same_v<InterestDim, continuous_dimension_type2>
1858 || std::is_same_v<InterestDim, continuous_dimension_type3>);
1859 if constexpr (std::is_same_v<
1860 InterestDim,
1861 typename evaluation_discrete_dimension_type1::
1862 continuous_dimension_type>) {
1863 return deriv_dim_1(spline_eval, spline_coef);
1864 } else if constexpr (std::is_same_v<
1865 InterestDim,
1866 typename evaluation_discrete_dimension_type2::
1867 continuous_dimension_type>) {
1868 return deriv_dim_2(spline_eval, spline_coef);
1869 } else if constexpr (std::is_same_v<
1870 InterestDim,
1871 typename evaluation_discrete_dimension_type3::
1872 continuous_dimension_type>) {
1873 return deriv_dim_3(spline_eval, spline_coef);
1874 }
1875 }
1876
1877 /**
1878 * @brief Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1879 *
1880 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1881 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1882 *
1883 * This is not a nD evaluation. This is a batched 3D differentiation.
1884 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1885 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1886 *
1887 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1888 *
1889 * @tparam InterestDim1 First dimension along which differentiation is performed.
1890 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1891 *
1892 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
1893 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1894 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1895 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1896 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1897 * the set of 3D spline coefficients retained to perform the evaluation).
1898 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1899 */
1900 template <
1901 class InterestDim1,
1902 class InterestDim2,
1903 class Layout1,
1904 class Layout2,
1905 class Layout3,
1906 class BatchedInterpolationDDom,
1907 class... CoordsDims>
1908 void deriv2(
1909 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1910 spline_eval,
1911 ddc::ChunkSpan<
1912 ddc::Coordinate<CoordsDims...> const,
1913 BatchedInterpolationDDom,
1914 Layout2,
1915 memory_space> const coords_eval,
1916 ddc::ChunkSpan<
1917 double const,
1918 batched_spline_domain_type<BatchedInterpolationDDom>,
1919 Layout3,
1920 memory_space> const spline_coef) const
1921 {
1922 static_assert(
1923 (std::is_same_v<
1924 InterestDim1,
1925 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1926 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1927 || (std::is_same_v<
1928 InterestDim2,
1929 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1930 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
1931 || (std::is_same_v<
1932 InterestDim1,
1933 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
1934 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1935 || (std::is_same_v<
1936 InterestDim2,
1937 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
1938 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
1939 || (std::is_same_v<
1940 InterestDim1,
1941 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1942 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1943 || (std::is_same_v<
1944 InterestDim2,
1945 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1946 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
1947
1948 if constexpr (
1949 (std::is_same_v<
1950 InterestDim1,
1951 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1952 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1953 || (std::is_same_v<
1954 InterestDim2,
1955 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1956 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
1957 return deriv_1_and_2(spline_eval, coords_eval, spline_coef);
1958 } else if constexpr (
1959 (std::is_same_v<
1960 InterestDim1,
1961 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
1962 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1963 || (std::is_same_v<
1964 InterestDim2,
1965 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
1966 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1967 return deriv_2_and_3(spline_eval, coords_eval, spline_coef);
1968 } else if constexpr (
1969 (std::is_same_v<
1970 InterestDim1,
1971 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1972 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
1973 || (std::is_same_v<
1974 InterestDim2,
1975 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1976 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
1977 return deriv_1_and_3(spline_eval, coords_eval, spline_coef);
1978 }
1979 }
1980
1981 /**
1982 * @brief Double-differentiate 3D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1983 *
1984 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
1985 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
1986 *
1987 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
1988 * This means that for each slice of spline_eval the evaluation is performed with
1989 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1990 *
1991 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1992 *
1993 * @tparam InterestDim1 First dimension along which differentiation is performed.
1994 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1995 *
1996 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
1997 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
1998 */
1999 template <
2000 class InterestDim1,
2001 class InterestDim2,
2002 class Layout1,
2003 class Layout2,
2004 class BatchedInterpolationDDom>
2005 void deriv2(
2006 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
2007 spline_eval,
2008 ddc::ChunkSpan<
2009 double const,
2010 batched_spline_domain_type<BatchedInterpolationDDom>,
2011 Layout2,
2012 memory_space> const spline_coef) const
2013 {
2014 static_assert(
2015 (std::is_same_v<
2016 InterestDim1,
2017 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2018 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
2019 || (std::is_same_v<
2020 InterestDim2,
2021 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2022 && std::is_same_v<InterestDim1, continuous_dimension_type2>)
2023 || (std::is_same_v<
2024 InterestDim1,
2025 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
2026 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2027 || (std::is_same_v<
2028 InterestDim2,
2029 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
2030 && std::is_same_v<InterestDim1, continuous_dimension_type3>)
2031 || (std::is_same_v<
2032 InterestDim1,
2033 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2034 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2035 || (std::is_same_v<
2036 InterestDim2,
2037 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2038 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
2039
2040 if constexpr (
2041 (std::is_same_v<
2042 InterestDim1,
2043 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2044 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
2045 || (std::is_same_v<
2046 InterestDim2,
2047 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2048 && std::is_same_v<InterestDim1, continuous_dimension_type2>)) {
2049 return deriv_1_and_2(spline_eval, spline_coef);
2050 } else if constexpr (
2051 (std::is_same_v<
2052 InterestDim1,
2053 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
2054 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2055 || (std::is_same_v<
2056 InterestDim2,
2057 typename evaluation_discrete_dimension_type2::continuous_dimension_type>
2058 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
2059 return deriv_2_and_3(spline_eval, spline_coef);
2060 } else if constexpr (
2061 (std::is_same_v<
2062 InterestDim1,
2063 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2064 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2065 || (std::is_same_v<
2066 InterestDim2,
2067 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2068 && std::is_same_v<InterestDim1, continuous_dimension_type3>)) {
2069 return deriv_1_and_3(spline_eval, spline_coef);
2070 }
2071 }
2072
2073 /**
2074 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimensions of interest.
2075 *
2076 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
2077 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
2078 *
2079 * This is not a nD evaluation. This is a batched 3D differentiation.
2080 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
2081 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
2082 *
2083 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
2084 *
2085 * @tparam InterestDim1 First dimension along which differentiation is performed.
2086 * @tparam InterestDim2 Second dimension along which differentiation is performed.
2087 * @tparam InterestDim3 Third dimension along which differentiation is performed.
2088 *
2089 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
2090 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
2091 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
2092 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
2093 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
2094 * the set of 3D spline coefficients retained to perform the evaluation).
2095 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
2096 */
2097 template <
2098 class InterestDim1,
2099 class InterestDim2,
2100 class InterestDim3,
2101 class Layout1,
2102 class Layout2,
2103 class Layout3,
2104 class BatchedInterpolationDDom,
2105 class... CoordsDims>
2106 void deriv3(
2107 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
2108 spline_eval,
2109 ddc::ChunkSpan<
2110 ddc::Coordinate<CoordsDims...> const,
2111 BatchedInterpolationDDom,
2112 Layout2,
2113 memory_space> const coords_eval,
2114 ddc::ChunkSpan<
2115 double const,
2116 batched_spline_domain_type<BatchedInterpolationDDom>,
2117 Layout3,
2118 memory_space> const spline_coef) const
2119 {
2120 static_assert(
2121 (std::is_same_v<
2122 InterestDim1,
2123 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2124 && std::is_same_v<InterestDim2, continuous_dimension_type2>
2125 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
2126 || (std::is_same_v<
2127 InterestDim3,
2128 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2129 && std::is_same_v<InterestDim1, continuous_dimension_type2>
2130 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2131 || (std::is_same_v<
2132 InterestDim2,
2133 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2134 && std::is_same_v<InterestDim3, continuous_dimension_type2>
2135 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
2136
2137 return deriv_1_2_3(spline_eval, coords_eval, spline_coef);
2138 }
2139
2140 /**
2141 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
2142 *
2143 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
2144 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
2145 *
2146 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
2147 * This means that for each slice of spline_eval the evaluation is performed with
2148 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
2149 *
2150 * Note: triple-differentiation other than cross-differentiation is not supported atm. See #440
2151 *
2152 * @tparam InterestDim1 First dimension along which differentiation is performed.
2153 * @tparam InterestDim2 Second dimension along which differentiation is performed.
2154 * @tparam InterestDim3 Third dimension along which differentiation is performed.
2155 *
2156 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
2157 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
2158 */
2159 template <
2160 class InterestDim1,
2161 class InterestDim2,
2162 class InterestDim3,
2163 class Layout1,
2164 class Layout2,
2165 class BatchedInterpolationDDom>
2166 void deriv3(
2167 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
2168 spline_eval,
2169 ddc::ChunkSpan<
2170 double const,
2171 batched_spline_domain_type<BatchedInterpolationDDom>,
2172 Layout2,
2173 memory_space> const spline_coef) const
2174 {
2175 static_assert(
2176 (std::is_same_v<
2177 InterestDim1,
2178 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2179 && std::is_same_v<InterestDim2, continuous_dimension_type2>
2180 && std::is_same_v<InterestDim3, continuous_dimension_type3>)
2181 || (std::is_same_v<
2182 InterestDim3,
2183 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2184 && std::is_same_v<InterestDim1, continuous_dimension_type2>
2185 && std::is_same_v<InterestDim2, continuous_dimension_type3>)
2186 || (std::is_same_v<
2187 InterestDim2,
2188 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
2189 && std::is_same_v<InterestDim3, continuous_dimension_type2>
2190 && std::is_same_v<InterestDim1, continuous_dimension_type3>));
2191
2192 return deriv_1_2_3(spline_eval, spline_coef);
2193 }
2194
2195 /** @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.
2196 *
2197 * 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.
2198 *
2199 * This is not a nD integration. This is a batched 3D integration.
2200 * This means that for each element of integrals, the integration is performed with the 3D set of
2201 * spline coefficients identified by the same DiscreteElement.
2202 *
2203 * @param[out] integrals The integrals of the 3D spline function on the subdomain of batch_domain. For practical reasons those are
2204 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
2205 * points represented by this domain are unused and irrelevant.
2206 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
2207 */
2208 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
2209 void integrate(
2210 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
2211 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
2212 spline_coef) const
2213 {
2214 static_assert(
2215 ddc::type_seq_contains_v<
2216 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>,
2217 to_type_seq_t<BatchedSplineDDom>>,
2218 "The spline coefficients domain must contain the bsplines dimensions");
2219 using batch_domain_type = ddc::
2220 remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2, bsplines_type3>;
2221 static_assert(
2222 std::is_same_v<batch_domain_type, BatchedDDom>,
2223 "The integrals domain must only contain the batch dimensions");
2224
2225 batch_domain_type batch_domain(integrals.domain());
2226 ddc::Chunk values1_alloc(
2227 ddc::DiscreteDomain<bsplines_type1>(spline_coef.domain()),
2228 ddc::KokkosAllocator<double, memory_space>());
2229 ddc::ChunkSpan values1 = values1_alloc.span_view();
2230 ddc::integrals(exec_space(), values1);
2231 ddc::Chunk values2_alloc(
2232 ddc::DiscreteDomain<bsplines_type2>(spline_coef.domain()),
2233 ddc::KokkosAllocator<double, memory_space>());
2234 ddc::ChunkSpan values2 = values2_alloc.span_view();
2235 ddc::integrals(exec_space(), values2);
2236 ddc::Chunk values3_alloc(
2237 ddc::DiscreteDomain<bsplines_type3>(spline_coef.domain()),
2238 ddc::KokkosAllocator<double, memory_space>());
2239 ddc::ChunkSpan values3 = values3_alloc.span_view();
2240 ddc::integrals(exec_space(), values3);
2241
2242 ddc::parallel_for_each(
2243 "ddc_splines_integrate_bsplines",
2244 exec_space(),
2245 batch_domain,
2246 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
2247 integrals(j) = 0;
2248 for (typename spline_domain_type1::discrete_element_type const i1 :
2249 values1.domain()) {
2250 for (typename spline_domain_type2::discrete_element_type const i2 :
2251 values2.domain()) {
2252 for (typename spline_domain_type3::discrete_element_type const i3 :
2253 values3.domain()) {
2254 integrals(j) += spline_coef(i1, i2, i3, j) * values1(i1)
2255 * values2(i2) * values3(i3);
2256 }
2257 }
2258 }
2259 });
2260 }
2261
2262private:
2263 /**
2264 * @brief Evaluate the function on B-splines at the coordinate given.
2265 *
2266 * This function firstly deals with the boundary conditions and calls the SplineEvaluator3D::eval_no_bc function
2267 * to evaluate.
2268 *
2269 * @param[in] coord_eval The 3D coordinate where we want to evaluate.
2270 * @param[in] spline_coef The B-splines coefficients of the function we want to evaluate.
2271 * @param[out] vals1 A ChunkSpan with the not-null values of each function of the spline in the first dimension.
2272 * @param[out] vals2 A ChunkSpan with the not-null values of each function of the spline in the second dimension.
2273 *
2274 * @return A double with the value of the function at the coordinate given.
2275 *
2276 * @see SplineBoundaryValue
2277 */
2278 template <class Layout, class... CoordsDims>
2279 KOKKOS_INLINE_FUNCTION double eval(
2280 ddc::Coordinate<CoordsDims...> coord_eval,
2281 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
2282 spline_coef) const
2283 {
2284 using Dim1 = continuous_dimension_type1;
2285 using Dim2 = continuous_dimension_type2;
2286 using Dim3 = continuous_dimension_type3;
2287 if constexpr (bsplines_type1::is_periodic()) {
2288 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()
2289 || ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
2290 ddc::get<Dim1>(coord_eval)
2291 -= Kokkos::floor(
2292 (ddc::get<Dim1>(coord_eval)
2293 - ddc::discrete_space<bsplines_type1>().rmin())
2294 / ddc::discrete_space<bsplines_type1>().length())
2295 * ddc::discrete_space<bsplines_type1>().length();
2296 }
2297 }
2298 if constexpr (bsplines_type2::is_periodic()) {
2299 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()
2300 || ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
2301 ddc::get<Dim2>(coord_eval)
2302 -= Kokkos::floor(
2303 (ddc::get<Dim2>(coord_eval)
2304 - ddc::discrete_space<bsplines_type2>().rmin())
2305 / ddc::discrete_space<bsplines_type2>().length())
2306 * ddc::discrete_space<bsplines_type2>().length();
2307 }
2308 }
2309 if constexpr (bsplines_type3::is_periodic()) {
2310 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()
2311 || ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
2312 ddc::get<Dim3>(coord_eval)
2313 -= Kokkos::floor(
2314 (ddc::get<Dim3>(coord_eval)
2315 - ddc::discrete_space<bsplines_type3>().rmin())
2316 / ddc::discrete_space<bsplines_type3>().length())
2317 * ddc::discrete_space<bsplines_type3>().length();
2318 }
2319 }
2320 if constexpr (!bsplines_type1::is_periodic()) {
2321 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()) {
2322 return m_lower_extrap_rule_1(coord_eval, spline_coef);
2323 }
2324 if (ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
2325 return m_upper_extrap_rule_1(coord_eval, spline_coef);
2326 }
2327 }
2328 if constexpr (!bsplines_type2::is_periodic()) {
2329 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()) {
2330 return m_lower_extrap_rule_2(coord_eval, spline_coef);
2331 }
2332 if (ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
2333 return m_upper_extrap_rule_2(coord_eval, spline_coef);
2334 }
2335 }
2336 if constexpr (!bsplines_type3::is_periodic()) {
2337 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()) {
2338 return m_lower_extrap_rule_3(coord_eval, spline_coef);
2339 }
2340 if (ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
2341 return m_upper_extrap_rule_3(coord_eval, spline_coef);
2342 }
2343 }
2344 return eval_no_bc<eval_type, eval_type, eval_type>(
2345 ddc::Coordinate<
2346 continuous_dimension_type1,
2347 continuous_dimension_type2,
2348 continuous_dimension_type3>(
2349 ddc::get<Dim1>(coord_eval),
2350 ddc::get<Dim2>(coord_eval),
2351 ddc::get<Dim3>(coord_eval)),
2352 spline_coef);
2353 }
2354
2355 /**
2356 * @brief Evaluate the function or its derivative at the coordinate given.
2357 *
2358 * @param[in] coord_eval The coordinate where we want to evaluate.
2359 * @param[in] splne_coef The B-splines coefficients of the function we want to evaluate.
2360 * @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`.
2361 * @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`.
2362 */
2363 template <class EvalType1, class EvalType2, class EvalType3, class Layout, class... CoordsDims>
2364 KOKKOS_INLINE_FUNCTION double eval_no_bc(
2365 ddc::Coordinate<CoordsDims...> const& coord_eval,
2366 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
2367 spline_coef) const
2368 {
2369 static_assert(
2370 std::is_same_v<EvalType1, eval_type> || std::is_same_v<EvalType1, eval_deriv_type>);
2371 static_assert(
2372 std::is_same_v<EvalType2, eval_type> || std::is_same_v<EvalType2, eval_deriv_type>);
2373 static_assert(
2374 std::is_same_v<EvalType3, eval_type> || std::is_same_v<EvalType3, eval_deriv_type>);
2375 ddc::DiscreteElement<bsplines_type1> jmin1;
2376 ddc::DiscreteElement<bsplines_type2> jmin2;
2377 ddc::DiscreteElement<bsplines_type3> jmin3;
2378
2379 std::array<double, bsplines_type1::degree() + 1> vals1_ptr;
2380 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>> const
2381 vals1(vals1_ptr.data());
2382 std::array<double, bsplines_type2::degree() + 1> vals2_ptr;
2383 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>> const
2384 vals2(vals2_ptr.data());
2385 std::array<double, bsplines_type3::degree() + 1> vals3_ptr;
2386 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type3::degree() + 1>> const
2387 vals3(vals3_ptr.data());
2388 ddc::Coordinate<continuous_dimension_type1> const coord_eval_interest1(coord_eval);
2389 ddc::Coordinate<continuous_dimension_type2> const coord_eval_interest2(coord_eval);
2390 ddc::Coordinate<continuous_dimension_type3> const coord_eval_interest3(coord_eval);
2391
2392 if constexpr (std::is_same_v<EvalType1, eval_type>) {
2393 jmin1 = ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
2394 } else if constexpr (std::is_same_v<EvalType1, eval_deriv_type>) {
2395 jmin1 = ddc::discrete_space<bsplines_type1>().eval_deriv(vals1, coord_eval_interest1);
2396 }
2397 if constexpr (std::is_same_v<EvalType2, eval_type>) {
2398 jmin2 = ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
2399 } else if constexpr (std::is_same_v<EvalType2, eval_deriv_type>) {
2400 jmin2 = ddc::discrete_space<bsplines_type2>().eval_deriv(vals2, coord_eval_interest2);
2401 }
2402 if constexpr (std::is_same_v<EvalType3, eval_type>) {
2403 jmin3 = ddc::discrete_space<bsplines_type3>().eval_basis(vals3, coord_eval_interest3);
2404 } else if constexpr (std::is_same_v<EvalType3, eval_deriv_type>) {
2405 jmin3 = ddc::discrete_space<bsplines_type3>().eval_deriv(vals3, coord_eval_interest3);
2406 }
2407
2408 double y = 0.0;
2409 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
2410 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
2411 for (std::size_t k = 0; k < bsplines_type3::degree() + 1; ++k) {
2412 y += spline_coef(
2413 ddc::DiscreteElement<
2414 bsplines_type1,
2415 bsplines_type2,
2416 bsplines_type3>(jmin1 + i, jmin2 + j, jmin3 + k))
2417 * vals1[i] * vals2[j] * vals3[k];
2418 }
2419 }
2420 }
2421 return y;
2422 }
2423};
2424
2425} // 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.