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