DDC 0.12.0
Loading...
Searching...
No Matches
spline_builder_2d.hpp
1// Copyright (C) The DDC development team, see COPYRIGHT.md file
2//
3// SPDX-License-Identifier: MIT
4
5#pragma once
6
7#include <cassert>
8#include <cstddef>
9#include <optional>
10
11#include <ddc/ddc.hpp>
12
15
16namespace ddc {
17
18/**
19 * @brief A class for creating a 2D spline approximation of a function.
20 *
21 * A class which contains an operator () which can be used to build a 2D spline approximation
22 * of a function. A 2D spline approximation uses a cross-product between two 1D SplineBuilder.
23 *
24 * @see SplineBuilder
25 */
26template <
27 class ExecSpace,
28 class MemorySpace,
29 class BSpline1,
30 class BSpline2,
31 class DDimI1,
32 class DDimI2,
33 ddc::BoundCond BcLower1,
34 ddc::BoundCond BcUpper1,
35 ddc::BoundCond BcLower2,
36 ddc::BoundCond BcUpper2,
37 ddc::SplineSolver Solver>
39{
40public:
41 /// @brief The type of the Kokkos execution space used by this class.
42 using exec_space = ExecSpace;
43
44 /// @brief The type of the Kokkos memory space used by this class.
45 using memory_space = MemorySpace;
46
47 /// @brief The type of the SplineBuilder used by this class to spline-approximate along first dimension.
48 using builder_type1 = ddc::
49 SplineBuilder<ExecSpace, MemorySpace, BSpline1, DDimI1, BcLower1, BcUpper1, Solver>;
50
51 /// @brief The type of the SplineBuilder used by this class to spline-approximate along second dimension.
52 using builder_type2 = ddc::
53 SplineBuilder<ExecSpace, MemorySpace, BSpline2, DDimI2, BcLower2, BcUpper2, Solver>;
54
55 /// @brief The type of the SplineBuilder used by this class to spline-approximate the second-dimension-derivatives along first dimension.
56 using builder_deriv_type1 = ddc::
57 SplineBuilder<ExecSpace, MemorySpace, BSpline1, DDimI1, BcLower1, BcUpper1, Solver>;
58
59 /// @brief The type of the first interpolation continuous dimension.
60 using continuous_dimension_type1 = builder_type1::continuous_dimension_type;
61
62 /// @brief The type of the second interpolation continuous dimension.
63 using continuous_dimension_type2 = builder_type2::continuous_dimension_type;
64
65 /// @brief The type of the first interpolation discrete dimension.
66 using interpolation_discrete_dimension_type1
67 = builder_type1::interpolation_discrete_dimension_type;
68
69 /// @brief The type of the second interpolation discrete dimension.
70 using interpolation_discrete_dimension_type2
71 = builder_type2::interpolation_discrete_dimension_type;
72
73 /// @brief The type of the B-splines in the first dimension.
74 using bsplines_type1 = builder_type1::bsplines_type;
75
76 /// @brief The type of the B-splines in the second dimension.
77 using bsplines_type2 = builder_type2::bsplines_type;
78
79 /// @brief The type of the Deriv domain on boundaries in the first dimension.
80 using deriv_type1 = builder_type1::deriv_type;
81
82 /// @brief The type of the Deriv domain on boundaries in the second dimension.
83 using deriv_type2 = builder_type2::deriv_type;
84
85 /// @brief The type of the domain for the interpolation mesh in the first dimension.
86 using interpolation_domain_type1 = builder_type1::interpolation_discrete_dimension_type;
87
88 /// @brief The type of the domain for the interpolation mesh in the second dimension.
89 using interpolation_domain_type2 = builder_type2::interpolation_discrete_dimension_type;
90
91 /// @brief The type of the domain for the interpolation mesh in the 2D dimension.
92 using interpolation_domain_type = ddc::DiscreteDomain<
93 interpolation_discrete_dimension_type1,
94 interpolation_discrete_dimension_type2>;
95
96 /**
97 * @brief The type of the whole domain representing interpolation points.
98 *
99 * @tparam The batched discrete domain on which the interpolation points are defined.
100 */
101 template <concepts::discrete_domain BatchedInterpolationDDom>
102 using batched_interpolation_domain_type = BatchedInterpolationDDom;
103
104 /**
105 * @brief The type of the batch domain (obtained by removing the dimensions of interest
106 * from the whole domain).
107 *
108 * @tparam The batched discrete domain on which the interpolation points are defined.
109 *
110 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y,
111 * this is DiscreteDomain<Z>.
112 */
113 template <concepts::discrete_domain BatchedInterpolationDDom>
114 using batch_domain_type = ddc::remove_dims_of_t<
115 BatchedInterpolationDDom,
116 interpolation_discrete_dimension_type1,
117 interpolation_discrete_dimension_type2>;
118
119 /**
120 * @brief The type of the whole spline domain (cartesian product of 2D spline domain
121 * and batch domain) preserving the order of dimensions.
122 *
123 * @tparam The batched discrete domain on which the interpolation points are defined.
124 *
125 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y
126 * (associated to B-splines tags BSplinesX and BSplinesY), this is DiscreteDomain<BSplinesX, BSplinesY, Z>
127 */
128 template <concepts::discrete_domain BatchedInterpolationDDom>
129 using batched_spline_domain_type
130 = ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
131 ddc::to_type_seq_t<BatchedInterpolationDDom>,
132 ddc::detail::TypeSeq<
133 interpolation_discrete_dimension_type1,
134 interpolation_discrete_dimension_type2>,
135 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2>>>;
136
137 /**
138 * @brief The type of the whole Derivs domain (cartesian product of the 1D Deriv domain
139 * and the associated batch domain) in the first dimension, preserving the order of dimensions.
140 *
141 * @tparam The batched discrete domain on which the interpolation points are defined.
142 *
143 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y,
144 * this is DiscreteDomain<Deriv<X>, Y, Z>.
145 */
146 template <concepts::discrete_domain BatchedInterpolationDDom>
147 using batched_derivs_domain_type1
148 = builder_type1::template batched_derivs_domain_type<BatchedInterpolationDDom>;
149
150 /**
151 * @brief The type of the whole Derivs domain (cartesian product of the 1D Deriv domain
152 * and the associated batch domain) in the second dimension, preserving the order of dimensions.
153 *
154 * @tparam The batched discrete domain on which the interpolation points are defined.
155 *
156 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y,
157 * this is DiscreteDomain<X, Deriv<Y>, Z>.
158 */
159 template <concepts::discrete_domain BatchedInterpolationDDom>
160 using batched_derivs_domain_type2 = ddc::replace_dim_of_t<
161 BatchedInterpolationDDom,
162 interpolation_discrete_dimension_type2,
163 deriv_type2>;
164
165 /**
166 * @brief The type of the whole Derivs domain (cartesian product of the 2D Deriv domain
167 * and the batch domain) in the second dimension, preserving the order of dimensions.
168 *
169 * @tparam The batched discrete domain on which the interpolation points are defined.
170 *
171 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and dimensions of interest X and Y,
172 * this is DiscreteDomain<Deriv<X>, Deriv<Y>, Z>.
173 */
174 template <concepts::discrete_domain BatchedInterpolationDDom>
175 using batched_derivs_domain_type
176 = ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_replace_t<
177 ddc::to_type_seq_t<BatchedInterpolationDDom>,
178 ddc::detail::TypeSeq<
179 interpolation_discrete_dimension_type1,
180 interpolation_discrete_dimension_type2>,
181 ddc::detail::TypeSeq<deriv_type1, deriv_type2>>>;
182
183private:
184 builder_type1 m_spline_builder1;
185 builder_deriv_type1 m_spline_builder_deriv1;
186 builder_type2 m_spline_builder2;
187
188public:
189 /**
190 * @brief Build a SplineBuilder2D acting on interpolation_domain.
191 *
192 * @param interpolation_domain The domain on which the interpolation points are defined, without the batch dimensions.
193 *
194 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
195 * of a chunk of right-hand-sides of the linear problem to be computed in parallel (chunks are treated
196 * by the linear solver one-after-the-other).
197 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
198 *
199 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
200 * define the size of a block used by the Block-Jacobi preconditioner.
201 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
202 *
203 * @see SplinesLinearProblemSparse
204 */
205 explicit SplineBuilder2D(
206 interpolation_domain_type const& interpolation_domain,
207 std::optional<std::size_t> cols_per_chunk = std::nullopt,
208 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
209 : m_spline_builder1(interpolation_domain, cols_per_chunk, preconditioner_max_block_size)
210 , m_spline_builder_deriv1(interpolation_domain)
211 , m_spline_builder2(interpolation_domain, cols_per_chunk, preconditioner_max_block_size)
212 {
213 }
214
215 /**
216 * @brief Build a SplineBuilder2D acting on the interpolation domain contained in batched_interpolation_domain.
217 *
218 * @param batched_interpolation_domain The domain on which the interpolation points are defined.
219 *
220 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
221 * of a chunk of right-hand-sides of the linear problem to be computed in parallel (chunks are treated
222 * by the linear solver one-after-the-other).
223 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
224 *
225 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
226 * define the size of a block used by the Block-Jacobi preconditioner.
227 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
228 *
229 * @see SplinesLinearProblemSparse
230 */
231 template <concepts::discrete_domain BatchedInterpolationDDom>
232 explicit SplineBuilder2D(
233 BatchedInterpolationDDom const& batched_interpolation_domain,
234 std::optional<std::size_t> cols_per_chunk = std::nullopt,
235 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
237 interpolation_domain_type(batched_interpolation_domain),
238 cols_per_chunk,
239 preconditioner_max_block_size)
240 {
241 }
242
243 /// @brief Copy-constructor is deleted.
244 SplineBuilder2D(SplineBuilder2D const& x) = delete;
245
246 /**
247 * @brief Move-constructs.
248 *
249 * @param x An rvalue to another SplineBuilder2D.
250 */
251 SplineBuilder2D(SplineBuilder2D&& x) = default;
252
253 /// @brief Destructs.
254 ~SplineBuilder2D() = default;
255
256 /// @brief Copy-assignment is deleted.
257 SplineBuilder2D& operator=(SplineBuilder2D const& x) = delete;
258
259 /** @brief Move-assigns.
260 *
261 * @param x An rvalue to another SplineBuilder.
262 * @return A reference to this object.
263 */
265
266 /**
267 * @brief Get the domain for the 2D interpolation mesh used by this class.
268 *
269 * This is 2D because it is defined along the dimensions of interest.
270 *
271 * @return The 2D domain for the interpolation mesh.
272 */
273 interpolation_domain_type interpolation_domain() const noexcept
274 {
275 return ddc::DiscreteDomain<interpolation_domain_type1, interpolation_domain_type2>(
276 m_spline_builder1.interpolation_domain(),
277 m_spline_builder2.interpolation_domain());
278 }
279
280 /**
281 * @brief Get the whole domain representing interpolation points.
282 *
283 * Values of the function must be provided on this domain in order
284 * to build a spline representation of the function (cartesian product of 2D interpolation_domain and batch_domain).
285 *
286 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
287 *
288 * @return The domain for the interpolation mesh.
289 */
290 template <concepts::discrete_domain BatchedInterpolationDDom>
291 BatchedInterpolationDDom batched_interpolation_domain(
292 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
293 {
294 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
295 return batched_interpolation_domain;
296 }
297
298 /**
299 * @brief Get the batch domain.
300 *
301 * Obtained by removing the dimensions of interest from the whole interpolation domain.
302 *
303 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
304 *
305 * @return The batch domain.
306 */
307 template <concepts::discrete_domain BatchedInterpolationDDom>
308 batch_domain_type<BatchedInterpolationDDom> batch_domain(
309 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
310 {
311 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
312 return ddc::remove_dims_of(batched_interpolation_domain, interpolation_domain());
313 }
314
315 /**
316 * @brief Get the 2D domain on which spline coefficients are defined.
317 *
318 * The 2D spline domain corresponding to the dimensions of interest.
319 *
320 * @return The 2D domain for the spline coefficients.
321 */
322 ddc::DiscreteDomain<bsplines_type1, bsplines_type2> spline_domain() const noexcept
323 {
324 return ddc::DiscreteDomain<bsplines_type1, bsplines_type2>(
325 ddc::discrete_space<bsplines_type1>().full_domain(),
326 ddc::discrete_space<bsplines_type2>().full_domain());
327 }
328
329 /**
330 * @brief Get the whole domain on which spline coefficients are defined.
331 *
332 * Spline approximations (spline-transformed functions) are computed on this domain.
333 *
334 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
335 *
336 * @return The domain for the spline coefficients.
337 */
338 template <concepts::discrete_domain BatchedInterpolationDDom>
339 batched_spline_domain_type<BatchedInterpolationDDom> batched_spline_domain(
340 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
341 {
342 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
343 return ddc::replace_dim_of<interpolation_discrete_dimension_type1, bsplines_type1>(
344 ddc::replace_dim_of<
345 interpolation_discrete_dimension_type2,
346 bsplines_type2>(batched_interpolation_domain, spline_domain()),
348 }
349
350 /**
351 * @brief Compute a 2D spline approximation of a function.
352 *
353 * Use the values of a function (defined on
354 * SplineBuilder2D::batched_interpolation_domain) and the derivatives of the
355 * function at the boundaries (in the case of BoundCond::HERMITE only)
356 * to calculate a 2D spline approximation of this function.
357 *
358 * The spline approximation is stored as a ChunkSpan of coefficients
359 * associated with B-splines.
360 *
361 * @param[out] spline
362 * The coefficients of the spline computed by this SplineBuilder.
363 * @param[in] vals
364 * The values of the function at the interpolation mesh.
365 * @param[in] derivs_min1
366 * The values of the derivatives at the lower boundary in the first dimension.
367 * @param[in] derivs_max1
368 * The values of the derivatives at the upper boundary in the first dimension.
369 * @param[in] derivs_min2
370 * The values of the derivatives at the lower boundary in the second dimension.
371 * @param[in] derivs_max2
372 * The values of the derivatives at the upper boundary in the second dimension.
373 * @param[in] mixed_derivs_min1_min2
374 * The values of the the cross-derivatives at the lower boundary in the first dimension
375 * and the lower boundary in the second dimension.
376 * @param[in] mixed_derivs_max1_min2
377 * The values of the the cross-derivatives at the upper boundary in the first dimension
378 * and the lower boundary in the second dimension.
379 * @param[in] mixed_derivs_min1_max2
380 * The values of the the cross-derivatives at the lower boundary in the first dimension
381 * and the upper boundary in the second dimension.
382 * @param[in] mixed_derivs_max1_max2
383 * The values of the the cross-derivatives at the upper boundary in the first dimension
384 * and the upper boundary in the second dimension.
385 */
386 template <class Layout, class BatchedInterpolationDDom>
387 void operator()(
388 ddc::ChunkSpan<
389 double,
390 batched_spline_domain_type<BatchedInterpolationDDom>,
391 Layout,
392 memory_space> spline,
393 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
394 std::optional<ddc::ChunkSpan<
395 double const,
396 batched_derivs_domain_type1<BatchedInterpolationDDom>,
397 Layout,
398 memory_space>> derivs_min1
399 = std::nullopt,
400 std::optional<ddc::ChunkSpan<
401 double const,
402 batched_derivs_domain_type1<BatchedInterpolationDDom>,
403 Layout,
404 memory_space>> derivs_max1
405 = std::nullopt,
406 std::optional<ddc::ChunkSpan<
407 double const,
408 batched_derivs_domain_type2<BatchedInterpolationDDom>,
409 Layout,
410 memory_space>> derivs_min2
411 = std::nullopt,
412 std::optional<ddc::ChunkSpan<
413 double const,
414 batched_derivs_domain_type2<BatchedInterpolationDDom>,
415 Layout,
416 memory_space>> derivs_max2
417 = std::nullopt,
418 std::optional<ddc::ChunkSpan<
419 double const,
420 batched_derivs_domain_type<BatchedInterpolationDDom>,
421 Layout,
422 memory_space>> mixed_derivs_min1_min2
423 = std::nullopt,
424 std::optional<ddc::ChunkSpan<
425 double const,
426 batched_derivs_domain_type<BatchedInterpolationDDom>,
427 Layout,
428 memory_space>> mixed_derivs_max1_min2
429 = std::nullopt,
430 std::optional<ddc::ChunkSpan<
431 double const,
432 batched_derivs_domain_type<BatchedInterpolationDDom>,
433 Layout,
434 memory_space>> mixed_derivs_min1_max2
435 = std::nullopt,
436 std::optional<ddc::ChunkSpan<
437 double const,
438 batched_derivs_domain_type<BatchedInterpolationDDom>,
439 Layout,
440 memory_space>> mixed_derivs_max1_max2
441 = std::nullopt) const;
442};
443
444
445template <
446 class ExecSpace,
447 class MemorySpace,
448 class BSpline1,
449 class BSpline2,
450 class DDimI1,
451 class DDimI2,
457template <class Layout, class BatchedInterpolationDDom>
458void SplineBuilder2D<
459 ExecSpace,
460 MemorySpace,
461 BSpline1,
462 BSpline2,
463 DDimI1,
464 DDimI2,
465 BcLower1,
466 BcUpper1,
467 BcLower2,
468 BcUpper2,
469 Solver>::
470operator()(
471 ddc::ChunkSpan<
472 double,
473 batched_spline_domain_type<BatchedInterpolationDDom>,
474 Layout,
475 memory_space> spline,
476 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
477 std::optional<ddc::ChunkSpan<
478 double const,
479 batched_derivs_domain_type1<BatchedInterpolationDDom>,
480 Layout,
481 memory_space>> const derivs_min1,
482 std::optional<ddc::ChunkSpan<
483 double const,
484 batched_derivs_domain_type1<BatchedInterpolationDDom>,
485 Layout,
486 memory_space>> const derivs_max1,
487 std::optional<ddc::ChunkSpan<
488 double const,
489 batched_derivs_domain_type2<BatchedInterpolationDDom>,
490 Layout,
491 memory_space>> const derivs_min2,
492 std::optional<ddc::ChunkSpan<
493 double const,
494 batched_derivs_domain_type2<BatchedInterpolationDDom>,
495 Layout,
496 memory_space>> const derivs_max2,
497 std::optional<ddc::ChunkSpan<
498 double const,
499 batched_derivs_domain_type<BatchedInterpolationDDom>,
500 Layout,
501 memory_space>> const mixed_derivs_min1_min2,
502 std::optional<ddc::ChunkSpan<
503 double const,
504 batched_derivs_domain_type<BatchedInterpolationDDom>,
505 Layout,
506 memory_space>> const mixed_derivs_max1_min2,
507 std::optional<ddc::ChunkSpan<
508 double const,
509 batched_derivs_domain_type<BatchedInterpolationDDom>,
510 Layout,
511 memory_space>> const mixed_derivs_min1_max2,
512 std::optional<ddc::ChunkSpan<
513 double const,
514 batched_derivs_domain_type<BatchedInterpolationDDom>,
515 Layout,
516 memory_space>> const mixed_derivs_max1_max2) const
517{
518 auto const batched_interpolation_domain = vals.domain();
519
520 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
521 assert(batch_domain_type<BatchedInterpolationDDom>(batched_interpolation_domain)
522 == batch_domain_type<BatchedInterpolationDDom>(spline.domain()));
523
524 if (batch_domain(batched_interpolation_domain).empty()) {
525 return;
526 }
527
528 // TODO: perform computations along dimension 1 on different streams ?
529 // Spline1-approximate derivs_min2 (to spline1_deriv_min)
530
531 auto const batched_interpolation_deriv_domain
532 = ddc::replace_dim_of<interpolation_discrete_dimension_type2, deriv_type2>(
533 batched_interpolation_domain,
534 ddc::DiscreteDomain<deriv_type2>(
535 ddc::DiscreteElement<deriv_type2>(builder_type2::s_odd),
536 ddc::DiscreteVector<deriv_type2>(bsplines_type2::degree() / 2)));
537
538 ddc::Chunk spline1_deriv_min_alloc(
539 m_spline_builder_deriv1.batched_spline_domain(batched_interpolation_deriv_domain),
540 ddc::KokkosAllocator<double, MemorySpace>());
541 auto spline1_deriv_min = spline1_deriv_min_alloc.span_view();
542 auto spline1_deriv_min_opt = std::optional(spline1_deriv_min.span_cview());
543 if constexpr (BcLower2 == ddc::BoundCond::HERMITE) {
544 m_spline_builder_deriv1(
545 spline1_deriv_min,
546 *derivs_min2,
547 mixed_derivs_min1_min2,
548 mixed_derivs_max1_min2);
549 } else {
550 spline1_deriv_min_opt = std::nullopt;
551 }
552
553 // Spline1-approximate vals (to spline1)
554 ddc::Chunk spline1_alloc(
555 m_spline_builder1.batched_spline_domain(batched_interpolation_domain),
556 ddc::KokkosAllocator<double, MemorySpace>());
557 ddc::ChunkSpan const spline1 = spline1_alloc.span_view();
558
559 m_spline_builder1(spline1, vals, derivs_min1, derivs_max1);
560
561 // Spline1-approximate derivs_max2 (to spline1_deriv_max)
562 ddc::Chunk spline1_deriv_max_alloc(
563 m_spline_builder_deriv1.batched_spline_domain(batched_interpolation_deriv_domain),
564 ddc::KokkosAllocator<double, MemorySpace>());
565 auto spline1_deriv_max = spline1_deriv_max_alloc.span_view();
566 auto spline1_deriv_max_opt = std::optional(spline1_deriv_max.span_cview());
567 if constexpr (BcUpper2 == ddc::BoundCond::HERMITE) {
568 m_spline_builder_deriv1(
569 spline1_deriv_max,
570 *derivs_max2,
571 mixed_derivs_min1_max2,
572 mixed_derivs_max1_max2);
573 } else {
574 spline1_deriv_max_opt = std::nullopt;
575 }
576
577 // Spline2-approximate spline1
578 m_spline_builder2(spline, spline1.span_cview(), spline1_deriv_min_opt, spline1_deriv_max_opt);
579}
580
581} // namespace ddc
friend class ChunkSpan
friend class Chunk
Definition chunk.hpp:81
friend class DiscreteDomain
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 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.
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.
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