DDC 0.12.0
Loading...
Searching...
No Matches
spline_evaluator_3d.hpp
1// Copyright (C) The DDC development team, see COPYRIGHT.md file
2//
3// SPDX-License-Identifier: MIT
4
5#pragma once
6
7#include <array>
8#include <cstddef>
9#include <type_traits>
10
11#include <ddc/ddc.hpp>
12
13#include <Kokkos_Core.hpp>
14
15#include "deriv.hpp"
16#include "integrals.hpp"
18
19namespace ddc {
20
21/**
22 * @brief A class to evaluate, differentiate or integrate a 3D spline function.
23 *
24 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a 3D 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 BSplines3 The discrete dimension representing the B-splines along the third dimension of interest.
31 * @tparam EvaluationDDim1 The first discrete dimension on which evaluation points are defined.
32 * @tparam EvaluationDDim2 The second discrete dimension on which evaluation points are defined.
33 * @tparam EvaluationDDim3 The third discrete dimension on which evaluation points are defined.
34 * @tparam LowerExtrapolationRule1 The lower extrapolation rule type along first dimension of interest.
35 * @tparam UpperExtrapolationRule1 The upper extrapolation rule type along first dimension of interest.
36 * @tparam LowerExtrapolationRule2 The lower extrapolation rule type along second dimension of interest.
37 * @tparam UpperExtrapolationRule2 The upper extrapolation rule type along second dimension of interest.
38 * @tparam LowerExtrapolationRule3 The lower extrapolation rule type along third dimension of interest.
39 * @tparam UpperExtrapolationRule3 The upper extrapolation rule type along third dimension of interest.
40 */
41template <
42 class ExecSpace,
43 class MemorySpace,
44 class BSplines1,
45 class BSplines2,
46 class BSplines3,
47 class EvaluationDDim1,
48 class EvaluationDDim2,
49 class EvaluationDDim3,
50 class LowerExtrapolationRule1,
51 class UpperExtrapolationRule1,
52 class LowerExtrapolationRule2,
53 class UpperExtrapolationRule2,
54 class LowerExtrapolationRule3,
55 class UpperExtrapolationRule3>
57{
58public:
59 /// @brief The type of the first evaluation continuous dimension used by this class.
60 using continuous_dimension_type1 = typename BSplines1::continuous_dimension_type;
61
62 /// @brief The type of the second evaluation continuous dimension used by this class.
63 using continuous_dimension_type2 = typename BSplines2::continuous_dimension_type;
64
65 /// @brief The type of the third evaluation continuous dimension used by this class.
66 using continuous_dimension_type3 = typename BSplines3::continuous_dimension_type;
67
68 /// @brief The type of the Kokkos execution space used by this class.
69 using exec_space = ExecSpace;
70
71 /// @brief The type of the Kokkos memory space used by this class.
72 using memory_space = MemorySpace;
73
74 /// @brief The type of the first discrete dimension of interest used by this class.
75 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
76
77 /// @brief The type of the second discrete dimension of interest used by this class.
78 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
79
80 /// @brief The type of the third discrete dimension of interest used by this class.
81 using evaluation_discrete_dimension_type3 = EvaluationDDim3;
82
83 /// @brief The discrete dimension representing the B-splines along first dimension.
84 using bsplines_type1 = BSplines1;
85
86 /// @brief The discrete dimension representing the B-splines along second dimension.
87 using bsplines_type2 = BSplines2;
88
89 /// @brief The discrete dimension representing the B-splines along third dimension.
90 using bsplines_type3 = BSplines3;
91
92 /// @brief The type of the domain for the 1D evaluation mesh along first dimension used by this class.
93 using evaluation_domain_type1 = ddc::DiscreteDomain<evaluation_discrete_dimension_type1>;
94
95 /// @brief The type of the domain for the 1D evaluation mesh along second dimension used by this class.
96 using evaluation_domain_type2 = ddc::DiscreteDomain<evaluation_discrete_dimension_type2>;
97
98 /// @brief The type of the domain for the 1D evaluation mesh along third dimension used by this class.
99 using evaluation_domain_type3 = ddc::DiscreteDomain<evaluation_discrete_dimension_type3>;
100
101 /// @brief The type of the domain for the 3D evaluation mesh used by this class.
102 using evaluation_domain_type = ddc::DiscreteDomain<
103 evaluation_discrete_dimension_type1,
104 evaluation_discrete_dimension_type2,
105 evaluation_discrete_dimension_type3>;
106
107 /**
108 * @brief The type of the whole domain representing evaluation points.
109 *
110 * @tparam The batched discrete domain on which the interpolation points are defined.
111 */
112 template <concepts::discrete_domain BatchedInterpolationDDom>
113 using batched_evaluation_domain_type = BatchedInterpolationDDom;
114
115 /// @brief The type of the 1D spline domain corresponding to the first dimension of interest.
116 using spline_domain_type1 = ddc::DiscreteDomain<bsplines_type1>;
117
118 /// @brief The type of the 1D spline domain corresponding to the second dimension of interest.
119 using spline_domain_type2 = ddc::DiscreteDomain<bsplines_type2>;
120
121 /// @brief The type of the 1D spline domain corresponding to the third dimension of interest.
122 using spline_domain_type3 = ddc::DiscreteDomain<bsplines_type3>;
123
124 /// @brief The type of the 3D spline domain corresponding to the dimensions of interest.
125 using spline_domain_type = ddc::DiscreteDomain<bsplines_type1, bsplines_type2, bsplines_type3>;
126
127 /**
128 * @brief The type of the batch domain (obtained by removing the dimensions of interest
129 * from the whole domain).
130 *
131 * @tparam The batched discrete domain on which the interpolation points are defined.
132 */
133 template <concepts::discrete_domain BatchedInterpolationDDom>
134 using batch_domain_type = typename ddc::remove_dims_of_t<
135 BatchedInterpolationDDom,
136 evaluation_discrete_dimension_type1,
137 evaluation_discrete_dimension_type2,
138 evaluation_discrete_dimension_type3>;
139
140 /**
141 * @brief The type of the whole spline domain (cartesian product of 3D spline domain
142 * and batch domain) preserving the underlying memory layout (order of dimensions).
143 *
144 * @tparam The batched discrete domain on which the interpolation points are defined.
145 */
146 template <concepts::discrete_domain BatchedInterpolationDDom>
147 using batched_spline_domain_type =
148 typename ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
149 ddc::to_type_seq_t<BatchedInterpolationDDom>,
150 ddc::detail::TypeSeq<
151 evaluation_discrete_dimension_type1,
152 evaluation_discrete_dimension_type2,
153 evaluation_discrete_dimension_type3>,
154 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>>>;
155
156 /// @brief The type of the extrapolation rule at the lower boundary along the first dimension.
157 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
158
159 /// @brief The type of the extrapolation rule at the upper boundary along the first dimension.
160 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
161
162 /// @brief The type of the extrapolation rule at the lower boundary along the second dimension.
163 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
164
165 /// @brief The type of the extrapolation rule at the upper boundary along the second dimension.
166 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
167
168 /// @brief The type of the extrapolation rule at the lower boundary along the third dimension.
169 using lower_extrapolation_rule_3_type = LowerExtrapolationRule3;
170
171 /// @brief The type of the extrapolation rule at the upper boundary along the third dimension.
172 using upper_extrapolation_rule_3_type = UpperExtrapolationRule3;
173
174private:
175 LowerExtrapolationRule1 m_lower_extrap_rule_1;
176
177 UpperExtrapolationRule1 m_upper_extrap_rule_1;
178
179 LowerExtrapolationRule2 m_lower_extrap_rule_2;
180
181 UpperExtrapolationRule2 m_upper_extrap_rule_2;
182
183 LowerExtrapolationRule3 m_lower_extrap_rule_3;
184
185 UpperExtrapolationRule3 m_upper_extrap_rule_3;
186
187public:
188 static_assert(
189 std::is_same_v<LowerExtrapolationRule1,
190 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
191 == bsplines_type1::is_periodic()
192 && std::is_same_v<
193 UpperExtrapolationRule1,
194 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type1>>
195 == bsplines_type1::is_periodic()
196 && std::is_same_v<
197 LowerExtrapolationRule2,
198 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
199 == bsplines_type2::is_periodic()
200 && std::is_same_v<
201 UpperExtrapolationRule2,
202 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type2>>
203 == bsplines_type2::is_periodic()
204 && std::is_same_v<
205 LowerExtrapolationRule3,
206 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
207 == bsplines_type3::is_periodic()
208 && std::is_same_v<
209 UpperExtrapolationRule3,
210 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type3>>
211 == bsplines_type3::is_periodic(),
212 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
213 static_assert(
214 std::is_invocable_r_v<
215 double,
216 LowerExtrapolationRule1,
217 ddc::Coordinate<continuous_dimension_type1>,
218 ddc::ChunkSpan<
219 double const,
220 spline_domain_type,
221 Kokkos::layout_right,
222 memory_space>>,
223 "LowerExtrapolationRule1::operator() has to be callable "
224 "with usual arguments.");
225 static_assert(
226 std::is_invocable_r_v<
227 double,
228 UpperExtrapolationRule1,
229 ddc::Coordinate<continuous_dimension_type1>,
230 ddc::ChunkSpan<
231 double const,
232 spline_domain_type,
233 Kokkos::layout_right,
234 memory_space>>,
235 "UpperExtrapolationRule1::operator() has to be callable "
236 "with usual arguments.");
237 static_assert(
238 std::is_invocable_r_v<
239 double,
240 LowerExtrapolationRule2,
241 ddc::Coordinate<continuous_dimension_type2>,
242 ddc::ChunkSpan<
243 double const,
244 spline_domain_type,
245 Kokkos::layout_right,
246 memory_space>>,
247 "LowerExtrapolationRule2::operator() has to be callable "
248 "with usual arguments.");
249 static_assert(
250 std::is_invocable_r_v<
251 double,
252 UpperExtrapolationRule2,
253 ddc::Coordinate<continuous_dimension_type2>,
254 ddc::ChunkSpan<
255 double const,
256 spline_domain_type,
257 Kokkos::layout_right,
258 memory_space>>,
259 "UpperExtrapolationRule2::operator() has to be callable "
260 "with usual arguments.");
261 static_assert(
262 std::is_invocable_r_v<
263 double,
264 LowerExtrapolationRule3,
265 ddc::Coordinate<continuous_dimension_type3>,
266 ddc::ChunkSpan<
267 double const,
268 spline_domain_type,
269 Kokkos::layout_right,
270 memory_space>>,
271 "LowerExtrapolationRule3::operator() has to be callable "
272 "with usual arguments.");
273 static_assert(
274 std::is_invocable_r_v<
275 double,
276 UpperExtrapolationRule3,
277 ddc::Coordinate<continuous_dimension_type3>,
278 ddc::ChunkSpan<
279 double const,
280 spline_domain_type,
281 Kokkos::layout_right,
282 memory_space>>,
283 "UpperExtrapolationRule3::operator() has to be callable "
284 "with usual arguments.");
285
286 /**
287 * @brief Build a SplineEvaluator3D acting on batched_spline_domain.
288 *
289 * @param lower_extrap_rule1 The extrapolation rule at the lower boundary along the first dimension.
290 * @param upper_extrap_rule1 The extrapolation rule at the upper boundary along the first dimension.
291 * @param lower_extrap_rule2 The extrapolation rule at the lower boundary along the second dimension.
292 * @param upper_extrap_rule2 The extrapolation rule at the upper boundary along the second dimension.
293 * @param lower_extrap_rule3 The extrapolation rule at the lower boundary along the third dimension.
294 * @param upper_extrap_rule3 The extrapolation rule at the upper boundary along the third dimension.
295 *
296 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
297 */
298 explicit SplineEvaluator3D(
299 LowerExtrapolationRule1 const& lower_extrap_rule1,
300 UpperExtrapolationRule1 const& upper_extrap_rule1,
301 LowerExtrapolationRule2 const& lower_extrap_rule2,
302 UpperExtrapolationRule2 const& upper_extrap_rule2,
303 LowerExtrapolationRule3 const& lower_extrap_rule3,
304 UpperExtrapolationRule3 const& upper_extrap_rule3)
305 : m_lower_extrap_rule_1(lower_extrap_rule1)
306 , m_upper_extrap_rule_1(upper_extrap_rule1)
307 , m_lower_extrap_rule_2(lower_extrap_rule2)
308 , m_upper_extrap_rule_2(upper_extrap_rule2)
309 , m_lower_extrap_rule_3(lower_extrap_rule3)
310 , m_upper_extrap_rule_3(upper_extrap_rule3)
311 {
312 }
313
314 /**
315 * @brief Copy-constructs.
316 *
317 * @param x A reference to another SplineEvaluator.
318 */
319 SplineEvaluator3D(SplineEvaluator3D const& x) = default;
320
321 /**
322 * @brief Move-constructs.
323 *
324 * @param x An rvalue to another SplineEvaluator.
325 */
327
328 /// @brief Destructs.
329 ~SplineEvaluator3D() = default;
330
331 /**
332 * @brief Copy-assigns.
333 *
334 * @param x A reference to another SplineEvaluator.
335 * @return A reference to this object.
336 */
337 SplineEvaluator3D& operator=(SplineEvaluator3D const& x) = default;
338
339 /**
340 * @brief Move-assigns.
341 *
342 * @param x An rvalue to another SplineEvaluator.
343 * @return A reference to this object.
344 */
346
347 /**
348 * @brief Get the lower extrapolation rule along the first dimension.
349 *
350 * 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.
351 *
352 * @return The lower extrapolation rule along the first dimension.
353 *
354 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
355 */
356 lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
357 {
358 return m_lower_extrap_rule_1;
359 }
360
361 /**
362 * @brief Get the upper extrapolation rule along the first dimension.
363 *
364 * 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.
365 *
366 * @return The upper extrapolation rule along the first dimension.
367 *
368 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
369 */
370 upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
371 {
372 return m_upper_extrap_rule_1;
373 }
374
375 /**
376 * @brief Get the lower extrapolation rule along the second dimension.
377 *
378 * 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.
379 *
380 * @return The lower extrapolation rule along the second dimension.
381 *
382 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
383 */
384 lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
385 {
386 return m_lower_extrap_rule_2;
387 }
388
389 /**
390 * @brief Get the upper extrapolation rule along the second dimension.
391 *
392 * 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.
393 *
394 * @return The upper extrapolation rule along the second dimension.
395 *
396 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
397 */
398 upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
399 {
400 return m_upper_extrap_rule_2;
401 }
402
403 /**
404 * @brief Get the lower extrapolation rule along the third dimension.
405 *
406 * 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.
407 *
408 * @return The lower extrapolation rule along the third dimension.
409 *
410 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
411 */
412 lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
413 {
414 return m_lower_extrap_rule_3;
415 }
416
417 /**
418 * @brief Get the upper extrapolation rule along the third dimension.
419 *
420 * 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.
421 *
422 * @return The upper extrapolation rule along the third dimension.
423 *
424 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
425 */
426 upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
427 {
428 return m_upper_extrap_rule_3;
429 }
430
431 /**
432 * @brief Evaluate 3D spline function (described by its spline coefficients) at a given coordinate.
433 *
434 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
435 *
436 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
437 *
438 * @param coord_eval The coordinate where the spline is evaluated. Note that only the components along the dimensions of interest are used.
439 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
440 *
441 * @return The value of the spline function at the desired coordinate.
442 */
443 template <class Layout, class... CoordsDims>
444 KOKKOS_FUNCTION double operator()(
445 ddc::Coordinate<CoordsDims...> const& coord_eval,
446 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
447 spline_coef) const
448 {
449 return eval(coord_eval, spline_coef);
450 }
451
452 /**
453 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
454 *
455 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
456 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
457 *
458 * This is not a nD evaluation. This is a batched 3D evaluation. This means that for each slice of coordinates
459 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 3D set of
460 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
461 *
462 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
463 *
464 * @param[out] spline_eval The values of the 3D spline function at the desired coordinates. For practical reasons those are
465 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
466 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
467 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
468 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
469 * the set of 3D spline coefficients retained to perform the evaluation).
470 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
471 */
472 template <
473 class Layout1,
474 class Layout2,
475 class Layout3,
476 class BatchedInterpolationDDom,
477 class... CoordsDims>
478 void operator()(
479 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
480 spline_eval,
481 ddc::ChunkSpan<
482 ddc::Coordinate<CoordsDims...> const,
483 BatchedInterpolationDDom,
484 Layout2,
485 memory_space> const coords_eval,
486 ddc::ChunkSpan<
487 double const,
488 batched_spline_domain_type<BatchedInterpolationDDom>,
489 Layout3,
490 memory_space> const spline_coef) const
491 {
492 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
493 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
494 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
495 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
496 ddc::parallel_for_each(
497 "ddc_splines_evaluate_3d",
498 exec_space(),
499 batch_domain,
500 KOKKOS_CLASS_LAMBDA(
501 typename batch_domain_type<
502 BatchedInterpolationDDom>::discrete_element_type const j) {
503 auto const spline_eval_3D = spline_eval[j];
504 auto const coords_eval_3D = coords_eval[j];
505 auto const spline_coef_3D = spline_coef[j];
506 for (auto const i1 : evaluation_domain1) {
507 for (auto const i2 : evaluation_domain2) {
508 for (auto const i3 : evaluation_domain3) {
509 spline_eval_3D(i1, i2, i3)
510 = eval(coords_eval_3D(i1, i2, i3), spline_coef_3D);
511 }
512 }
513 }
514 });
515 }
516
517 /**
518 * @brief Evaluate 3D spline function (described by its spline coefficients) on a mesh.
519 *
520 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
521 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
522 *
523 * This is not a multidimensional evaluation. This is a batched 3D evaluation.
524 * This means that for each slice of spline_eval the evaluation is performed with
525 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
526 *
527 * Remark: calling SplineBuilder3D then SplineEvaluator3D corresponds to a 3D spline interpolation.
528 *
529 * @param[out] spline_eval The values of the 3D spline function at their coordinates.
530 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
531 */
532 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
533 void operator()(
534 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
535 spline_eval,
536 ddc::ChunkSpan<
537 double const,
538 batched_spline_domain_type<BatchedInterpolationDDom>,
539 Layout2,
540 memory_space> const spline_coef) const
541 {
542 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
543 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
544 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
545 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
546 ddc::parallel_for_each(
547 "ddc_splines_evaluate_3d",
548 exec_space(),
549 batch_domain,
550 KOKKOS_CLASS_LAMBDA(
551 typename batch_domain_type<
552 BatchedInterpolationDDom>::discrete_element_type const j) {
553 auto const spline_eval_3D = spline_eval[j];
554 auto const spline_coef_3D = spline_coef[j];
555 for (auto const i1 : evaluation_domain1) {
556 for (auto const i2 : evaluation_domain2) {
557 for (auto const i3 : evaluation_domain3) {
558 ddc::Coordinate<
559 continuous_dimension_type1,
560 continuous_dimension_type2,
561 continuous_dimension_type3>
562 coord_eval_3D(
563 ddc::coordinate(i1),
564 ddc::coordinate(i2),
565 ddc::coordinate(i3));
566 spline_eval_3D(i1, i2, i3)
567 = eval(coord_eval_3D(i1, i2, i3), spline_coef_3D);
568 }
569 }
570 }
571 });
572 }
573
574 /**
575 * @brief Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along the dimensions of interest.
576 *
577 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be
578 * obtained via various methods, such as using a SplineBuilder3D.
579 *
580 * @param deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
581 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
582 * @param coord_eval The coordinate where the spline is differentiated. Note that only the components along the dimensions of interest are used.
583 * @param spline_coef A ChunkSpan storing the 3D spline coefficients.
584 *
585 * @return The derivative of the spline function at the desired coordinate.
586 */
587 template <class DElem, class Layout, class... CoordsDims>
588 KOKKOS_FUNCTION double deriv(
589 DElem const& deriv_order,
590 ddc::Coordinate<CoordsDims...> const& coord_eval,
591 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
592 spline_coef) const
593 {
594 static_assert(is_discrete_element_v<DElem>);
595 return eval_no_bc(deriv_order, coord_eval, spline_coef);
596 }
597
598 /**
599 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
600 *
601 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
602 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
603 *
604 * This is not a nD differentiation. This is a batched 3D differentiation.
605 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
606 * the differentiation is performed with the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
607 *
608 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
609 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
610 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates. For practical reasons those are
611 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
612 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
613 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
614 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
615 * the set of 3D spline coefficients retained to perform the evaluation).
616 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
617 */
618 template <
619 class DElem,
620 class Layout1,
621 class Layout2,
622 class Layout3,
623 class BatchedInterpolationDDom,
624 class... CoordsDims>
625 void deriv(
626 DElem const& deriv_order,
627 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
628 spline_eval,
629 ddc::ChunkSpan<
630 ddc::Coordinate<CoordsDims...> const,
631 BatchedInterpolationDDom,
632 Layout2,
633 memory_space> const coords_eval,
634 ddc::ChunkSpan<
635 double const,
636 batched_spline_domain_type<BatchedInterpolationDDom>,
637 Layout3,
638 memory_space> const spline_coef) const
639 {
640 static_assert(is_discrete_element_v<DElem>);
641
642 batch_domain_type<BatchedInterpolationDDom> const batch_domain(coords_eval.domain());
643 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
644 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
645 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
646 ddc::parallel_for_each(
647 "ddc_splines_cross_differentiate_3d",
648 exec_space(),
649 batch_domain,
650 KOKKOS_CLASS_LAMBDA(
651 typename batch_domain_type<
652 BatchedInterpolationDDom>::discrete_element_type const j) {
653 auto const spline_eval_3D = spline_eval[j];
654 auto const coords_eval_3D = coords_eval[j];
655 auto const spline_coef_3D = spline_coef[j];
656 for (auto const i1 : evaluation_domain1) {
657 for (auto const i2 : evaluation_domain2) {
658 for (auto const i3 : evaluation_domain3) {
659 spline_eval_3D(i1, i2, i3) = eval_no_bc(
660 deriv_order,
661 coords_eval_3D(i1, i2, i3),
662 spline_coef_3D);
663 }
664 }
665 }
666 });
667 }
668
669 /**
670 * @brief Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimensions of interest.
671 *
672 * The spline coefficients represent a 3D spline function defined on a cartesian product of batch_domain and B-splines
673 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
674 *
675 * This is not a multidimensional differentiation. This is a batched 3D differentiation.
676 * This means that for each slice of spline_eval the differentiation is performed with
677 * the 3D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
678 *
679 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
680 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
681 * @param[out] spline_eval The derivatives of the 3D spline function at the desired coordinates.
682 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
683 */
684 template <class DElem, class Layout1, class Layout2, class BatchedInterpolationDDom>
685 void deriv(
686 DElem const& deriv_order,
687 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
688 spline_eval,
689 ddc::ChunkSpan<
690 double const,
691 batched_spline_domain_type<BatchedInterpolationDDom>,
692 Layout2,
693 memory_space> const spline_coef) const
694 {
695 static_assert(is_discrete_element_v<DElem>);
696
697 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
698 evaluation_domain_type1 const evaluation_domain1(spline_eval.domain());
699 evaluation_domain_type2 const evaluation_domain2(spline_eval.domain());
700 evaluation_domain_type3 const evaluation_domain3(spline_eval.domain());
701 ddc::parallel_for_each(
702 "ddc_splines_cross_differentiate_3d",
703 exec_space(),
704 batch_domain,
705 KOKKOS_CLASS_LAMBDA(
706 typename batch_domain_type<
707 BatchedInterpolationDDom>::discrete_element_type const j) {
708 auto const spline_eval_3D = spline_eval[j];
709 auto const spline_coef_3D = spline_coef[j];
710 for (auto const i1 : evaluation_domain1) {
711 for (auto const i2 : evaluation_domain2) {
712 for (auto const i3 : evaluation_domain3) {
713 ddc::Coordinate<
714 continuous_dimension_type1,
715 continuous_dimension_type2,
716 continuous_dimension_type3>
717 coord_eval_3D(
718 ddc::coordinate(i1),
719 ddc::coordinate(i2),
720 ddc::coordinate(i3));
721 spline_eval_3D(i1, i2, i3)
722 = eval_no_bc(deriv_order, coord_eval_3D, spline_coef_3D);
723 }
724 }
725 }
726 });
727 }
728
729 /** @brief Perform batched 3D integrations of a spline function (described by its spline coefficients) along the dimensions of interest and store results on a subdomain of batch_domain.
730 *
731 * The spline coefficients represent a 3D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder3D.
732 *
733 * This is not a nD integration. This is a batched 3D integration.
734 * This means that for each element of integrals, the integration is performed with the 3D set of
735 * spline coefficients identified by the same DiscreteElement.
736 *
737 * @param[out] integrals The integrals of the 3D spline function on the subdomain of batch_domain. For practical reasons those are
738 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
739 * points represented by this domain are unused and irrelevant.
740 * @param[in] spline_coef A ChunkSpan storing the 3D spline coefficients.
741 */
742 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
743 void integrate(
744 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
745 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
746 spline_coef) const
747 {
748 static_assert(
749 ddc::type_seq_contains_v<
750 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>,
751 to_type_seq_t<BatchedSplineDDom>>,
752 "The spline coefficients domain must contain the bsplines dimensions");
753 using batch_domain_type = ddc::
754 remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2, bsplines_type3>;
755 static_assert(
756 std::is_same_v<batch_domain_type, BatchedDDom>,
757 "The integrals domain must only contain the batch dimensions");
758
759 batch_domain_type batch_domain(integrals.domain());
760 ddc::Chunk values1_alloc(
761 ddc::DiscreteDomain<bsplines_type1>(spline_coef.domain()),
762 ddc::KokkosAllocator<double, memory_space>());
763 ddc::ChunkSpan values1 = values1_alloc.span_view();
764 ddc::integrals(exec_space(), values1);
765 ddc::Chunk values2_alloc(
766 ddc::DiscreteDomain<bsplines_type2>(spline_coef.domain()),
767 ddc::KokkosAllocator<double, memory_space>());
768 ddc::ChunkSpan values2 = values2_alloc.span_view();
769 ddc::integrals(exec_space(), values2);
770 ddc::Chunk values3_alloc(
771 ddc::DiscreteDomain<bsplines_type3>(spline_coef.domain()),
772 ddc::KokkosAllocator<double, memory_space>());
773 ddc::ChunkSpan values3 = values3_alloc.span_view();
774 ddc::integrals(exec_space(), values3);
775
776 ddc::parallel_for_each(
777 "ddc_splines_integrate_bsplines",
778 exec_space(),
779 batch_domain,
780 KOKKOS_LAMBDA(typename batch_domain_type::discrete_element_type const j) {
781 integrals(j) = 0;
782 for (typename spline_domain_type1::discrete_element_type const i1 :
783 values1.domain()) {
784 for (typename spline_domain_type2::discrete_element_type const i2 :
785 values2.domain()) {
786 for (typename spline_domain_type3::discrete_element_type const i3 :
787 values3.domain()) {
788 integrals(j) += spline_coef(i1, i2, i3, j) * values1(i1)
789 * values2(i2) * values3(i3);
790 }
791 }
792 }
793 });
794 }
795
796private:
797 /**
798 * @brief Evaluate the function on B-splines at the coordinate given.
799 *
800 * This function firstly deals with the boundary conditions and calls the SplineEvaluator3D::eval_no_bc function
801 * to evaluate.
802 *
803 * @param[in] coord_eval The 3D coordinate where we want to evaluate.
804 * @param[in] spline_coef The B-splines coefficients of the function we want to evaluate.
805 * @param[out] vals1 A ChunkSpan with the not-null values of each function of the spline in the first dimension.
806 * @param[out] vals2 A ChunkSpan with the not-null values of each function of the spline in the second dimension.
807 *
808 * @return A double with the value of the function at the coordinate given.
809 *
810 * @see SplineBoundaryValue
811 */
812 template <class Layout, class... CoordsDims>
813 KOKKOS_INLINE_FUNCTION double eval(
814 ddc::Coordinate<CoordsDims...> coord_eval,
815 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
816 spline_coef) const
817 {
818 using Dim1 = continuous_dimension_type1;
819 using Dim2 = continuous_dimension_type2;
820 using Dim3 = continuous_dimension_type3;
821 if constexpr (bsplines_type1::is_periodic()) {
822 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()
823 || ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
824 ddc::get<Dim1>(coord_eval)
825 -= Kokkos::floor(
826 (ddc::get<Dim1>(coord_eval)
827 - ddc::discrete_space<bsplines_type1>().rmin())
828 / ddc::discrete_space<bsplines_type1>().length())
829 * ddc::discrete_space<bsplines_type1>().length();
830 }
831 }
832 if constexpr (bsplines_type2::is_periodic()) {
833 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()
834 || ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
835 ddc::get<Dim2>(coord_eval)
836 -= Kokkos::floor(
837 (ddc::get<Dim2>(coord_eval)
838 - ddc::discrete_space<bsplines_type2>().rmin())
839 / ddc::discrete_space<bsplines_type2>().length())
840 * ddc::discrete_space<bsplines_type2>().length();
841 }
842 }
843 if constexpr (bsplines_type3::is_periodic()) {
844 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()
845 || ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
846 ddc::get<Dim3>(coord_eval)
847 -= Kokkos::floor(
848 (ddc::get<Dim3>(coord_eval)
849 - ddc::discrete_space<bsplines_type3>().rmin())
850 / ddc::discrete_space<bsplines_type3>().length())
851 * ddc::discrete_space<bsplines_type3>().length();
852 }
853 }
854 if constexpr (!bsplines_type1::is_periodic()) {
855 if (ddc::get<Dim1>(coord_eval) < ddc::discrete_space<bsplines_type1>().rmin()) {
856 return m_lower_extrap_rule_1(coord_eval, spline_coef);
857 }
858 if (ddc::get<Dim1>(coord_eval) > ddc::discrete_space<bsplines_type1>().rmax()) {
859 return m_upper_extrap_rule_1(coord_eval, spline_coef);
860 }
861 }
862 if constexpr (!bsplines_type2::is_periodic()) {
863 if (ddc::get<Dim2>(coord_eval) < ddc::discrete_space<bsplines_type2>().rmin()) {
864 return m_lower_extrap_rule_2(coord_eval, spline_coef);
865 }
866 if (ddc::get<Dim2>(coord_eval) > ddc::discrete_space<bsplines_type2>().rmax()) {
867 return m_upper_extrap_rule_2(coord_eval, spline_coef);
868 }
869 }
870 if constexpr (!bsplines_type3::is_periodic()) {
871 if (ddc::get<Dim3>(coord_eval) < ddc::discrete_space<bsplines_type3>().rmin()) {
872 return m_lower_extrap_rule_3(coord_eval, spline_coef);
873 }
874 if (ddc::get<Dim3>(coord_eval) > ddc::discrete_space<bsplines_type3>().rmax()) {
875 return m_upper_extrap_rule_3(coord_eval, spline_coef);
876 }
877 }
878 return eval_no_bc(
879 ddc::DiscreteElement<>(),
880 ddc::Coordinate<
881 continuous_dimension_type1,
882 continuous_dimension_type2,
883 continuous_dimension_type3>(
884 ddc::get<Dim1>(coord_eval),
885 ddc::get<Dim2>(coord_eval),
886 ddc::get<Dim3>(coord_eval)),
887 spline_coef);
888 }
889
890 /**
891 * @brief Evaluate the function or its derivative at the coordinate given.
892 *
893 * @param[in] deriv_order A DiscreteElement containing the orders of derivation for each of the dimensions of interest.
894 * If one of the dimensions is not present, its corresponding order of derivation is considered to be 0.
895 * @param[in] coord_eval The coordinate where we want to evaluate.
896 * @param[in] splne_coef The B-splines coefficients of the function we want to evaluate.
897 */
898 template <class... DerivDims, class Layout, class... CoordsDims>
899 KOKKOS_INLINE_FUNCTION double eval_no_bc(
900 ddc::DiscreteElement<DerivDims...> const& deriv_order,
901 ddc::Coordinate<CoordsDims...> const& coord_eval,
902 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
903 spline_coef) const
904 {
905 using deriv_dim1 = Deriv<continuous_dimension_type1>;
906 using deriv_dim2 = Deriv<continuous_dimension_type2>;
907 using deriv_dim3 = Deriv<continuous_dimension_type3>;
908 using deriv_dims = detail::TypeSeq<DerivDims...>;
909
910 // Check that the tags are valid
911 static_assert(
912 (in_tags_v<DerivDims, ddc::detail::TypeSeq<deriv_dim1, deriv_dim2, deriv_dim3>>
913 && ...),
914 "The only valid dimensions for deriv_order are Deriv<Dim1>, Deriv<Dim2> and "
915 "Deriv<Dim3>");
916
917 ddc::DiscreteElement<bsplines_type1> jmin1;
918 ddc::DiscreteElement<bsplines_type2> jmin2;
919 ddc::DiscreteElement<bsplines_type3> jmin3;
920
921 std::array<double, bsplines_type1::degree() + 1> vals1_ptr;
922 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>> const
923 vals1(vals1_ptr.data());
924 std::array<double, bsplines_type2::degree() + 1> vals2_ptr;
925 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>> const
926 vals2(vals2_ptr.data());
927 std::array<double, bsplines_type3::degree() + 1> vals3_ptr;
928 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type3::degree() + 1>> const
929 vals3(vals3_ptr.data());
930 ddc::Coordinate<continuous_dimension_type1> const coord_eval_interest1(coord_eval);
931 ddc::Coordinate<continuous_dimension_type2> const coord_eval_interest2(coord_eval);
932 ddc::Coordinate<continuous_dimension_type3> const coord_eval_interest3(coord_eval);
933
934 if constexpr (!in_tags_v<deriv_dim1, deriv_dims>) {
935 jmin1 = ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
936 } else {
937 auto const order1 = deriv_order.template uid<deriv_dim1>();
938 KOKKOS_ASSERT(order1 > 0 && order1 <= bsplines_type1::degree())
939
940 std::array<double, (bsplines_type1::degree() + 1) * (bsplines_type1::degree() + 1)>
941 derivs1_ptr;
942 Kokkos::mdspan<
943 double,
944 Kokkos::extents<
945 std::size_t,
946 bsplines_type1::degree() + 1,
947 Kokkos::dynamic_extent>> const derivs1(derivs1_ptr.data(), order1 + 1);
948
949 jmin1 = ddc::discrete_space<bsplines_type1>()
950 .eval_basis_and_n_derivs(derivs1, coord_eval_interest1, order1);
951
952 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
953 vals1[i] = DDC_MDSPAN_ACCESS_OP(derivs1, i, order1);
954 }
955 }
956
957 if constexpr (!in_tags_v<deriv_dim2, deriv_dims>) {
958 jmin2 = ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
959 } else {
960 auto const order2 = deriv_order.template uid<deriv_dim2>();
961 KOKKOS_ASSERT(order2 > 0 && order2 <= bsplines_type2::degree())
962
963 std::array<double, (bsplines_type2::degree() + 1) * (bsplines_type2::degree() + 1)>
964 derivs2_ptr;
965 Kokkos::mdspan<
966 double,
967 Kokkos::extents<
968 std::size_t,
969 bsplines_type2::degree() + 1,
970 Kokkos::dynamic_extent>> const derivs2(derivs2_ptr.data(), order2 + 1);
971
972 jmin2 = ddc::discrete_space<bsplines_type2>()
973 .eval_basis_and_n_derivs(derivs2, coord_eval_interest2, order2);
974
975 for (std::size_t i = 0; i < bsplines_type2::degree() + 1; ++i) {
976 vals2[i] = DDC_MDSPAN_ACCESS_OP(derivs2, i, order2);
977 }
978 }
979
980 if constexpr (!in_tags_v<deriv_dim3, deriv_dims>) {
981 jmin3 = ddc::discrete_space<bsplines_type3>().eval_basis(vals3, coord_eval_interest3);
982 } else {
983 auto const order3 = deriv_order.template uid<deriv_dim3>();
984 KOKKOS_ASSERT(order3 > 0 && order3 <= bsplines_type3::degree())
985
986 std::array<double, (bsplines_type3::degree() + 1) * (bsplines_type3::degree() + 1)>
987 derivs3_ptr;
988 Kokkos::mdspan<
989 double,
990 Kokkos::extents<
991 std::size_t,
992 bsplines_type3::degree() + 1,
993 Kokkos::dynamic_extent>> const derivs3(derivs3_ptr.data(), order3 + 1);
994
995 jmin3 = ddc::discrete_space<bsplines_type3>()
996 .eval_basis_and_n_derivs(derivs3, coord_eval_interest3, order3);
997
998 for (std::size_t i = 0; i < bsplines_type3::degree() + 1; ++i) {
999 vals3[i] = DDC_MDSPAN_ACCESS_OP(derivs3, i, order3);
1000 }
1001 }
1002
1003 double y = 0.0;
1004 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
1005 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
1006 for (std::size_t k = 0; k < bsplines_type3::degree() + 1; ++k) {
1007 y += spline_coef(
1008 ddc::DiscreteElement<
1009 bsplines_type1,
1010 bsplines_type2,
1011 bsplines_type3>(jmin1 + i, jmin2 + j, jmin3 + k))
1012 * vals1[i] * vals2[j] * vals3[k];
1013 }
1014 }
1015 }
1016 return y;
1017 }
1018};
1019
1020} // 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 3D spline function.
upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
Get the upper extrapolation rule along the first dimension.
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 3D spline function (described by its spline coefficients) on a mesh along the dimension...
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator3D & operator=(SplineEvaluator3D &&x)=default
Move-assigns.
upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
Get the upper extrapolation rule along the third dimension.
SplineEvaluator3D(LowerExtrapolationRule1 const &lower_extrap_rule1, UpperExtrapolationRule1 const &upper_extrap_rule1, LowerExtrapolationRule2 const &lower_extrap_rule2, UpperExtrapolationRule2 const &upper_extrap_rule2, LowerExtrapolationRule3 const &lower_extrap_rule3, UpperExtrapolationRule3 const &upper_extrap_rule3)
Build a SplineEvaluator3D acting on batched_spline_domain.
~SplineEvaluator3D()=default
Destructs.
lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
Get the lower extrapolation rule along the first dimension.
upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
Get the upper extrapolation rule along the second dimension.
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 3D spline function (described by its spline coefficients) on a mesh along the dimension...
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 3D spline function (described by its spline coefficients) at a given coordinate along t...
SplineEvaluator3D(SplineEvaluator3D const &x)=default
Copy-constructs.
KOKKOS_FUNCTION double operator()(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) at a given coordinate.
lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
Get the lower extrapolation rule along the third dimension.
SplineEvaluator3D(SplineEvaluator3D &&x)=default
Move-constructs.
void integrate(ddc::ChunkSpan< double, BatchedDDom, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, BatchedSplineDDom, Layout2, memory_space > const spline_coef) const
Perform batched 3D integrations of a spline function (described by its spline coefficients) along the...
lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
Get the lower extrapolation rule along the second dimension.
SplineEvaluator3D & operator=(SplineEvaluator3D const &x)=default
Copy-assigns.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) on a mesh.
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