DDC 0.12.0
Loading...
Searching...
No Matches
spline_evaluator.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 spline function.
23 *
24 * A class which contains an operator () which can be used to evaluate, differentiate or integrate a 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 BSplines The discrete dimension representing the B-splines.
29 * @tparam EvaluationDDim The discrete dimension on which evaluation points are defined.
30 * @tparam LowerExtrapolationRule The lower extrapolation rule type.
31 * @tparam UpperExtrapolationRule The upper extrapolation rule type.
32 */
33template <
34 class ExecSpace,
35 class MemorySpace,
36 class BSplines,
37 class EvaluationDDim,
38 class LowerExtrapolationRule,
39 class UpperExtrapolationRule>
41{
42public:
43 /// @brief The type of the Kokkos execution space used by this class.
44 using exec_space = ExecSpace;
45
46 /// @brief The type of the Kokkos memory space used by this class.
47 using memory_space = MemorySpace;
48
49 /// @brief The type of the evaluation continuous dimension (continuous dimension of interest) used by this class.
50 using continuous_dimension_type = typename BSplines::continuous_dimension_type;
51
52 /// @brief The type of the evaluation discrete dimension (discrete dimension of interest) used by this class.
53 using evaluation_discrete_dimension_type = EvaluationDDim;
54
55 /// @brief The discrete dimension representing the B-splines.
56 using bsplines_type = BSplines;
57
58 /// @brief The type of the domain for the 1D evaluation mesh used by this class.
59 using evaluation_domain_type = ddc::DiscreteDomain<evaluation_discrete_dimension_type>;
60
61 /**
62 * @brief The type of the whole domain representing evaluation points.
63 *
64 * @tparam The batched discrete domain on which the interpolation points are defined.
65 */
66 template <concepts::discrete_domain BatchedInterpolationDDom>
67 using batched_evaluation_domain_type = BatchedInterpolationDDom;
68
69 /// @brief The type of the 1D spline domain corresponding to the dimension of interest.
70 using spline_domain_type = ddc::DiscreteDomain<bsplines_type>;
71
72 /**
73 * @brief The type of the batch domain (obtained by removing the dimension of interest
74 * from the whole domain).
75 *
76 * @tparam The batched discrete domain on which the interpolation points are defined.
77 */
78 template <concepts::discrete_domain BatchedInterpolationDDom>
79 using batch_domain_type
80 = ddc::remove_dims_of_t<BatchedInterpolationDDom, evaluation_discrete_dimension_type>;
81
82 /**
83 * @brief The type of the whole spline domain (cartesian product of 1D spline domain
84 * and batch domain) preserving the order of dimensions.
85 *
86 * @tparam The batched discrete domain on which the interpolation points are defined.
87 */
88 template <concepts::discrete_domain BatchedInterpolationDDom>
89 using batched_spline_domain_type = ddc::replace_dim_of_t<
90 BatchedInterpolationDDom,
91 evaluation_discrete_dimension_type,
92 bsplines_type>;
93
94 /// @brief The type of the extrapolation rule at the lower boundary.
95 using lower_extrapolation_rule_type = LowerExtrapolationRule;
96
97 /// @brief The type of the extrapolation rule at the upper boundary.
98 using upper_extrapolation_rule_type = UpperExtrapolationRule;
99
100
101private:
102 LowerExtrapolationRule m_lower_extrap_rule;
103
104 UpperExtrapolationRule m_upper_extrap_rule;
105
106public:
107 static_assert(
108 std::is_same_v<LowerExtrapolationRule,
109 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
110 == bsplines_type::is_periodic()
111 && std::is_same_v<
112 UpperExtrapolationRule,
113 typename ddc::PeriodicExtrapolationRule<continuous_dimension_type>>
114 == bsplines_type::is_periodic(),
115 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
116 static_assert(
117 std::is_invocable_r_v<
118 double,
119 LowerExtrapolationRule,
120 ddc::Coordinate<continuous_dimension_type>,
121 ddc::ChunkSpan<
122 double const,
123 spline_domain_type,
124 Kokkos::layout_right,
125 memory_space>>,
126 "LowerExtrapolationRule::operator() has to be callable with usual arguments.");
127 static_assert(
128 std::is_invocable_r_v<
129 double,
130 UpperExtrapolationRule,
131 ddc::Coordinate<continuous_dimension_type>,
132 ddc::ChunkSpan<
133 double const,
134 spline_domain_type,
135 Kokkos::layout_right,
136 memory_space>>,
137 "UpperExtrapolationRule::operator() has to be callable with usual arguments.");
138
139 /**
140 * @brief Build a SplineEvaluator acting on batched_spline_domain.
141 *
142 * @param lower_extrap_rule The extrapolation rule at the lower boundary.
143 * @param upper_extrap_rule The extrapolation rule at the upper boundary.
144 *
145 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
146 */
147 explicit SplineEvaluator(
148 LowerExtrapolationRule const& lower_extrap_rule,
149 UpperExtrapolationRule const& upper_extrap_rule)
150 : m_lower_extrap_rule(lower_extrap_rule)
151 , m_upper_extrap_rule(upper_extrap_rule)
152 {
153 }
154
155 /**
156 * @brief Copy-constructs.
157 *
158 * @param x A reference to another SplineEvaluator.
159 */
160 SplineEvaluator(SplineEvaluator const& x) = default;
161
162 /**
163 * @brief Move-constructs.
164 *
165 * @param x An rvalue to another SplineEvaluator.
166 */
167 SplineEvaluator(SplineEvaluator&& x) = default;
168
169 /// @brief Destructs
170 ~SplineEvaluator() = default;
171
172 /**
173 * @brief Copy-assigns.
174 *
175 * @param x A reference to another SplineEvaluator.
176 * @return A reference to this object.
177 */
178 SplineEvaluator& operator=(SplineEvaluator const& x) = default;
179
180 /**
181 * @brief Move-assigns.
182 *
183 * @param x An rvalue to another SplineEvaluator.
184 * @return A reference to this object.
185 */
187
188 /**
189 * @brief Get the lower extrapolation rule.
190 *
191 * 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.
192 *
193 * @return The lower extrapolation rule.
194 *
195 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
196 */
197 lower_extrapolation_rule_type lower_extrapolation_rule() const
198 {
199 return m_lower_extrap_rule;
200 }
201
202 /**
203 * @brief Get the upper extrapolation rule.
204 *
205 * 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.
206 *
207 * @return The upper extrapolation rule.
208 *
209 * @see NullExtrapolationRule ConstantExtrapolationRule PeriodicExtrapolationRule
210 */
211 upper_extrapolation_rule_type upper_extrapolation_rule() const
212 {
213 return m_upper_extrap_rule;
214 }
215
216 /**
217 * @brief Evaluate 1D spline function (described by its spline coefficients) at a given coordinate.
218 *
219 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
220 *
221 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
222 *
223 * @param coord_eval The coordinate where the spline is evaluated. Note that only the component along the dimension of interest is used.
224 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
225 *
226 * @return The value of the spline function at the desired coordinate.
227 */
228 template <class Layout, class... CoordsDims>
229 KOKKOS_FUNCTION double operator()(
230 ddc::Coordinate<CoordsDims...> const& coord_eval,
231 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
232 spline_coef) const
233 {
234 return eval(coord_eval, spline_coef);
235 }
236
237 /**
238 * @brief Evaluate spline function (described by its spline coefficients) on a mesh.
239 *
240 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
241 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
242 *
243 * This is not a multidimensional evaluation. This is a batched 1D evaluation. This means that for each slice of coordinates
244 * identified by a batch_domain_type::discrete_element_type, the evaluation is performed with the 1D set of
245 * spline coefficients identified by the same batch_domain_type::discrete_element_type.
246 *
247 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline interpolation.
248 *
249 * @param[out] spline_eval The values of the spline function at the desired coordinates. For practical reasons those are
250 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
251 * @param[in] coords_eval The coordinates where the spline is evaluated. Those are
252 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
253 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
254 * the set of 1D spline coefficients retained to perform the evaluation).
255 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
256 */
257 template <
258 class Layout1,
259 class Layout2,
260 class Layout3,
261 class BatchedInterpolationDDom,
262 class... CoordsDims>
263 void operator()(
264 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
265 spline_eval,
266 ddc::ChunkSpan<
267 ddc::Coordinate<CoordsDims...> const,
268 BatchedInterpolationDDom,
269 Layout2,
270 memory_space> const coords_eval,
271 ddc::ChunkSpan<
272 double const,
273 batched_spline_domain_type<BatchedInterpolationDDom>,
274 Layout3,
275 memory_space> const spline_coef) const
276 {
277 evaluation_domain_type const evaluation_domain(spline_eval.domain());
278 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
279
280 ddc::parallel_for_each(
281 "ddc_splines_evaluate",
282 exec_space(),
283 batch_domain,
284 KOKKOS_CLASS_LAMBDA(
285 typename batch_domain_type<
286 BatchedInterpolationDDom>::discrete_element_type const j) {
287 auto const spline_eval_1D = spline_eval[j];
288 auto const coords_eval_1D = coords_eval[j];
289 auto const spline_coef_1D = spline_coef[j];
290 for (auto const i : evaluation_domain) {
291 spline_eval_1D(i) = eval(coords_eval_1D(i), spline_coef_1D);
292 }
293 });
294 }
295
296 /**
297 * @brief Evaluate a spline function (described by its spline coefficients) on a mesh.
298 *
299 * The spline coefficients represent a spline function defined on a cartesian
300 * product of batch_domain and B-splines (basis splines). They can be obtained
301 * via various methods, such as using a SplineBuilder.
302 *
303 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
304 * This means that for each slice of spline_eval the evaluation is performed with
305 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
306 *
307 * Remark: calling SplineBuilder then SplineEvaluator corresponds to a spline
308 * interpolation.
309 *
310 * @param[out] spline_eval The values of the spline function at the coordinates
311 * of the mesh.
312 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
313 */
314 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
315 void operator()(
316 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
317 spline_eval,
318 ddc::ChunkSpan<
319 double const,
320 batched_spline_domain_type<BatchedInterpolationDDom>,
321 Layout2,
322 memory_space> const spline_coef) const
323 {
324 evaluation_domain_type const evaluation_domain(spline_eval.domain());
325 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
326
327 ddc::parallel_for_each(
328 "ddc_splines_evaluate",
329 exec_space(),
330 batch_domain,
331 KOKKOS_CLASS_LAMBDA(
332 typename batch_domain_type<
333 BatchedInterpolationDDom>::discrete_element_type const j) {
334 auto const spline_eval_1D = spline_eval[j];
335 auto const spline_coef_1D = spline_coef[j];
336 for (auto const i : evaluation_domain) {
337 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
338 = ddc::coordinate(i);
339 spline_eval_1D(i) = eval(coord_eval_1D, spline_coef_1D);
340 }
341 });
342 }
343
344#if defined(DDC_BUILD_DEPRECATED_CODE)
345 /**
346 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
347 *
348 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
349 * obtained via various methods, such as using a SplineBuilder.
350 *
351 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
352 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
353 *
354 * @return The derivative of the spline function at the desired coordinate.
355 */
356 template <class Layout, class... CoordsDims>
357 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] KOKKOS_FUNCTION double
358 deriv(ddc::Coordinate<CoordsDims...> const& coord_eval,
359 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const spline_coef)
360 const
361 {
362 return deriv(
363 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
364 coord_eval,
365 spline_coef);
366 }
367
368 /**
369 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
370 *
371 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
372 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
373 *
374 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
375 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
376 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
377 *
378 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
379 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
380 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
381 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
382 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
383 * the set of 1D spline coefficients retained to perform the evaluation).
384 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
385 */
386 template <
387 class Layout1,
388 class Layout2,
389 class Layout3,
390 class BatchedInterpolationDDom,
391 class... CoordsDims>
392 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
393 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
394 spline_eval,
395 ddc::ChunkSpan<
396 ddc::Coordinate<CoordsDims...> const,
397 BatchedInterpolationDDom,
398 Layout2,
399 memory_space> const coords_eval,
400 ddc::ChunkSpan<
401 double const,
402 batched_spline_domain_type<BatchedInterpolationDDom>,
403 Layout3,
404 memory_space> const spline_coef) const
405 {
406 return deriv(
407 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
408 spline_eval,
409 coords_eval,
410 spline_coef);
411 }
412
413 /**
414 * @brief Differentiate spline function (described by its spline coefficients) on a mesh.
415 *
416 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
417 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
418 *
419 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
420 * This is not a multidimensional evaluation. This is a batched 1D evaluation.
421 * This means that for each slice of spline_eval the evaluation is performed with
422 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
423 *
424 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
425 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
426 */
427 template <class Layout1, class Layout2, class BatchedInterpolationDDom>
428 [[deprecated("Use deriv(DiscreteElement<Deriv<Dim>>(1), ...) instead")]] void deriv(
429 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
430 spline_eval,
431 ddc::ChunkSpan<
432 double const,
433 batched_spline_domain_type<BatchedInterpolationDDom>,
434 Layout2,
435 memory_space> const spline_coef) const
436 {
437 return deriv(
438 ddc::DiscreteElement<Deriv<continuous_dimension_type>>(1),
439 spline_eval,
440 spline_coef);
441 }
442#endif
443
444 /**
445 * @brief Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
446 *
447 * The spline coefficients represent a 1D spline function defined on a B-splines (basis splines). They can be
448 * obtained via various methods, such as using a SplineBuilder.
449 *
450 * @param deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
451 * If the dimension is not present, the order of derivation is considered to be 0.
452 * @param coord_eval The coordinate where the spline is differentiated. Note that only the component along the dimension of interest is used.
453 * @param spline_coef A ChunkSpan storing the 1D spline coefficients.
454 *
455 * @return The derivative of the spline function at the desired coordinate.
456 */
457 template <class DElem, class Layout, class... CoordsDims>
458 KOKKOS_FUNCTION double deriv(
459 DElem const& deriv_order,
460 ddc::Coordinate<CoordsDims...> const& coord_eval,
461 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
462 spline_coef) const
463 {
464 static_assert(is_discrete_element_v<DElem>);
465 return eval_no_bc(deriv_order, coord_eval, spline_coef);
466 }
467
468 /**
469 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
470 *
471 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
472 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
473 *
474 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
475 * This means that for each slice of coordinates identified by a batch_domain_type::discrete_element_type,
476 * the derivation is performed with the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
477 *
478 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
479 * If the dimension is not present, the order of derivation is considered to be 0.
480 * @param[out] spline_eval The derivatives of the spline function at the desired coordinates. For practical reasons those are
481 * stored in a ChunkSpan defined on a batched_evaluation_domain_type.
482 * @param[in] coords_eval The coordinates where the spline is differentiated. Those are
483 * stored in a ChunkSpan defined on a batched_evaluation_domain_type. Note that the coordinates of the
484 * points represented by this domain are unused and irrelevant (but the points themselves (DiscreteElement) are used to select
485 * the set of 1D spline coefficients retained to perform the evaluation).
486 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
487 */
488 template <
489 class DElem,
490 class Layout1,
491 class Layout2,
492 class Layout3,
493 class BatchedInterpolationDDom,
494 class... CoordsDims>
495 void deriv(
496 DElem const& deriv_order,
497 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
498 spline_eval,
499 ddc::ChunkSpan<
500 ddc::Coordinate<CoordsDims...> const,
501 BatchedInterpolationDDom,
502 Layout2,
503 memory_space> const coords_eval,
504 ddc::ChunkSpan<
505 double const,
506 batched_spline_domain_type<BatchedInterpolationDDom>,
507 Layout3,
508 memory_space> const spline_coef) const
509 {
510 static_assert(is_discrete_element_v<DElem>);
511
512 evaluation_domain_type const evaluation_domain(spline_eval.domain());
513 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
514
515 ddc::parallel_for_each(
516 "ddc_splines_differentiate",
517 exec_space(),
518 batch_domain,
519 KOKKOS_CLASS_LAMBDA(
520 typename batch_domain_type<
521 BatchedInterpolationDDom>::discrete_element_type const j) {
522 auto const spline_eval_1D = spline_eval[j];
523 auto const coords_eval_1D = coords_eval[j];
524 auto const spline_coef_1D = spline_coef[j];
525 for (auto const i : evaluation_domain) {
526 spline_eval_1D(i)
527 = eval_no_bc(deriv_order, coords_eval_1D(i), spline_coef_1D);
528 }
529 });
530 }
531
532 /**
533 * @brief Differentiate 1D spline function (described by its spline coefficients) on a mesh.
534 *
535 * The spline coefficients represent a spline function defined on a cartesian product of batch_domain and B-splines
536 * (basis splines). They can be obtained via various methods, such as using a SplineBuilder.
537 *
538 * The derivation is not performed in a multidimensional way (in any sense). This is a batched 1D derivation.
539 * This means that for each slice of spline_eval the evaluation is performed with
540 * the 1D set of spline coefficients identified by the same batch_domain_type::discrete_element_type.
541 *
542 * @param[in] deriv_order A DiscreteElement containing the order of derivation for the dimension of interest.
543 * If the dimension is not present, the order of derivation is considered to be 0.
544 * @param[out] spline_eval The derivatives of the spline function at the coordinates.
545 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
546 */
547 template <class DElem, class Layout1, class Layout2, class BatchedInterpolationDDom>
548 void deriv(
549 DElem const& deriv_order,
550 ddc::ChunkSpan<double, BatchedInterpolationDDom, Layout1, memory_space> const
551 spline_eval,
552 ddc::ChunkSpan<
553 double const,
554 batched_spline_domain_type<BatchedInterpolationDDom>,
555 Layout2,
556 memory_space> const spline_coef) const
557 {
558 static_assert(is_discrete_element_v<DElem>);
559
560 evaluation_domain_type const evaluation_domain(spline_eval.domain());
561 batch_domain_type<BatchedInterpolationDDom> const batch_domain(spline_eval.domain());
562
563 ddc::parallel_for_each(
564 "ddc_splines_differentiate",
565 exec_space(),
566 batch_domain,
567 KOKKOS_CLASS_LAMBDA(
568 typename batch_domain_type<
569 BatchedInterpolationDDom>::discrete_element_type const j) {
570 auto const spline_eval_1D = spline_eval[j];
571 auto const spline_coef_1D = spline_coef[j];
572 for (auto const i : evaluation_domain) {
573 ddc::Coordinate<continuous_dimension_type> coord_eval_1D
574 = ddc::coordinate(i);
575 spline_eval_1D(i) = eval_no_bc(deriv_order, coord_eval_1D, spline_coef_1D);
576 }
577 });
578 }
579
580 /** @brief Perform batched 1D integrations of a spline function (described by its spline coefficients) along the dimension of interest and store results on a subdomain of batch_domain.
581 *
582 * The spline coefficients represent a spline function defined on a B-splines (basis splines). They can be obtained via the SplineBuilder.
583 *
584 * The integration is not performed in a multidimensional way (in any sense). This is a batched 1D integration.
585 * This means that for each element of integrals, the integration is performed with the 1D set of
586 * spline coefficients identified by the same DiscreteElement.
587 *
588 * @param[out] integrals The integrals of the spline function on the subdomain of batch_domain. For practical reasons those are
589 * stored in a ChunkSpan defined on a batch_domain_type. Note that the coordinates of the
590 * points represented by this domain are unused and irrelevant.
591 * @param[in] spline_coef A ChunkSpan storing the spline coefficients.
592 */
593 template <class Layout1, class Layout2, class BatchedDDom, class BatchedSplineDDom>
594 void integrate(
595 ddc::ChunkSpan<double, BatchedDDom, Layout1, memory_space> const integrals,
596 ddc::ChunkSpan<double const, BatchedSplineDDom, Layout2, memory_space> const
597 spline_coef) const
598 {
599 static_assert(
600 ddc::in_tags_v<bsplines_type, to_type_seq_t<BatchedSplineDDom>>,
601 "The spline coefficients domain must contain the bsplines dimension");
602 using batch_domain_type = ddc::remove_dims_of_t<BatchedSplineDDom, bsplines_type>;
603 static_assert(
604 std::is_same_v<batch_domain_type, BatchedDDom>,
605 "The integrals domain must only contain the batch dimensions");
606
607 BatchedDDom const batch_domain(integrals.domain());
608 ddc::Chunk values_alloc(
609 ddc::DiscreteDomain<bsplines_type>(spline_coef.domain()),
610 ddc::KokkosAllocator<double, memory_space>());
611 ddc::ChunkSpan const values = values_alloc.span_view();
612 ddc::integrals(exec_space(), values);
613
614 ddc::parallel_for_each(
615 "ddc_splines_integrate",
616 exec_space(),
617 batch_domain,
618 KOKKOS_LAMBDA(typename BatchedDDom::discrete_element_type const j) {
619 integrals(j) = 0;
620 for (typename spline_domain_type::discrete_element_type const i :
621 values.domain()) {
622 integrals(j) += spline_coef(i, j) * values(i);
623 }
624 });
625 }
626
627private:
628 template <class Layout, class... CoordsDims>
629 KOKKOS_INLINE_FUNCTION double eval(
630 ddc::Coordinate<CoordsDims...> const& coord_eval,
631 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
632 spline_coef) const
633 {
634 ddc::Coordinate<continuous_dimension_type> coord_eval_interest(coord_eval);
635 if constexpr (bsplines_type::is_periodic()) {
636 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()
637 || coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
638 coord_eval_interest -= Kokkos::floor(
639 (coord_eval_interest
640 - ddc::discrete_space<bsplines_type>().rmin())
641 / ddc::discrete_space<bsplines_type>().length())
642 * ddc::discrete_space<bsplines_type>().length();
643 }
644 } else {
645 if (coord_eval_interest < ddc::discrete_space<bsplines_type>().rmin()) {
646 return m_lower_extrap_rule(coord_eval_interest, spline_coef);
647 }
648 if (coord_eval_interest > ddc::discrete_space<bsplines_type>().rmax()) {
649 return m_upper_extrap_rule(coord_eval_interest, spline_coef);
650 }
651 }
652 return eval_no_bc(ddc::DiscreteElement<>(), coord_eval_interest, spline_coef);
653 }
654
655 template <class... DerivDims, class Layout, class... CoordsDims>
656 KOKKOS_INLINE_FUNCTION double eval_no_bc(
657 ddc::DiscreteElement<DerivDims...> const& deriv_order,
658 ddc::Coordinate<CoordsDims...> const& coord_eval,
659 ddc::ChunkSpan<double const, spline_domain_type, Layout, memory_space> const
660 spline_coef) const
661 {
662 using deriv_dim = ddc::Deriv<continuous_dimension_type>;
663 static_assert(
664 sizeof...(DerivDims) == 0
665 || type_seq_same_v<
666 detail::TypeSeq<DerivDims...>,
667 detail::TypeSeq<deriv_dim>>,
668 "The only valid dimension for deriv_order is Deriv<Dim>");
669
670 ddc::DiscreteElement<bsplines_type> jmin;
671 std::array<double, bsplines_type::degree() + 1> vals_ptr;
672 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
673 vals(vals_ptr.data());
674 ddc::Coordinate<continuous_dimension_type> const coord_eval_interest(coord_eval);
675
676 if constexpr (sizeof...(DerivDims) == 0) {
677 jmin = ddc::discrete_space<bsplines_type>().eval_basis(vals, coord_eval_interest);
678 } else {
679 auto const order = deriv_order.uid();
680 KOKKOS_ASSERT(order > 0 && order <= bsplines_type::degree())
681
682 std::array<double, (bsplines_type::degree() + 1) * (bsplines_type::degree() + 1)>
683 derivs_ptr;
684 Kokkos::mdspan<
685 double,
686 Kokkos::extents<
687 std::size_t,
688 bsplines_type::degree() + 1,
689 Kokkos::dynamic_extent>> const derivs(derivs_ptr.data(), order + 1);
690
691 jmin = ddc::discrete_space<bsplines_type>()
692 .eval_basis_and_n_derivs(derivs, coord_eval_interest, order);
693
694 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
695 vals[i] = DDC_MDSPAN_ACCESS_OP(derivs, i, order);
696 }
697 }
698
699 double y = 0.0;
700 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
701 y += spline_coef(ddc::DiscreteElement<bsplines_type>(jmin + i)) * vals[i];
702 }
703 return y;
704 }
705};
706
707} // 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 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