DDC 0.12.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 = typename 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 typename 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 /// Calculate offset so that the matrix is diagonally dominant
204 void compute_offset(interpolation_domain_type const& interpolation_domain, int& offset);
205
206public:
207 /**
208 * @brief Build a SplineBuilder acting on interpolation_domain.
209 *
210 * @param interpolation_domain The domain on which the interpolation points are defined.
211 *
212 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
213 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
214 * by the linear solver one-after-the-other).
215 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
216 *
217 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
218 * define the size of a block used by the Block-Jacobi preconditioner.
219 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
220 *
221 * @see MatrixSparse
222 */
223 explicit SplineBuilder(
224 interpolation_domain_type const& interpolation_domain,
225 std::optional<std::size_t> cols_per_chunk = std::nullopt,
226 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
227 : m_interpolation_domain(interpolation_domain)
228 , m_dx((ddc::discrete_space<BSplines>().rmax() - ddc::discrete_space<BSplines>().rmin())
229 / ddc::discrete_space<BSplines>().ncells())
230 {
231 static_assert(
232 ((BcLower == BoundCond::PERIODIC) == (BcUpper == BoundCond::PERIODIC)),
233 "Incompatible boundary conditions");
234 check_valid_grid();
235
236 compute_offset(this->interpolation_domain(), m_offset);
237
238 // Calculate block sizes
239 int lower_block_size;
240 int upper_block_size;
241 if constexpr (bsplines_type::is_uniform()) {
242 upper_block_size = compute_block_sizes_uniform(BcLower, s_nbe_xmin);
243 lower_block_size = compute_block_sizes_uniform(BcUpper, s_nbe_xmax);
244 } else {
245 upper_block_size = compute_block_sizes_non_uniform(BcLower, s_nbe_xmin);
246 lower_block_size = compute_block_sizes_non_uniform(BcUpper, s_nbe_xmax);
247 }
248 allocate_matrix(
249 lower_block_size,
250 upper_block_size,
251 cols_per_chunk,
252 preconditioner_max_block_size);
253 }
254
255 /**
256 * @brief Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
257 *
258 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
259 *
260 * @param cols_per_chunk A parameter used by the slicer (internal to the solver) to define the size
261 * of a chunk of right-hand sides of the linear problem to be computed in parallel (chunks are treated
262 * by the linear solver one-after-the-other).
263 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
264 *
265 * @param preconditioner_max_block_size A parameter used by the slicer (internal to the solver) to
266 * define the size of a block used by the Block-Jacobi preconditioner.
267 * This value is optional. If no value is provided then the default value is chosen by the requested solver.
268 *
269 * @see MatrixSparse
270 */
271 template <concepts::discrete_domain BatchedInterpolationDDom>
272 explicit SplineBuilder(
273 BatchedInterpolationDDom const& batched_interpolation_domain,
274 std::optional<std::size_t> cols_per_chunk = std::nullopt,
275 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt)
277 interpolation_domain_type(batched_interpolation_domain),
278 cols_per_chunk,
279 preconditioner_max_block_size)
280 {
281 }
282
283
284 /// @brief Copy-constructor is deleted.
285 SplineBuilder(SplineBuilder const& x) = delete;
286
287 /** @brief Move-constructs.
288 *
289 * @param x An rvalue to another SplineBuilder.
290 */
291 SplineBuilder(SplineBuilder&& x) = default;
292
293 /// @brief Destructs.
294 ~SplineBuilder() = default;
295
296 /// @brief Copy-assignment is deleted.
297 SplineBuilder& operator=(SplineBuilder const& x) = delete;
298
299 /** @brief Move-assigns.
300 *
301 * @param x An rvalue to another SplineBuilder.
302 * @return A reference to this object.
303 */
304 SplineBuilder& operator=(SplineBuilder&& x) = default;
305
306 /**
307 * @brief Get the domain for the 1D interpolation mesh used by this class.
308 *
309 * This is 1D because it is defined along the dimension of interest.
310 *
311 * @return The 1D domain for the interpolation mesh.
312 */
313 interpolation_domain_type interpolation_domain() const noexcept
314 {
315 return m_interpolation_domain;
316 }
317
318 /**
319 * @brief Get the whole domain representing interpolation points.
320 *
321 * Values of the function must be provided on this domain in order
322 * to build a spline representation of the function (cartesian product of 1D interpolation_domain and batch_domain).
323 *
324 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
325 *
326 * @return The domain for the interpolation mesh.
327 */
328 template <concepts::discrete_domain BatchedInterpolationDDom>
329 BatchedInterpolationDDom batched_interpolation_domain(
330 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
331 {
332 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
333 return batched_interpolation_domain;
334 }
335
336 /**
337 * @brief Get the batch domain.
338 *
339 * Obtained by removing the dimension of interest from the whole interpolation domain.
340 *
341 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
342 *
343 * @return The batch domain.
344 */
345 template <class BatchedInterpolationDDom>
346 batch_domain_type<BatchedInterpolationDDom> batch_domain(
347 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
348 {
349 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
350 return ddc::remove_dims_of(batched_interpolation_domain, interpolation_domain());
351 }
352
353 /**
354 * @brief Get the 1D domain on which spline coefficients are defined.
355 *
356 * The 1D spline domain corresponding to the dimension of interest.
357 *
358 * @return The 1D domain for the spline coefficients.
359 */
360 ddc::DiscreteDomain<bsplines_type> spline_domain() const noexcept
361 {
362 return ddc::discrete_space<bsplines_type>().full_domain();
363 }
364
365 /**
366 * @brief Get the whole domain on which spline coefficients are defined.
367 *
368 * Spline approximations (spline-transformed functions) are computed on this domain.
369 *
370 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
371 *
372 * @return The domain for the spline coefficients.
373 */
374 template <class BatchedInterpolationDDom>
375 batched_spline_domain_type<BatchedInterpolationDDom> batched_spline_domain(
376 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
377 {
378 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
379 return ddc::replace_dim_of<
380 interpolation_discrete_dimension_type,
381 bsplines_type>(batched_interpolation_domain, spline_domain());
382 }
383
384private:
385 /**
386 * @brief Get the whole domain on which spline coefficients are defined, with the dimension of interest being the leading dimension.
387 *
388 * 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.
389 *
390 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
391 *
392 * @return The (transposed) domain for the spline coefficients.
393 */
394 template <class BatchedInterpolationDDom>
395 batched_spline_tr_domain_type<BatchedInterpolationDDom> batched_spline_tr_domain(
396 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
397 {
398 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
399 return batched_spline_tr_domain_type<BatchedInterpolationDDom>(
400 ddc::replace_dim_of<bsplines_type, bsplines_type>(
401 batched_spline_domain(batched_interpolation_domain),
402 ddc::DiscreteDomain<bsplines_type>(
403 ddc::DiscreteElement<bsplines_type>(0),
404 ddc::DiscreteVector<bsplines_type>(
405 m_matrix->required_number_of_rhs_rows()))));
406 }
407
408public:
409 /**
410 * @brief Get the whole domain on which derivatives on lower boundary are defined.
411 *
412 * This is only used with BoundCond::HERMITE boundary conditions.
413 *
414 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
415 *
416 * @return The domain for the Derivs values.
417 */
418 template <class BatchedInterpolationDDom>
419 batched_derivs_domain_type<BatchedInterpolationDDom> batched_derivs_xmin_domain(
420 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
421 {
422 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
423 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
424 batched_interpolation_domain,
425 ddc::DiscreteDomain<deriv_type>(
426 ddc::DiscreteElement<deriv_type>(1),
427 ddc::DiscreteVector<deriv_type>(s_nbv_xmin)));
428 }
429
430 /**
431 * @brief Get the whole domain on which derivatives on upper boundary are defined.
432 *
433 * This is only used with BoundCond::HERMITE boundary conditions.
434 *
435 * @param batched_interpolation_domain The whole domain on which the interpolation points are defined.
436 *
437 * @return The domain for the Derivs values.
438 */
439 template <class BatchedInterpolationDDom>
440 batched_derivs_domain_type<BatchedInterpolationDDom> batched_derivs_xmax_domain(
441 BatchedInterpolationDDom const& batched_interpolation_domain) const noexcept
442 {
443 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
444 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
445 batched_interpolation_domain,
446 ddc::DiscreteDomain<deriv_type>(
447 ddc::DiscreteElement<deriv_type>(1),
448 ddc::DiscreteVector<deriv_type>(s_nbv_xmax)));
449 }
450
451 /**
452 * @brief Compute a spline approximation of a function.
453 *
454 * Use the values of a function (defined on
455 * SplineBuilder::batched_interpolation_domain) and the derivatives of the
456 * function at the boundaries (in the case of BoundCond::HERMITE only, defined
457 * on SplineBuilder::batched_derivs_xmin_domain and SplineBuilder::batched_derivs_xmax_domain)
458 * to calculate a spline approximation of this function.
459 *
460 * The spline approximation is stored as a ChunkSpan of coefficients
461 * associated with B-splines.
462 *
463 * @param[out] spline The coefficients of the spline computed by this SplineBuilder.
464 * @param[in] vals The values of the function on the interpolation mesh.
465 * @param[in] derivs_xmin The values of the derivatives at the lower boundary
466 * (used only with BoundCond::HERMITE lower boundary condition).
467 * @param[in] derivs_xmax The values of the derivatives at the upper boundary
468 * (used only with BoundCond::HERMITE upper boundary condition).
469 */
470 template <class Layout, class BatchedInterpolationDDom>
471 void operator()(
472 ddc::ChunkSpan<
473 double,
474 batched_spline_domain_type<BatchedInterpolationDDom>,
475 Layout,
476 memory_space> spline,
477 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
478 std::optional<ddc::ChunkSpan<
479 double const,
480 batched_derivs_domain_type<BatchedInterpolationDDom>,
481 Layout,
482 memory_space>> derivs_xmin
483 = std::nullopt,
484 std::optional<ddc::ChunkSpan<
485 double const,
486 batched_derivs_domain_type<BatchedInterpolationDDom>,
487 Layout,
488 memory_space>> derivs_xmax
489 = std::nullopt) const;
490
491 /**
492 * @brief Compute the quadrature coefficients associated to the b-splines used by this SplineBuilder.
493 *
494 * Those coefficients can be used to perform integration way faster than SplineEvaluator::integrate().
495 *
496 * This function solves matrix equation A^t*Q=integral_bsplines. In case of HERMITE boundary conditions,
497 * integral_bsplines contains the integral coefficients at the boundaries, and Q thus has to
498 * be split in three parts (quadrature coefficients for the derivatives at lower boundary,
499 * for the values inside the domain and for the derivatives at upper boundary).
500 *
501 * A discrete function f can then be integrated using sum_j Q_j*f_j for j in interpolation_domain.
502 * If boundary condition is HERMITE, sum_j Qderiv_j*(d^j f/dx^j) for j in derivs_domain
503 * must be added at the boundary.
504 *
505 * Please refer to section 2.8.1 of Emily's Bourne phd (https://theses.fr/2022AIXM0412) for more information and to
506 * the (Non)PeriodicSplineBuilderTest for example usage to compute integrals.
507 *
508 * @tparam OutMemorySpace The Kokkos::MemorySpace on which the quadrature coefficients are be returned
509 * (but they are computed on ExecSpace then copied).
510 *
511 * @return A tuple containing the three Chunks containing the quadrature coefficients (if HERMITE
512 * is not used, first and third are empty).
513 */
514 template <class OutMemorySpace = MemorySpace>
515 std::tuple<
516 ddc::Chunk<
517 double,
519 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
520 ddc::KokkosAllocator<double, OutMemorySpace>>,
521 ddc::Chunk<
522 double,
523 ddc::DiscreteDomain<InterpolationDDim>,
524 ddc::KokkosAllocator<double, OutMemorySpace>>,
525 ddc::Chunk<
526 double,
528 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
529 ddc::KokkosAllocator<double, OutMemorySpace>>>
531
532private:
533 static int compute_block_sizes_uniform(ddc::BoundCond bound_cond, int nbc);
534
535 static int compute_block_sizes_non_uniform(ddc::BoundCond bound_cond, int nbc);
536
537 void allocate_matrix(
538 int lower_block_size,
539 int upper_block_size,
540 std::optional<std::size_t> cols_per_chunk = std::nullopt,
541 std::optional<unsigned int> preconditioner_max_block_size = std::nullopt);
542
543 void build_matrix_system();
544
545 void check_valid_grid();
546
547 template <class KnotElement>
548 static void check_n_points_in_cell(int n_points_in_cell, KnotElement current_cell_end_idx);
549};
550
551template <
552 class ExecSpace,
553 class MemorySpace,
554 class BSplines,
555 class InterpolationDDim,
556 ddc::BoundCond BcLower,
557 ddc::BoundCond BcUpper,
558 SplineSolver Solver>
559void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
560 compute_offset(interpolation_domain_type const& interpolation_domain, int& offset)
561{
562 if constexpr (bsplines_type::is_periodic()) {
563 // Calculate offset so that the matrix is diagonally dominant
564 std::array<double, bsplines_type::degree() + 1> values_ptr;
565 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const
566 values(values_ptr.data());
567 ddc::DiscreteElement<interpolation_discrete_dimension_type> start(
568 interpolation_domain.front());
569 auto jmin = ddc::discrete_space<BSplines>()
570 .eval_basis(values, ddc::coordinate(start + BSplines::degree()));
571 if constexpr (bsplines_type::degree() % 2 == 0) {
572 offset = jmin.uid() - start.uid() + bsplines_type::degree() / 2 - BSplines::degree();
573 } else {
574 int const mid = bsplines_type::degree() / 2;
575 offset = jmin.uid() - start.uid()
576 + (DDC_MDSPAN_ACCESS_OP(values, mid) > DDC_MDSPAN_ACCESS_OP(values, mid + 1)
577 ? mid
578 : mid + 1)
579 - BSplines::degree();
580 }
581 } else {
582 offset = 0;
583 }
584}
585
586template <
587 class ExecSpace,
588 class MemorySpace,
589 class BSplines,
590 class InterpolationDDim,
591 ddc::BoundCond BcLower,
592 ddc::BoundCond BcUpper,
593 SplineSolver Solver>
594int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
595 compute_block_sizes_uniform(ddc::BoundCond const bound_cond, int const nbc)
596{
597 if (bound_cond == ddc::BoundCond::PERIODIC) {
598 return static_cast<int>(bsplines_type::degree()) / 2;
599 }
600
601 if (bound_cond == ddc::BoundCond::HERMITE
602 || bound_cond == ddc::BoundCond::HOMOGENEOUS_HERMITE) {
603 return nbc;
604 }
605
606 if (bound_cond == ddc::BoundCond::GREVILLE) {
607 return static_cast<int>(bsplines_type::degree()) - 1;
608 }
609
610 throw std::runtime_error("ddc::BoundCond not handled");
611}
612
613template <
614 class ExecSpace,
615 class MemorySpace,
616 class BSplines,
617 class InterpolationDDim,
618 ddc::BoundCond BcLower,
619 ddc::BoundCond BcUpper,
620 SplineSolver Solver>
621int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
622 compute_block_sizes_non_uniform(ddc::BoundCond const bound_cond, int const nbc)
623{
624 if (bound_cond == ddc::BoundCond::PERIODIC || bound_cond == ddc::BoundCond::GREVILLE) {
625 return static_cast<int>(bsplines_type::degree()) - 1;
626 }
627
628 if (bound_cond == ddc::BoundCond::HERMITE
629 || bound_cond == ddc::BoundCond::HOMOGENEOUS_HERMITE) {
630 return nbc + 1;
631 }
632
633 throw std::runtime_error("ddc::BoundCond not handled");
634}
635
636template <
637 class ExecSpace,
638 class MemorySpace,
639 class BSplines,
640 class InterpolationDDim,
641 ddc::BoundCond BcLower,
642 ddc::BoundCond BcUpper,
643 SplineSolver Solver>
644void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
645 allocate_matrix(
646 [[maybe_unused]] int lower_block_size,
647 [[maybe_unused]] int upper_block_size,
648 std::optional<std::size_t> cols_per_chunk,
649 std::optional<unsigned int> preconditioner_max_block_size)
650{
651 // Special case: linear spline
652 // No need for matrix assembly
653 // (disabled)
654 // if constexpr (bsplines_type::degree() == 1) {
655 // return;
656 // }
657
658 if constexpr (Solver == ddc::SplineSolver::LAPACK) {
659 int upper_band_width;
660 if (bsplines_type::is_uniform()) {
661 upper_band_width = bsplines_type::degree() / 2;
662 } else {
663 upper_band_width = bsplines_type::degree() - 1;
664 }
665 if constexpr (bsplines_type::is_periodic()) {
666 m_matrix = ddc::detail::SplinesLinearProblemMaker::make_new_periodic_band_matrix<
667 ExecSpace>(
668 ddc::discrete_space<BSplines>().nbasis(),
669 upper_band_width,
670 upper_band_width,
671 bsplines_type::is_uniform());
672 } else {
673 m_matrix = ddc::detail::SplinesLinearProblemMaker::
674 make_new_block_matrix_with_band_main_block<ExecSpace>(
675 ddc::discrete_space<BSplines>().nbasis(),
676 upper_band_width,
677 upper_band_width,
678 bsplines_type::is_uniform(),
679 lower_block_size,
680 upper_block_size);
681 }
682 } else if constexpr (Solver == ddc::SplineSolver::GINKGO) {
683 m_matrix = ddc::detail::SplinesLinearProblemMaker::make_new_sparse<ExecSpace>(
684 ddc::discrete_space<BSplines>().nbasis(),
685 cols_per_chunk,
686 preconditioner_max_block_size);
687 }
688
689 build_matrix_system();
690
691 m_matrix->setup_solver();
692}
693
694template <
695 class ExecSpace,
696 class MemorySpace,
697 class BSplines,
698 class InterpolationDDim,
699 ddc::BoundCond BcLower,
700 ddc::BoundCond BcUpper,
701 SplineSolver Solver>
702void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
703 build_matrix_system()
704{
705 // Hermite boundary conditions at xmin, if any
706 if constexpr (
708 std::array<double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
709 derivs_ptr;
710 ddc::DSpan2D const
711 derivs(derivs_ptr.data(),
712 bsplines_type::degree() + 1,
713 bsplines_type::degree() / 2 + 1);
714 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
715 derivs,
716 ddc::discrete_space<BSplines>().rmin(),
717 s_nbe_xmin);
718
719 // In order to improve the condition number of the matrix, we normalize
720 // all derivatives by multiplying the i-th derivative by dx^i
721 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
722 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
723 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *= ddc::detail::ipow(m_dx, j);
724 }
725 }
726
727 if constexpr (s_nbe_xmin > 0) {
728 // iterate only to deg as last bspline is 0
729 for (std::size_t i = 0; i < s_nbe_xmin; ++i) {
730 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
731 m_matrix->set_element(i, j, DDC_MDSPAN_ACCESS_OP(derivs, j, i + s_odd));
732 }
733 }
734 }
735 }
736
737 // Interpolation points
738 std::array<double, bsplines_type::degree() + 1> values_ptr;
739 Kokkos::mdspan<double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>> const values(
740 values_ptr.data());
741
742 int start = interpolation_domain().front().uid();
743 ddc::host_for_each(interpolation_domain(), [&](auto ix) {
744 auto jmin = ddc::discrete_space<BSplines>().eval_basis(
745 values,
746 ddc::coordinate(ddc::DiscreteElement<interpolation_discrete_dimension_type>(ix)));
747 for (std::size_t s = 0; s < bsplines_type::degree() + 1; ++s) {
748 int const j = ddc::detail::
749 modulo(int(jmin.uid() - m_offset + s),
750 static_cast<int>(ddc::discrete_space<BSplines>().nbasis()));
751 m_matrix->set_element(
752 ix.uid() - start + s_nbe_xmin,
753 j,
754 DDC_MDSPAN_ACCESS_OP(values, s));
755 }
756 });
757
758 // Hermite boundary conditions at xmax, if any
759 if constexpr (
761 std::array<double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
762 derivs_ptr;
763 Kokkos::mdspan<
764 double,
765 Kokkos::extents<
766 std::size_t,
767 bsplines_type::degree() + 1,
768 bsplines_type::degree() / 2 + 1>> const derivs(derivs_ptr.data());
769
770 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
771 derivs,
772 ddc::discrete_space<BSplines>().rmax(),
773 s_nbe_xmax);
774
775 // In order to improve the condition number of the matrix, we normalize
776 // all derivatives by multiplying the i-th derivative by dx^i
777 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
778 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
779 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *= ddc::detail::ipow(m_dx, j);
780 }
781 }
782
783 if constexpr (s_nbe_xmax > 0) {
784 int const i0 = ddc::discrete_space<BSplines>().nbasis() - s_nbe_xmax;
785 int const j0 = ddc::discrete_space<BSplines>().nbasis() - bsplines_type::degree();
786 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
787 for (std::size_t i = 0; i < s_nbe_xmax; ++i) {
788 m_matrix->set_element(
789 i0 + i,
790 j0 + j,
791 DDC_MDSPAN_ACCESS_OP(derivs, j + 1, i + s_odd));
792 }
793 }
794 }
795 }
796}
797
798template <
799 class ExecSpace,
800 class MemorySpace,
801 class BSplines,
802 class InterpolationDDim,
806template <class Layout, class BatchedInterpolationDDom>
807void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
808operator()(
809 ddc::ChunkSpan<
810 double,
811 batched_spline_domain_type<BatchedInterpolationDDom>,
812 Layout,
813 memory_space> spline,
814 ddc::ChunkSpan<double const, BatchedInterpolationDDom, Layout, memory_space> vals,
815 std::optional<ddc::ChunkSpan<
816 double const,
817 batched_derivs_domain_type<BatchedInterpolationDDom>,
818 Layout,
819 memory_space>> const derivs_xmin,
820 std::optional<ddc::ChunkSpan<
821 double const,
822 batched_derivs_domain_type<BatchedInterpolationDDom>,
823 Layout,
824 memory_space>> const derivs_xmax) const
825{
826 auto const batched_interpolation_domain = vals.domain();
827
828 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
829 assert(batch_domain_type<BatchedInterpolationDDom>(batched_interpolation_domain)
830 == batch_domain_type<BatchedInterpolationDDom>(spline.domain()));
831
832 if (batch_domain(batched_interpolation_domain).empty()) {
833 return;
834 }
835
836 assert(vals.template extent<interpolation_discrete_dimension_type>()
837 == ddc::discrete_space<bsplines_type>().nbasis() - s_nbe_xmin - s_nbe_xmax);
838
839 if constexpr (BcLower == BoundCond::HERMITE) {
840 assert(ddc::DiscreteElement<deriv_type>(derivs_xmin->domain().front()).uid() == s_odd);
841 assert(derivs_xmin.has_value() || s_nbe_xmin == 0);
842 } else {
843 assert(!derivs_xmin.has_value() || derivs_xmin->template extent<deriv_type>() == 0);
844 }
845 if constexpr (BcUpper == BoundCond::HERMITE) {
846 assert(ddc::DiscreteElement<deriv_type>(derivs_xmax->domain().front()).uid() == s_odd);
847 assert(derivs_xmax.has_value() || s_nbe_xmax == 0);
848 } else {
849 assert(!derivs_xmax.has_value() || derivs_xmax->template extent<deriv_type>() == 0);
850 }
851
852 // Hermite boundary conditions at xmin, if any
853 // NOTE: For consistency with the linear system, the i-th derivative
854 // provided by the user must be multiplied by dx^i
855 if constexpr (BcLower == BoundCond::HERMITE) {
856 assert(derivs_xmin->template extent<deriv_type>() == s_nbe_xmin);
857 auto derivs_xmin_values = *derivs_xmin;
858 auto const dx_proxy = m_dx;
859 auto const odd_proxy = s_odd;
860 ddc::parallel_for_each(
861 "ddc_splines_hermite_compute_lower_coefficients",
862 exec_space(),
863 batch_domain(batched_interpolation_domain),
864 KOKKOS_LAMBDA(
865 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
866 j) {
867 for (int i = 0; i < s_nbe_xmin; ++i) {
868 spline(ddc::DiscreteElement<bsplines_type>(i), j)
869 = derivs_xmin_values(
870 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
871 j)
872 * ddc::detail::ipow(dx_proxy, i + odd_proxy);
873 }
874 });
875 } else if constexpr (BcLower == BoundCond::HOMOGENEOUS_HERMITE) {
876 ddc::DiscreteDomain<bsplines_type> const dx_splines(
877 ddc::DiscreteElement<bsplines_type>(0),
878 ddc::DiscreteVector<bsplines_type>(s_nbe_xmin));
879 batched_spline_domain_type<BatchedInterpolationDDom> const
880 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
881 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
882 }
883
884 // Fill spline with vals (to work in spline afterward and preserve vals)
885 ddc::parallel_fill(
886 exec_space(),
887 spline[ddc::DiscreteDomain<bsplines_type>(
888 ddc::DiscreteElement<bsplines_type>(s_nbe_xmin),
889 ddc::DiscreteVector<bsplines_type>(m_offset))],
890 0.);
891 // NOTE: We rely on Kokkos::deep_copy because ddc::parallel_deepcopy do not support
892 // different domain-typed Chunks.
893 Kokkos::deep_copy(
894 exec_space(),
895 spline[ddc::DiscreteDomain<bsplines_type>(
896 ddc::DiscreteElement<bsplines_type>(s_nbe_xmin + m_offset),
897 ddc::DiscreteVector<bsplines_type>(static_cast<std::size_t>(
898 vals.domain()
899 .template extent<
900 interpolation_discrete_dimension_type>())))]
901 .allocation_kokkos_view(),
902 vals.allocation_kokkos_view());
903
904
905
906 // Hermite boundary conditions at xmax, if any
907 // NOTE: For consistency with the linear system, the i-th derivative
908 // provided by the user must be multiplied by dx^i
909 auto const& nbasis_proxy = ddc::discrete_space<bsplines_type>().nbasis();
910 if constexpr (BcUpper == BoundCond::HERMITE) {
911 assert(derivs_xmax->template extent<deriv_type>() == s_nbe_xmax);
912 auto derivs_xmax_values = *derivs_xmax;
913 auto const dx_proxy = m_dx;
914 auto const odd_proxy = s_odd;
915 ddc::parallel_for_each(
916 "ddc_splines_hermite_compute_upper_coefficients",
917 exec_space(),
918 batch_domain(batched_interpolation_domain),
919 KOKKOS_LAMBDA(
920 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
921 j) {
922 for (int i = 0; i < s_nbe_xmax; ++i) {
923 spline(ddc::DiscreteElement<bsplines_type>(nbasis_proxy - s_nbe_xmax + i),
924 j)
925 = derivs_xmax_values(
926 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
927 j)
928 * ddc::detail::ipow(dx_proxy, i + odd_proxy);
929 }
930 });
931 } else if constexpr (BcUpper == BoundCond::HOMOGENEOUS_HERMITE) {
932 ddc::DiscreteDomain<bsplines_type> const dx_splines(
933 ddc::DiscreteElement<bsplines_type>(nbasis_proxy - s_nbe_xmax),
934 ddc::DiscreteVector<bsplines_type>(s_nbe_xmax));
935 batched_spline_domain_type<BatchedInterpolationDDom> const
936 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
937 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
938 }
939
940 // 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
941 auto const& offset_proxy = m_offset;
942 ddc::Chunk spline_tr_alloc(
943 batched_spline_tr_domain(batched_interpolation_domain),
944 ddc::KokkosAllocator<double, memory_space>());
945 ddc::ChunkSpan const spline_tr = spline_tr_alloc.span_view();
946 ddc::parallel_for_each(
947 "ddc_splines_transpose_rhs",
948 exec_space(),
949 batch_domain(batched_interpolation_domain),
950 KOKKOS_LAMBDA(
951 typename batch_domain_type<
952 BatchedInterpolationDDom>::discrete_element_type const j) {
953 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
954 spline_tr(ddc::DiscreteElement<bsplines_type>(i), j)
955 = spline(ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j);
956 }
957 });
958 // Create a 2D Kokkos::View to manage spline_tr as a matrix
959 Kokkos::View<double**, Kokkos::LayoutRight, exec_space> const bcoef_section(
960 spline_tr.data_handle(),
961 static_cast<std::size_t>(spline_tr.template extent<bsplines_type>()),
962 batch_domain(batched_interpolation_domain).size());
963 // Compute spline coef
964 m_matrix->solve(bcoef_section, false);
965 // Transpose back spline_tr into spline.
966 ddc::parallel_for_each(
967 "ddc_splines_transpose_back_rhs",
968 exec_space(),
969 batch_domain(batched_interpolation_domain),
970 KOKKOS_LAMBDA(
971 typename batch_domain_type<
972 BatchedInterpolationDDom>::discrete_element_type const j) {
973 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
974 spline(ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j)
975 = spline_tr(ddc::DiscreteElement<bsplines_type>(i), j);
976 }
977 });
978
979 // Duplicate the lower spline coefficients to the upper side in case of periodic boundaries
980 if (bsplines_type::is_periodic()) {
981 ddc::parallel_for_each(
982 "ddc_splines_periodic_rows_duplicate_rhs",
983 exec_space(),
984 batch_domain(batched_interpolation_domain),
985 KOKKOS_LAMBDA(
986 typename batch_domain_type<
987 BatchedInterpolationDDom>::discrete_element_type const j) {
988 if (offset_proxy != 0) {
989 for (int i = 0; i < offset_proxy; ++i) {
990 spline(ddc::DiscreteElement<bsplines_type>(i), j) = spline(
991 ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i),
992 j);
993 }
994 for (std::size_t i = offset_proxy; i < bsplines_type::degree(); ++i) {
995 spline(ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i), j)
996 = spline(ddc::DiscreteElement<bsplines_type>(i), j);
997 }
998 }
999 for (std::size_t i(0); i < bsplines_type::degree(); ++i) {
1000 ddc::DiscreteElement<bsplines_type> const i_start(i);
1001 ddc::DiscreteElement<bsplines_type> const i_end(nbasis_proxy + i);
1002
1003 spline(i_end, j) = spline(i_start, j);
1004 }
1005 });
1006 }
1007}
1008
1009template <
1010 class ExecSpace,
1011 class MemorySpace,
1012 class BSplines,
1013 class InterpolationDDim,
1017template <class OutMemorySpace>
1018std::tuple<
1019 ddc::Chunk<
1020 double,
1022 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
1023 ddc::KokkosAllocator<double, OutMemorySpace>>,
1024 ddc::Chunk<
1025 double,
1026 ddc::DiscreteDomain<InterpolationDDim>,
1027 ddc::KokkosAllocator<double, OutMemorySpace>>,
1028 ddc::Chunk<
1029 double,
1031 ddc::Deriv<typename InterpolationDDim::continuous_dimension_type>>,
1032 ddc::KokkosAllocator<double, OutMemorySpace>>>
1033SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1035{
1036 // Compute integrals of bsplines
1037 ddc::Chunk integral_bsplines(spline_domain(), ddc::KokkosAllocator<double, MemorySpace>());
1038 ddc::integrals(ExecSpace(), integral_bsplines.span_view());
1039
1040 // Remove additional B-splines in the periodic case (cf. UniformBSplines::full_domain() documentation)
1041 ddc::ChunkSpan const integral_bsplines_without_periodic_additional_bsplines
1042 = integral_bsplines[spline_domain().take_first(
1043 ddc::DiscreteVector<bsplines_type>(m_matrix->size()))];
1044
1045 // Allocate mirror with additional rows (cf. SplinesLinearProblem3x3Blocks documentation)
1046 Kokkos::View<double**, Kokkos::LayoutRight, MemorySpace> const
1047 integral_bsplines_mirror_with_additional_allocation(
1048 "integral_bsplines_mirror_with_additional_allocation",
1049 m_matrix->required_number_of_rhs_rows(),
1050 1);
1051
1052 // Extract relevant subview
1053 Kokkos::View<double*, Kokkos::LayoutRight, MemorySpace> const integral_bsplines_mirror
1054 = Kokkos::
1055 subview(integral_bsplines_mirror_with_additional_allocation,
1056 std::
1057 pair {static_cast<std::size_t>(0),
1058 integral_bsplines_without_periodic_additional_bsplines
1059 .size()},
1060 0);
1061
1062 // Solve matrix equation A^t*X=integral_bsplines
1063 Kokkos::deep_copy(
1064 integral_bsplines_mirror,
1065 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view());
1066 m_matrix->solve(integral_bsplines_mirror_with_additional_allocation, true);
1067 Kokkos::deep_copy(
1068 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view(),
1069 integral_bsplines_mirror);
1070
1071 // Slice into three ChunkSpan corresponding to lower derivatives, function values and upper derivatives
1072 ddc::ChunkSpan const coefficients_derivs_xmin
1073 = integral_bsplines_without_periodic_additional_bsplines[spline_domain().take_first(
1074 ddc::DiscreteVector<bsplines_type>(s_nbv_xmin))];
1075 ddc::ChunkSpan const coefficients = integral_bsplines_without_periodic_additional_bsplines
1077 .remove_first(ddc::DiscreteVector<bsplines_type>(s_nbv_xmin))
1078 .take_first(
1079 ddc::DiscreteVector<bsplines_type>(
1080 ddc::discrete_space<bsplines_type>().nbasis() - s_nbv_xmin
1081 - s_nbv_xmax))];
1082 ddc::ChunkSpan const coefficients_derivs_xmax
1083 = integral_bsplines_without_periodic_additional_bsplines
1085 .remove_first(
1086 ddc::DiscreteVector<bsplines_type>(
1087 s_nbv_xmin + coefficients.size()))
1088 .take_first(ddc::DiscreteVector<bsplines_type>(s_nbv_xmax))];
1089
1090 // Multiply derivatives coefficients by dx^n
1091 auto const dx_proxy = m_dx;
1092 auto const odd_proxy = s_odd;
1093 ddc::parallel_for_each(
1094 exec_space(),
1095 coefficients_derivs_xmin.domain(),
1096 KOKKOS_LAMBDA(ddc::DiscreteElement<bsplines_type> i) {
1097 coefficients_derivs_xmin(i) *= ddc::detail::
1098 ipow(dx_proxy,
1099 static_cast<std::size_t>(get<bsplines_type>(
1100 (i - coefficients_derivs_xmin.domain().front()) + odd_proxy)));
1101 });
1102 ddc::parallel_for_each(
1103 exec_space(),
1104 coefficients_derivs_xmax.domain(),
1105 KOKKOS_LAMBDA(ddc::DiscreteElement<bsplines_type> i) {
1106 coefficients_derivs_xmax(i) *= ddc::detail::
1107 ipow(dx_proxy,
1108 static_cast<std::size_t>(get<bsplines_type>(
1109 (i - coefficients_derivs_xmax.domain().front()) + odd_proxy)));
1110 });
1111
1112 ddc::DiscreteElement<deriv_type> const first_deriv(s_odd);
1113 // Allocate Chunk on deriv_type and interpolation_discrete_dimension_type and copy quadrature coefficients into it
1114 ddc::Chunk coefficients_derivs_xmin_out(
1116 deriv_type>(first_deriv, ddc::DiscreteVector<deriv_type>(s_nbv_xmin)),
1117 ddc::KokkosAllocator<double, OutMemorySpace>());
1118 ddc::Chunk coefficients_out(
1119 interpolation_domain().take_first(
1120 ddc::DiscreteVector<interpolation_discrete_dimension_type>(
1121 coefficients.size())),
1122 ddc::KokkosAllocator<double, OutMemorySpace>());
1123 ddc::Chunk coefficients_derivs_xmax_out(
1125 deriv_type>(first_deriv, ddc::DiscreteVector<deriv_type>(s_nbv_xmax)),
1126 ddc::KokkosAllocator<double, OutMemorySpace>());
1127 Kokkos::deep_copy(
1128 coefficients_derivs_xmin_out.allocation_kokkos_view(),
1129 coefficients_derivs_xmin.allocation_kokkos_view());
1130 Kokkos::deep_copy(
1131 coefficients_out.allocation_kokkos_view(),
1132 coefficients.allocation_kokkos_view());
1133 Kokkos::deep_copy(
1134 coefficients_derivs_xmax_out.allocation_kokkos_view(),
1135 coefficients_derivs_xmax.allocation_kokkos_view());
1136 return std::make_tuple(
1137 std::move(coefficients_derivs_xmin_out),
1138 std::move(coefficients_out),
1139 std::move(coefficients_derivs_xmax_out));
1140}
1141
1142template <
1143 class ExecSpace,
1144 class MemorySpace,
1145 class BSplines,
1146 class InterpolationDDim,
1150template <class KnotElement>
1151void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1152 check_n_points_in_cell(int const n_points_in_cell, KnotElement const current_cell_end_idx)
1153{
1154 if (n_points_in_cell > BSplines::degree() + 1) {
1155 KnotElement const rmin_idx = ddc::discrete_space<BSplines>().break_point_domain().front();
1156 int const failed_cell = (current_cell_end_idx - rmin_idx).value();
1157 throw std::runtime_error(
1158 "The spline problem is overconstrained. There are "
1159 + std::to_string(n_points_in_cell) + " points in the " + std::to_string(failed_cell)
1160 + "-th cell.");
1161 }
1162}
1163
1164template <
1165 class ExecSpace,
1166 class MemorySpace,
1167 class BSplines,
1168 class InterpolationDDim,
1169 ddc::BoundCond BcLower,
1170 ddc::BoundCond BcUpper,
1171 SplineSolver Solver>
1172void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1173 check_valid_grid()
1174{
1175 std::size_t const n_interp_points = interpolation_domain().size();
1176 std::size_t const expected_npoints
1177 = ddc::discrete_space<BSplines>().nbasis() - s_nbe_xmin - s_nbe_xmax;
1178 if (n_interp_points != expected_npoints) {
1179 throw std::runtime_error(
1180 "Incorrect number of points supplied to NonUniformInterpolationPoints. "
1181 "(Received : "
1182 + std::to_string(n_interp_points)
1183 + ", expected : " + std::to_string(expected_npoints));
1184 }
1185 int n_points_in_cell = 0;
1186 auto current_cell_end_idx = ddc::discrete_space<BSplines>().break_point_domain().front() + 1;
1187 ddc::host_for_each(interpolation_domain(), [&](auto idx) {
1188 ddc::Coordinate<continuous_dimension_type> const point = ddc::coordinate(idx);
1189 if (point > ddc::coordinate(current_cell_end_idx)) {
1190 // Check the points found in the previous cell
1191 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1192 // Initialise the number of points in the subsequent cell, including the new point
1193 n_points_in_cell = 1;
1194 // Move to the next cell
1195 current_cell_end_idx += 1;
1196 } else if (point == ddc::coordinate(current_cell_end_idx)) {
1197 // Check the points found in the previous cell including the point on the boundary
1198 check_n_points_in_cell(n_points_in_cell + 1, current_cell_end_idx);
1199 // Initialise the number of points in the subsequent cell, including the point on the boundary
1200 n_points_in_cell = 1;
1201 // Move to the next cell
1202 current_cell_end_idx += 1;
1203 } else {
1204 // Indicate that the point is in the cell
1205 n_points_in_cell += 1;
1206 }
1207 });
1208 // Check the number of points in the final cell
1209 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1210}
1211
1212} // 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(SplineBuilder const &x)=delete
Copy-constructor is deleted.
SplineBuilder(BatchedInterpolationDDom const &batched_interpolation_domain, std::optional< std::size_t > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditioner_max_block_size=std::nullopt)
Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
static constexpr int s_nbv_xmax
The number of input values defining the boundary condition at the upper bound.
static constexpr bool s_odd
Indicates if the degree of the splines is odd or even.
static constexpr int s_nbv_xmin
The number of input values defining the boundary condition at the lower bound.
interpolation_domain_type interpolation_domain() const noexcept
Get the domain for the 1D interpolation mesh used by this class.
batched_derivs_domain_type< BatchedInterpolationDDom > batched_derivs_xmin_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain on which derivatives on lower boundary are defined.
static constexpr ddc::BoundCond s_bc_xmax
The boundary condition implemented at the upper bound.
void operator()(ddc::ChunkSpan< double, batched_spline_domain_type< BatchedInterpolationDDom >, Layout, memory_space > spline, ddc::ChunkSpan< double const, BatchedInterpolationDDom, Layout, memory_space > vals, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > derivs_xmin=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > derivs_xmax=std::nullopt) const
Compute a spline approximation of a function.
batched_spline_domain_type< BatchedInterpolationDDom > batched_spline_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain on which spline coefficients are defined.
ddc::DiscreteDomain< bsplines_type > spline_domain() const noexcept
Get the 1D domain on which spline coefficients are defined.
~SplineBuilder()=default
Destructs.
static constexpr int s_nbe_xmin
The number of equations defining the boundary condition at the lower bound.
Storage class of the static attributes of the discrete dimension.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_last_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the last support knot associated to a DiscreteElement identifying a B-splin...
Impl(ddc::Coordinate< CDim > rmin, ddc::Coordinate< CDim > rmax, std::size_t ncells)
Constructs a spline basis (B-splines) with n equidistant knots over .
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmax() const noexcept
Returns the coordinate of the upper bound of the domain on which the B-splines are defined.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis(DSpan1D values, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-splines at a given coordinate.
KOKKOS_INLINE_FUNCTION ddc::DiscreteDomain< knot_discrete_dimension_type > break_point_domain() const
Returns the discrete domain which describes the break points.
KOKKOS_INLINE_FUNCTION ddc::Coordinate< CDim > rmin() const noexcept
Returns the coordinate of the lower bound of the domain on which the B-splines are defined.
~Impl()=default
Destructs.
KOKKOS_INLINE_FUNCTION std::size_t nbasis() const noexcept
Returns the number of basis functions.
Impl(Impl const &x)=default
Copy-constructs.
KOKKOS_INLINE_FUNCTION std::size_t size() const noexcept
Returns the number of elements necessary to construct a spline representation of a function.
Impl(Impl &&x)=default
Move-constructs.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_basis_and_n_derivs(ddc::DSpan2D derivs, ddc::Coordinate< CDim > const &x, std::size_t n) const
Evaluates non-zero B-spline values and derivatives at a given coordinate.
KOKKOS_INLINE_FUNCTION ddc::DiscreteElement< knot_discrete_dimension_type > get_first_support_knot(discrete_element_type const &ix) const
Returns the coordinate of the first support knot associated to a DiscreteElement identifying a B-spli...
KOKKOS_INLINE_FUNCTION double length() const noexcept
Returns the length of the domain.
Impl(Impl< DDim, OriginMemorySpace > const &impl)
Copy-constructs from another Impl with a different Kokkos memory space.
KOKKOS_INLINE_FUNCTION std::size_t ncells() const noexcept
Returns the number of cells over which the B-splines are defined.
KOKKOS_INLINE_FUNCTION discrete_element_type eval_deriv(DSpan1D derivs, ddc::Coordinate< CDim > const &x) const
Evaluates non-zero B-spline derivatives at a given coordinate.
Impl & operator=(Impl &&x)=default
Move-assigns.
KOKKOS_INLINE_FUNCTION discrete_domain_type full_domain() const
Returns the discrete domain including eventual additional B-splines in the periodic case.
Impl & operator=(Impl const &x)=default
Copy-assigns.
The type of a uniform 1D spline basis (B-spline).
static constexpr bool is_uniform() noexcept
Indicates if the B-splines are uniform or not (this is the case here).
static constexpr bool is_periodic() noexcept
Indicates if the B-splines are periodic or not.
static constexpr std::size_t degree() noexcept
The degree of B-splines.
UniformPointSampling models a uniform discretization of the provided continuous dimension.
The top-level namespace of DDC.
constexpr int n_boundary_equations(ddc::BoundCond const bc, std::size_t const degree)
Return the number of equations needed to describe a given boundary condition.
constexpr bool is_uniform_bsplines_v
Indicates if a tag corresponds to uniform B-splines or not.
BoundCond
An enum representing a spline boundary condition.
@ HOMOGENEOUS_HERMITE
Homogeneous Hermite boundary condition (derivatives are 0)
@ GREVILLE
Use Greville points instead of conditions on derivative for B-Spline interpolation.
@ HERMITE
Hermite boundary condition.
@ PERIODIC
Periodic boundary condition u(1)=u(n)
ddc::ChunkSpan< double, ddc::DiscreteDomain< DDim >, Layout, MemorySpace > integrals(ExecSpace const &execution_space, ddc::ChunkSpan< double, ddc::DiscreteDomain< DDim >, Layout, MemorySpace > int_vals)
Compute the integrals of the B-splines.
SplineSolver
An enum determining the backend solver of a SplineBuilder or SplineBuilder2d.
@ LAPACK
Enum member to identify the LAPACK-based solver (direct method)
@ GINKGO
Enum member to identify the Ginkgo-based solver (iterative method)
constexpr bool is_non_uniform_bsplines_v
Indicates if a tag corresponds to non-uniform B-splines or not.
A templated struct representing a discrete dimension storing the derivatives of a function along a co...
Definition deriv.hpp:15
If the type DDim is a B-spline, defines type to the discrete dimension of the associated knots.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos, ddc::Coordinate< DimNI > eval_pos_not_interest_min, ddc::Coordinate< DimNI > eval_pos_not_interest_max)
Instantiate a ConstantExtrapolationRule.
KOKKOS_FUNCTION double operator()(CoordType coord_extrap, ddc::ChunkSpan< double const, ddc::DiscreteDomain< BSplines1, BSplines2 >, Layout, MemorySpace > const spline_coef) const
Get the value of the function on B-splines at a coordinate outside the domain.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos)
Instantiate a ConstantExtrapolationRule.
KOKKOS_FUNCTION double operator()(CoordType pos, ddc::ChunkSpan< double const, ddc::DiscreteDomain< BSplines >, Layout, MemorySpace > const spline_coef) const
Get the value of the function on B-splines at a coordinate outside the domain.
ConstantExtrapolationRule(ddc::Coordinate< DimI > eval_pos)
Instantiate a ConstantExtrapolationRule.
A functor describing a null extrapolation boundary value for 1D spline evaluator.
KOKKOS_FUNCTION double operator()(CoordType, ChunkSpan) const
Evaluates the spline at a coordinate outside of the domain.
KOKKOS_FUNCTION double operator()(CoordType, ChunkSpan) const