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 = 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 =
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;
206 void compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset);
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
229 interpolation_domain_type
const& interpolation_domain,
230 std::optional<std::size_t> cols_per_chunk = std::nullopt,
231 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
232 : m_interpolation_domain(interpolation_domain)
233 , m_dx((
ddc::discrete_space<BSplines>().rmax() -
ddc::discrete_space<BSplines>().rmin())
234 /
ddc::discrete_space<BSplines>().ncells())
235 , m_label(std::move(label))
239 "Incompatible boundary conditions");
242 compute_offset(
this->interpolation_domain(), m_offset);
245 int lower_block_size;
246 int upper_block_size;
247 if constexpr (bsplines_type::is_uniform()) {
248 upper_block_size = compute_block_sizes_uniform(BcLower,
s_nbe_xmin);
249 lower_block_size = compute_block_sizes_uniform(BcUpper,
s_nbe_xmax);
251 upper_block_size = compute_block_sizes_non_uniform(BcLower,
s_nbe_xmin);
252 lower_block_size = compute_block_sizes_non_uniform(BcUpper,
s_nbe_xmax);
258 preconditioner_max_block_size);
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
278 interpolation_domain_type
const& interpolation_domain,
279 std::optional<std::size_t> cols_per_chunk = std::nullopt,
280 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
283 interpolation_domain,
285 preconditioner_max_block_size)
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307 template <
concepts::discrete_domain BatchedInterpolationDDom>
310 BatchedInterpolationDDom
const& batched_interpolation_domain,
311 std::optional<std::size_t> cols_per_chunk = std::nullopt,
312 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
315 interpolation_domain_type(batched_interpolation_domain),
317 preconditioner_max_block_size)
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337 template <
concepts::discrete_domain BatchedInterpolationDDom>
339 BatchedInterpolationDDom
const& batched_interpolation_domain,
340 std::optional<std::size_t> cols_per_chunk = std::nullopt,
341 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt)
344 interpolation_domain_type(batched_interpolation_domain),
346 preconditioner_max_block_size)
355
356
357
367
368
369
370
374
375
376
377
378
379
382 return m_interpolation_domain;
386
387
388
389
390
391
392
393
394
395 template <
concepts::discrete_domain BatchedInterpolationDDom>
397 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
400 return batched_interpolation_domain;
404
405
406
407
408
409
410
411
412 template <
class BatchedInterpolationDDom>
413 batch_domain_type<BatchedInterpolationDDom>
batch_domain(
414 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
421
422
423
424
425
426
429 return ddc::discrete_space<bsplines_type>().full_domain();
433
434
435
436
437
438
439
440
441 template <
class BatchedInterpolationDDom>
443 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
446 return ddc::replace_dim_of<
447 interpolation_discrete_dimension_type,
453
454
455
456
457
458
459
460
461 template <
class BatchedInterpolationDDom>
462 batched_spline_tr_domain_type<BatchedInterpolationDDom> batched_spline_tr_domain(
463 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
466 return batched_spline_tr_domain_type<BatchedInterpolationDDom>(
467 ddc::replace_dim_of<bsplines_type, bsplines_type>(
468 batched_spline_domain(batched_interpolation_domain),
470 ddc::DiscreteElement<bsplines_type>(0),
472 m_matrix->required_number_of_rhs_rows()))));
477
478
479
480
481
482
483
484
485 template <
class BatchedInterpolationDDom>
487 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
490 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
491 batched_interpolation_domain,
493 ddc::DiscreteElement<deriv_type>(1),
498
499
500
501
502
503
504
505
506 template <
class BatchedInterpolationDDom>
508 BatchedInterpolationDDom
const& batched_interpolation_domain)
const noexcept
511 return ddc::replace_dim_of<interpolation_discrete_dimension_type, deriv_type>(
512 batched_interpolation_domain,
514 ddc::DiscreteElement<deriv_type>(1),
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537 template <
class Layout,
class BatchedInterpolationDDom>
541 batched_spline_domain_type<BatchedInterpolationDDom>,
543 memory_space> spline,
544 ddc::
ChunkSpan<
double const, BatchedInterpolationDDom, Layout, memory_space> vals,
547 batched_derivs_domain_type<BatchedInterpolationDDom>,
549 memory_space>> derivs_xmin
553 batched_derivs_domain_type<BatchedInterpolationDDom>,
555 memory_space>> derivs_xmax
556 = std::nullopt)
const;
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581 template <
class OutMemorySpace = MemorySpace>
600 static int compute_block_sizes_uniform(
ddc::
BoundCond bound_cond,
int nbc);
602 static int compute_block_sizes_non_uniform(
ddc::
BoundCond bound_cond,
int nbc);
604 void allocate_matrix(
605 int lower_block_size,
606 int upper_block_size,
607 std::optional<std::size_t> cols_per_chunk = std::nullopt,
608 std::optional<
unsigned int> preconditioner_max_block_size = std::nullopt);
610 void build_matrix_system();
612 void check_valid_grid();
614 template <
class KnotElement>
615 static void check_n_points_in_cell(
int n_points_in_cell, KnotElement current_cell_end_idx);
622 class InterpolationDDim,
626void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
627 compute_offset(interpolation_domain_type
const& interpolation_domain,
int& offset)
629 if constexpr (bsplines_type::is_periodic()) {
631 std::array<
double, bsplines_type::degree() + 1> values_ptr;
632 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>>
const
633 values(values_ptr.data());
634 ddc::DiscreteElement<interpolation_discrete_dimension_type> start(
635 interpolation_domain.front());
636 auto jmin =
ddc::discrete_space<BSplines>()
637 .eval_basis(values,
ddc::coordinate(start + BSplines::degree()));
638 if constexpr (bsplines_type::degree() % 2 == 0) {
639 offset = jmin.uid() - start.uid() + bsplines_type::degree() / 2 - BSplines::degree();
641 int const mid = bsplines_type::degree() / 2;
642 offset = jmin.uid() - start.uid()
643 + (DDC_MDSPAN_ACCESS_OP(values, mid) > DDC_MDSPAN_ACCESS_OP(values, mid + 1)
646 - BSplines::degree();
657 class InterpolationDDim,
661int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
662 compute_block_sizes_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
665 return static_cast<
int>(bsplines_type::degree()) / 2;
674 return static_cast<
int>(bsplines_type::degree()) - 1;
677 throw std::runtime_error(
"ddc::BoundCond not handled");
684 class InterpolationDDim,
688int SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
689 compute_block_sizes_non_uniform(
ddc::
BoundCond const bound_cond,
int const nbc)
692 return static_cast<
int>(bsplines_type::degree()) - 1;
700 throw std::runtime_error(
"ddc::BoundCond not handled");
707 class InterpolationDDim,
711void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
713 [[maybe_unused]]
int lower_block_size,
714 [[maybe_unused]]
int upper_block_size,
715 std::optional<std::size_t> cols_per_chunk,
716 std::optional<
unsigned int> preconditioner_max_block_size)
726 int upper_band_width;
727 if (bsplines_type::is_uniform()) {
728 upper_band_width = bsplines_type::degree() / 2;
730 upper_band_width = bsplines_type::degree() - 1;
732 if constexpr (bsplines_type::is_periodic()) {
733 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_periodic_band_matrix<
735 ddc::discrete_space<BSplines>().nbasis(),
738 bsplines_type::is_uniform());
740 m_matrix =
ddc::detail::SplinesLinearProblemMaker::
741 make_new_block_matrix_with_band_main_block<ExecSpace>(
742 ddc::discrete_space<BSplines>().nbasis(),
745 bsplines_type::is_uniform(),
750 m_matrix =
ddc::detail::SplinesLinearProblemMaker::make_new_sparse<ExecSpace>(
751 ddc::discrete_space<BSplines>().nbasis(),
753 preconditioner_max_block_size);
756 build_matrix_system();
758 m_matrix->setup_solver();
765 class InterpolationDDim,
769void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
770 build_matrix_system()
775 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
778 derivs(derivs_ptr.data(),
779 bsplines_type::degree() + 1,
780 bsplines_type::degree() / 2 + 1);
781 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
783 ddc::discrete_space<BSplines>().rmin(),
788 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
789 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
790 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *=
ddc::detail::ipow(m_dx, j);
796 for (std::size_t i = 0; i <
s_nbe_xmin; ++i) {
797 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
798 m_matrix->set_element(i, j, DDC_MDSPAN_ACCESS_OP(derivs, j, i + s_odd));
805 std::array<
double, bsplines_type::degree() + 1> values_ptr;
806 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type::degree() + 1>>
const values(
811 auto jmin =
ddc::discrete_space<BSplines>().eval_basis(
813 ddc::coordinate(
ddc::DiscreteElement<interpolation_discrete_dimension_type>(ix)));
814 for (std::size_t s = 0; s < bsplines_type::degree() + 1; ++s) {
815 int const j =
ddc::detail::
816 modulo(
int(jmin.uid() - m_offset + s),
817 static_cast<
int>(
ddc::discrete_space<BSplines>().nbasis()));
818 m_matrix->set_element(
821 DDC_MDSPAN_ACCESS_OP(values, s));
828 std::array<
double, (bsplines_type::degree() / 2 + 1) * (bsplines_type::degree() + 1)>
834 bsplines_type::degree() + 1,
835 bsplines_type::degree() / 2 + 1>>
const derivs(derivs_ptr.data());
837 ddc::discrete_space<BSplines>().eval_basis_and_n_derivs(
839 ddc::discrete_space<BSplines>().rmax(),
844 for (std::size_t i = 0; i < bsplines_type::degree() + 1; ++i) {
845 for (std::size_t j = 1; j < bsplines_type::degree() / 2 + 1; ++j) {
846 DDC_MDSPAN_ACCESS_OP(derivs, i, j) *=
ddc::detail::ipow(m_dx, j);
851 int const i0 =
ddc::discrete_space<BSplines>().nbasis() -
s_nbe_xmax;
852 int const j0 =
ddc::discrete_space<BSplines>().nbasis() - bsplines_type::degree();
853 for (std::size_t j = 0; j < bsplines_type::degree(); ++j) {
854 for (std::size_t i = 0; i <
s_nbe_xmax; ++i) {
855 m_matrix->set_element(
858 DDC_MDSPAN_ACCESS_OP(derivs, j + 1, i +
s_odd));
873template <
class Layout,
class BatchedInterpolationDDom>
874void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
878 batched_spline_domain_type<BatchedInterpolationDDom>,
880 memory_space> spline,
881 ddc::
ChunkSpan<
double const, BatchedInterpolationDDom, Layout, memory_space> vals,
884 batched_derivs_domain_type<BatchedInterpolationDDom>,
886 memory_space>>
const derivs_xmin,
889 batched_derivs_domain_type<BatchedInterpolationDDom>,
891 memory_space>>
const derivs_xmax)
const
893 auto const batched_interpolation_domain = vals.domain();
896 assert(batch_domain_type<BatchedInterpolationDDom>(batched_interpolation_domain)
897 == batch_domain_type<BatchedInterpolationDDom>(spline.domain()));
899 if (batch_domain(batched_interpolation_domain).empty()) {
903 assert(vals.
template extent<interpolation_discrete_dimension_type>()
907 assert(
ddc::DiscreteElement<deriv_type>(derivs_xmin->domain().front()).uid() ==
s_odd);
908 assert(derivs_xmin.has_value() ||
s_nbe_xmin == 0);
910 assert(!derivs_xmin.has_value() || derivs_xmin->
template extent<deriv_type>() == 0);
913 assert(
ddc::DiscreteElement<deriv_type>(derivs_xmax->domain().front()).uid() ==
s_odd);
914 assert(derivs_xmax.has_value() ||
s_nbe_xmax == 0);
916 assert(!derivs_xmax.has_value() || derivs_xmax->
template extent<deriv_type>() == 0);
923 assert(derivs_xmin->
template extent<deriv_type>() ==
s_nbe_xmin);
924 auto derivs_xmin_values = *derivs_xmin;
925 auto const dx_proxy = m_dx;
926 auto const odd_proxy =
s_odd;
927 ddc::parallel_for_each(
928 "ddc_splines_hermite_compute_lower_coefficients",
930 batch_domain(batched_interpolation_domain),
932 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type j) {
934 spline(
ddc::DiscreteElement<bsplines_type>(i), j)
935 = derivs_xmin_values(
936 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
938 *
ddc::detail::ipow(dx_proxy, i + odd_proxy);
943 ddc::DiscreteElement<bsplines_type>(0),
945 batched_spline_domain_type<BatchedInterpolationDDom>
const
946 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
947 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
966 interpolation_discrete_dimension_type>())))]
967 .allocation_kokkos_view(),
968 vals.allocation_kokkos_view());
975 auto const& nbasis_proxy =
ddc::discrete_space<bsplines_type>().nbasis();
977 assert(derivs_xmax->
template extent<deriv_type>() ==
s_nbe_xmax);
978 auto derivs_xmax_values = *derivs_xmax;
979 auto const dx_proxy = m_dx;
980 auto const odd_proxy =
s_odd;
981 ddc::parallel_for_each(
982 "ddc_splines_hermite_compute_upper_coefficients",
984 batch_domain(batched_interpolation_domain),
986 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type j) {
988 spline(
ddc::DiscreteElement<bsplines_type>(nbasis_proxy -
s_nbe_xmax + i),
990 = derivs_xmax_values(
991 ddc::DiscreteElement<deriv_type>(i + odd_proxy),
993 *
ddc::detail::ipow(dx_proxy, i + odd_proxy);
998 ddc::DiscreteElement<bsplines_type>(nbasis_proxy -
s_nbe_xmax),
1000 batched_spline_domain_type<BatchedInterpolationDDom>
const
1001 dx_spline_domain(dx_splines, batch_domain(batched_interpolation_domain));
1002 ddc::parallel_fill(exec_space(), spline[dx_spline_domain], 0.0);
1006 auto const& offset_proxy = m_offset;
1008 m_label +
" > spline_tr (ddc::SplineBuilder::operator())",
1009 batched_spline_tr_domain(batched_interpolation_domain),
1011 ddc::
ChunkSpan const spline_tr = spline_tr_alloc.span_view();
1012 ddc::parallel_for_each(
1013 m_label +
" > ddc_splines_transpose_rhs",
1015 batch_domain(batched_interpolation_domain),
1017 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
const j) {
1018 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
1019 spline_tr(
ddc::DiscreteElement<bsplines_type>(i), j)
1020 = spline(
ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j);
1024 Kokkos::View<
double**, Kokkos::LayoutRight, exec_space>
const bcoef_section(
1025 spline_tr.data_handle(),
1026 static_cast<std::size_t>(spline_tr.
template extent<bsplines_type>()),
1027 batch_domain(batched_interpolation_domain).size());
1029 m_matrix->solve(bcoef_section,
false);
1031 ddc::parallel_for_each(
1032 m_label +
" > ddc_splines_transpose_back_rhs",
1034 batch_domain(batched_interpolation_domain),
1036 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
const j) {
1037 for (std::size_t i = 0; i < nbasis_proxy; ++i) {
1038 spline(
ddc::DiscreteElement<bsplines_type>(i + offset_proxy), j)
1039 = spline_tr(
ddc::DiscreteElement<bsplines_type>(i), j);
1044 if (bsplines_type::is_periodic()) {
1045 ddc::parallel_for_each(
1046 m_label +
" > ddc_splines_periodic_rows_duplicate_rhs",
1048 batch_domain(batched_interpolation_domain),
1050 batch_domain_type<BatchedInterpolationDDom>::discrete_element_type
const
1052 if (offset_proxy != 0) {
1053 for (
int i = 0; i < offset_proxy; ++i) {
1054 spline(
ddc::DiscreteElement<bsplines_type>(i), j) = spline(
1055 ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i),
1058 for (std::size_t i = offset_proxy; i < bsplines_type::degree(); ++i) {
1059 spline(
ddc::DiscreteElement<bsplines_type>(nbasis_proxy + i), j)
1060 = spline(
ddc::DiscreteElement<bsplines_type>(i), j);
1063 for (std::size_t i(0); i < bsplines_type::degree(); ++i) {
1064 ddc::DiscreteElement<bsplines_type>
const i_start(i);
1065 ddc::DiscreteElement<bsplines_type>
const i_end(nbasis_proxy + i);
1067 spline(i_end, j) = spline(i_start, j);
1081template <
class OutMemorySpace>
1097SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1102 ddc::integrals(ExecSpace(), integral_bsplines.span_view());
1105 ddc::
ChunkSpan const integral_bsplines_without_periodic_additional_bsplines
1110 Kokkos::View<
double**, Kokkos::LayoutRight, MemorySpace>
const
1111 integral_bsplines_mirror_with_additional_allocation(
1112 m_label +
" > integral_bsplines_mirror_with_additional_allocation",
1113 m_matrix->required_number_of_rhs_rows(),
1117 Kokkos::View<
double*, Kokkos::LayoutRight, MemorySpace>
const integral_bsplines_mirror
1119 subview(integral_bsplines_mirror_with_additional_allocation,
1121 pair {
static_cast<std::size_t>(0),
1122 integral_bsplines_without_periodic_additional_bsplines
1128 integral_bsplines_mirror,
1129 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view());
1130 m_matrix->solve(integral_bsplines_mirror_with_additional_allocation,
true);
1132 integral_bsplines_without_periodic_additional_bsplines.allocation_kokkos_view(),
1133 integral_bsplines_mirror);
1137 = integral_bsplines_without_periodic_additional_bsplines[
spline_domain().take_first(
1139 ddc::
ChunkSpan const coefficients = integral_bsplines_without_periodic_additional_bsplines
1147 = integral_bsplines_without_periodic_additional_bsplines
1155 auto const dx_proxy = m_dx;
1156 auto const odd_proxy =
s_odd;
1157 ddc::parallel_for_each(
1159 coefficients_derivs_xmin.domain(),
1160 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1161 coefficients_derivs_xmin(i) *=
ddc::detail::
1163 static_cast<std::size_t>(get<bsplines_type>(
1164 (i - coefficients_derivs_xmin.domain().front()) + odd_proxy)));
1166 ddc::parallel_for_each(
1168 coefficients_derivs_xmax.domain(),
1169 KOKKOS_LAMBDA(
ddc::DiscreteElement<bsplines_type> i) {
1170 coefficients_derivs_xmax(i) *=
ddc::detail::
1172 static_cast<std::size_t>(get<bsplines_type>(
1173 (i - coefficients_derivs_xmax.domain().front()) + odd_proxy)));
1176 ddc::DiscreteElement<deriv_type>
const first_deriv(
s_odd);
1178 ddc::
Chunk coefficients_derivs_xmin_out(
1185 coefficients.size())),
1187 ddc::
Chunk coefficients_derivs_xmax_out(
1192 coefficients_derivs_xmin_out.allocation_kokkos_view(),
1193 coefficients_derivs_xmin.allocation_kokkos_view());
1195 coefficients_out.allocation_kokkos_view(),
1196 coefficients.allocation_kokkos_view());
1198 coefficients_derivs_xmax_out.allocation_kokkos_view(),
1199 coefficients_derivs_xmax.allocation_kokkos_view());
1200 return std::make_tuple(
1201 std::move(coefficients_derivs_xmin_out),
1202 std::move(coefficients_out),
1203 std::move(coefficients_derivs_xmax_out));
1214template <
class KnotElement>
1215void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1216 check_n_points_in_cell(
int const n_points_in_cell, KnotElement
const current_cell_end_idx)
1218 if (n_points_in_cell > BSplines::degree() + 1) {
1219 KnotElement
const rmin_idx =
ddc::discrete_space<BSplines>().break_point_domain().front();
1220 int const failed_cell = (current_cell_end_idx - rmin_idx).value();
1221 throw std::runtime_error(
1222 "The spline problem is overconstrained. There are "
1223 + std::to_string(n_points_in_cell) +
" points in the " + std::to_string(failed_cell)
1232 class InterpolationDDim,
1236void SplineBuilder<ExecSpace, MemorySpace, BSplines, InterpolationDDim, BcLower, BcUpper, Solver>::
1240 std::size_t
const expected_npoints
1242 if (n_interp_points != expected_npoints) {
1243 throw std::runtime_error(
1244 "Incorrect number of points supplied to NonUniformInterpolationPoints. "
1246 + std::to_string(n_interp_points)
1247 +
", expected : " + std::to_string(expected_npoints));
1249 int n_points_in_cell = 0;
1250 auto current_cell_end_idx =
ddc::discrete_space<BSplines>().break_point_domain().front() + 1;
1252 ddc::Coordinate<continuous_dimension_type>
const point =
ddc::coordinate(idx);
1253 if (point >
ddc::coordinate(current_cell_end_idx)) {
1255 check_n_points_in_cell(n_points_in_cell, current_cell_end_idx);
1257 n_points_in_cell = 1;
1259 current_cell_end_idx += 1;
1260 }
else if (point ==
ddc::coordinate(current_cell_end_idx)) {
1262 check_n_points_in_cell(n_points_in_cell + 1, current_cell_end_idx);
1264 n_points_in_cell = 1;
1266 current_cell_end_idx += 1;
1269 n_points_in_cell += 1;
1273 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(std::string label, BatchedInterpolationDDom const &batched_interpolation_domain, std::optional< std::size_t > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditioner_max_block_size=std::nullopt)
Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
SplineBuilder(SplineBuilder const &x)=delete
Copy-constructor is deleted.
SplineBuilder(BatchedInterpolationDDom const &batched_interpolation_domain, std::optional< std::size_t > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditioner_max_block_size=std::nullopt)
Build a SplineBuilder acting on the interpolation domain contained by batched_interpolation_domain.
static constexpr int s_nbv_xmax
The number of input values defining the boundary condition at the upper bound.
static constexpr bool s_odd
Indicates if the degree of the splines is odd or even.
static constexpr int s_nbv_xmin
The number of input values defining the boundary condition at the lower bound.
SplineBuilder(std::string label, interpolation_domain_type const &interpolation_domain, std::optional< std::size_t > cols_per_chunk=std::nullopt, std::optional< unsigned int > preconditioner_max_block_size=std::nullopt)
Build a SplineBuilder acting on interpolation_domain.
interpolation_domain_type interpolation_domain() const noexcept
Get the domain for the 1D interpolation mesh used by this class.
batched_derivs_domain_type< BatchedInterpolationDDom > batched_derivs_xmin_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain on which derivatives on lower boundary are defined.
static constexpr ddc::BoundCond s_bc_xmax
The boundary condition implemented at the upper bound.
void operator()(ddc::ChunkSpan< double, batched_spline_domain_type< BatchedInterpolationDDom >, Layout, memory_space > spline, ddc::ChunkSpan< double const, BatchedInterpolationDDom, Layout, memory_space > vals, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > derivs_xmin=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > derivs_xmax=std::nullopt) const
Compute a spline approximation of a function.
batched_spline_domain_type< BatchedInterpolationDDom > batched_spline_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain on which spline coefficients are defined.
ddc::DiscreteDomain< bsplines_type > spline_domain() const noexcept
Get the 1D domain on which spline coefficients are defined.
~SplineBuilder()=default
Destructs.
static constexpr int s_nbe_xmin
The number of equations defining the boundary condition at the lower bound.
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.