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