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