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