13#include <Kokkos_Core.hpp>
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
47 class EvaluationDDim1,
48 class EvaluationDDim2,
49 class EvaluationDDim3,
50 class LowerExtrapolationRule1,
51 class UpperExtrapolationRule1,
52 class LowerExtrapolationRule2,
53 class UpperExtrapolationRule2,
54 class LowerExtrapolationRule3,
55 class UpperExtrapolationRule3>
60 using continuous_dimension_type1 =
typename BSplines1::continuous_dimension_type;
63 using continuous_dimension_type2 =
typename BSplines2::continuous_dimension_type;
66 using continuous_dimension_type3 =
typename BSplines3::continuous_dimension_type;
69 using exec_space = ExecSpace;
72 using memory_space = MemorySpace;
75 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
78 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
81 using evaluation_discrete_dimension_type3 = EvaluationDDim3;
84 using bsplines_type1 = BSplines1;
87 using bsplines_type2 = BSplines2;
90 using bsplines_type3 = BSplines3;
93 using evaluation_domain_type1 =
ddc::
DiscreteDomain<evaluation_discrete_dimension_type1>;
96 using evaluation_domain_type2 =
ddc::
DiscreteDomain<evaluation_discrete_dimension_type2>;
99 using evaluation_domain_type3 =
ddc::
DiscreteDomain<evaluation_discrete_dimension_type3>;
103 evaluation_discrete_dimension_type1,
104 evaluation_discrete_dimension_type2,
105 evaluation_discrete_dimension_type3>;
108
109
110
111
112 template <concepts::discrete_domain BatchedInterpolationDDom>
113 using batched_evaluation_domain_type = BatchedInterpolationDDom;
125 using spline_domain_type =
ddc::
DiscreteDomain<bsplines_type1, bsplines_type2, bsplines_type3>;
128
129
130
131
132
133 template <concepts::discrete_domain BatchedInterpolationDDom>
134 using batch_domain_type =
typename ddc::remove_dims_of_t<
135 BatchedInterpolationDDom,
136 evaluation_discrete_dimension_type1,
137 evaluation_discrete_dimension_type2,
138 evaluation_discrete_dimension_type3>;
141
142
143
144
145
146 template <concepts::discrete_domain BatchedInterpolationDDom>
147 using batched_spline_domain_type =
148 typename ddc::detail::convert_type_seq_to_discrete_domain_t<
ddc::type_seq_replace_t<
149 ddc::to_type_seq_t<BatchedInterpolationDDom>,
150 ddc::detail::TypeSeq<
151 evaluation_discrete_dimension_type1,
152 evaluation_discrete_dimension_type2,
153 evaluation_discrete_dimension_type3>,
154 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>>>;
157 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
160 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
163 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
166 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
169 using lower_extrapolation_rule_3_type = LowerExtrapolationRule3;
172 using upper_extrapolation_rule_3_type = UpperExtrapolationRule3;
175 LowerExtrapolationRule1 m_lower_extrap_rule_1;
177 UpperExtrapolationRule1 m_upper_extrap_rule_1;
179 LowerExtrapolationRule2 m_lower_extrap_rule_2;
181 UpperExtrapolationRule2 m_upper_extrap_rule_2;
183 LowerExtrapolationRule3 m_lower_extrap_rule_3;
185 UpperExtrapolationRule3 m_upper_extrap_rule_3;
189 std::is_same_v<LowerExtrapolationRule1,
191 == bsplines_type1::is_periodic()
193 UpperExtrapolationRule1,
195 == bsplines_type1::is_periodic()
197 LowerExtrapolationRule2,
199 == bsplines_type2::is_periodic()
201 UpperExtrapolationRule2,
203 == bsplines_type2::is_periodic()
205 LowerExtrapolationRule3,
207 == bsplines_type3::is_periodic()
209 UpperExtrapolationRule3,
211 == bsplines_type3::is_periodic(),
212 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
214 std::is_invocable_r_v<
216 LowerExtrapolationRule1,
217 ddc::Coordinate<continuous_dimension_type1>,
221 Kokkos::layout_right,
223 "LowerExtrapolationRule1::operator() has to be callable "
224 "with usual arguments.");
226 std::is_invocable_r_v<
228 UpperExtrapolationRule1,
229 ddc::Coordinate<continuous_dimension_type1>,
233 Kokkos::layout_right,
235 "UpperExtrapolationRule1::operator() has to be callable "
236 "with usual arguments.");
238 std::is_invocable_r_v<
240 LowerExtrapolationRule2,
241 ddc::Coordinate<continuous_dimension_type2>,
245 Kokkos::layout_right,
247 "LowerExtrapolationRule2::operator() has to be callable "
248 "with usual arguments.");
250 std::is_invocable_r_v<
252 UpperExtrapolationRule2,
253 ddc::Coordinate<continuous_dimension_type2>,
257 Kokkos::layout_right,
259 "UpperExtrapolationRule2::operator() has to be callable "
260 "with usual arguments.");
262 std::is_invocable_r_v<
264 LowerExtrapolationRule3,
265 ddc::Coordinate<continuous_dimension_type3>,
269 Kokkos::layout_right,
271 "LowerExtrapolationRule3::operator() has to be callable "
272 "with usual arguments.");
274 std::is_invocable_r_v<
276 UpperExtrapolationRule3,
277 ddc::Coordinate<continuous_dimension_type3>,
281 Kokkos::layout_right,
283 "UpperExtrapolationRule3::operator() has to be callable "
284 "with usual arguments.");
287
288
289
290
291
292
293
294
295
296
297
299 LowerExtrapolationRule1
const& lower_extrap_rule1,
300 UpperExtrapolationRule1
const& upper_extrap_rule1,
301 LowerExtrapolationRule2
const& lower_extrap_rule2,
302 UpperExtrapolationRule2
const& upper_extrap_rule2,
303 LowerExtrapolationRule3
const& lower_extrap_rule3,
304 UpperExtrapolationRule3
const& upper_extrap_rule3)
305 : m_lower_extrap_rule_1(lower_extrap_rule1)
306 , m_upper_extrap_rule_1(upper_extrap_rule1)
307 , m_lower_extrap_rule_2(lower_extrap_rule2)
308 , m_upper_extrap_rule_2(upper_extrap_rule2)
309 , m_lower_extrap_rule_3(lower_extrap_rule3)
310 , m_upper_extrap_rule_3(upper_extrap_rule3)
315
316
317
318
322
323
324
325
332
333
334
335
336
340
341
342
343
344
348
349
350
351
352
353
354
355
358 return m_lower_extrap_rule_1;
362
363
364
365
366
367
368
369
372 return m_upper_extrap_rule_1;
376
377
378
379
380
381
382
383
386 return m_lower_extrap_rule_2;
390
391
392
393
394
395
396
397
400 return m_upper_extrap_rule_2;
404
405
406
407
408
409
410
411
414 return m_lower_extrap_rule_3;
418
419
420
421
422
423
424
425
428 return m_upper_extrap_rule_3;
432
433
434
435
436
437
438
439
440
441
442
443 template <
class Layout,
class... CoordsDims>
445 ddc::Coordinate<CoordsDims...>
const& coord_eval,
446 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
449 return eval(coord_eval, spline_coef);
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
476 class BatchedInterpolationDDom,
479 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
482 ddc::Coordinate<CoordsDims...>
const,
483 BatchedInterpolationDDom,
485 memory_space>
const coords_eval,
488 batched_spline_domain_type<BatchedInterpolationDDom>,
490 memory_space>
const spline_coef)
const
492 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
493 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
494 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
495 evaluation_domain_type3
const evaluation_domain3(spline_eval.domain());
496 ddc::parallel_for_each(
497 "ddc_splines_evaluate_3d",
501 typename batch_domain_type<
502 BatchedInterpolationDDom>::discrete_element_type
const j) {
503 auto const spline_eval_3D = spline_eval[j];
504 auto const coords_eval_3D = coords_eval[j];
505 auto const spline_coef_3D = spline_coef[j];
506 for (
auto const i1 : evaluation_domain1) {
507 for (
auto const i2 : evaluation_domain2) {
508 for (
auto const i3 : evaluation_domain3) {
509 spline_eval_3D(i1, i2, i3)
510 = eval(coords_eval_3D(i1, i2, i3), spline_coef_3D);
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532 template <
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
534 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
538 batched_spline_domain_type<BatchedInterpolationDDom>,
540 memory_space>
const spline_coef)
const
542 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
543 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
544 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
545 evaluation_domain_type3
const evaluation_domain3(spline_eval.domain());
546 ddc::parallel_for_each(
547 "ddc_splines_evaluate_3d",
551 typename batch_domain_type<
552 BatchedInterpolationDDom>::discrete_element_type
const j) {
553 auto const spline_eval_3D = spline_eval[j];
554 auto const spline_coef_3D = spline_coef[j];
555 for (
auto const i1 : evaluation_domain1) {
556 for (
auto const i2 : evaluation_domain2) {
557 for (
auto const i3 : evaluation_domain3) {
559 continuous_dimension_type1,
560 continuous_dimension_type2,
561 continuous_dimension_type3>
565 ddc::coordinate(i3));
566 spline_eval_3D(i1, i2, i3)
567 = eval(coord_eval_3D(i1, i2, i3), spline_coef_3D);
575
576
577
578
579
580
581
582
583
584
585
586
587 template <
class DElem,
class Layout,
class... CoordsDims>
588 KOKKOS_FUNCTION
double deriv(
589 DElem
const& deriv_order,
590 ddc::Coordinate<CoordsDims...>
const& coord_eval,
591 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
594 static_assert(is_discrete_element_v<DElem>);
595 return eval_no_bc(deriv_order, coord_eval, spline_coef);
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
623 class BatchedInterpolationDDom,
626 DElem
const& deriv_order,
627 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
630 ddc::Coordinate<CoordsDims...>
const,
631 BatchedInterpolationDDom,
633 memory_space>
const coords_eval,
636 batched_spline_domain_type<BatchedInterpolationDDom>,
638 memory_space>
const spline_coef)
const
640 static_assert(is_discrete_element_v<DElem>);
642 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
643 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
644 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
645 evaluation_domain_type3
const evaluation_domain3(spline_eval.domain());
646 ddc::parallel_for_each(
647 "ddc_splines_cross_differentiate_3d",
651 typename batch_domain_type<
652 BatchedInterpolationDDom>::discrete_element_type
const j) {
653 auto const spline_eval_3D = spline_eval[j];
654 auto const coords_eval_3D = coords_eval[j];
655 auto const spline_coef_3D = spline_coef[j];
656 for (
auto const i1 : evaluation_domain1) {
657 for (
auto const i2 : evaluation_domain2) {
658 for (
auto const i3 : evaluation_domain3) {
659 spline_eval_3D(i1, i2, i3) = eval_no_bc(
661 coords_eval_3D(i1, i2, i3),
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684 template <
class DElem,
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
686 DElem
const& deriv_order,
687 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
691 batched_spline_domain_type<BatchedInterpolationDDom>,
693 memory_space>
const spline_coef)
const
695 static_assert(is_discrete_element_v<DElem>);
697 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
698 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
699 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
700 evaluation_domain_type3
const evaluation_domain3(spline_eval.domain());
701 ddc::parallel_for_each(
702 "ddc_splines_cross_differentiate_3d",
706 typename batch_domain_type<
707 BatchedInterpolationDDom>::discrete_element_type
const j) {
708 auto const spline_eval_3D = spline_eval[j];
709 auto const spline_coef_3D = spline_coef[j];
710 for (
auto const i1 : evaluation_domain1) {
711 for (
auto const i2 : evaluation_domain2) {
712 for (
auto const i3 : evaluation_domain3) {
714 continuous_dimension_type1,
715 continuous_dimension_type2,
716 continuous_dimension_type3>
720 ddc::coordinate(i3));
721 spline_eval_3D(i1, i2, i3)
722 = eval_no_bc(deriv_order, coord_eval_3D, spline_coef_3D);
730
731
732
733
734
735
736
737
738
739
740
741
742 template <
class Layout1,
class Layout2,
class BatchedDDom,
class BatchedSplineDDom>
744 ddc::
ChunkSpan<
double, BatchedDDom, Layout1, memory_space>
const integrals,
745 ddc::
ChunkSpan<
double const, BatchedSplineDDom, Layout2, memory_space>
const
749 ddc::type_seq_contains_v<
750 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2, bsplines_type3>,
751 to_type_seq_t<BatchedSplineDDom>>,
752 "The spline coefficients domain must contain the bsplines dimensions");
753 using batch_domain_type =
ddc::
754 remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2, bsplines_type3>;
756 std::is_same_v<batch_domain_type, BatchedDDom>,
757 "The integrals domain must only contain the batch dimensions");
759 batch_domain_type batch_domain(integrals.domain());
764 ddc::integrals(exec_space(), values1);
769 ddc::integrals(exec_space(), values2);
774 ddc::integrals(exec_space(), values3);
776 ddc::parallel_for_each(
777 "ddc_splines_integrate_bsplines",
780 KOKKOS_LAMBDA(
typename batch_domain_type::discrete_element_type
const j) {
782 for (
typename spline_domain_type1::discrete_element_type
const i1 :
784 for (
typename spline_domain_type2::discrete_element_type
const i2 :
786 for (
typename spline_domain_type3::discrete_element_type
const i3 :
788 integrals(j) += spline_coef(i1, i2, i3, j) * values1(i1)
789 * values2(i2) * values3(i3);
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812 template <
class Layout,
class... CoordsDims>
813 KOKKOS_INLINE_FUNCTION
double eval(
814 ddc::Coordinate<CoordsDims...> coord_eval,
815 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
818 using Dim1 = continuous_dimension_type1;
819 using Dim2 = continuous_dimension_type2;
820 using Dim3 = continuous_dimension_type3;
821 if constexpr (bsplines_type1::is_periodic()) {
822 if (
ddc::get<Dim1>(coord_eval) <
ddc::discrete_space<bsplines_type1>().rmin()
823 ||
ddc::get<Dim1>(coord_eval) >
ddc::discrete_space<bsplines_type1>().rmax()) {
824 ddc::get<Dim1>(coord_eval)
826 (
ddc::get<Dim1>(coord_eval)
827 -
ddc::discrete_space<bsplines_type1>().rmin())
828 /
ddc::discrete_space<bsplines_type1>().length())
829 *
ddc::discrete_space<bsplines_type1>().length();
832 if constexpr (bsplines_type2::is_periodic()) {
833 if (
ddc::get<Dim2>(coord_eval) <
ddc::discrete_space<bsplines_type2>().rmin()
834 ||
ddc::get<Dim2>(coord_eval) >
ddc::discrete_space<bsplines_type2>().rmax()) {
835 ddc::get<Dim2>(coord_eval)
837 (
ddc::get<Dim2>(coord_eval)
838 -
ddc::discrete_space<bsplines_type2>().rmin())
839 /
ddc::discrete_space<bsplines_type2>().length())
840 *
ddc::discrete_space<bsplines_type2>().length();
843 if constexpr (bsplines_type3::is_periodic()) {
844 if (
ddc::get<Dim3>(coord_eval) <
ddc::discrete_space<bsplines_type3>().rmin()
845 ||
ddc::get<Dim3>(coord_eval) >
ddc::discrete_space<bsplines_type3>().rmax()) {
846 ddc::get<Dim3>(coord_eval)
848 (
ddc::get<Dim3>(coord_eval)
849 -
ddc::discrete_space<bsplines_type3>().rmin())
850 /
ddc::discrete_space<bsplines_type3>().length())
851 *
ddc::discrete_space<bsplines_type3>().length();
854 if constexpr (!bsplines_type1::is_periodic()) {
855 if (
ddc::get<Dim1>(coord_eval) <
ddc::discrete_space<bsplines_type1>().rmin()) {
856 return m_lower_extrap_rule_1(coord_eval, spline_coef);
858 if (
ddc::get<Dim1>(coord_eval) >
ddc::discrete_space<bsplines_type1>().rmax()) {
859 return m_upper_extrap_rule_1(coord_eval, spline_coef);
862 if constexpr (!bsplines_type2::is_periodic()) {
863 if (
ddc::get<Dim2>(coord_eval) <
ddc::discrete_space<bsplines_type2>().rmin()) {
864 return m_lower_extrap_rule_2(coord_eval, spline_coef);
866 if (
ddc::get<Dim2>(coord_eval) >
ddc::discrete_space<bsplines_type2>().rmax()) {
867 return m_upper_extrap_rule_2(coord_eval, spline_coef);
870 if constexpr (!bsplines_type3::is_periodic()) {
871 if (
ddc::get<Dim3>(coord_eval) <
ddc::discrete_space<bsplines_type3>().rmin()) {
872 return m_lower_extrap_rule_3(coord_eval, spline_coef);
874 if (
ddc::get<Dim3>(coord_eval) >
ddc::discrete_space<bsplines_type3>().rmax()) {
875 return m_upper_extrap_rule_3(coord_eval, spline_coef);
881 continuous_dimension_type1,
882 continuous_dimension_type2,
883 continuous_dimension_type3>(
884 ddc::get<Dim1>(coord_eval),
885 ddc::get<Dim2>(coord_eval),
886 ddc::get<Dim3>(coord_eval)),
891
892
893
894
895
896
897
898 template <
class... DerivDims,
class Layout,
class... CoordsDims>
899 KOKKOS_INLINE_FUNCTION
double eval_no_bc(
900 ddc::DiscreteElement<DerivDims...>
const& deriv_order,
901 ddc::Coordinate<CoordsDims...>
const& coord_eval,
902 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
905 using deriv_dim1 =
Deriv<continuous_dimension_type1>;
906 using deriv_dim2 =
Deriv<continuous_dimension_type2>;
907 using deriv_dim3 =
Deriv<continuous_dimension_type3>;
908 using deriv_dims = detail::TypeSeq<DerivDims...>;
912 (in_tags_v<DerivDims,
ddc::detail::TypeSeq<deriv_dim1, deriv_dim2, deriv_dim3>>
914 "The only valid dimensions for deriv_order are Deriv<Dim1>, Deriv<Dim2> and "
917 ddc::DiscreteElement<bsplines_type1> jmin1;
918 ddc::DiscreteElement<bsplines_type2> jmin2;
919 ddc::DiscreteElement<bsplines_type3> jmin3;
921 std::array<
double, bsplines_type1::degree() + 1> vals1_ptr;
922 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>>
const
923 vals1(vals1_ptr.data());
924 std::array<
double, bsplines_type2::degree() + 1> vals2_ptr;
925 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>>
const
926 vals2(vals2_ptr.data());
927 std::array<
double, bsplines_type3::degree() + 1> vals3_ptr;
928 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type3::degree() + 1>>
const
929 vals3(vals3_ptr.data());
930 ddc::Coordinate<continuous_dimension_type1>
const coord_eval_interest1(coord_eval);
931 ddc::Coordinate<continuous_dimension_type2>
const coord_eval_interest2(coord_eval);
932 ddc::Coordinate<continuous_dimension_type3>
const coord_eval_interest3(coord_eval);
934 if constexpr (!in_tags_v<deriv_dim1, deriv_dims>) {
935 jmin1 =
ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
937 auto const order1 = deriv_order.
template uid<deriv_dim1>();
938 KOKKOS_ASSERT(order1 > 0 && order1 <= bsplines_type1::degree())
940 std::array<
double, (bsplines_type1::degree() + 1) * (bsplines_type1::degree() + 1)>
946 bsplines_type1::degree() + 1,
947 Kokkos::dynamic_extent>>
const derivs1(derivs1_ptr.data(), order1 + 1);
949 jmin1 =
ddc::discrete_space<bsplines_type1>()
950 .eval_basis_and_n_derivs(derivs1, coord_eval_interest1, order1);
952 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
953 vals1[i] = DDC_MDSPAN_ACCESS_OP(derivs1, i, order1);
957 if constexpr (!in_tags_v<deriv_dim2, deriv_dims>) {
958 jmin2 =
ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
960 auto const order2 = deriv_order.
template uid<deriv_dim2>();
961 KOKKOS_ASSERT(order2 > 0 && order2 <= bsplines_type2::degree())
963 std::array<
double, (bsplines_type2::degree() + 1) * (bsplines_type2::degree() + 1)>
969 bsplines_type2::degree() + 1,
970 Kokkos::dynamic_extent>>
const derivs2(derivs2_ptr.data(), order2 + 1);
972 jmin2 =
ddc::discrete_space<bsplines_type2>()
973 .eval_basis_and_n_derivs(derivs2, coord_eval_interest2, order2);
975 for (std::size_t i = 0; i < bsplines_type2::degree() + 1; ++i) {
976 vals2[i] = DDC_MDSPAN_ACCESS_OP(derivs2, i, order2);
980 if constexpr (!in_tags_v<deriv_dim3, deriv_dims>) {
981 jmin3 =
ddc::discrete_space<bsplines_type3>().eval_basis(vals3, coord_eval_interest3);
983 auto const order3 = deriv_order.
template uid<deriv_dim3>();
984 KOKKOS_ASSERT(order3 > 0 && order3 <= bsplines_type3::degree())
986 std::array<
double, (bsplines_type3::degree() + 1) * (bsplines_type3::degree() + 1)>
992 bsplines_type3::degree() + 1,
993 Kokkos::dynamic_extent>>
const derivs3(derivs3_ptr.data(), order3 + 1);
995 jmin3 =
ddc::discrete_space<bsplines_type3>()
996 .eval_basis_and_n_derivs(derivs3, coord_eval_interest3, order3);
998 for (std::size_t i = 0; i < bsplines_type3::degree() + 1; ++i) {
999 vals3[i] = DDC_MDSPAN_ACCESS_OP(derivs3, i, order3);
1004 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
1005 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
1006 for (std::size_t k = 0; k < bsplines_type3::degree() + 1; ++k) {
1008 ddc::DiscreteElement<
1011 bsplines_type3>(jmin1 + i, jmin2 + j, jmin3 + k))
1012 * vals1[i] * vals2[j] * vals3[k];
friend class DiscreteDomain
KOKKOS_DEFAULTED_FUNCTION constexpr DiscreteElement()=default
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 2D spline approximation of a function.
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_type1< BatchedInterpolationDDom >, Layout, memory_space > > derivs_min1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1< BatchedInterpolationDDom >, Layout, memory_space > > derivs_max1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2< BatchedInterpolationDDom >, Layout, memory_space > > derivs_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2< BatchedInterpolationDDom >, Layout, memory_space > > derivs_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_max2=std::nullopt) const
Compute a 2D spline approximation of a function.
BatchedInterpolationDDom batched_interpolation_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain representing interpolation points.
ddc::DiscreteDomain< bsplines_type1, bsplines_type2 > spline_domain() const noexcept
Get the 2D domain on which spline coefficients are defined.
SplineBuilder2D(SplineBuilder2D &&x)=default
Move-constructs.
SplineBuilder2D & operator=(SplineBuilder2D &&x)=default
Move-assigns.
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.
SplineBuilder2D(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 SplineBuilder2D acting on the interpolation domain contained in batched_interpolation_domain.
SplineBuilder2D(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 SplineBuilder2D acting on interpolation_domain.
~SplineBuilder2D()=default
Destructs.
SplineBuilder2D(SplineBuilder2D const &x)=delete
Copy-constructor is deleted.
interpolation_domain_type interpolation_domain() const noexcept
Get the domain for the 2D interpolation mesh used by this class.
batch_domain_type< BatchedInterpolationDDom > batch_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the batch domain.
SplineBuilder2D & operator=(SplineBuilder2D const &x)=delete
Copy-assignment is deleted.
A class for creating a 3D spline approximation of a function.
SplineBuilder3D(SplineBuilder3D const &x)=delete
Copy-constructor is deleted.
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_type1< BatchedInterpolationDDom >, Layout, memory_space > > derivs_min1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1< BatchedInterpolationDDom >, Layout, memory_space > > derivs_max1=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2< BatchedInterpolationDDom >, Layout, memory_space > > derivs_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2< BatchedInterpolationDDom >, Layout, memory_space > > derivs_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type3< BatchedInterpolationDDom >, Layout, memory_space > > derivs_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type3< BatchedInterpolationDDom >, Layout, memory_space > > derivs_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_2< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_2< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_min2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_2< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_2< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_max2=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min2_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type2_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max2_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type1_3< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_min2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_min2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_max2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_max2_min3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_min2_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_min2_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_min1_max2_max3=std::nullopt, std::optional< ddc::ChunkSpan< double const, batched_derivs_domain_type< BatchedInterpolationDDom >, Layout, memory_space > > mixed_derivs_max1_max2_max3=std::nullopt) const
Compute a 3D spline approximation of a function.
BatchedInterpolationDDom batched_interpolation_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the whole domain representing interpolation points.
~SplineBuilder3D()=default
Destructs.
SplineBuilder3D(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 SplineBuilder3D acting on the interpolation domain contained in batched_interpolation_domain.
ddc::DiscreteDomain< bsplines_type1, bsplines_type2, bsplines_type3 > spline_domain() const noexcept
Get the 3D domain on which spline coefficients are defined.
interpolation_domain_type interpolation_domain() const noexcept
Get the domain for the 3D interpolation mesh used by this class.
SplineBuilder3D & operator=(SplineBuilder3D const &x)=delete
Copy-assignment is deleted.
SplineBuilder3D(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 SplineBuilder3D acting on interpolation_domain.
SplineBuilder3D(SplineBuilder3D &&x)=default
Move-constructs.
SplineBuilder3D & operator=(SplineBuilder3D &&x)=default
Move-assigns.
batch_domain_type< BatchedInterpolationDDom > batch_domain(BatchedInterpolationDDom const &batched_interpolation_domain) const noexcept
Get the batch domain.
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.
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.
A class to evaluate, differentiate or integrate a 2D spline function.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Evaluate 2D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator2D(SplineEvaluator2D &&x)=default
Move-constructs.
lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
Get the lower extrapolation rule along the first dimension.
SplineEvaluator2D & operator=(SplineEvaluator2D const &x)=default
Copy-assigns.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Differentiate 2D spline function (described by its spline coefficients) on a mesh along the dimension...
SplineEvaluator2D(SplineEvaluator2D const &x)=default
Copy-constructs.
~SplineEvaluator2D()=default
Destructs.
KOKKOS_FUNCTION double deriv(DElem const &deriv_order, ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Differentiate 2D spline function (described by its spline coefficients) at a given coordinate along t...
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Differentiate 2D spline function (described by its spline coefficients) on a mesh along the dimension...
void integrate(ddc::ChunkSpan< double, BatchedDDom, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, BatchedSplineDDom, Layout2, memory_space > const spline_coef) const
Perform batched 2D integrations of a spline function (described by its spline coefficients) along the...
upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
Get the upper extrapolation rule along the second dimension.
upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
Get the upper extrapolation rule along the first dimension.
lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
Get the lower extrapolation rule along the second dimension.
SplineEvaluator2D & operator=(SplineEvaluator2D &&x)=default
Move-assigns.
KOKKOS_FUNCTION double operator()(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Evaluate 2D spline function (described by its spline coefficients) at a given coordinate.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Evaluate 2D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator2D(LowerExtrapolationRule1 const &lower_extrap_rule1, UpperExtrapolationRule1 const &upper_extrap_rule1, LowerExtrapolationRule2 const &lower_extrap_rule2, UpperExtrapolationRule2 const &upper_extrap_rule2)
Build a SplineEvaluator2D acting on batched_spline_domain.
A class to evaluate, differentiate or integrate a 3D spline function.
upper_extrapolation_rule_1_type upper_extrapolation_rule_dim_1() const
Get the upper extrapolation rule along the first dimension.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimension...
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) on a mesh.
SplineEvaluator3D & operator=(SplineEvaluator3D &&x)=default
Move-assigns.
upper_extrapolation_rule_3_type upper_extrapolation_rule_dim_3() const
Get the upper extrapolation rule along the third dimension.
SplineEvaluator3D(LowerExtrapolationRule1 const &lower_extrap_rule1, UpperExtrapolationRule1 const &upper_extrap_rule1, LowerExtrapolationRule2 const &lower_extrap_rule2, UpperExtrapolationRule2 const &upper_extrap_rule2, LowerExtrapolationRule3 const &lower_extrap_rule3, UpperExtrapolationRule3 const &upper_extrap_rule3)
Build a SplineEvaluator3D acting on batched_spline_domain.
~SplineEvaluator3D()=default
Destructs.
lower_extrapolation_rule_1_type lower_extrapolation_rule_dim_1() const
Get the lower extrapolation rule along the first dimension.
upper_extrapolation_rule_2_type upper_extrapolation_rule_dim_2() const
Get the upper extrapolation rule along the second dimension.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Differentiate 3D spline function (described by its spline coefficients) on a mesh along the dimension...
KOKKOS_FUNCTION double deriv(DElem const &deriv_order, ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Differentiate 3D spline function (described by its spline coefficients) at a given coordinate along t...
SplineEvaluator3D(SplineEvaluator3D const &x)=default
Copy-constructs.
KOKKOS_FUNCTION double operator()(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) at a given coordinate.
lower_extrapolation_rule_3_type lower_extrapolation_rule_dim_3() const
Get the lower extrapolation rule along the third dimension.
SplineEvaluator3D(SplineEvaluator3D &&x)=default
Move-constructs.
void integrate(ddc::ChunkSpan< double, BatchedDDom, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, BatchedSplineDDom, Layout2, memory_space > const spline_coef) const
Perform batched 3D integrations of a spline function (described by its spline coefficients) along the...
lower_extrapolation_rule_2_type lower_extrapolation_rule_dim_2() const
Get the lower extrapolation rule along the second dimension.
SplineEvaluator3D & operator=(SplineEvaluator3D const &x)=default
Copy-assigns.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Evaluate 3D spline function (described by its spline coefficients) on a mesh.
A class to evaluate, differentiate or integrate a spline function.
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Evaluate spline function (described by its spline coefficients) on a mesh.
upper_extrapolation_rule_type upper_extrapolation_rule() const
Get the upper extrapolation rule.
SplineEvaluator & operator=(SplineEvaluator const &x)=default
Copy-assigns.
SplineEvaluator & operator=(SplineEvaluator &&x)=default
Move-assigns.
SplineEvaluator(LowerExtrapolationRule const &lower_extrap_rule, UpperExtrapolationRule const &upper_extrap_rule)
Build a SplineEvaluator acting on batched_spline_domain.
lower_extrapolation_rule_type lower_extrapolation_rule() const
Get the lower extrapolation rule.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) on a mesh.
KOKKOS_FUNCTION double operator()(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Evaluate 1D spline function (described by its spline coefficients) at a given coordinate.
KOKKOS_FUNCTION double deriv(DElem const &deriv_order, ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) at a given coordinate.
SplineEvaluator(SplineEvaluator const &x)=default
Copy-constructs.
SplineEvaluator(SplineEvaluator &&x)=default
Move-constructs.
void deriv(DElem const &deriv_order, ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< ddc::Coordinate< CoordsDims... > const, BatchedInterpolationDDom, Layout2, memory_space > const coords_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout3, memory_space > const spline_coef) const
Differentiate 1D spline function (described by its spline coefficients) on a mesh.
void integrate(ddc::ChunkSpan< double, BatchedDDom, Layout1, memory_space > const integrals, ddc::ChunkSpan< double const, BatchedSplineDDom, Layout2, memory_space > const spline_coef) const
Perform batched 1D integrations of a spline function (described by its spline coefficients) along the...
void operator()(ddc::ChunkSpan< double, BatchedInterpolationDDom, Layout1, memory_space > const spline_eval, ddc::ChunkSpan< double const, batched_spline_domain_type< BatchedInterpolationDDom >, Layout2, memory_space > const spline_coef) const
Evaluate a spline function (described by its spline coefficients) on a mesh.
~SplineEvaluator()=default
Destructs.
#define DDC_BUILD_DEPRECATED_CODE
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.