DDC 0.14.0
Loading...
Searching...
No Matches
spline_builder.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 <cassert>
9#include <cstddef>
10#include <memory>
11#include <optional>
12#include <stdexcept>
13#include <string>
14#include <tuple>
15#include <utility>
16
17#include <ddc/ddc.hpp>
18
19#include <Kokkos_Core.hpp>
20
21#include "deriv.hpp"
22#include "integrals.hpp"
23#include "math_tools.hpp"
27#include "view.hpp"
28
29namespace ddc {
30
31/**
32 * @brief An enum determining the backend solver of a SplineBuilder or SplineBuilder2d.
33 *
34 * An enum determining the backend solver of a SplineBuilder or SplineBuilder2d.
35 */
36enum class SplineSolver {
37 GINKGO, ///< Enum member to identify the Ginkgo-based solver (iterative method)
38 LAPACK ///< Enum member to identify the LAPACK-based solver (direct method)
39};
40
41/**
42 * @brief A class for creating a spline approximation of a function.
43 *
44 * A class which contains an operator () which can be used to build a spline approximation
45 * of a function. A spline approximation is represented by coefficients stored in a Chunk
46 * of B-splines. The spline is constructed such that it respects the boundary conditions
47 * BcLower and BcUpper, and it interpolates the function at the points on the interpolation_discrete_dimension
48 * associated with interpolation_discrete_dimension_type.
49 * @tparam ExecSpace The Kokkos execution space on which the spline approximation is performed.
50 * @tparam MemorySpace The Kokkos memory space on which the data (interpolation function and splines coefficients) is stored.
51 * @tparam BSplines The discrete dimension representing the B-splines.
52 * @tparam InterpolationDDim The discrete dimension on which interpolation points are defined.
53 * @tparam BcLower The lower boundary condition.
54 * @tparam BcUpper The upper boundary condition.
55 * @tparam Solver The SplineSolver giving the backend used to perform the spline approximation.
56 */
57template <
58 class ExecSpace,
59 class MemorySpace,
60 class BSplines,
61 class InterpolationDDim,
62 ddc::BoundCond BcLower,
63 ddc::BoundCond BcUpper,
64 SplineSolver Solver>
65class SplineBuilder
66{
67 static_assert(
68 (BSplines::is_periodic() && (BcLower == ddc::BoundCond::PERIODIC)
69 && (BcUpper == ddc::BoundCond::PERIODIC))
70 || (!BSplines::is_periodic() && (BcLower != ddc::BoundCond::PERIODIC)
71 && (BcUpper != ddc::BoundCond::PERIODIC)));
72
73public:
74 /// @brief The type of the Kokkos execution space used by this class.
75 using exec_space = ExecSpace;
76
77 /// @brief The type of the Kokkos memory space used by this class.
78 using memory_space = MemorySpace;
79
80 /// @brief The type of the interpolation continuous dimension (continuous dimension of interest) used by this class.
81 using continuous_dimension_type = InterpolationDDim::continuous_dimension_type;
82
83 /// @brief The type of the interpolation discrete dimension (discrete dimension of interest) used by this class.
84 using interpolation_discrete_dimension_type = InterpolationDDim;
85
86 /// @brief The discrete dimension representing the B-splines.
87 using bsplines_type = BSplines;
88
89 /// @brief The type of the Deriv dimension at the boundaries.
90 using deriv_type = ddc::Deriv<continuous_dimension_type>;
91
92 /// @brief The type of the domain for the 1D interpolation mesh used by this class.
93 using interpolation_domain_type = ddc::DiscreteDomain<interpolation_discrete_dimension_type>;
94
95 /**
96 * @brief The type of the whole domain representing interpolation points.
97 *
98 * @tparam The batched discrete domain on which the interpolation points are defined.
99 */
100 template <concepts::discrete_domain BatchedInterpolationDDom>
101 using batched_interpolation_domain_type = BatchedInterpolationDDom;
102
103 /**
104 * @brief The type of the batch domain (obtained by removing the dimension of interest
105 * from the whole domain).
106 *
107 * @tparam The batched discrete domain on which the interpolation points are defined.
108 *
109 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and a dimension of interest Y,
110 * this is DiscreteDomain<X,Z>
111 */
112 template <concepts::discrete_domain BatchedInterpolationDDom>
113 using batch_domain_type = ddc::
114 remove_dims_of_t<BatchedInterpolationDDom, interpolation_discrete_dimension_type>;
115
116 /**
117 * @brief The type of the whole spline domain (cartesian product of 1D spline domain
118 * and batch domain) preserving the underlying memory layout (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 a dimension of interest Y
123 * (associated to a B-splines tag BSplinesY), this is DiscreteDomain<X,BSplinesY,Z>.
124 */
125 template <concepts::discrete_domain BatchedInterpolationDDom>
126 using batched_spline_domain_type = ddc::replace_dim_of_t<
127 BatchedInterpolationDDom,
128 interpolation_discrete_dimension_type,
129 bsplines_type>;
130
131private:
132 /**
133 * @brief The type of the whole spline domain (cartesian product of the 1D spline domain
134 * and the batch domain) with 1D spline dimension being the leading dimension.
135 *
136 * @tparam The batched discrete domain on which the interpolation points are defined.
137 *
138 * Example: For batched_interpolation_domain_type = DiscreteDomain<X,Y,Z> and a dimension of interest Y
139 * (associated to a B-splines tag BSplinesY), this is DiscreteDomain<BSplinesY,X,Z>.
140 */
141 template <concepts::discrete_domain BatchedInterpolationDDom>
142 using batched_spline_tr_domain_type
143 = ddc::detail::convert_type_seq_to_discrete_domain_t<ddc::type_seq_merge_t<
144 ddc::detail::TypeSeq<bsplines_type>,
145 ddc::type_seq_remove_t<
146 ddc::to_type_seq_t<BatchedInterpolationDDom>,
147 ddc::detail::TypeSeq<interpolation_discrete_dimension_type>>>>;
148
149public:
150 /**
151 * @brief The type of the whole Deriv domain (cartesian product of 1D Deriv domain
152 * and batch domain) preserving the underlying memory layout (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 a dimension of interest Y,
157 * this is DiscreteDomain<X,Deriv<Y>,Z>
158 */
159 template <concepts::discrete_domain BatchedInterpolationDDom>
160 using batched_derivs_domain_type = ddc::replace_dim_of_t<
161 BatchedInterpolationDDom,
162 interpolation_discrete_dimension_type,
163 deriv_type>;
164
165 /// @brief Indicates if the degree of the splines is odd or even.
166 static constexpr bool s_odd = BSplines::degree() % 2;
167
168 /// @brief The number of equations defining the boundary condition at the lower bound.
169 static constexpr int s_nbe_xmin = n_boundary_equations(BcLower, BSplines::degree());
170
171 /// @brief The number of equations defining the boundary condition at the upper bound.
172 static constexpr int s_nbe_xmax = n_boundary_equations(BcUpper, BSplines::degree());
173
174 /// @brief The number of input values defining the boundary condition at the lower bound.
175 static constexpr int s_nbv_xmin = BcLower == BoundCond::HOMOGENEOUS_HERMITE
176 ? 0
177 : n_boundary_equations(BcLower, BSplines::degree());
178
179 /// @brief The number of input values defining the boundary condition at the upper bound.
180 static constexpr int s_nbv_xmax = BcUpper == BoundCond::HOMOGENEOUS_HERMITE
181 ? 0
182 : n_boundary_equations(BcUpper, BSplines::degree());
183
184 /// @brief The boundary condition implemented at the lower bound.
185 static constexpr ddc::BoundCond s_bc_xmin = BcLower;
186
187 /// @brief The boundary condition implemented at the upper bound.
188 static constexpr ddc::BoundCond s_bc_xmax = BcUpper;
189
190 /// @brief The SplineSolver giving the backend used to perform the spline approximation.
191 static constexpr SplineSolver s_spline_solver = Solver;
192
193private:
194 interpolation_domain_type m_interpolation_domain;
195
196 int m_offset = 0;
197
198 double m_dx; // average cell size for normalization of derivatives
199
200 // interpolator specific
201 std::unique_ptr<ddc::detail::SplinesLinearProblem<exec_space>> m_matrix;
202
203 std::string m_label;
204
205 /// Calculate offset so that the matrix is diagonally dominant
206 void compute_offset(interpolation_domain_type const& interpolation_domain, int& offset);
207
208public:
209 /**
210 * @brief Build a SplineBuilder acting on interpolation_domain.
211 *
212 * @param label A label used to tag parallel regions and memory allocations for profiling.
213 *
214 * @param interpolation_domain The domain on which the interpolation points are defined.
215 *
216 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
217 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
218 * by the linear solver one-after-the-other).
219 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
220 *
221 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
222 * define the size of a block used by the Block-Jacobi preconditioner.
223 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
224 *
225 * @see MatrixSparse
226 */
227 explicit SplineBuilder(
228 std::string label,
229 interpolation_domain_type const& interpolation_domain,
230 std::optional<std::size_t> cols_per_chunk = std::nullopt,
231 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
232 : m_interpolation_domain(interpolation_domain)
233 , m_dx((ddc::discrete_space<BSplines>().rmax() - ddc::discrete_space<BSplines>().rmin())
234 / ddc::discrete_space<BSplines>().ncells())
235 , m_label(std::move(label))
236 {
237 static_assert(
238 ((BcLower == BoundCond::PERIODIC) == (BcUpper == BoundCond::PERIODIC)),
239 "Incompatible boundary conditions");
240 check_valid_grid();
241
242 compute_offset(this->interpolation_domain(), m_offset);
243
244 // Calculate block sizes
245 int lower_block_size;
246 int upper_block_size;
247 if constexpr (bsplines_type::is_uniform()) {
248 upper_block_size = compute_block_sizes_uniform(BcLower, s_nbe_xmin);
249 lower_block_size = compute_block_sizes_uniform(BcUpper, s_nbe_xmax);
250 } else {
251 upper_block_size = compute_block_sizes_non_uniform(BcLower, s_nbe_xmin);
252 lower_block_size = compute_block_sizes_non_uniform(BcUpper, s_nbe_xmax);
253 }
254 allocate_matrix(
255 lower_block_size,
256 upper_block_size,
257 cols_per_chunk,
258 preconditioner_max_block_size);
259 }
260
261 /**
262 * @brief Build a SplineBuilder acting on interpolation_domain.
263 *
264 * @param interpolation_domain The domain on which the interpolation points are defined.
265 *
266 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
267 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
268 * by the linear solver one-after-the-other).
269 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
270 *
271 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
272 * define the size of a block used by the Block-Jacobi preconditioner.
273 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
274 *
275 * @see MatrixSparse
276 */
277 explicit SplineBuilder(
278 interpolation_domain_type const& interpolation_domain,
279 std::optional<std::size_t> cols_per_chunk = std::nullopt,
280 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
282 "no-label",
283 interpolation_domain,
284 cols_per_chunk,
285 preconditioner_max_block_size)
286 {
287 }
288
289 /**
290 * @brief Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
291 *
292 * @param label A label used to tag parallel regions and memory allocations for profiling.
293 *
294 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
295 *
296 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
297 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
298 * by the linear solver one-after-the-other).
299 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
300 *
301 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
302 * define the size of a block used by the Block-Jacobi preconditioner.
303 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
304 *
305 * @see MatrixSparse
306 */
307 template <concepts::discrete_domain BatchedInterpolationDDom>
308 explicit SplineBuilder(
309 std::string label,
310 BatchedInterpolationDDom const& batched_interpolation_domain,
311 std::optional<std::size_t> cols_per_chunk = std::nullopt,
312 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
314 std::move(label),
315 interpolation_domain_type(batched_interpolation_domain),
316 cols_per_chunk,
317 preconditioner_max_block_size)
318 {
319 }
320
321 /**
322 * @brief Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
323 *
324 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
325 *
326 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
327 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
328 * by the linear solver one-after-the-other).
329 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
330 *
331 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
332 * define the size of a block used by the Block-Jacobi preconditioner.
333 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
334 *
335 * @see MatrixSparse
336 */
337 template <concepts::discrete_domain BatchedInterpolationDDom>
338 explicit SplineBuilder(
339 BatchedInterpolationDDom const& batched_interpolation_domain,
340 std::optional<std::size_t> cols_per_chunk = std::nullopt,
341 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
343 "no-label",
344 interpolation_domain_type(batched_interpolation_domain),
345 cols_per_chunk,
346 preconditioner_max_block_size)
347 {
348 }
349
350
351 /// @brief Copy-constructor is deleted.
352 SplineBuilder(SplineBuilder const& x) = delete;
353
354 /** @brief Move-constructs.
355 *
356 * @param x An rvalue to another SplineBuilder.
357 */
358 SplineBuilder(SplineBuilder&& x) = default;
359
360 /// @brief Destructs.
361 ~SplineBuilder() = default;
362
363 /// @brief Copy-assignment is deleted.
364 SplineBuilder& operator=(SplineBuilder const& x) = delete;
365
366 /** @brief Move-assigns.
367 *
368 * @param x An rvalue to another SplineBuilder.
369 * @return A reference to this object.
370 */
371 SplineBuilder& operator=(SplineBuilder&& x) = default;
372
373 /**
374 * @brief Get the domain for the 1D interpolation mesh used by this class.
375 *
376 * This is 1D because it is defined along the dimension of interest.
377 *
378 * @return The 1D domain for the interpolation mesh.
379 */
380 interpolation_domain_type interpolation_domain() const noexcept
381 {
382 return m_interpolation_domain;
383 }
384
385 /**
386 * @brief Get the whole domain representing interpolation points.
387 *
388 * Values of the function must be provided on this domain in order
389 * to build a spline representation of the function (cartesian product of 1D interpolation_domain and batch_domain).
390 *
391 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
392 *
393 * @return The domain for the interpolation mesh.
394 */
395 template <concepts::discrete_domain BatchedInterpolationDDom>
396 BatchedInterpolationDDom batched_interpolation_domain(
397 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
398 {
399 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
400 return batched_interpolation_domain;
401 }
402
403 /**
404 * @brief Get the batch domain.
405 *
406 * Obtained by removing the dimension of interest from the whole interpolation domain.
407 *
408 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
409 *
410 * @return The batch domain.
411 */
412 template <class BatchedInterpolationDDom>
413 batch_domain_type<BatchedInterpolationDDom> batch_domain(
414 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
415 {
416 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
417 return ddc::remove_dims_of(batched_interpolation_domain, interpolation_domain());
418 }
419
420 /**
421 * @brief Get the 1D domain on which spline coefficients are defined.
422 *
423 * The 1D spline domain corresponding to the dimension of interest.
424 *
425 * @return The 1D domain for the spline coefficients.
426 */
427 ddc::DiscreteDomain<bsplines_type> spline_domain() const noexcept
428 {
429 return ddc::discrete_space<bsplines_type>().full_domain();
430 }
431
432 /**
433 * @brief Get the whole domain on which spline coefficients are defined.
434 *
435 * Spline approximations (spline-transformed functions) are computed on this domain.
436 *
437 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
438 *
439 * @return The domain for the spline coefficients.
440 */
441 template <class BatchedInterpolationDDom>
442 batched_spline_domain_type<BatchedInterpolationDDom> batched_spline_domain(
443 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
444 {
445 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
446 return ddc::replace_dim_of<
447 interpolation_discrete_dimension_type,
448 bsplines_type>(batched_interpolation_domain, spline_domain());
449 }
450
451private:
452 /**
453 * @brief Get the whole domain on which spline coefficients are defined, with the dimension of interest being the leading dimension.
454 *
455 * This is used internally due to solver limitation and because it may be beneficial to computation performance. For LAPACK backend and non-periodic boundary condition, we are using SplinesLinearSolver3x3Blocks which requires upper_block_size additional rows for internal operations.
456 *
457 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
458 *
459 * @return The (transposed) domain for the spline coefficients.
460 */
461 template <class BatchedInterpolationDDom>
462 batched_spline_tr_domain_type<BatchedInterpolationDDom> batched_spline_tr_domain(
463 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
464 {
465 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
466 return batched_spline_tr_domain_type<BatchedInterpolationDDom>(
467 ddc::replace_dim_of<bsplines_type, bsplines_type>(
468 batched_spline_domain(batched_interpolation_domain),
469 ddc::DiscreteDomain<bsplines_type>(
470 ddc::DiscreteElement<bsplines_type>(0),
471 ddc::DiscreteVector<bsplines_type>(
472 m_matrix->required_number_of_rhs_rows()))));
473 }
474
475public:
476 /**
477 * @brief Get the whole domain on which derivatives on lower boundary are defined.
478 *
479 * This is only used with BoundCond::HERMITE boundary conditions.
480 *
481 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
482 *
483 * @return The domain for the Derivs values.
484 */
485 template <class BatchedInterpolationDDom>
486 batched_derivs_domain_type<BatchedInterpolationDDom> batched_derivs_xmin_domain(
487 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
488 {
489 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
490 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
491 batched_interpolation_domain,
492 ddc::DiscreteDomain<deriv_type>(
493 ddc::DiscreteElement<deriv_type>(1),
494 ddc::DiscreteVector<deriv_type>(s_nbv_xmin)));
495 }
496
497 /**
498 * @brief Get the whole domain on which derivatives on upper boundary are defined.
499 *
500 * This is only used with BoundCond::HERMITE boundary conditions.
501 *
502 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
503 *
504 * @return The domain for the Derivs values.
505 */
506 template <class BatchedInterpolationDDom>
507 batched_derivs_domain_type<BatchedInterpolationDDom> batched_derivs_xmax_domain(
508 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
509 {
510 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
511 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
512 batched_interpolation_domain,
513 ddc::DiscreteDomain<deriv_type>(
514 ddc::DiscreteElement<deriv_type>(1),
515 ddc::DiscreteVector<deriv_type>(s_nbv_xmax)));
516 }
517
518 /**
519 * @brief Compute a spline approximation of a function.
520 *
521 * Use the values of a function (defined on
522 * SplineBuilder::batched_interpolation_domain) and the derivatives of the
523 * function at the boundaries (in the case of BoundCond::HERMITE only, defined
524 * on SplineBuilder::batched_derivs_xmin_domain and SplineBuilder::batched_derivs_xmax_domain)
525 * to calculate a spline approximation of this function.
526 *
527 * The spline approximation is stored as a ChunkSpan of coefficients
528 * associated with B-splines.
529 *
530 * @param[out] spline The coefficients of the spline computed by this SplineBuilder.
531 * @param[in] vals The values of the function on the interpolation mesh.
532 * @param[in] derivs_xmin The values of the derivatives at the lower boundary
533 * (used only with BoundCond::HERMITE lower boundary condition).
534 * @param[in] derivs_xmax The values of the derivatives at the upper boundary
535 * (used only with BoundCond::HERMITE upper boundary condition).
536 */
537 template <class Layout, class BatchedInterpolationDDom>
538 void operator()(
539 ddc::ChunkSpan<
540 double,
541 batched_spline_domain_type<BatchedInterpolationDDom>,
542 Layout,
543 memory_space> spline,
544 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
545 std::optional<ddc::ChunkSpan<
546 double const,
547 batched_derivs_domain_type<BatchedInterpolationDDom>,
548 Layout,
549 memory_space>> derivs_xmin
550 = std::nullopt,
551 std::optional<ddc::ChunkSpan<
552 double const,
553 batched_derivs_domain_type<BatchedInterpolationDDom>,
554 Layout,
555 memory_space>> derivs_xmax
556 = std::nullopt) const;
557
558 /**
559 * @brief Compute the quadrature coefficients associated to the b-splines used by this SplineBuilder.
560 *
561 * Those coefficients can be used to perform integration way faster than SplineEvaluator::integrate().
562 *
563 * This function solves matrix equation A^t*Q=integral_bsplines. In case of HERMITE boundary conditions,
564 * integral_bsplines contains the integral coefficients at the boundaries, and Q thus has to
565 * be split in three parts (quadrature coefficients for the derivatives at lower boundary,
566 * for the values inside the domain and for the derivatives at upper boundary).
567 *
568 * A discrete function f can then be integrated using sum_j Q_j*f_j for j in interpolation_domain.
569 * If boundary condition is HERMITE, sum_j Qderiv_j*(d^j f/dx^j) for j in derivs_domain
570 * must be added at the boundary.
571 *
572 * Please refer to section 2.8.1 of Emily's Bourne phd (https://theses.fr/2022AIXM0412) for more information and to
573 * the (Non)PeriodicSplineBuilderTest for example usage to compute integrals.
574 *
575 * @tparam OutMemorySpace The Kokkos::MemorySpace on which the quadrature coefficients are be returned
576 * (but they are computed on ExecSpace then copied).
577 *
578 * @return A tuple containing the three Chunks containing the quadrature coefficients (if HERMITE
579 * is not used, first and third are empty).
580 */
581 template <class OutMemorySpace = MemorySpace>
582 std::tuple<
583 ddc::Chunk<
584 double,
586 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
587 ddc::KokkosAllocator<double, OutMemorySpace>>,
588 ddc::Chunk<
589 double,
590 ddc::DiscreteDomain<InterpolationDDim>,
591 ddc::KokkosAllocator<double, OutMemorySpace>>,
592 ddc::Chunk<
593 double,
595 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
596 ddc::KokkosAllocator<double, OutMemorySpace>>>
598
599private:
600 static int compute_block_sizes_uniform(ddc::BoundCond bound_cond, int nbc);
601
602 static int compute_block_sizes_non_uniform(ddc::BoundCond bound_cond, int nbc);
603
604 void allocate_matrix(
605 int lower_block_size,
606 int upper_block_size,
607 std::optional<std::size_t> cols_per_chunk = std::nullopt,
608 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt);
609
610 void build_matrix_system();
611
612 void check_valid_grid();
613
614 template <class KnotElement>
615 static void check_n_points_in_cell(int n_points_in_cell, KnotElement current_cell_end_idx);
616};
617
618template <
619 class ExecSpace,
620 class MemorySpace,
621 class BSplines,
622 class InterpolationDDim,
623 ddc::BoundCond BcLower,
624 ddc::BoundCond BcUpper,
625 SplineSolver Solver>
626void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
627 compute_offset(interpolation_domain_type const& interpolation_domain, int& offset)
628{
629 if constexpr (bsplines_type::is_periodic()) {
630 // Calculate offset so that the matrix is diagonally dominant
631 std::array<double, bsplines_type::degree() + 1> values_ptr;
632 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
633 values(values_ptr.data());
634 ddc::DiscreteElement<interpolation_discrete_dimension_type> start(
635 interpolation_domain.front());
636 auto jmin = ddc::discrete_space<BSplines>()
637 .eval_basis(values, ddc::coordinate(start + BSplines::degree()));
638 if constexpr (bsplines_type::degree() % 2 == 0) {
639 offset = jmin.uid() - start.uid() + bsplines_type::degree() / 2 - BSplines::degree();
640 } else {
641 int const mid = bsplines_type::degree() / 2;
642 offset = jmin.uid() - start.uid()
643 + (DDC_MDSPAN_ACCESS_OP(values, mid) > DDC_MDSPAN_ACCESS_OP(values, mid + 1)
644 ? mid
645 : mid + 1)
646 - BSplines::degree();
647 }
648 } else {
649 offset = 0;
650 }
651}
652
653template <
654 class ExecSpace,
655 class MemorySpace,
656 class BSplines,
657 class InterpolationDDim,
658 ddc::BoundCond BcLower,
659 ddc::BoundCond BcUpper,
660 SplineSolver Solver>
661int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
662 compute_block_sizes_uniform(ddc::BoundCond const bound_cond, int const nbc)
663{
664 if (bound_cond == ddc::BoundCond::PERIODIC) {
665 return static_cast<int>(bsplines_type::degree()) / 2;
666 }
667
668 if (bound_cond == ddc::BoundCond::HERMITE
669 || bound_cond == ddc::BoundCond::HOMOGENEOUS_HERMITE) {
670 return nbc;
671 }
672
673 if (bound_cond == ddc::BoundCond::GREVILLE) {
674 return static_cast<int>(bsplines_type::degree()) - 1;
675 }
676
677 throw std::runtime_error("ddc::BoundCond not handled");
678}
679
680template <
681 class ExecSpace,
682 class MemorySpace,
683 class BSplines,
684 class InterpolationDDim,
685 ddc::BoundCond BcLower,
686 ddc::BoundCond BcUpper,
687 SplineSolver Solver>
688int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
689 compute_block_sizes_non_uniform(ddc::BoundCond const bound_cond, int const nbc)
690{
691 if (bound_cond == ddc::BoundCond::PERIODIC || bound_cond == ddc::BoundCond::GREVILLE) {
692 return static_cast<int>(bsplines_type::degree()) - 1;
693 }
694
695 if (bound_cond == ddc::BoundCond::HERMITE
696 || bound_cond == ddc::BoundCond::HOMOGENEOUS_HERMITE) {
697 return nbc + 1;
698 }
699
700 throw std::runtime_error("ddc::BoundCond not handled");
701}
702
703template <
704 class ExecSpace,
705 class MemorySpace,
706 class BSplines,
707 class InterpolationDDim,
708 ddc::BoundCond BcLower,
709 ddc::BoundCond BcUpper,
710 SplineSolver Solver>
711void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
712 allocate_matrix(
713 [[maybe_unused]] int lower_block_size,
714 [[maybe_unused]] int upper_block_size,
715 std::optional<std::size_t> cols_per_chunk,
716 std::optional<unsigned int> preconditioner_max_block_size)
717{
718 // Special case: linear spline
719 // No need for matrix assembly
720 // (disabled)
721 // if constexpr (bsplines_type::degree() == 1) {
722 // return;
723 // }
724
725 if constexpr (Solver == ddc::SplineSolver::LAPACK) {
726 int upper_band_width;
727 if (bsplines_type::is_uniform()) {
728 upper_band_width = bsplines_type::degree() / 2;
729 } else {
730 upper_band_width = bsplines_type::degree() - 1;
731 }
732 if constexpr (bsplines_type::is_periodic()) {
733 m_matrix = ddc::detail::SplinesLinearProblemMaker::make_new_periodic_band_matrix<
734 ExecSpace>(
735 ddc::discrete_space<BSplines>().nbasis(),
736 upper_band_width,
737 upper_band_width,
738 bsplines_type::is_uniform());
739 } else {
740 m_matrix = ddc::detail::SplinesLinearProblemMaker::
741 make_new_block_matrix_with_band_main_block<ExecSpace>(
742 ddc::discrete_space<BSplines>().nbasis(),
743 upper_band_width,
744 upper_band_width,
745 bsplines_type::is_uniform(),
746 lower_block_size,
747 upper_block_size);
748 }
749 } else if constexpr (Solver == ddc::SplineSolver::GINKGO) {
750 m_matrix = ddc::detail::SplinesLinearProblemMaker::make_new_sparse<ExecSpace>(
751 ddc::discrete_space<BSplines>().nbasis(),
752 cols_per_chunk,
753 preconditioner_max_block_size);
754 }
755
756 build_matrix_system();
757
758 m_matrix->setup_solver();
759}
760
761template <
762 class ExecSpace,
763 class MemorySpace,
764 class BSplines,
765 class InterpolationDDim,
766 ddc::BoundCond BcLower,
767 ddc::BoundCond BcUpper,
768 SplineSolver Solver>
769void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
770 build_matrix_system()
771{
772 // Hermite boundary conditions at xmin, if any
773 if constexpr (
775 std::array<double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
776 derivs_ptr;
777 ddc::DSpan2D const
778 derivs(derivs_ptr.data(),
779 bsplines_type::degree() + 1,
780 bsplines_type::degree() / 2 + 1);
781 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
782 derivs,
783 ddc::discrete_space<BSplines>().rmin(),
784 s_nbe_xmin);
785
786 // In order to improve the condition number of the matrix, we normalize
787 // all derivatives by multiplying the i-th derivative by dx^i
788 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
789 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
790 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *= ddc::detail::ipow(m_dx, j);
791 }
792 }
793
794 if constexpr (s_nbe_xmin > 0) {
795 // iterate only to deg as last bspline is 0
796 for (std::size_t i = 0; i < s_nbe_xmin; ++i) {
797 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
798 m_matrix->set_element(i, j, DDC_MDSPAN_ACCESS_OP(derivs, j, i + s_odd));
799 }
800 }
801 }
802 }
803
804 // Interpolation points
805 std::array<double, bsplines_type::degree() + 1> values_ptr;
806 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const values(
807 values_ptr.data());
808
809 int start = interpolation_domain().front().uid();
810 ddc::host_for_each(interpolation_domain(), [&](auto ix) {
811 auto jmin = ddc::discrete_space<BSplines>().eval_basis(
812 values,
813 ddc::coordinate(ddc::DiscreteElement<interpolation_discrete_dimension_type>(ix)));
814 for (std::size_t s = 0; s < bsplines_type::degree() + 1; ++s) {
815 int const j = ddc::detail::
816 modulo(int(jmin.uid() - m_offset + s),
817 static_cast<int>(ddc::discrete_space<BSplines>().nbasis()));
818 m_matrix->set_element(
819 ix.uid() - start + s_nbe_xmin,
820 j,
821 DDC_MDSPAN_ACCESS_OP(values, s));
822 }
823 });
824
825 // Hermite boundary conditions at xmax, if any
826 if constexpr (
828 std::array<double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
829 derivs_ptr;
830 Kokkos::mdspan<
831 double,
832 Kokkos::extents<
833 std::size_t,
834 bsplines_type::degree() + 1,
835 bsplines_type::degree() / 2 + 1>> const derivs(derivs_ptr.data());
836
837 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
838 derivs,
839 ddc::discrete_space<BSplines>().rmax(),
840 s_nbe_xmax);
841
842 // In order to improve the condition number of the matrix, we normalize
843 // all derivatives by multiplying the i-th derivative by dx^i
844 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
845 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
846 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *= ddc::detail::ipow(m_dx, j);
847 }
848 }
849
850 if constexpr (s_nbe_xmax > 0) {
851 int const i0 = ddc::discrete_space<BSplines>().nbasis() - s_nbe_xmax;
852 int const j0 = ddc::discrete_space<BSplines>().nbasis() - bsplines_type::degree();
853 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
854 for (std::size_t i = 0; i < s_nbe_xmax; ++i) {
855 m_matrix->set_element(
856 i0 + i,
857 j0 + j,
858 DDC_MDSPAN_ACCESS_OP(derivs, j + 1, i + s_odd));
859 }
860 }
861 }
862 }
863}
864
865template <
866 class ExecSpace,
867 class MemorySpace,
868 class BSplines,
869 class InterpolationDDim,
873template <class Layout, class BatchedInterpolationDDom>
874void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
875operator()(
876 ddc::ChunkSpan<
877 double,
878 batched_spline_domain_type<BatchedInterpolationDDom>,
879 Layout,
880 memory_space> spline,
881 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
882 std::optional<ddc::ChunkSpan<
883 double const,
884 batched_derivs_domain_type<BatchedInterpolationDDom>,
885 Layout,
886 memory_space>> const derivs_xmin,
887 std::optional<ddc::ChunkSpan<
888 double const,
889 batched_derivs_domain_type<BatchedInterpolationDDom>,
890 Layout,
891 memory_space>> const derivs_xmax) const
892{
893 auto const batched_interpolation_domain = vals.domain();
894
895 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
896 assert(batch_domain_type<BatchedInterpolationDDom>(batched_interpolation_domain)
897 == batch_domain_type<BatchedInterpolationDDom>(spline.domain()));
898
899 if (batch_domain(batched_interpolation_domain).empty()) {
900 return;
901 }
902
903 assert(vals.template extent<interpolation_discrete_dimension_type>()
904 == ddc::discrete_space<bsplines_type>().nbasis() - s_nbe_xmin - s_nbe_xmax);
905
906 if constexpr (BcLower == BoundCond::HERMITE) {
907 assert(ddc::DiscreteElement<deriv_type>(derivs_xmin->domain().front()).uid() == s_odd);
908 assert(derivs_xmin.has_value() || s_nbe_xmin == 0);
909 } else {
910 assert(!derivs_xmin.has_value() || derivs_xmin->template extent<deriv_type>() == 0);
911 }
912 if constexpr (BcUpper == BoundCond::HERMITE) {
913 assert(ddc::DiscreteElement<deriv_type>(derivs_xmax->domain().front()).uid() == s_odd);
914 assert(derivs_xmax.has_value() || s_nbe_xmax == 0);
915 } else {
916 assert(!derivs_xmax.has_value() || derivs_xmax->template extent<deriv_type>() == 0);
917 }
918
919 // Hermite boundary conditions at xmin, if any
920 // NOTE: For consistency with the linear system, the i-th derivative
921 // provided by the user must be multiplied by dx^i
922 if constexpr (BcLower == BoundCond::HERMITE) {
923 assert(derivs_xmin->template extent<deriv_type>() == s_nbe_xmin);
924 auto derivs_xmin_values = *derivs_xmin;
925 auto const dx_proxy = m_dx;
926 auto const odd_proxy = s_odd;
927 ddc::parallel_for_each(
928 "ddc_splines_hermite_compute_lower_coefficients",
929 exec_space(),
930 batch_domain(batched_interpolation_domain),
931 KOKKOS_LAMBDA(
932 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type j) {
933 for (int i = 0; i < s_nbe_xmin; ++i) {
934 spline(ddc::DiscreteElement<bsplines_type>(i), j)
935 = derivs_xmin_values(
936 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
937 j)
938 * ddc::detail::ipow(dx_proxy, i + odd_proxy);
939 }
940 });
941 } else if constexpr (BcLower == BoundCond::HOMOGENEOUS_HERMITE) {
942 ddc::DiscreteDomain<bsplines_type> const dx_splines(
943 ddc::DiscreteElement<bsplines_type>(0),
944 ddc::DiscreteVector<bsplines_type>(s_nbe_xmin));
945 batched_spline_domain_type<BatchedInterpolationDDom> const
946 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
947 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
948 }
949
950 // Fill spline with vals (to work in spline afterward and preserve vals)
951 ddc::parallel_fill(
952 exec_space(),
953 spline[ddc::DiscreteDomain<bsplines_type>(
954 ddc::DiscreteElement<bsplines_type>(s_nbe_xmin),
955 ddc::DiscreteVector<bsplines_type>(m_offset))],
956 0.);
957 // NOTE: We rely on Kokkos::deep_copy because ddc::parallel_deepcopy do not support
958 // different domain-typed Chunks.
959 Kokkos::deep_copy(
960 exec_space(),
961 spline[ddc::DiscreteDomain<bsplines_type>(
962 ddc::DiscreteElement<bsplines_type>(s_nbe_xmin + m_offset),
963 ddc::DiscreteVector<bsplines_type>(static_cast<std::size_t>(
964 vals.domain()
965 .template extent<
966 interpolation_discrete_dimension_type>())))]
967 .allocation_kokkos_view(),
968 vals.allocation_kokkos_view());
969
970
971
972 // Hermite boundary conditions at xmax, if any
973 // NOTE: For consistency with the linear system, the i-th derivative
974 // provided by the user must be multiplied by dx^i
975 auto const& nbasis_proxy = ddc::discrete_space<bsplines_type>().nbasis();
976 if constexpr (BcUpper == BoundCond::HERMITE) {
977 assert(derivs_xmax->template extent<deriv_type>() == s_nbe_xmax);
978 auto derivs_xmax_values = *derivs_xmax;
979 auto const dx_proxy = m_dx;
980 auto const odd_proxy = s_odd;
981 ddc::parallel_for_each(
982 "ddc_splines_hermite_compute_upper_coefficients",
983 exec_space(),
984 batch_domain(batched_interpolation_domain),
985 KOKKOS_LAMBDA(
986 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type j) {
987 for (int i = 0; i < s_nbe_xmax; ++i) {
988 spline(ddc::DiscreteElement<bsplines_type>(nbasis_proxy - s_nbe_xmax + i),
989 j)
990 = derivs_xmax_values(
991 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
992 j)
993 * ddc::detail::ipow(dx_proxy, i + odd_proxy);
994 }
995 });
996 } else if constexpr (BcUpper == BoundCond::HOMOGENEOUS_HERMITE) {
997 ddc::DiscreteDomain<bsplines_type> const dx_splines(
998 ddc::DiscreteElement<bsplines_type>(nbasis_proxy - s_nbe_xmax),
999 ddc::DiscreteVector<bsplines_type>(s_nbe_xmax));
1000 batched_spline_domain_type<BatchedInterpolationDDom> const
1001 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
1002 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
1003 }
1004
1005 // Allocate and fill a transposed version of spline in order to get dimension of interest as last dimension (optimal for GPU, necessary for Ginkgo). Also select only relevant rows in case of periodic boundaries
1006 auto const& offset_proxy = m_offset;
1007 ddc::Chunk spline_tr_alloc(
1008 m_label + " > spline_tr (ddc::SplineBuilder::operator())",
1009 batched_spline_tr_domain(batched_interpolation_domain),
1010 ddc::KokkosAllocator<double, memory_space>());
1011 ddc::ChunkSpan const spline_tr = spline_tr_alloc.span_view();
1012 ddc::parallel_for_each(
1013 m_label + " > ddc_splines_transpose_rhs",
1014 exec_space(),
1015 batch_domain(batched_interpolation_domain),
1016 KOKKOS_LAMBDA(
1017 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type const j) {
1018 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
1019 spline_tr(ddc::DiscreteElement<bsplines_type>(i), j)
1020 = spline(ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j);
1021 }
1022 });
1023 // Create a 2D Kokkos::View to manage spline_tr as a matrix
1024 Kokkos::View<double**, Kokkos::LayoutRight, exec_space> const bcoef_section(
1025 spline_tr.data_handle(),
1026 static_cast<std::size_t>(spline_tr.template extent<bsplines_type>()),
1027 batch_domain(batched_interpolation_domain).size());
1028 // Compute spline coef
1029 m_matrix->solve(bcoef_section, false);
1030 // Transpose back spline_tr into spline.
1031 ddc::parallel_for_each(
1032 m_label + " > ddc_splines_transpose_back_rhs",
1033 exec_space(),
1034 batch_domain(batched_interpolation_domain),
1035 KOKKOS_LAMBDA(
1036 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type const j) {
1037 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
1038 spline(ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j)
1039 = spline_tr(ddc::DiscreteElement<bsplines_type>(i), j);
1040 }
1041 });
1042
1043 // Duplicate the lower spline coefficients to the upper side in case of periodic boundaries
1044 if (bsplines_type::is_periodic()) {
1045 ddc::parallel_for_each(
1046 m_label + " > ddc_splines_periodic_rows_duplicate_rhs",
1047 exec_space(),
1048 batch_domain(batched_interpolation_domain),
1049 KOKKOS_LAMBDA(
1050 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type const
1051 j) {
1052 if (offset_proxy != 0) {
1053 for (int i = 0; i < offset_proxy; ++i) {
1054 spline(ddc::DiscreteElement<bsplines_type>(i), j) = spline(
1055 ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i),
1056 j);
1057 }
1058 for (std::size_t i = offset_proxy; i < bsplines_type::degree(); ++i) {
1059 spline(ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i), j)
1060 = spline(ddc::DiscreteElement<bsplines_type>(i), j);
1061 }
1062 }
1063 for (std::size_t i(0); i < bsplines_type::degree(); ++i) {
1064 ddc::DiscreteElement<bsplines_type> const i_start(i);
1065 ddc::DiscreteElement<bsplines_type> const i_end(nbasis_proxy + i);
1066
1067 spline(i_end, j) = spline(i_start, j);
1068 }
1069 });
1070 }
1071}
1072
1073template <
1074 class ExecSpace,
1075 class MemorySpace,
1076 class BSplines,
1077 class InterpolationDDim,
1081template <class OutMemorySpace>
1082std::tuple<
1083 ddc::Chunk<
1084 double,
1086 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
1087 ddc::KokkosAllocator<double, OutMemorySpace>>,
1088 ddc::Chunk<
1089 double,
1090 ddc::DiscreteDomain<InterpolationDDim>,
1091 ddc::KokkosAllocator<double, OutMemorySpace>>,
1092 ddc::Chunk<
1093 double,
1095 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
1096 ddc::KokkosAllocator<double, OutMemorySpace>>>
1097SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1099{
1100 // Compute integrals of bsplines
1101 ddc::Chunk integral_bsplines(spline_domain(), ddc::KokkosAllocator<double, MemorySpace>());
1102 ddc::integrals(ExecSpace(), integral_bsplines.span_view());
1103
1104 // Remove additional B-splines in the periodic case (cf. UniformBSplines::full_domain() documentation)
1105 ddc::ChunkSpan const integral_bsplines_without_periodic_additional_bsplines
1106 = integral_bsplines[spline_domain().take_first(
1107 ddc::DiscreteVector<bsplines_type>(m_matrix->size()))];
1108
1109 // Allocate mirror with additional rows (cf. SplinesLinearProblem3x3Blocks documentation)
1110 Kokkos::View<double**, Kokkos::LayoutRight, MemorySpace> const
1111 integral_bsplines_mirror_with_additional_allocation(
1112 m_label + " > integral_bsplines_mirror_with_additional_allocation",
1113 m_matrix->required_number_of_rhs_rows(),
1114 1);
1115
1116 // Extract relevant subview
1117 Kokkos::View<double*, Kokkos::LayoutRight, MemorySpace> const integral_bsplines_mirror
1118 = Kokkos::
1119 subview(integral_bsplines_mirror_with_additional_allocation,
1120 std::
1121 pair {static_cast<std::size_t>(0),
1122 integral_bsplines_without_periodic_additional_bsplines
1123 .size()},
1124 0);
1125
1126 // Solve matrix equation A^t*X=integral_bsplines
1127 Kokkos::deep_copy(
1128 integral_bsplines_mirror,
1129 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view());
1130 m_matrix->solve(integral_bsplines_mirror_with_additional_allocation, true);
1131 Kokkos::deep_copy(
1132 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view(),
1133 integral_bsplines_mirror);
1134
1135 // Slice into three ChunkSpan corresponding to lower derivatives, function values and upper derivatives
1136 ddc::ChunkSpan const coefficients_derivs_xmin
1137 = integral_bsplines_without_periodic_additional_bsplines[spline_domain().take_first(
1138 ddc::DiscreteVector<bsplines_type>(s_nbv_xmin))];
1139 ddc::ChunkSpan const coefficients = integral_bsplines_without_periodic_additional_bsplines
1141 .remove_first(ddc::DiscreteVector<bsplines_type>(s_nbv_xmin))
1142 .take_first(
1143 ddc::DiscreteVector<bsplines_type>(
1144 ddc::discrete_space<bsplines_type>().nbasis() - s_nbv_xmin
1145 - s_nbv_xmax))];
1146 ddc::ChunkSpan const coefficients_derivs_xmax
1147 = integral_bsplines_without_periodic_additional_bsplines
1149 .remove_first(
1150 ddc::DiscreteVector<bsplines_type>(
1151 s_nbv_xmin + coefficients.size()))
1152 .take_first(ddc::DiscreteVector<bsplines_type>(s_nbv_xmax))];
1153
1154 // Multiply derivatives coefficients by dx^n
1155 auto const dx_proxy = m_dx;
1156 auto const odd_proxy = s_odd;
1157 ddc::parallel_for_each(
1158 exec_space(),
1159 coefficients_derivs_xmin.domain(),
1160 KOKKOS_LAMBDA(ddc::DiscreteElement<bsplines_type> i) {
1161 coefficients_derivs_xmin(i) *= ddc::detail::
1162 ipow(dx_proxy,
1163 static_cast<std::size_t>(get<bsplines_type>(
1164 (i - coefficients_derivs_xmin.domain().front()) + odd_proxy)));
1165 });
1166 ddc::parallel_for_each(
1167 exec_space(),
1168 coefficients_derivs_xmax.domain(),
1169 KOKKOS_LAMBDA(ddc::DiscreteElement<bsplines_type> i) {
1170 coefficients_derivs_xmax(i) *= ddc::detail::
1171 ipow(dx_proxy,
1172 static_cast<std::size_t>(get<bsplines_type>(
1173 (i - coefficients_derivs_xmax.domain().front()) + odd_proxy)));
1174 });
1175
1176 ddc::DiscreteElement<deriv_type> const first_deriv(s_odd);
1177 // Allocate Chunk on deriv_type and interpolation_discrete_dimension_type and copy quadrature coefficients into it
1178 ddc::Chunk coefficients_derivs_xmin_out(
1180 deriv_type>(first_deriv, ddc::DiscreteVector<deriv_type>(s_nbv_xmin)),
1181 ddc::KokkosAllocator<double, OutMemorySpace>());
1182 ddc::Chunk coefficients_out(
1183 interpolation_domain().take_first(
1184 ddc::DiscreteVector<interpolation_discrete_dimension_type>(
1185 coefficients.size())),
1186 ddc::KokkosAllocator<double, OutMemorySpace>());
1187 ddc::Chunk coefficients_derivs_xmax_out(
1189 deriv_type>(first_deriv, ddc::DiscreteVector<deriv_type>(s_nbv_xmax)),
1190 ddc::KokkosAllocator<double, OutMemorySpace>());
1191 Kokkos::deep_copy(
1192 coefficients_derivs_xmin_out.allocation_kokkos_view(),
1193 coefficients_derivs_xmin.allocation_kokkos_view());
1194 Kokkos::deep_copy(
1195 coefficients_out.allocation_kokkos_view(),
1196 coefficients.allocation_kokkos_view());
1197 Kokkos::deep_copy(
1198 coefficients_derivs_xmax_out.allocation_kokkos_view(),
1199 coefficients_derivs_xmax.allocation_kokkos_view());
1200 return std::make_tuple(
1201 std::move(coefficients_derivs_xmin_out),
1202 std::move(coefficients_out),
1203 std::move(coefficients_derivs_xmax_out));
1204}
1205
1206template <
1207 class ExecSpace,
1208 class MemorySpace,
1209 class BSplines,
1210 class InterpolationDDim,
1214template <class KnotElement>
1215void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1216 check_n_points_in_cell(int const n_points_in_cell, KnotElement const current_cell_end_idx)
1217{
1218 if (n_points_in_cell > BSplines::degree() + 1) {
1219 KnotElement const rmin_idx = ddc::discrete_space<BSplines>().break_point_domain().front();
1220 int const failed_cell = (current_cell_end_idx - rmin_idx).value();
1221 throw std::runtime_error(
1222 "The spline problem is overconstrained. There are "
1223 + std::to_string(n_points_in_cell) + " points in the " + std::to_string(failed_cell)
1224 + "-th cell.");
1225 }
1226}
1227
1228template <
1229 class ExecSpace,
1230 class MemorySpace,
1231 class BSplines,
1232 class InterpolationDDim,
1233 ddc::BoundCond BcLower,
1234 ddc::BoundCond BcUpper,
1235 SplineSolver Solver>
1236void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1237 check_valid_grid()
1238{
1239 std::size_t const n_interp_points = interpolation_domain().size();
1240 std::size_t const expected_npoints
1241 = ddc::discrete_space<BSplines>().nbasis() - s_nbe_xmin - s_nbe_xmax;
1242 if (n_interp_points != expected_npoints) {
1243 throw std::runtime_error(
1244 "Incorrect number of points supplied to NonUniformInterpolationPoints. "
1245 "(Received : "
1246 + std::to_string(n_interp_points)
1247 + ", expected : " + std::to_string(expected_npoints));
1248 }
1249 int n_points_in_cell = 0;
1250 auto current_cell_end_idx = ddc::discrete_space<BSplines>().break_point_domain().front() + 1;
1251 ddc::host_for_each(interpolation_domain(), [&](auto idx) {
1252 ddc::Coordinate<continuous_dimension_type> const point = ddc::coordinate(idx);
1253 if (point > ddc::coordinate(current_cell_end_idx)) {
1254 // Check the points found in the previous cell
1255 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1256 // Initialise the number of points in the subsequent cell, including the new point
1257 n_points_in_cell = 1;
1258 // Move to the next cell
1259 current_cell_end_idx += 1;
1260 } else if (point == ddc::coordinate(current_cell_end_idx)) {
1261 // Check the points found in the previous cell including the point on the boundary
1262 check_n_points_in_cell(n_points_in_cell + 1, current_cell_end_idx);
1263 // Initialise the number of points in the subsequent cell, including the point on the boundary
1264 n_points_in_cell = 1;
1265 // Move to the next cell
1266 current_cell_end_idx += 1;
1267 } else {
1268 // Indicate that the point is in the cell
1269 n_points_in_cell += 1;
1270 }
1271 });
1272 // Check the number of points in the final cell
1273 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1274}
1275
1276} // 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 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