20#include <Kokkos_Core.hpp>
23#include "integrals.hpp"
24#include "math_tools.hpp"
25#include "spline_boundary_conditions.hpp"
26#include "splines_linear_problem.hpp"
27#include "splines_linear_problem_maker.hpp"
33
34
35
36
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
62 class InterpolationDDim,
76 using exec_space = ExecSpace;
79 using memory_space = MemorySpace;
82 using continuous_dimension_type =
typename InterpolationDDim::continuous_dimension_type;
85 using interpolation_discrete_dimension_type = InterpolationDDim;
88 using bsplines_type = BSplines;
91 using deriv_type =
ddc::
Deriv<continuous_dimension_type>;
94 using interpolation_domain_type =
ddc::
DiscreteDomain<interpolation_discrete_dimension_type>;
97
98
99
100
102 class BatchedInterpolationDDom,
103 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
104 using batched_interpolation_domain_type = BatchedInterpolationDDom;
107
108
109
110
111
112
113
114
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>;
122
123
124
125
126
127
128
129
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,
140
141
142
143
144
145
146
147
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>>>>;
160
161
162
163
164
165
166
167
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,
177 static constexpr bool s_odd = BSplines::degree() % 2;
180 static constexpr int s_nbc_xmin = n_boundary_equations(BcLower, BSplines::degree());
183 static constexpr int s_nbc_xmax = n_boundary_equations(BcUpper, BSplines::degree());
195 interpolation_domain_type m_interpolation_domain;
202 std::unique_ptr<
ddc::detail::SplinesLinearProblem<exec_space>> m_matrix;
205 void compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset);
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
225 interpolation_domain_type
const& interpolation_domain,
226 std::optional<std::size_t> cols_per_chunk = std::nullopt,
227 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
228 : m_interpolation_domain(interpolation_domain)
229 , m_dx((
ddc::discrete_space<BSplines>().rmax() -
ddc::discrete_space<BSplines>().rmin())
230 /
ddc::discrete_space<BSplines>().ncells())
234 "Incompatible boundary conditions");
237 compute_offset(
this->interpolation_domain(), m_offset);
240 int lower_block_size;
241 int upper_block_size;
242 if constexpr (bsplines_type::is_uniform()) {
243 upper_block_size = compute_block_sizes_uniform(BcLower,
s_nbc_xmin);
244 lower_block_size = compute_block_sizes_uniform(BcUpper,
s_nbc_xmax);
246 upper_block_size = compute_block_sizes_non_uniform(BcLower,
s_nbc_xmin);
247 lower_block_size = compute_block_sizes_non_uniform(BcUpper,
s_nbc_xmax);
253 preconditioner_max_block_size);
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
273 class BatchedInterpolationDDom,
274 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
276 BatchedInterpolationDDom
const& batched_interpolation_domain,
277 std::optional<std::size_t> cols_per_chunk = std::nullopt,
278 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
280 interpolation_domain_type(batched_interpolation_domain),
282 preconditioner_max_block_size)
291
292
293
303
304
305
306
310
311
312
313
314
315
318 return m_interpolation_domain;
322
323
324
325
326
327
328
329
330
332 class BatchedInterpolationDDom,
333 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
335 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
337 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
338 return batched_interpolation_domain;
342
343
344
345
346
347
348
349
350 template <
class BatchedInterpolationDDom>
351 batch_domain_type<BatchedInterpolationDDom>
batch_domain(
352 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
354 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
359
360
361
362
363
364
367 return ddc::discrete_space<bsplines_type>().full_domain();
371
372
373
374
375
376
377
378
379 template <
class BatchedInterpolationDDom>
381 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
383 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
384 return ddc::replace_dim_of<
385 interpolation_discrete_dimension_type,
391
392
393
394
395
396
397
398
399 template <
class BatchedInterpolationDDom>
400 batched_spline_tr_domain_type<BatchedInterpolationDDom> batched_spline_tr_domain(
401 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
403 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
404 return batched_spline_tr_domain_type<BatchedInterpolationDDom>(
405 ddc::replace_dim_of<bsplines_type, bsplines_type>(
406 batched_spline_domain(batched_interpolation_domain),
408 ddc::DiscreteElement<bsplines_type>(0),
410 m_matrix->required_number_of_rhs_rows()))));
415
416
417
418
419
420
421
422
423 template <
class BatchedInterpolationDDom>
425 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
427 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
428 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
429 batched_interpolation_domain,
431 ddc::DiscreteElement<deriv_type>(1),
436
437
438
439
440
441
442
443
444 template <
class BatchedInterpolationDDom>
446 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
448 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
449 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
450 batched_interpolation_domain,
452 ddc::DiscreteElement<deriv_type>(1),
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475 template <
class Layout,
class BatchedInterpolationDDom>
479 batched_spline_domain_type<BatchedInterpolationDDom>,
481 memory_space> spline,
482 ddc::
ChunkSpan<
double const, BatchedInterpolationDDom, Layout, memory_space> vals,
485 batched_derivs_domain_type<BatchedInterpolationDDom>,
487 memory_space>> derivs_xmin
491 batched_derivs_domain_type<BatchedInterpolationDDom>,
493 memory_space>> derivs_xmax
494 = std::nullopt)
const;
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519 template <
class OutMemorySpace = MemorySpace>
538 static int compute_block_sizes_uniform(
ddc::
BoundCond bound_cond,
int nbc);
540 static int compute_block_sizes_non_uniform(
ddc::
BoundCond bound_cond,
int nbc);
542 void allocate_matrix(
543 int lower_block_size,
544 int upper_block_size,
545 std::optional<std::size_t> cols_per_chunk = std::nullopt,
546 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt);
548 void build_matrix_system();
550 void check_valid_grid();
552 template <
class KnotElement>
553 static void check_n_points_in_cell(
int n_points_in_cell, KnotElement current_cell_end_idx);
560 class InterpolationDDim,
564void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
565 compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset)
567 if constexpr (bsplines_type::is_periodic()) {
569 std::array<
double, bsplines_type::degree() + 1> values_ptr;
570 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>>
const
571 values(values_ptr.data());
572 ddc::DiscreteElement<interpolation_discrete_dimension_type> start(
573 interpolation_domain.front());
574 auto jmin =
ddc::discrete_space<BSplines>()
575 .eval_basis(values,
ddc::coordinate(start + BSplines::degree()));
576 if constexpr (bsplines_type::degree() % 2 == 0) {
577 offset = jmin.uid() - start.uid() + bsplines_type::degree() / 2 - BSplines::degree();
579 int const mid = bsplines_type::degree() / 2;
580 offset = jmin.uid() - start.uid()
581 + (DDC_MDSPAN_ACCESS_OP(values, mid) > DDC_MDSPAN_ACCESS_OP(values, mid + 1)
584 - BSplines::degree();
595 class InterpolationDDim,
599int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
600 compute_block_sizes_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
603 return static_cast<
int>(bsplines_type::degree()) / 2;
611 return static_cast<
int>(bsplines_type::degree()) - 1;
614 throw std::runtime_error(
"ddc::BoundCond not handled");
621 class InterpolationDDim,
625int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
626 compute_block_sizes_non_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
629 return static_cast<
int>(bsplines_type::degree()) - 1;
636 throw std::runtime_error(
"ddc::BoundCond not handled");
643 class InterpolationDDim,
647void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
649 [[maybe_unused]]
int lower_block_size,
650 [[maybe_unused]]
int upper_block_size,
651 std::optional<std::size_t> cols_per_chunk,
652 std::optional<
unsigned int> preconditioner_max_block_size)
662 int upper_band_width;
663 if (bsplines_type::is_uniform()) {
664 upper_band_width = bsplines_type::degree() / 2;
666 upper_band_width = bsplines_type::degree() - 1;
668 if constexpr (bsplines_type::is_periodic()) {
669 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_periodic_band_matrix<
671 ddc::discrete_space<BSplines>().nbasis(),
674 bsplines_type::is_uniform());
676 m_matrix =
ddc::detail::SplinesLinearProblemMaker::
677 make_new_block_matrix_with_band_main_block<ExecSpace>(
678 ddc::discrete_space<BSplines>().nbasis(),
681 bsplines_type::is_uniform(),
686 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_sparse<ExecSpace>(
687 ddc::discrete_space<BSplines>().nbasis(),
689 preconditioner_max_block_size);
692 build_matrix_system();
694 m_matrix->setup_solver();
701 class InterpolationDDim,
705void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
706 build_matrix_system()
710 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
713 derivs(derivs_ptr.data(),
714 bsplines_type::degree() + 1,
715 bsplines_type::degree() / 2 + 1);
716 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
718 ddc::discrete_space<BSplines>().rmin(),
723 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
724 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
725 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *=
ddc::detail::ipow(m_dx, j);
730 for (std::size_t i = 0; i <
s_nbc_xmin; ++i) {
731 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
732 m_matrix->set_element(
735 DDC_MDSPAN_ACCESS_OP(derivs, j, s_nbc_xmin - i - 1 + s_odd));
741 std::array<
double, bsplines_type::degree() + 1> values_ptr;
742 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>>
const values(
747 auto jmin =
ddc::discrete_space<BSplines>().eval_basis(
749 ddc::coordinate(
ddc::DiscreteElement<interpolation_discrete_dimension_type>(ix)));
750 for (std::size_t s = 0; s < bsplines_type::degree() + 1; ++s) {
751 int const j =
ddc::detail::
752 modulo(
int(jmin.uid() - m_offset + s),
753 static_cast<
int>(
ddc::discrete_space<BSplines>().nbasis()));
754 m_matrix->set_element(
757 DDC_MDSPAN_ACCESS_OP(values, s));
763 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
769 bsplines_type::degree() + 1,
770 bsplines_type::degree() / 2 + 1>>
const derivs(derivs_ptr.data());
772 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
774 ddc::discrete_space<BSplines>().rmax(),
779 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
780 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
781 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *=
ddc::detail::ipow(m_dx, j);
785 int const i0 =
ddc::discrete_space<BSplines>().nbasis() -
s_nbc_xmax;
786 int const j0 =
ddc::discrete_space<BSplines>().nbasis() - bsplines_type::degree();
787 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
788 for (std::size_t i = 0; i <
s_nbc_xmax; ++i) {
789 m_matrix->set_element(
792 DDC_MDSPAN_ACCESS_OP(derivs, j + 1, i + s_odd));
806template <
class Layout,
class BatchedInterpolationDDom>
807void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
811 batched_spline_domain_type<BatchedInterpolationDDom>,
813 memory_space> spline,
814 ddc::
ChunkSpan<
double const, BatchedInterpolationDDom, Layout, memory_space> vals,
817 batched_derivs_domain_type<BatchedInterpolationDDom>,
819 memory_space>>
const derivs_xmin,
822 batched_derivs_domain_type<BatchedInterpolationDDom>,
824 memory_space>>
const derivs_xmax)
const
826 auto const batched_interpolation_domain = vals.domain();
828 assert(interpolation_domain() == interpolation_domain_type(batched_interpolation_domain));
830 assert(vals.
template extent<interpolation_discrete_dimension_type>()
831 == ddc::discrete_space<bsplines_type>().nbasis() - s_nbc_xmin - s_nbc_xmax);
833 assert((BcLower == ddc::BoundCond::HERMITE)
834 != (!derivs_xmin.has_value() || derivs_xmin->
template extent<deriv_type>() == 0));
835 assert((BcUpper == ddc::BoundCond::HERMITE)
836 != (!derivs_xmax.has_value() || derivs_xmax->
template extent<deriv_type>() == 0));
838 assert(ddc::DiscreteElement<deriv_type>(derivs_xmin->domain().front()).uid() == 1);
841 assert(ddc::DiscreteElement<deriv_type>(derivs_xmax->domain().front()).uid() == 1);
848 assert(derivs_xmin->
template extent<deriv_type>() == s_nbc_xmin);
849 auto derivs_xmin_values = *derivs_xmin;
850 auto const dx_proxy = m_dx;
851 ddc::parallel_for_each(
852 "ddc_splines_hermite_compute_lower_coefficients",
854 batch_domain(batched_interpolation_domain),
856 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
859 spline(
ddc::DiscreteElement<bsplines_type>(
s_nbc_xmin - i), j)
860 = derivs_xmin_values(
ddc::DiscreteElement<deriv_type>(i), j)
861 *
ddc::detail::ipow(dx_proxy, i +
s_odd - 1);
882 interpolation_discrete_dimension_type>())))]
883 .allocation_kokkos_view(),
884 vals.allocation_kokkos_view());
891 auto const& nbasis_proxy =
ddc::discrete_space<bsplines_type>().nbasis();
893 assert(derivs_xmax->
template extent<deriv_type>() == s_nbc_xmax);
894 auto derivs_xmax_values = *derivs_xmax;
895 auto const dx_proxy = m_dx;
896 ddc::parallel_for_each(
897 "ddc_splines_hermite_compute_upper_coefficients",
899 batch_domain(batched_interpolation_domain),
901 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
904 spline(
ddc::DiscreteElement<bsplines_type>(nbasis_proxy -
s_nbc_xmax - i),
906 = derivs_xmax_values(
ddc::DiscreteElement<deriv_type>(i + 1), j)
907 *
ddc::detail::ipow(dx_proxy, i +
s_odd);
913 auto const& offset_proxy = m_offset;
915 batched_spline_tr_domain(batched_interpolation_domain),
917 ddc::
ChunkSpan const spline_tr = spline_tr_alloc.span_view();
918 ddc::parallel_for_each(
919 "ddc_splines_transpose_rhs",
921 batch_domain(batched_interpolation_domain),
923 typename batch_domain_type<
924 BatchedInterpolationDDom>::discrete_element_type
const j) {
925 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
926 spline_tr(
ddc::DiscreteElement<bsplines_type>(i), j)
927 = spline(
ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j);
931 Kokkos::View<
double**, Kokkos::LayoutRight, exec_space>
const bcoef_section(
932 spline_tr.data_handle(),
933 static_cast<std::size_t>(spline_tr.
template extent<bsplines_type>()),
934 batch_domain(batched_interpolation_domain).size());
936 m_matrix->solve(bcoef_section,
false);
938 ddc::parallel_for_each(
939 "ddc_splines_transpose_back_rhs",
941 batch_domain(batched_interpolation_domain),
943 typename batch_domain_type<
944 BatchedInterpolationDDom>::discrete_element_type
const j) {
945 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
946 spline(
ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j)
947 = spline_tr(
ddc::DiscreteElement<bsplines_type>(i), j);
952 if (bsplines_type::is_periodic()) {
953 ddc::parallel_for_each(
954 "ddc_splines_periodic_rows_duplicate_rhs",
956 batch_domain(batched_interpolation_domain),
958 typename batch_domain_type<
959 BatchedInterpolationDDom>::discrete_element_type
const j) {
960 if (offset_proxy != 0) {
961 for (
int i = 0; i < offset_proxy; ++i) {
962 spline(
ddc::DiscreteElement<bsplines_type>(i), j) = spline(
963 ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i),
966 for (std::size_t i = offset_proxy; i < bsplines_type::degree(); ++i) {
967 spline(
ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i), j)
968 = spline(
ddc::DiscreteElement<bsplines_type>(i), j);
971 for (std::size_t i(0); i < bsplines_type::degree(); ++i) {
972 ddc::DiscreteElement<bsplines_type>
const i_start(i);
973 ddc::DiscreteElement<bsplines_type>
const i_end(nbasis_proxy + i);
975 spline(i_end, j) = spline(i_start, j);
989template <
class OutMemorySpace>
1005SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1010 ddc::integrals(ExecSpace(), integral_bsplines.span_view());
1013 ddc::
ChunkSpan const integral_bsplines_without_periodic_additional_bsplines
1018 Kokkos::View<
double**, Kokkos::LayoutRight, MemorySpace>
const
1019 integral_bsplines_mirror_with_additional_allocation(
1020 "integral_bsplines_mirror_with_additional_allocation",
1021 m_matrix->required_number_of_rhs_rows(),
1025 Kokkos::View<
double*, Kokkos::LayoutRight, MemorySpace>
const integral_bsplines_mirror
1027 subview(integral_bsplines_mirror_with_additional_allocation,
1029 pair {
static_cast<std::size_t>(0),
1030 integral_bsplines_without_periodic_additional_bsplines
1036 integral_bsplines_mirror,
1037 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view());
1038 m_matrix->solve(integral_bsplines_mirror_with_additional_allocation,
true);
1040 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view(),
1041 integral_bsplines_mirror);
1045 = integral_bsplines_without_periodic_additional_bsplines[
spline_domain().take_first(
1047 ddc::
ChunkSpan const coefficients = integral_bsplines_without_periodic_additional_bsplines
1055 = integral_bsplines_without_periodic_additional_bsplines
1064 ddc::parallel_for_each(
1066 coefficients_derivs_xmin.domain(),
1067 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1068 ddc::Coordinate<continuous_dimension_type>
const dx
1069 =
ddc::distance_at_right(interpolation_domain_proxy.front() + 1);
1070 coefficients_derivs_xmin(i) *=
ddc::detail::
1072 static_cast<std::size_t>(get<bsplines_type>(
1073 i - coefficients_derivs_xmin.domain().front() + 1)));
1075 ddc::parallel_for_each(
1077 coefficients_derivs_xmax.domain(),
1078 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1079 ddc::Coordinate<continuous_dimension_type>
const dx
1080 =
ddc::distance_at_left(interpolation_domain_proxy.back() - 1);
1081 coefficients_derivs_xmax(i) *=
ddc::detail::
1083 static_cast<std::size_t>(get<bsplines_type>(
1084 i - coefficients_derivs_xmax.domain().front() + 1)));
1088 ddc::
Chunk coefficients_derivs_xmin_out(
1090 ddc::DiscreteElement<deriv_type>(1),
1096 coefficients.size())),
1098 ddc::
Chunk coefficients_derivs_xmax_out(
1100 ddc::DiscreteElement<deriv_type>(1),
1104 coefficients_derivs_xmin_out.allocation_kokkos_view(),
1105 coefficients_derivs_xmin.allocation_kokkos_view());
1107 coefficients_out.allocation_kokkos_view(),
1108 coefficients.allocation_kokkos_view());
1110 coefficients_derivs_xmax_out.allocation_kokkos_view(),
1111 coefficients_derivs_xmax.allocation_kokkos_view());
1112 return std::make_tuple(
1113 std::move(coefficients_derivs_xmin_out),
1114 std::move(coefficients_out),
1115 std::move(coefficients_derivs_xmax_out));
1126template <
class KnotElement>
1127void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1128 check_n_points_in_cell(
int const n_points_in_cell, KnotElement
const current_cell_end_idx)
1130 if (n_points_in_cell > BSplines::degree() + 1) {
1131 KnotElement
const rmin_idx =
ddc::discrete_space<BSplines>().break_point_domain().front();
1132 int const failed_cell = (current_cell_end_idx - rmin_idx).value();
1133 throw std::runtime_error(
1134 "The spline problem is overconstrained. There are "
1135 + std::to_string(n_points_in_cell) +
" points in the " + std::to_string(failed_cell)
1144 class InterpolationDDim,
1148void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1152 std::size_t
const expected_npoints
1154 if (n_interp_points != expected_npoints) {
1155 throw std::runtime_error(
1156 "Incorrect number of points supplied to NonUniformInterpolationPoints. "
1158 + std::to_string(n_interp_points)
1159 +
", expected : " + std::to_string(expected_npoints));
1161 int n_points_in_cell = 0;
1162 auto current_cell_end_idx =
ddc::discrete_space<BSplines>().break_point_domain().front() + 1;
1164 ddc::Coordinate<continuous_dimension_type>
const point =
ddc::coordinate(idx);
1165 if (point >
ddc::coordinate(current_cell_end_idx)) {
1167 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1169 n_points_in_cell = 1;
1171 current_cell_end_idx += 1;
1172 }
else if (point ==
ddc::coordinate(current_cell_end_idx)) {
1174 check_n_points_in_cell(n_points_in_cell + 1, current_cell_end_idx);
1176 n_points_in_cell = 1;
1178 current_cell_end_idx += 1;
1181 n_points_in_cell += 1;
1185 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
friend class DiscreteDomain
KOKKOS_FUNCTION constexpr bool operator!=(DiscreteVector< OTags... > const &rhs) const noexcept
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.
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.
static constexpr int s_nbc_xmin
The number of equations defining the boundary condition 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 bool s_odd
Indicates if the degree of the splines is odd or even.
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.
static constexpr int s_nbc_xmax
The number of equations defining the boundary condition 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.
The top-level namespace of DDC.
BoundCond
An enum representing a spline boundary condition.
@ 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)
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)
A templated struct representing a discrete dimension storing the derivatives of a function along a co...