19#include <Kokkos_Core.hpp>
32
33
34
35
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
61 class InterpolationDDim,
75 using exec_space = ExecSpace;
78 using memory_space = MemorySpace;
81 using continuous_dimension_type =
typename InterpolationDDim::continuous_dimension_type;
84 using interpolation_discrete_dimension_type = InterpolationDDim;
87 using bsplines_type = BSplines;
90 using deriv_type =
ddc::
Deriv<continuous_dimension_type>;
93 using interpolation_domain_type =
ddc::
DiscreteDomain<interpolation_discrete_dimension_type>;
96
97
98
99
100 template <concepts::discrete_domain BatchedInterpolationDDom>
101 using batched_interpolation_domain_type = BatchedInterpolationDDom;
104
105
106
107
108
109
110
111
112 template <concepts::discrete_domain BatchedInterpolationDDom>
113 using batch_domain_type =
ddc::
114 remove_dims_of_t<BatchedInterpolationDDom, interpolation_discrete_dimension_type>;
117
118
119
120
121
122
123
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,
133
134
135
136
137
138
139
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>>>>;
151
152
153
154
155
156
157
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,
166 static constexpr bool s_odd = BSplines::degree() % 2;
169 static constexpr int s_nbe_xmin = n_boundary_equations(BcLower, BSplines::degree());
172 static constexpr int s_nbe_xmax = n_boundary_equations(BcUpper, BSplines::degree());
177 : n_boundary_equations(BcLower, BSplines::degree());
182 : n_boundary_equations(BcUpper, BSplines::degree());
194 interpolation_domain_type m_interpolation_domain;
201 std::unique_ptr<
ddc::detail::SplinesLinearProblem<exec_space>> m_matrix;
204 void compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset);
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
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())
233 "Incompatible boundary conditions");
236 compute_offset(
this->interpolation_domain(), m_offset);
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);
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);
252 preconditioner_max_block_size);
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271 template <
concepts::discrete_domain BatchedInterpolationDDom>
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),
279 preconditioner_max_block_size)
288
289
290
300
301
302
303
307
308
309
310
311
312
315 return m_interpolation_domain;
319
320
321
322
323
324
325
326
327
328 template <
concepts::discrete_domain BatchedInterpolationDDom>
330 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
333 return batched_interpolation_domain;
337
338
339
340
341
342
343
344
345 template <
class BatchedInterpolationDDom>
346 batch_domain_type<BatchedInterpolationDDom>
batch_domain(
347 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
354
355
356
357
358
359
362 return ddc::discrete_space<bsplines_type>().full_domain();
366
367
368
369
370
371
372
373
374 template <
class BatchedInterpolationDDom>
376 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
379 return ddc::replace_dim_of<
380 interpolation_discrete_dimension_type,
386
387
388
389
390
391
392
393
394 template <
class BatchedInterpolationDDom>
395 batched_spline_tr_domain_type<BatchedInterpolationDDom> batched_spline_tr_domain(
396 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
399 return batched_spline_tr_domain_type<BatchedInterpolationDDom>(
400 ddc::replace_dim_of<bsplines_type, bsplines_type>(
401 batched_spline_domain(batched_interpolation_domain),
403 ddc::DiscreteElement<bsplines_type>(0),
405 m_matrix->required_number_of_rhs_rows()))));
410
411
412
413
414
415
416
417
418 template <
class BatchedInterpolationDDom>
420 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
423 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
424 batched_interpolation_domain,
426 ddc::DiscreteElement<deriv_type>(1),
431
432
433
434
435
436
437
438
439 template <
class BatchedInterpolationDDom>
441 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
444 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
445 batched_interpolation_domain,
447 ddc::DiscreteElement<deriv_type>(1),
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470 template <
class Layout,
class BatchedInterpolationDDom>
474 batched_spline_domain_type<BatchedInterpolationDDom>,
476 memory_space> spline,
477 ddc::
ChunkSpan<
double const, BatchedInterpolationDDom, Layout, memory_space> vals,
480 batched_derivs_domain_type<BatchedInterpolationDDom>,
482 memory_space>> derivs_xmin
486 batched_derivs_domain_type<BatchedInterpolationDDom>,
488 memory_space>> derivs_xmax
489 = std::nullopt)
const;
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514 template <
class OutMemorySpace = MemorySpace>
533 static int compute_block_sizes_uniform(
ddc::
BoundCond bound_cond,
int nbc);
535 static int compute_block_sizes_non_uniform(
ddc::
BoundCond bound_cond,
int nbc);
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);
543 void build_matrix_system();
545 void check_valid_grid();
547 template <
class KnotElement>
548 static void check_n_points_in_cell(
int n_points_in_cell, KnotElement current_cell_end_idx);
555 class InterpolationDDim,
559void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
560 compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset)
562 if constexpr (bsplines_type::is_periodic()) {
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();
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)
579 - BSplines::degree();
590 class InterpolationDDim,
594int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
595 compute_block_sizes_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
598 return static_cast<
int>(bsplines_type::degree()) / 2;
607 return static_cast<
int>(bsplines_type::degree()) - 1;
610 throw std::runtime_error(
"ddc::BoundCond not handled");
617 class InterpolationDDim,
621int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
622 compute_block_sizes_non_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
625 return static_cast<
int>(bsplines_type::degree()) - 1;
633 throw std::runtime_error(
"ddc::BoundCond not handled");
640 class InterpolationDDim,
644void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
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)
659 int upper_band_width;
660 if (bsplines_type::is_uniform()) {
661 upper_band_width = bsplines_type::degree() / 2;
663 upper_band_width = bsplines_type::degree() - 1;
665 if constexpr (bsplines_type::is_periodic()) {
666 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_periodic_band_matrix<
668 ddc::discrete_space<BSplines>().nbasis(),
671 bsplines_type::is_uniform());
673 m_matrix =
ddc::detail::SplinesLinearProblemMaker::
674 make_new_block_matrix_with_band_main_block<ExecSpace>(
675 ddc::discrete_space<BSplines>().nbasis(),
678 bsplines_type::is_uniform(),
683 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_sparse<ExecSpace>(
684 ddc::discrete_space<BSplines>().nbasis(),
686 preconditioner_max_block_size);
689 build_matrix_system();
691 m_matrix->setup_solver();
698 class InterpolationDDim,
702void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
703 build_matrix_system()
708 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
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(
716 ddc::discrete_space<BSplines>().rmin(),
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);
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));
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(
744 auto jmin =
ddc::discrete_space<BSplines>().eval_basis(
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(
754 DDC_MDSPAN_ACCESS_OP(values, s));
761 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
767 bsplines_type::degree() + 1,
768 bsplines_type::degree() / 2 + 1>>
const derivs(derivs_ptr.data());
770 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
772 ddc::discrete_space<BSplines>().rmax(),
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);
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(
791 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();
829 assert(batch_domain_type<BatchedInterpolationDDom>(batched_interpolation_domain)
830 == batch_domain_type<BatchedInterpolationDDom>(spline.domain()));
832 if (batch_domain(batched_interpolation_domain).empty()) {
836 assert(vals.
template extent<interpolation_discrete_dimension_type>()
840 assert(
ddc::DiscreteElement<deriv_type>(derivs_xmin->domain().front()).uid() ==
s_odd);
841 assert(derivs_xmin.has_value() ||
s_nbe_xmin == 0);
843 assert(!derivs_xmin.has_value() || derivs_xmin->
template extent<deriv_type>() == 0);
846 assert(
ddc::DiscreteElement<deriv_type>(derivs_xmax->domain().front()).uid() ==
s_odd);
847 assert(derivs_xmax.has_value() ||
s_nbe_xmax == 0);
849 assert(!derivs_xmax.has_value() || derivs_xmax->
template extent<deriv_type>() == 0);
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",
863 batch_domain(batched_interpolation_domain),
865 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
868 spline(
ddc::DiscreteElement<bsplines_type>(i), j)
869 = derivs_xmin_values(
870 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
872 *
ddc::detail::ipow(dx_proxy, i + odd_proxy);
877 ddc::DiscreteElement<bsplines_type>(0),
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);
900 interpolation_discrete_dimension_type>())))]
901 .allocation_kokkos_view(),
902 vals.allocation_kokkos_view());
909 auto const& nbasis_proxy =
ddc::discrete_space<bsplines_type>().nbasis();
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",
918 batch_domain(batched_interpolation_domain),
920 typename batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
923 spline(
ddc::DiscreteElement<bsplines_type>(nbasis_proxy -
s_nbe_xmax + i),
925 = derivs_xmax_values(
926 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
928 *
ddc::detail::ipow(dx_proxy, i + odd_proxy);
933 ddc::DiscreteElement<bsplines_type>(nbasis_proxy -
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);
941 auto const& offset_proxy = m_offset;
943 batched_spline_tr_domain(batched_interpolation_domain),
945 ddc::
ChunkSpan const spline_tr = spline_tr_alloc.span_view();
946 ddc::parallel_for_each(
947 "ddc_splines_transpose_rhs",
949 batch_domain(batched_interpolation_domain),
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);
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());
964 m_matrix->solve(bcoef_section,
false);
966 ddc::parallel_for_each(
967 "ddc_splines_transpose_back_rhs",
969 batch_domain(batched_interpolation_domain),
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);
980 if (bsplines_type::is_periodic()) {
981 ddc::parallel_for_each(
982 "ddc_splines_periodic_rows_duplicate_rhs",
984 batch_domain(batched_interpolation_domain),
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),
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);
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);
1003 spline(i_end, j) = spline(i_start, j);
1017template <
class OutMemorySpace>
1033SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1038 ddc::integrals(ExecSpace(), integral_bsplines.span_view());
1041 ddc::
ChunkSpan const integral_bsplines_without_periodic_additional_bsplines
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(),
1053 Kokkos::View<
double*, Kokkos::LayoutRight, MemorySpace>
const integral_bsplines_mirror
1055 subview(integral_bsplines_mirror_with_additional_allocation,
1057 pair {
static_cast<std::size_t>(0),
1058 integral_bsplines_without_periodic_additional_bsplines
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);
1068 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view(),
1069 integral_bsplines_mirror);
1073 = integral_bsplines_without_periodic_additional_bsplines[
spline_domain().take_first(
1075 ddc::
ChunkSpan const coefficients = integral_bsplines_without_periodic_additional_bsplines
1083 = integral_bsplines_without_periodic_additional_bsplines
1091 auto const dx_proxy = m_dx;
1092 auto const odd_proxy =
s_odd;
1093 ddc::parallel_for_each(
1095 coefficients_derivs_xmin.domain(),
1096 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1097 coefficients_derivs_xmin(i) *=
ddc::detail::
1099 static_cast<std::size_t>(get<bsplines_type>(
1100 (i - coefficients_derivs_xmin.domain().front()) + odd_proxy)));
1102 ddc::parallel_for_each(
1104 coefficients_derivs_xmax.domain(),
1105 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1106 coefficients_derivs_xmax(i) *=
ddc::detail::
1108 static_cast<std::size_t>(get<bsplines_type>(
1109 (i - coefficients_derivs_xmax.domain().front()) + odd_proxy)));
1112 ddc::DiscreteElement<deriv_type>
const first_deriv(
s_odd);
1114 ddc::
Chunk coefficients_derivs_xmin_out(
1121 coefficients.size())),
1123 ddc::
Chunk coefficients_derivs_xmax_out(
1128 coefficients_derivs_xmin_out.allocation_kokkos_view(),
1129 coefficients_derivs_xmin.allocation_kokkos_view());
1131 coefficients_out.allocation_kokkos_view(),
1132 coefficients.allocation_kokkos_view());
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));
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)
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)
1168 class InterpolationDDim,
1172void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1176 std::size_t
const expected_npoints
1178 if (n_interp_points != expected_npoints) {
1179 throw std::runtime_error(
1180 "Incorrect number of points supplied to NonUniformInterpolationPoints. "
1182 + std::to_string(n_interp_points)
1183 +
", expected : " + std::to_string(expected_npoints));
1185 int n_points_in_cell = 0;
1186 auto current_cell_end_idx =
ddc::discrete_space<BSplines>().break_point_domain().front() + 1;
1188 ddc::Coordinate<continuous_dimension_type>
const point =
ddc::coordinate(idx);
1189 if (point >
ddc::coordinate(current_cell_end_idx)) {
1191 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1193 n_points_in_cell = 1;
1195 current_cell_end_idx += 1;
1196 }
else if (point ==
ddc::coordinate(current_cell_end_idx)) {
1198 check_n_points_in_cell(n_points_in_cell + 1, current_cell_end_idx);
1200 n_points_in_cell = 1;
1202 current_cell_end_idx += 1;
1205 n_points_in_cell += 1;
1209 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 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.
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.
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...
If the type DDim is a B-spline, defines type to the discrete dimension of the associated knots.