DDC 0.5.2
Loading...
Searching...
No Matches
spline_evaluator_2d.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 2D spline function.
22 *
23 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a 2D 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 EvaluationDDim1 The first discrete dimension on which evaluation points are defined.
30 * @tparam EvaluationDDim2 The second discrete dimension on which evaluation points are defined.
31 * @tparam LowerExtrapolationRule1 The lower extrapolation rule type along first dimension of interest.
32 * @tparam UpperExtrapolationRule1 The upper extrapolation rule type along first dimension of interest.
33 * @tparam LowerExtrapolationRule2 The lower extrapolation rule type along second dimension of interest.
34 * @tparam UpperExtrapolationRule2 The upper extrapolation rule type along second dimension of interest.
35 */
36template <
37 class ExecSpace,
38 class MemorySpace,
39 class BSplines1,
40 class BSplines2,
41 class EvaluationDDim1,
42 class EvaluationDDim2,
43 class LowerExtrapolationRule1,
44 class UpperExtrapolationRule1,
45 class LowerExtrapolationRule2,
46 class UpperExtrapolationRule2>
48{
49private:
50 /**
51 * @brief Tag to indicate that the value of the spline should be evaluated.
52 */
53 struct eval_type
54 {
55 };
56
57 /**
58 * @brief Tag to indicate that derivative of the spline should be evaluated.
59 */
60 struct eval_deriv_type
61 {
62 };
63
64public:
65 /// @brief The type of the first evaluation continuous dimension used by this class.
66 using continuous_dimension_type1 = typename BSplines1::continuous_dimension_type;
67
68 /// @brief The type of the second evaluation continuous dimension used by this class.
69 using continuous_dimension_type2 = typename BSplines2::continuous_dimension_type;
70
71 /// @brief The type of the Kokkos execution space used by this class.
72 using exec_space = ExecSpace;
73
74 /// @brief The type of the Kokkos memory space used by this class.
75 using memory_space = MemorySpace;
76
77 /// @brief The type of the first discrete dimension of interest used by this class.
78 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
79
80 /// @brief The type of the second discrete dimension of interest used by this class.
81 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
82
83 /// @brief The discrete dimension representing the B-splines along first dimension.
84 using bsplines_type1 = BSplines1;
85
86 /// @brief The discrete dimension representing the B-splines along second dimension.
87 using bsplines_type2 = BSplines2;
88
89 /// @brief The type of the domain for the 1D evaluation mesh along first dimension used by this class.
90 using evaluation_domain_type1 = ddc::DiscreteDomain<evaluation_discrete_dimension_type1>;
91
92 /// @brief The type of the domain for the 1D evaluation mesh along second dimension used by this class.
93 using evaluation_domain_type2 = ddc::DiscreteDomain<evaluation_discrete_dimension_type2>;
94
95 /// @brief The type of the domain for the 2D evaluation mesh used by this class.
96 using evaluation_domain_type = ddc::DiscreteDomain<
97 evaluation_discrete_dimension_type1,
98 evaluation_discrete_dimension_type2>;
99
100 /**
101 * @brief The type of the whole domain representing evaluation points.
102 *
103 * @tparam The batched discrete domain on which the interpolation points are defined.
104 */
105 template <
106 class BatchedInterpolationDDom,
107 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
108 using batched_evaluation_domain_type = BatchedInterpolationDDom;
109
110 /// @brief The type of the 1D spline domain corresponding to the first dimension of interest.
111 using spline_domain_type1 = ddc::DiscreteDomain<bsplines_type1>;
112
113 /// @brief The type of the 1D spline domain corresponding to the second dimension of interest.
114 using spline_domain_type2 = ddc::DiscreteDomain<bsplines_type2>;
115
116 /// @brief The type of the 2D spline domain corresponding to the dimensions of interest.
117 using spline_domain_type = ddc::DiscreteDomain<bsplines_type1, bsplines_type2>;
118
119 /**
120 * @brief The type of the batch domain (obtained by removing the dimensions of interest
121 * from the whole domain).
122 *
123 * @tparam The batched discrete domain on which the interpolation points are defined.
124 */
125 template <
126 class BatchedInterpolationDDom,
127 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
128 using batch_domain_type = typename ddc::remove_dims_of_t<
129 BatchedInterpolationDDom,
130 evaluation_discrete_dimension_type1,
131 evaluation_discrete_dimension_type2>;
132
133 /**
134 * @brief The type of the whole spline domain (cartesian product of 2D spline domain
135 * and batch domain) preserving the underlying memory layout (order of dimensions).
136 *
137 * @tparam The batched discrete domain on which the interpolation points are defined.
138 */
139 template <
140 class BatchedInterpolationDDom,
141 class = std::enable_if_t<ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
142 using batched_spline_domain_type =
143 typename ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
144 ddc::to_type_seq_t<BatchedInterpolationDDom>,
145 ddc::detail::TypeSeq<
146 evaluation_discrete_dimension_type1,
147 evaluation_discrete_dimension_type2>,
148 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2>>>;
149
150 /// @brief The type of the extrapolation rule at the lower boundary along the first dimension.
151 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
152
153 /// @brief The type of the extrapolation rule at the upper boundary along the first dimension.
154 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
155
156 /// @brief The type of the extrapolation rule at the lower boundary along the second dimension.
157 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
158
159 /// @brief The type of the extrapolation rule at the upper boundary along the second dimension.
160 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
161
162private:
163 LowerExtrapolationRule1 m_lower_extrap_rule_1;
164
165 UpperExtrapolationRule1 m_upper_extrap_rule_1;
166
167 LowerExtrapolationRule2 m_lower_extrap_rule_2;
168
169 UpperExtrapolationRule2 m_upper_extrap_rule_2;
170
171public:
172 static_assert(
173 std::is_same_v<LowerExtrapolationRule1,
174 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
175 == bsplines_type1::is_periodic()
176 && std::is_same_v<
177 UpperExtrapolationRule1,
178 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
179 == bsplines_type1::is_periodic()
180 && std::is_same_v<
181 LowerExtrapolationRule2,
182 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
183 == bsplines_type2::is_periodic()
184 && std::is_same_v<
185 UpperExtrapolationRule2,
186 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
187 == bsplines_type2::is_periodic(),
188 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
189 static_assert(
190 std::is_invocable_r_v<
191 double,
192 LowerExtrapolationRule1,
193 ddc::Coordinate<continuous_dimension_type1>,
194 ddc::ChunkSpan<
195 double const,
196 spline_domain_type,
197 Kokkos::layout_right,
198 memory_space>>,
199 "LowerExtrapolationRule1::operator() has to be callable "
200 "with usual arguments.");
201 static_assert(
202 std::is_invocable_r_v<
203 double,
204 UpperExtrapolationRule1,
205 ddc::Coordinate<continuous_dimension_type1>,
206 ddc::ChunkSpan<
207 double const,
208 spline_domain_type,
209 Kokkos::layout_right,
210 memory_space>>,
211 "UpperExtrapolationRule1::operator() has to be callable "
212 "with usual arguments.");
213 static_assert(
214 std::is_invocable_r_v<
215 double,
216 LowerExtrapolationRule2,
217 ddc::Coordinate<continuous_dimension_type2>,
218 ddc::ChunkSpan<
219 double const,
220 spline_domain_type,
221 Kokkos::layout_right,
222 memory_space>>,
223 "LowerExtrapolationRule2::operator() has to be callable "
224 "with usual arguments.");
225 static_assert(
226 std::is_invocable_r_v<
227 double,
228 UpperExtrapolationRule2,
229 ddc::Coordinate<continuous_dimension_type2>,
230 ddc::ChunkSpan<
231 double const,
232 spline_domain_type,
233 Kokkos::layout_right,
234 memory_space>>,
235 "UpperExtrapolationRule2::operator() has to be callable "
236 "with usual arguments.");
237
238 /**
239 * @brief Build a SplineEvaluator2D acting on batched_spline_domain.
240 *
241 * @param lower_extrap_rule1 The extrapolation rule at the lower boundary along the first dimension.
242 * @param upper_extrap_rule1 The extrapolation rule at the upper boundary along the first dimension.
243 * @param lower_extrap_rule2 The extrapolation rule at the lower boundary along the second dimension.
244 * @param upper_extrap_rule2 The extrapolation rule at the upper boundary along the second dimension.
245 *
246 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
247 */
248 explicit SplineEvaluator2D(
249 LowerExtrapolationRule1 const& lower_extrap_rule1,
250 UpperExtrapolationRule1 const& upper_extrap_rule1,
251 LowerExtrapolationRule2 const& lower_extrap_rule2,
252 UpperExtrapolationRule2 const& upper_extrap_rule2)
253 : m_lower_extrap_rule_1(lower_extrap_rule1)
254 , m_upper_extrap_rule_1(upper_extrap_rule1)
255 , m_lower_extrap_rule_2(lower_extrap_rule2)
256 , m_upper_extrap_rule_2(upper_extrap_rule2)
257 {
258 }
259
260 /**
261 * @brief Copy-constructs.
262 *
263 * @param x A reference to another SplineEvaluator.
264 */
265 SplineEvaluator2D(SplineEvaluator2D const& x) = default;
266
267 /**
268 * @brief Move-constructs.
269 *
270 * @param x An rvalue to another SplineEvaluator.
271 */
273
274 /// @brief Destructs.
275 ~SplineEvaluator2D() = default;
276
277 /**
278 * @brief Copy-assigns.
279 *
280 * @param x A reference to another SplineEvaluator.
281 * @return A reference to this object.
282 */
283 SplineEvaluator2D& operator=(SplineEvaluator2D const& x) = default;
284
285 /**
286 * @brief Move-assigns.
287 *
288 * @param x An rvalue to another SplineEvaluator.
289 * @return A reference to this object.
290 */
292
293 /**
294 * @brief Get the lower extrapolation rule along the first dimension.
295 *
296 * 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.
297 *
298 * @return The lower extrapolation rule along the first dimension.
299 *
300 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
301 */
302 lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
303 {
304 return m_lower_extrap_rule_1;
305 }
306
307 /**
308 * @brief Get the upper extrapolation rule along the first dimension.
309 *
310 * 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.
311 *
312 * @return The upper extrapolation rule along the first dimension.
313 *
314 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
315 */
316 upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
317 {
318 return m_upper_extrap_rule_1;
319 }
320
321 /**
322 * @brief Get the lower extrapolation rule along the second dimension.
323 *
324 * 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.
325 *
326 * @return The lower extrapolation rule along the second dimension.
327 *
328 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
329 */
330 lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
331 {
332 return m_lower_extrap_rule_2;
333 }
334
335 /**
336 * @brief Get the upper extrapolation rule along the second dimension.
337 *
338 * 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.
339 *
340 * @return The upper extrapolation rule along the second dimension.
341 *
342 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
343 */
344 upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
345 {
346 return m_upper_extrap_rule_2;
347 }
348
349 /**
350 * @brief Evaluate 2D spline function (described by its spline coefficients) at a given coordinate.
351 *
352 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
353 *
354 * Remark: calling SplineBuilder2D then SplineEvaluator2D corresponds to a 2D spline interpolation.
355 *
356 * @param coord_eval The coordinate where the spline is evaluated. Note that only the components along the dimensions of interest are used.
357 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
358 *
359 * @return The value of the spline function at the desired coordinate.
360 */
361 template <class Layout, class... CoordsDims>
362 KOKKOS_FUNCTION double operator()(
363 ddc::Coordinate<CoordsDims...> const& coord_eval,
364 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
365 spline_coef) const
366 {
367 return eval(coord_eval, spline_coef);
368 }
369
370 /**
371 * @brief Evaluate 2D spline function (described by its spline coefficients) on a mesh.
372 *
373 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
374 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
375 *
376 * This is not a nD evaluation. This is a batched 2D evaluation. This means that for each slice of coordinates
377 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 2D set of
378 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
379 *
380 * Remark: calling SplineBuilder2D then SplineEvaluator2D corresponds to a 2D spline interpolation.
381 *
382 * @param[out] spline_eval The values of the 2D spline function at the desired coordinates. For practical reasons those are
383 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
384 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
385 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
386 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
387 * the set of 2D spline coefficients retained to perform the evaluation).
388 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
389 */
390 template <
391 class Layout1,
392 class Layout2,
393 class Layout3,
394 class BatchedInterpolationDDom,
395 class... CoordsDims>
396 void operator()(
397 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
398 spline_eval,
399 ddc::ChunkSpan<
400 ddc::Coordinate<CoordsDims...> const,
401 BatchedInterpolationDDom,
402 Layout2,
403 memory_space> const coords_eval,
404 ddc::ChunkSpan<
405 double const,
406 batched_spline_domain_type<BatchedInterpolationDDom>,
407 Layout3,
408 memory_space> const spline_coef) const
409 {
410 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
411 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
412 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
413 ddc::parallel_for_each(
414 "ddc_splines_evaluate_2d",
415 exec_space(),
416 batch_domain,
417 KOKKOS_CLASS_LAMBDA(
418 typename batch_domain_type<
419 BatchedInterpolationDDom>::discrete_element_type const j) {
420 const auto spline_eval_2D = spline_eval[j];
421 const auto coords_eval_2D = coords_eval[j];
422 const auto spline_coef_2D = spline_coef[j];
423 for (auto const i1 : evaluation_domain1) {
424 for (auto const i2 : evaluation_domain2) {
425 spline_eval_2D(i1, i2) = eval(coords_eval_2D(i1, i2), spline_coef_2D);
426 }
427 }
428 });
429 }
430
431 /**
432 * @brief Evaluate 2D spline function (described by its spline coefficients) on a mesh.
433 *
434 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
435 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
436 *
437 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
438 * This means that for each slice of spline_eval the evaluation is performed with
439 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
440 *
441 * Remark: calling SplineBuilder2D then SplineEvaluator2D corresponds to a 2D spline interpolation.
442 *
443 * @param[out] spline_eval The values of the 2D spline function at their coordinates.
444 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
445 */
446 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
447 void operator()(
448 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
449 spline_eval,
450 ddc::ChunkSpan<
451 double const,
452 batched_spline_domain_type<BatchedInterpolationDDom>,
453 Layout2,
454 memory_space> const spline_coef) const
455 {
456 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
457 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
458 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
459 ddc::parallel_for_each(
460 "ddc_splines_evaluate_2d",
461 exec_space(),
462 batch_domain,
463 KOKKOS_CLASS_LAMBDA(
464 typename batch_domain_type<
465 BatchedInterpolationDDom>::discrete_element_type const j) {
466 const auto spline_eval_2D = spline_eval[j];
467 const auto spline_coef_2D = spline_coef[j];
468 for (auto const i1 : evaluation_domain1) {
469 for (auto const i2 : evaluation_domain2) {
470 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
471 coord_eval_2D(ddc::coordinate(i1), ddc::coordinate(i2));
472 spline_eval_2D(i1, i2) = eval(coord_eval_2D(i1, i2), spline_coef_2D);
473 }
474 }
475 });
476 }
477
478 /**
479 * @brief Differentiate 2D spline function (described by its spline coefficients) at a given coordinate along first dimension of interest.
480 *
481 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be
482 * obtained via various methods, such as using a SplineBuilder2D.
483 *
484 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
485 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
486 *
487 * @return The derivative of the spline function at the desired coordinate.
488 */
489 template <class Layout, class... CoordsDims>
491 ddc::Coordinate<CoordsDims...> const& coord_eval,
492 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
493 spline_coef) const
494 {
495 return eval_no_bc<eval_deriv_type, eval_type>(coord_eval, spline_coef);
496 }
497
498 /**
499 * @brief Differentiate 2D spline function (described by its spline coefficients) at a given coordinate along second dimension of interest.
500 *
501 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be
502 * obtained via various methods, such as using a SplineBuilder2D.
503 *
504 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
505 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
506 *
507 * @return The derivative of the spline function at the desired coordinate.
508 */
509 template <class Layout, class... CoordsDims>
511 ddc::Coordinate<CoordsDims...> const& coord_eval,
512 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
513 spline_coef) const
514 {
515 return eval_no_bc<eval_type, eval_deriv_type>(coord_eval, spline_coef);
516 }
517
518 /**
519 * @brief Cross-differentiate 2D spline function (described by its spline coefficients) at a given coordinate.
520 *
521 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be
522 * obtained via various methods, such as using a SplineBuilder2D.
523 *
524 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
525 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
526 *
527 * @return The derivative of the spline function at the desired coordinate.
528 */
529 template <class Layout, class... CoordsDims>
531 ddc::Coordinate<CoordsDims...> const& coord_eval,
532 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
533 spline_coef) const
534 {
535 return eval_no_bc<eval_deriv_type, eval_deriv_type>(coord_eval, spline_coef);
536 }
537
538 /**
539 * @brief Differentiate 2D spline function (described by its spline coefficients) at a given coordinate along a specified dimension of interest.
540 *
541 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be
542 * obtained via various methods, such as using a SplineBuilder2D.
543 *
544 * @tparam InterestDim Dimension along which differentiation is performed.
545 *
546 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
547 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
548 *
549 * @return The derivative of the spline function at the desired coordinate.
550 */
551 template <class InterestDim, class Layout, class... CoordsDims>
552 KOKKOS_FUNCTION double deriv(
553 ddc::Coordinate<CoordsDims...> const& coord_eval,
554 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
555 spline_coef) const
556 {
557 static_assert(
558 std::is_same_v<InterestDim, continuous_dimension_type1>
559 || std::is_same_v<InterestDim, continuous_dimension_type2>);
560 if constexpr (std::is_same_v<
561 InterestDim,
562 typename evaluation_discrete_dimension_type1::
563 continuous_dimension_type>) {
564 return deriv_dim_1(coord_eval, spline_coef);
565 } else if constexpr (std::is_same_v<
566 InterestDim,
567 typename evaluation_discrete_dimension_type2::
568 continuous_dimension_type>) {
569 return deriv_dim_2(coord_eval, spline_coef);
570 }
571 }
572
573 /**
574 * @brief Double-differentiate 2D spline function (described by its spline coefficients) at a given coordinate along specified dimensions of interest.
575 *
576 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be
577 * obtained via various methods, such as using a SplineBuilder2D.
578 *
579 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
580 *
581 * @tparam InterestDim1 First dimension along which differentiation is performed.
582 * @tparam InterestDim2 Second dimension along which differentiation is performed.
583 *
584 * @param coord_eval The coordinate where the spline is double-differentiated. Note that only the components along the dimensions of interest are used.
585 * @param spline_coef A ChunkSpan storing the 2D spline coefficients.
586 *
587 * @return The derivative of the spline function at the desired coordinate.
588 */
589 template <class InterestDim1, class InterestDim2, class Layout, class... CoordsDims>
590 KOKKOS_FUNCTION double deriv2(
591 ddc::Coordinate<CoordsDims...> const& coord_eval,
592 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
593 spline_coef) const
594 {
595 static_assert(
596 (std::is_same_v<
597 InterestDim1,
598 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
599 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
600 || (std::is_same_v<
601 InterestDim2,
602 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
603 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
604 return deriv_1_and_2(coord_eval, spline_coef);
605 }
606
607 /**
608 * @brief Differentiate 2D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
609 *
610 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
611 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
612 *
613 * This is not a nD evaluation. This is a batched 2D differentiation.
614 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
615 * the differentiation is performed with the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
616 *
617 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates. For practical reasons those are
618 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
619 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
620 * the set of 2D spline coefficients retained to perform the evaluation).
621 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
622 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
623 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
624 * the set of 2D spline coefficients retained to perform the evaluation).
625 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
626 */
627 template <
628 class Layout1,
629 class Layout2,
630 class Layout3,
631 class BatchedInterpolationDDom,
632 class... CoordsDims>
633 void deriv_dim_1(
634 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
635 spline_eval,
636 ddc::ChunkSpan<
637 ddc::Coordinate<CoordsDims...> const,
638 BatchedInterpolationDDom,
639 Layout2,
640 memory_space> const coords_eval,
641 ddc::ChunkSpan<
642 double const,
643 batched_spline_domain_type<BatchedInterpolationDDom>,
644 Layout3,
645 memory_space> const spline_coef) const
646 {
647 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
648 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
649 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
650 ddc::parallel_for_each(
651 "ddc_splines_differentiate_2d_dim_1",
652 exec_space(),
653 batch_domain,
654 KOKKOS_CLASS_LAMBDA(
655 typename batch_domain_type<
656 BatchedInterpolationDDom>::discrete_element_type const j) {
657 const auto spline_eval_2D = spline_eval[j];
658 const auto coords_eval_2D = coords_eval[j];
659 const auto spline_coef_2D = spline_coef[j];
660 for (auto const i1 : evaluation_domain1) {
661 for (auto const i2 : evaluation_domain2) {
662 spline_eval_2D(i1, i2) = eval_no_bc<
663 eval_deriv_type,
664 eval_type>(coords_eval_2D(i1, i2), spline_coef_2D);
665 }
666 }
667 });
668 }
669
670 /**
671 * @brief Differentiate 2D spline function (described by its spline coefficients) on a mesh along first dimension of interest.
672 *
673 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
674 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
675 *
676 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
677 * This means that for each slice of spline_eval the evaluation is performed with
678 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
679 *
680 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates.
681 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
682 */
683 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
684 void deriv_dim_1(
685 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
686 spline_eval,
687 ddc::ChunkSpan<
688 double const,
689 batched_spline_domain_type<BatchedInterpolationDDom>,
690 Layout2,
691 memory_space> const spline_coef) const
692 {
693 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
694 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
695 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
696 ddc::parallel_for_each(
697 "ddc_splines_differentiate_2d_dim_1",
698 exec_space(),
699 batch_domain,
700 KOKKOS_CLASS_LAMBDA(
701 typename batch_domain_type<
702 BatchedInterpolationDDom>::discrete_element_type const j) {
703 const auto spline_eval_2D = spline_eval[j];
704 const auto spline_coef_2D = spline_coef[j];
705 for (auto const i1 : evaluation_domain1) {
706 for (auto const i2 : evaluation_domain2) {
707 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
708 coord_eval_2D(ddc::coordinate(i1), ddc::coordinate(i2));
709 spline_eval_2D(i1, i2) = eval_no_bc<
710 eval_deriv_type,
711 eval_type>(coord_eval_2D, spline_coef_2D);
712 }
713 }
714 });
715 }
716
717 /**
718 * @brief Differentiate 2D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
719 *
720 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
721 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
722 *
723 * This is not a nD differentiation. This is a batched 2D differentiation.
724 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
725 * the differentiation is performed with the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
726 *
727 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates. For practical reasons those are
728 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
729 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
730 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
731 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
732 * the set of 2D spline coefficients retained to perform the evaluation).
733 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
734 */
735 template <
736 class Layout1,
737 class Layout2,
738 class Layout3,
739 class BatchedInterpolationDDom,
740 class... CoordsDims>
741 void deriv_dim_2(
742 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
743 spline_eval,
744 ddc::ChunkSpan<
745 ddc::Coordinate<CoordsDims...> const,
746 BatchedInterpolationDDom,
747 Layout2,
748 memory_space> const coords_eval,
749 ddc::ChunkSpan<
750 double const,
751 batched_spline_domain_type<BatchedInterpolationDDom>,
752 Layout3,
753 memory_space> const spline_coef) const
754 {
755 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
756 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
757 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
758 ddc::parallel_for_each(
759 "ddc_splines_differentiate_2d_dim_2",
760 exec_space(),
761 batch_domain,
762 KOKKOS_CLASS_LAMBDA(
763 typename batch_domain_type<
764 BatchedInterpolationDDom>::discrete_element_type const j) {
765 const auto spline_eval_2D = spline_eval[j];
766 const auto coords_eval_2D = coords_eval[j];
767 const auto spline_coef_2D = spline_coef[j];
768 for (auto const i1 : evaluation_domain1) {
769 for (auto const i2 : evaluation_domain2) {
770 spline_eval_2D(i1, i2) = eval_no_bc<
771 eval_type,
772 eval_deriv_type>(coords_eval_2D(i1, i2), spline_coef_2D);
773 }
774 }
775 });
776 }
777
778 /**
779 * @brief Differentiate 2D spline function (described by its spline coefficients) on a mesh along second dimension of interest.
780 *
781 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
782 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
783 *
784 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
785 * This means that for each slice of spline_eval the evaluation is performed with
786 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
787 *
788 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates.
789 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
790 */
791 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
792 void deriv_dim_2(
793 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
794 spline_eval,
795 ddc::ChunkSpan<
796 double const,
797 batched_spline_domain_type<BatchedInterpolationDDom>,
798 Layout2,
799 memory_space> const spline_coef) const
800 {
801 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
802 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
803 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
804 ddc::parallel_for_each(
805 "ddc_splines_differentiate_2d_dim_2",
806 exec_space(),
807 batch_domain,
808 KOKKOS_CLASS_LAMBDA(
809 typename batch_domain_type<
810 BatchedInterpolationDDom>::discrete_element_type const j) {
811 const auto spline_eval_2D = spline_eval[j];
812 const auto spline_coef_2D = spline_coef[j];
813 for (auto const i1 : evaluation_domain1) {
814 for (auto const i2 : evaluation_domain2) {
815 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
816 coord_eval_2D(ddc::coordinate(i1), ddc::coordinate(i2));
817 spline_eval_2D(i1, i2) = eval_no_bc<
818 eval_type,
819 eval_deriv_type>(coord_eval_2D, spline_coef_2D);
820 }
821 }
822 });
823 }
824
825 /**
826 * @brief Cross-differentiate 2D spline function (described by its spline coefficients) on a mesh along dimensions of interest.
827 *
828 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
829 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
830 *
831 * This is not a nD cross-differentiation. This is a batched 2D cross-differentiation.
832 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
833 * the cross-differentiation is performed with the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
834 *
835 * @param[out] spline_eval The cross-derivatives of the 2D spline function at the desired coordinates. For practical reasons those are
836 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
837 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
838 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
839 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
840 * the set of 2D spline coefficients retained to perform the evaluation).
841 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
842 */
843 template <
844 class Layout1,
845 class Layout2,
846 class Layout3,
847 class BatchedInterpolationDDom,
848 class... CoordsDims>
849 void deriv_1_and_2(
850 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
851 spline_eval,
852 ddc::ChunkSpan<
853 ddc::Coordinate<CoordsDims...> const,
854 BatchedInterpolationDDom,
855 Layout2,
856 memory_space> const coords_eval,
857 ddc::ChunkSpan<
858 double const,
859 batched_spline_domain_type<BatchedInterpolationDDom>,
860 Layout3,
861 memory_space> const spline_coef) const
862 {
863 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
864 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
865 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
866 ddc::parallel_for_each(
867 "ddc_splines_cross_differentiate",
868 exec_space(),
869 batch_domain,
870 KOKKOS_CLASS_LAMBDA(
871 typename batch_domain_type<
872 BatchedInterpolationDDom>::discrete_element_type const j) {
873 const auto spline_eval_2D = spline_eval[j];
874 const auto coords_eval_2D = coords_eval[j];
875 const auto spline_coef_2D = spline_coef[j];
876 for (auto const i1 : evaluation_domain1) {
877 for (auto const i2 : evaluation_domain2) {
878 spline_eval_2D(i1, i2) = eval_no_bc<
879 eval_deriv_type,
880 eval_deriv_type>(coords_eval_2D(i1, i2), spline_coef_2D);
881 }
882 }
883 });
884 }
885
886 /**
887 * @brief Cross-differentiate 2D spline function (described by its spline coefficients) on a mesh along dimensions of interest.
888 *
889 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
890 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
891 *
892 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
893 * This means that for each slice of spline_eval the evaluation is performed with
894 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
895 *
896 * @param[out] spline_eval The cross-derivatives of the 2D spline function at the desired coordinates.
897 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
898 */
899 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
900 void deriv_1_and_2(
901 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
902 spline_eval,
903 ddc::ChunkSpan<
904 double const,
905 batched_spline_domain_type<BatchedInterpolationDDom>,
906 Layout2,
907 memory_space> const spline_coef) const
908 {
909 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
910 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
911 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
912 ddc::parallel_for_each(
913 "ddc_splines_cross_differentiate",
914 exec_space(),
915 batch_domain,
916 KOKKOS_CLASS_LAMBDA(
917 typename batch_domain_type<
918 BatchedInterpolationDDom>::discrete_element_type const j) {
919 const auto spline_eval_2D = spline_eval[j];
920 const auto spline_coef_2D = spline_coef[j];
921 for (auto const i1 : evaluation_domain1) {
922 for (auto const i2 : evaluation_domain2) {
923 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
924 coord_eval_2D(ddc::coordinate(i1), ddc::coordinate(i2));
925 spline_eval_2D(i1, i2) = eval_no_bc<
926 eval_deriv_type,
927 eval_deriv_type>(coord_eval_2D, spline_coef_2D);
928 }
929 }
930 });
931 }
932
933 /**
934 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
935 *
936 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
937 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
938 *
939 * This is not a nD evaluation. This is a batched 2D differentiation.
940 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
941 * the differentiation is performed with the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
942 *
943 * @tparam InterestDim Dimension along which differentiation is performed.
944 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates. For practical reasons those are
945 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
946 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
947 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
948 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
949 * the set of 2D spline coefficients retained to perform the evaluation).
950 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
951 */
952 template <
953 class InterestDim,
954 class Layout1,
955 class Layout2,
956 class Layout3,
957 class BatchedInterpolationDDom,
958 class... CoordsDims>
959 void deriv(
960 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
961 spline_eval,
962 ddc::ChunkSpan<
963 ddc::Coordinate<CoordsDims...> const,
964 BatchedInterpolationDDom,
965 Layout2,
966 memory_space> const coords_eval,
967 ddc::ChunkSpan<
968 double const,
969 batched_spline_domain_type<BatchedInterpolationDDom>,
970 Layout3,
971 memory_space> const spline_coef) const
972 {
973 static_assert(
974 std::is_same_v<
975 InterestDim,
976 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
977 || std::is_same_v<InterestDim, continuous_dimension_type2>);
978 if constexpr (std::is_same_v<
979 InterestDim,
980 typename evaluation_discrete_dimension_type1::
981 continuous_dimension_type>) {
982 return deriv_dim_1(spline_eval, coords_eval, spline_coef);
983 } else if constexpr (std::is_same_v<
984 InterestDim,
985 typename evaluation_discrete_dimension_type2::
986 continuous_dimension_type>) {
987 return deriv_dim_2(spline_eval, coords_eval, spline_coef);
988 }
989 }
990
991 /**
992 * @brief Differentiate spline function (described by its spline coefficients) on a mesh along a specified dimension of interest.
993 *
994 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
995 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
996 *
997 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
998 * This means that for each slice of spline_eval the evaluation is performed with
999 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1000 *
1001 * @tparam InterestDim Dimension along which differentiation is performed.
1002 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates.
1003 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
1004 */
1005 template <class InterestDim, class Layout1, class Layout2, class BatchedInterpolationDDom>
1006 void deriv(
1007 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1008 spline_eval,
1009 ddc::ChunkSpan<
1010 double const,
1011 batched_spline_domain_type<BatchedInterpolationDDom>,
1012 Layout2,
1013 memory_space> const spline_coef) const
1014 {
1015 static_assert(
1016 std::is_same_v<
1017 InterestDim,
1018 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1019 || std::is_same_v<InterestDim, continuous_dimension_type2>);
1020 if constexpr (std::is_same_v<
1021 InterestDim,
1022 typename evaluation_discrete_dimension_type1::
1023 continuous_dimension_type>) {
1024 return deriv_dim_1(spline_eval, spline_coef);
1025 } else if constexpr (std::is_same_v<
1026 InterestDim,
1027 typename evaluation_discrete_dimension_type2::
1028 continuous_dimension_type>) {
1029 return deriv_dim_2(spline_eval, spline_coef);
1030 }
1031 }
1032
1033 /**
1034 * @brief Double-differentiate 2D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1035 *
1036 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
1037 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
1038 *
1039 * This is not a nD evaluation. This is a batched 2D differentiation.
1040 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
1041 * the differentiation is performed with the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1042 *
1043 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1044 *
1045 * @tparam InterestDim1 First dimension along which differentiation is performed.
1046 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1047 *
1048 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates. For practical reasons those are
1049 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
1050 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
1051 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
1052 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
1053 * the set of 2D spline coefficients retained to perform the evaluation).
1054 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
1055 */
1056 template <
1057 class InterestDim1,
1058 class InterestDim2,
1059 class Layout1,
1060 class Layout2,
1061 class Layout3,
1062 class BatchedInterpolationDDom,
1063 class... CoordsDims>
1064 void deriv2(
1065 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1066 spline_eval,
1067 ddc::ChunkSpan<
1068 ddc::Coordinate<CoordsDims...> const,
1069 BatchedInterpolationDDom,
1070 Layout2,
1071 memory_space> const coords_eval,
1072 ddc::ChunkSpan<
1073 double const,
1074 batched_spline_domain_type<BatchedInterpolationDDom>,
1075 Layout3,
1076 memory_space> const spline_coef) const
1077 {
1078 static_assert(
1079 (std::is_same_v<
1080 InterestDim1,
1081 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1082 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1083 || (std::is_same_v<
1084 InterestDim2,
1085 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1086 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
1087 return deriv_1_and_2(spline_eval, coords_eval, spline_coef);
1088 }
1089
1090 /**
1091 * @brief Double-differentiate 2D spline function (described by its spline coefficients) on a mesh along specified dimensions of interest.
1092 *
1093 * The spline coefficients represent a 2D spline function defined on a cartesian product of batch_domain and B-splines
1094 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
1095 *
1096 * This is not a multidimensional evaluation. This is a batched 2D evaluation.
1097 * This means that for each slice of spline_eval the evaluation is performed with
1098 * the 2D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
1099 *
1100 * Note: double-differentiation other than cross-differentiation is not supported atm. See #440
1101 *
1102 * @tparam InterestDim1 First dimension along which differentiation is performed.
1103 * @tparam InterestDim2 Second dimension along which differentiation is performed.
1104 *
1105 * @param[out] spline_eval The derivatives of the 2D spline function at the desired coordinates.
1106 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
1107 */
1108 template <
1109 class InterestDim1,
1110 class InterestDim2,
1111 class Layout1,
1112 class Layout2,
1113 class BatchedInterpolationDDom>
1114 void deriv2(
1115 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
1116 spline_eval,
1117 ddc::ChunkSpan<
1118 double const,
1119 batched_spline_domain_type<BatchedInterpolationDDom>,
1120 Layout2,
1121 memory_space> const spline_coef) const
1122 {
1123 static_assert(
1124 (std::is_same_v<
1125 InterestDim1,
1126 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1127 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1128 || (std::is_same_v<
1129 InterestDim2,
1130 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1131 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
1132 return deriv_1_and_2(spline_eval, spline_coef);
1133 }
1134
1135 /** @brief Perform batched 2D integrations of a spline function (described by its spline coefficients) along the dimensions of interest and store results on a subdomain of batch_domain.
1136 *
1137 * The spline coefficients represent a 2D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder2D.
1138 *
1139 * This is not a nD integration. This is a batched 2D integration.
1140 * This means that for each element of integrals, the integration is performed with the 2D set of
1141 * spline coefficients identified by the same DiscreteElement.
1142 *
1143 * @param[out] integrals The integrals of the 2D spline function on the subdomain of batch_domain. For practical reasons those are
1144 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
1145 * points represented by this domain are unused and irrelevant.
1146 * @param[in] spline_coef A ChunkSpan storing the 2D spline coefficients.
1147 */
1148 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
1149 void integrate(
1150 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
1151 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
1152 spline_coef) const
1153 {
1154 static_assert(
1155 ddc::type_seq_contains_v<
1156 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2>,
1157 to_type_seq_t<BatchedSplineDDom>>,
1158 "The spline coefficients domain must contain the bsplines dimensions");
1159 using batch_domain_type
1160 = ddc::remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2>;
1161 static_assert(
1162 std::is_same_v<batch_domain_type, BatchedDDom>,
1163 "The integrals domain must only contain the batch dimensions");
1164
1165 batch_domain_type batch_domain(integrals.domain());
1166 ddc::Chunk values1_alloc(
1167 ddc::DiscreteDomain<bsplines_type1>(spline_coef.domain()),
1168 ddc::KokkosAllocator<double, memory_space>());
1169 ddc::ChunkSpan values1 = values1_alloc.span_view();
1170 ddc::integrals(exec_space(), values1);
1171 ddc::Chunk values2_alloc(
1172 ddc::DiscreteDomain<bsplines_type2>(spline_coef.domain()),
1173 ddc::KokkosAllocator<double, memory_space>());
1174 ddc::ChunkSpan values2 = values2_alloc.span_view();
1175 ddc::integrals(exec_space(), values2);
1176
1177 ddc::parallel_for_each(
1178 "ddc_splines_integrate_bsplines",
1179 exec_space(),
1180 batch_domain,
1181 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
1182 integrals(j) = 0;
1183 for (typename spline_domain_type1::discrete_element_type const i1 :
1184 values1.domain()) {
1185 for (typename spline_domain_type2::discrete_element_type const i2 :
1186 values2.domain()) {
1187 integrals(j) += spline_coef(i1, i2, j) * values1(i1) * values2(i2);
1188 }
1189 }
1190 });
1191 }
1192
1193private:
1194 /**
1195 * @brief Evaluate the function on B-splines at the coordinate given.
1196 *
1197 * This function firstly deals with the boundary conditions and calls the SplineEvaluator2D::eval_no_bc function
1198 * to evaluate.
1199 *
1200 * @param[in] coord_eval The 2D coordinate where we want to evaluate.
1201 * @param[in] spline_coef The B-splines coefficients of the function we want to evaluate.
1202 * @param[out] vals1 A ChunkSpan with the not-null values of each function of the spline in the first dimension.
1203 * @param[out] vals2 A ChunkSpan with the not-null values of each function of the spline in the second dimension.
1204 *
1205 * @return A double with the value of the function at the coordinate given.
1206 *
1207 * @see SplineBoundaryValue
1208 */
1209 template <class Layout, class... CoordsDims>
1210 KOKKOS_INLINE_FUNCTION double eval(
1211 ddc::Coordinate<CoordsDims...> coord_eval,
1212 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
1213 spline_coef) const
1214 {
1215 using Dim1 = continuous_dimension_type1;
1216 using Dim2 = continuous_dimension_type2;
1217 if constexpr (bsplines_type1::is_periodic()) {
1218 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()
1219 || ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
1220 ddc::get<Dim1>(coord_eval)
1221 -= Kokkos::floor(
1222 (ddc::get<Dim1>(coord_eval)
1223 - ddc::discrete_space<bsplines_type1>().rmin())
1224 / ddc::discrete_space<bsplines_type1>().length())
1225 * ddc::discrete_space<bsplines_type1>().length();
1226 }
1227 }
1228 if constexpr (bsplines_type2::is_periodic()) {
1229 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()
1230 || ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
1231 ddc::get<Dim2>(coord_eval)
1232 -= Kokkos::floor(
1233 (ddc::get<Dim2>(coord_eval)
1234 - ddc::discrete_space<bsplines_type2>().rmin())
1235 / ddc::discrete_space<bsplines_type2>().length())
1236 * ddc::discrete_space<bsplines_type2>().length();
1237 }
1238 }
1239 if constexpr (!bsplines_type1::is_periodic()) {
1240 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()) {
1241 return m_lower_extrap_rule_1(coord_eval, spline_coef);
1242 }
1243 if (ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
1244 return m_upper_extrap_rule_1(coord_eval, spline_coef);
1245 }
1246 }
1247 if constexpr (!bsplines_type2::is_periodic()) {
1248 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()) {
1249 return m_lower_extrap_rule_2(coord_eval, spline_coef);
1250 }
1251 if (ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
1252 return m_upper_extrap_rule_2(coord_eval, spline_coef);
1253 }
1254 }
1255 return eval_no_bc<eval_type, eval_type>(
1256 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>(
1257 ddc::get<Dim1>(coord_eval),
1258 ddc::get<Dim2>(coord_eval)),
1259 spline_coef);
1260 }
1261
1262 /**
1263 * @brief Evaluate the function or its derivative at the coordinate given.
1264 *
1265 * @param[in] coord_eval The coordinate where we want to evaluate.
1266 * @param[in] splne_coef The B-splines coefficients of the function we want to evaluate.
1267 * @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`.
1268 * @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`.
1269 */
1270 template <class EvalType1, class EvalType2, class Layout, class... CoordsDims>
1271 KOKKOS_INLINE_FUNCTION double eval_no_bc(
1272 ddc::Coordinate<CoordsDims...> const& coord_eval,
1273 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
1274 spline_coef) const
1275 {
1276 static_assert(
1277 std::is_same_v<EvalType1, eval_type> || std::is_same_v<EvalType1, eval_deriv_type>);
1278 static_assert(
1279 std::is_same_v<EvalType2, eval_type> || std::is_same_v<EvalType2, eval_deriv_type>);
1280 ddc::DiscreteElement<bsplines_type1> jmin1;
1281 ddc::DiscreteElement<bsplines_type2> jmin2;
1282
1283 std::array<double, bsplines_type1::degree() + 1> vals1_ptr;
1284 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>> const
1285 vals1(vals1_ptr.data());
1286 std::array<double, bsplines_type2::degree() + 1> vals2_ptr;
1287 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>> const
1288 vals2(vals2_ptr.data());
1289 ddc::Coordinate<continuous_dimension_type1> const coord_eval_interest1(coord_eval);
1290 ddc::Coordinate<continuous_dimension_type2> const coord_eval_interest2(coord_eval);
1291
1292 if constexpr (std::is_same_v<EvalType1, eval_type>) {
1293 jmin1 = ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
1294 } else if constexpr (std::is_same_v<EvalType1, eval_deriv_type>) {
1295 jmin1 = ddc::discrete_space<bsplines_type1>().eval_deriv(vals1, coord_eval_interest1);
1296 }
1297 if constexpr (std::is_same_v<EvalType2, eval_type>) {
1298 jmin2 = ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
1299 } else if constexpr (std::is_same_v<EvalType2, eval_deriv_type>) {
1300 jmin2 = ddc::discrete_space<bsplines_type2>().eval_deriv(vals2, coord_eval_interest2);
1301 }
1302
1303 double y = 0.0;
1304 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
1305 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
1306 y += spline_coef(
1307 ddc::DiscreteElement<
1308 bsplines_type1,
1309 bsplines_type2>(jmin1 + i, jmin2 + j))
1310 * vals1[i] * vals2[j];
1311 }
1312 }
1313 return y;
1314 }
1315};
1316
1317} // namespace ddc
friend class DiscreteDomain
A class to evaluate, differentiate or integrate a 2D spline function.
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 2D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator2D(SplineEvaluator2D &&x)=default
Move-constructs.
lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
Get the lower extrapolation rule along the first dimension.
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...
SplineEvaluator2D & operator=(SplineEvaluator2D const &x)=default
Copy-assigns.
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 2D spline function (described by its spline coefficients) at a given coordinate ...
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 2D spline function (described by its spline coefficients) at a given coordinate along s...
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 2D spline function (described by its spline coefficients) on a mesh along first dimensi...
SplineEvaluator2D(SplineEvaluator2D const &x)=default
Copy-constructs.
~SplineEvaluator2D()=default
Destructs.
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 2D spline function (described by its spline coefficients) on a mesh along specif...
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 2D spline function (described by its spline coefficients) at a given coordinate along a...
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 2D spline function (described by its spline coefficients) on a mesh along first dimensi...
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 2D spline function (described by its spline coefficients) on a mesh along specif...
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 2D spline function (described by its spline coefficients) on a mesh along dimensi...
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 2D integrations of a spline function (described by its spline coefficients) along the...
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...
upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
Get the upper extrapolation rule along the second dimension.
upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
Get the upper extrapolation rule along the first dimension.
lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
Get the lower extrapolation rule along the second dimension.
SplineEvaluator2D & operator=(SplineEvaluator2D &&x)=default
Move-assigns.
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 2D spline function (described by its spline coefficients) at a given coordinate.
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 2D spline function (described by its spline coefficients) on a mesh along 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 2D spline function (described by its spline coefficients) at a given coordinate.
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 2D 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 2D spline function (described by its spline coefficients) on a mesh.
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 2D spline function (described by its spline coefficients) at a given coordinate along f...
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 2D spline function (described by its spline coefficients) on a mesh along second dimens...
SplineEvaluator2D(LowerExtrapolationRule1 const &lower_extrap_rule1, UpperExtrapolationRule1 const &upper_extrap_rule1, LowerExtrapolationRule2 const &lower_extrap_rule2, UpperExtrapolationRule2 const &upper_extrap_rule2)
Build a SplineEvaluator2D acting on batched_spline_domain.
The top-level namespace of DDC.