13#include <Kokkos_Core.hpp>
15#include "integrals.hpp"
16#include "periodic_extrapolation_rule.hpp"
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
41 class EvaluationDDim1,
42 class EvaluationDDim2,
43 class LowerExtrapolationRule1,
44 class UpperExtrapolationRule1,
45 class LowerExtrapolationRule2,
46 class UpperExtrapolationRule2>
51
52
58
59
60 struct eval_deriv_type
66 using continuous_dimension_type1 =
typename BSplines1::continuous_dimension_type;
69 using continuous_dimension_type2 =
typename BSplines2::continuous_dimension_type;
72 using exec_space = ExecSpace;
75 using memory_space = MemorySpace;
78 using evaluation_discrete_dimension_type1 = EvaluationDDim1;
81 using evaluation_discrete_dimension_type2 = EvaluationDDim2;
84 using bsplines_type1 = BSplines1;
87 using bsplines_type2 = BSplines2;
90 using evaluation_domain_type1 =
ddc::
DiscreteDomain<evaluation_discrete_dimension_type1>;
93 using evaluation_domain_type2 =
ddc::
DiscreteDomain<evaluation_discrete_dimension_type2>;
97 evaluation_discrete_dimension_type1,
98 evaluation_discrete_dimension_type2>;
101
102
103
104
106 class BatchedInterpolationDDom,
107 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
108 using batched_evaluation_domain_type = BatchedInterpolationDDom;
120
121
122
123
124
126 class BatchedInterpolationDDom,
127 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
128 using batch_domain_type =
typename ddc::remove_dims_of_t<
129 BatchedInterpolationDDom,
130 evaluation_discrete_dimension_type1,
131 evaluation_discrete_dimension_type2>;
134
135
136
137
138
140 class BatchedInterpolationDDom,
141 class = std::enable_if_t<
ddc::is_discrete_domain_v<BatchedInterpolationDDom>>>
142 using batched_spline_domain_type =
143 typename ddc::detail::convert_type_seq_to_discrete_domain_t<
ddc::type_seq_replace_t<
144 ddc::to_type_seq_t<BatchedInterpolationDDom>,
145 ddc::detail::TypeSeq<
146 evaluation_discrete_dimension_type1,
147 evaluation_discrete_dimension_type2>,
148 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2>>>;
151 using lower_extrapolation_rule_1_type = LowerExtrapolationRule1;
154 using upper_extrapolation_rule_1_type = UpperExtrapolationRule1;
157 using lower_extrapolation_rule_2_type = LowerExtrapolationRule2;
160 using upper_extrapolation_rule_2_type = UpperExtrapolationRule2;
163 LowerExtrapolationRule1 m_lower_extrap_rule_1;
165 UpperExtrapolationRule1 m_upper_extrap_rule_1;
167 LowerExtrapolationRule2 m_lower_extrap_rule_2;
169 UpperExtrapolationRule2 m_upper_extrap_rule_2;
173 std::is_same_v<LowerExtrapolationRule1,
175 == bsplines_type1::is_periodic()
177 UpperExtrapolationRule1,
179 == bsplines_type1::is_periodic()
181 LowerExtrapolationRule2,
183 == bsplines_type2::is_periodic()
185 UpperExtrapolationRule2,
187 == bsplines_type2::is_periodic(),
188 "PeriodicExtrapolationRule has to be used if and only if dimension is periodic");
190 std::is_invocable_r_v<
192 LowerExtrapolationRule1,
193 ddc::Coordinate<continuous_dimension_type1>,
197 Kokkos::layout_right,
199 "LowerExtrapolationRule1::operator() has to be callable "
200 "with usual arguments.");
202 std::is_invocable_r_v<
204 UpperExtrapolationRule1,
205 ddc::Coordinate<continuous_dimension_type1>,
209 Kokkos::layout_right,
211 "UpperExtrapolationRule1::operator() has to be callable "
212 "with usual arguments.");
214 std::is_invocable_r_v<
216 LowerExtrapolationRule2,
217 ddc::Coordinate<continuous_dimension_type2>,
221 Kokkos::layout_right,
223 "LowerExtrapolationRule2::operator() has to be callable "
224 "with usual arguments.");
226 std::is_invocable_r_v<
228 UpperExtrapolationRule2,
229 ddc::Coordinate<continuous_dimension_type2>,
233 Kokkos::layout_right,
235 "UpperExtrapolationRule2::operator() has to be callable "
236 "with usual arguments.");
239
240
241
242
243
244
245
246
247
249 LowerExtrapolationRule1
const& lower_extrap_rule1,
250 UpperExtrapolationRule1
const& upper_extrap_rule1,
251 LowerExtrapolationRule2
const& lower_extrap_rule2,
252 UpperExtrapolationRule2
const& upper_extrap_rule2)
253 : m_lower_extrap_rule_1(lower_extrap_rule1)
254 , m_upper_extrap_rule_1(upper_extrap_rule1)
255 , m_lower_extrap_rule_2(lower_extrap_rule2)
256 , m_upper_extrap_rule_2(upper_extrap_rule2)
261
262
263
264
268
269
270
271
278
279
280
281
282
286
287
288
289
290
294
295
296
297
298
299
300
301
304 return m_lower_extrap_rule_1;
308
309
310
311
312
313
314
315
318 return m_upper_extrap_rule_1;
322
323
324
325
326
327
328
329
332 return m_lower_extrap_rule_2;
336
337
338
339
340
341
342
343
346 return m_upper_extrap_rule_2;
350
351
352
353
354
355
356
357
358
359
360
361 template <
class Layout,
class... CoordsDims>
363 ddc::Coordinate<CoordsDims...>
const& coord_eval,
364 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
367 return eval(coord_eval, spline_coef);
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
394 class BatchedInterpolationDDom,
397 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
400 ddc::Coordinate<CoordsDims...>
const,
401 BatchedInterpolationDDom,
403 memory_space>
const coords_eval,
406 batched_spline_domain_type<BatchedInterpolationDDom>,
408 memory_space>
const spline_coef)
const
410 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
411 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
412 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
413 ddc::parallel_for_each(
414 "ddc_splines_evaluate_2d",
418 typename batch_domain_type<
419 BatchedInterpolationDDom>::discrete_element_type
const j) {
420 const auto spline_eval_2D = spline_eval[j];
421 const auto coords_eval_2D = coords_eval[j];
422 const auto spline_coef_2D = spline_coef[j];
423 for (
auto const i1 : evaluation_domain1) {
424 for (
auto const i2 : evaluation_domain2) {
425 spline_eval_2D(i1, i2) = eval(coords_eval_2D(i1, i2), spline_coef_2D);
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446 template <
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
448 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
452 batched_spline_domain_type<BatchedInterpolationDDom>,
454 memory_space>
const spline_coef)
const
456 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
457 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
458 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
459 ddc::parallel_for_each(
460 "ddc_splines_evaluate_2d",
464 typename batch_domain_type<
465 BatchedInterpolationDDom>::discrete_element_type
const j) {
466 const auto spline_eval_2D = spline_eval[j];
467 const auto spline_coef_2D = spline_coef[j];
468 for (
auto const i1 : evaluation_domain1) {
469 for (
auto const i2 : evaluation_domain2) {
470 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
471 coord_eval_2D(
ddc::coordinate(i1),
ddc::coordinate(i2));
472 spline_eval_2D(i1, i2) = eval(coord_eval_2D(i1, i2), spline_coef_2D);
479
480
481
482
483
484
485
486
487
488
489 template <
class Layout,
class... CoordsDims>
491 ddc::Coordinate<CoordsDims...>
const& coord_eval,
492 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
495 return eval_no_bc<eval_deriv_type, eval_type>(coord_eval, spline_coef);
499
500
501
502
503
504
505
506
507
508
509 template <
class Layout,
class... CoordsDims>
511 ddc::Coordinate<CoordsDims...>
const& coord_eval,
512 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
515 return eval_no_bc<eval_type, eval_deriv_type>(coord_eval, spline_coef);
519
520
521
522
523
524
525
526
527
528
529 template <
class Layout,
class... CoordsDims>
531 ddc::Coordinate<CoordsDims...>
const& coord_eval,
532 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
535 return eval_no_bc<eval_deriv_type, eval_deriv_type>(coord_eval, spline_coef);
539
540
541
542
543
544
545
546
547
548
549
550
551 template <
class InterestDim,
class Layout,
class... CoordsDims>
553 ddc::Coordinate<CoordsDims...>
const& coord_eval,
554 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
558 std::is_same_v<InterestDim, continuous_dimension_type1>
559 || std::is_same_v<InterestDim, continuous_dimension_type2>);
560 if constexpr (std::is_same_v<
562 typename evaluation_discrete_dimension_type1::
563 continuous_dimension_type>) {
564 return deriv_dim_1(coord_eval, spline_coef);
565 }
else if constexpr (std::is_same_v<
567 typename evaluation_discrete_dimension_type2::
568 continuous_dimension_type>) {
569 return deriv_dim_2(coord_eval, spline_coef);
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589 template <
class InterestDim1,
class InterestDim2,
class Layout,
class... CoordsDims>
591 ddc::Coordinate<CoordsDims...>
const& coord_eval,
592 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
598 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
599 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
602 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
603 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
604 return deriv_1_and_2(coord_eval, spline_coef);
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
631 class BatchedInterpolationDDom,
634 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
637 ddc::Coordinate<CoordsDims...>
const,
638 BatchedInterpolationDDom,
640 memory_space>
const coords_eval,
643 batched_spline_domain_type<BatchedInterpolationDDom>,
645 memory_space>
const spline_coef)
const
647 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
648 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
649 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
650 ddc::parallel_for_each(
651 "ddc_splines_differentiate_2d_dim_1",
655 typename batch_domain_type<
656 BatchedInterpolationDDom>::discrete_element_type
const j) {
657 const auto spline_eval_2D = spline_eval[j];
658 const auto coords_eval_2D = coords_eval[j];
659 const auto spline_coef_2D = spline_coef[j];
660 for (
auto const i1 : evaluation_domain1) {
661 for (
auto const i2 : evaluation_domain2) {
662 spline_eval_2D(i1, i2) = eval_no_bc<
664 eval_type>(coords_eval_2D(i1, i2), spline_coef_2D);
671
672
673
674
675
676
677
678
679
680
681
682
683 template <
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
685 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
689 batched_spline_domain_type<BatchedInterpolationDDom>,
691 memory_space>
const spline_coef)
const
693 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
694 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
695 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
696 ddc::parallel_for_each(
697 "ddc_splines_differentiate_2d_dim_1",
701 typename batch_domain_type<
702 BatchedInterpolationDDom>::discrete_element_type
const j) {
703 const auto spline_eval_2D = spline_eval[j];
704 const auto spline_coef_2D = spline_coef[j];
705 for (
auto const i1 : evaluation_domain1) {
706 for (
auto const i2 : evaluation_domain2) {
707 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
708 coord_eval_2D(
ddc::coordinate(i1),
ddc::coordinate(i2));
709 spline_eval_2D(i1, i2) = eval_no_bc<
711 eval_type>(coord_eval_2D, spline_coef_2D);
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
739 class BatchedInterpolationDDom,
742 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
745 ddc::Coordinate<CoordsDims...>
const,
746 BatchedInterpolationDDom,
748 memory_space>
const coords_eval,
751 batched_spline_domain_type<BatchedInterpolationDDom>,
753 memory_space>
const spline_coef)
const
755 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
756 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
757 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
758 ddc::parallel_for_each(
759 "ddc_splines_differentiate_2d_dim_2",
763 typename batch_domain_type<
764 BatchedInterpolationDDom>::discrete_element_type
const j) {
765 const auto spline_eval_2D = spline_eval[j];
766 const auto coords_eval_2D = coords_eval[j];
767 const auto spline_coef_2D = spline_coef[j];
768 for (
auto const i1 : evaluation_domain1) {
769 for (
auto const i2 : evaluation_domain2) {
770 spline_eval_2D(i1, i2) = eval_no_bc<
772 eval_deriv_type>(coords_eval_2D(i1, i2), spline_coef_2D);
779
780
781
782
783
784
785
786
787
788
789
790
791 template <
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
793 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
797 batched_spline_domain_type<BatchedInterpolationDDom>,
799 memory_space>
const spline_coef)
const
801 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
802 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
803 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
804 ddc::parallel_for_each(
805 "ddc_splines_differentiate_2d_dim_2",
809 typename batch_domain_type<
810 BatchedInterpolationDDom>::discrete_element_type
const j) {
811 const auto spline_eval_2D = spline_eval[j];
812 const auto spline_coef_2D = spline_coef[j];
813 for (
auto const i1 : evaluation_domain1) {
814 for (
auto const i2 : evaluation_domain2) {
815 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
816 coord_eval_2D(
ddc::coordinate(i1),
ddc::coordinate(i2));
817 spline_eval_2D(i1, i2) = eval_no_bc<
819 eval_deriv_type>(coord_eval_2D, spline_coef_2D);
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
847 class BatchedInterpolationDDom,
850 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
853 ddc::Coordinate<CoordsDims...>
const,
854 BatchedInterpolationDDom,
856 memory_space>
const coords_eval,
859 batched_spline_domain_type<BatchedInterpolationDDom>,
861 memory_space>
const spline_coef)
const
863 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(coords_eval.domain());
864 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
865 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
866 ddc::parallel_for_each(
867 "ddc_splines_cross_differentiate",
871 typename batch_domain_type<
872 BatchedInterpolationDDom>::discrete_element_type
const j) {
873 const auto spline_eval_2D = spline_eval[j];
874 const auto coords_eval_2D = coords_eval[j];
875 const auto spline_coef_2D = spline_coef[j];
876 for (
auto const i1 : evaluation_domain1) {
877 for (
auto const i2 : evaluation_domain2) {
878 spline_eval_2D(i1, i2) = eval_no_bc<
880 eval_deriv_type>(coords_eval_2D(i1, i2), spline_coef_2D);
887
888
889
890
891
892
893
894
895
896
897
898
899 template <
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
901 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
905 batched_spline_domain_type<BatchedInterpolationDDom>,
907 memory_space>
const spline_coef)
const
909 batch_domain_type<BatchedInterpolationDDom>
const batch_domain(spline_eval.domain());
910 evaluation_domain_type1
const evaluation_domain1(spline_eval.domain());
911 evaluation_domain_type2
const evaluation_domain2(spline_eval.domain());
912 ddc::parallel_for_each(
913 "ddc_splines_cross_differentiate",
917 typename batch_domain_type<
918 BatchedInterpolationDDom>::discrete_element_type
const j) {
919 const auto spline_eval_2D = spline_eval[j];
920 const auto spline_coef_2D = spline_coef[j];
921 for (
auto const i1 : evaluation_domain1) {
922 for (
auto const i2 : evaluation_domain2) {
923 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>
924 coord_eval_2D(
ddc::coordinate(i1),
ddc::coordinate(i2));
925 spline_eval_2D(i1, i2) = eval_no_bc<
927 eval_deriv_type>(coord_eval_2D, spline_coef_2D);
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
957 class BatchedInterpolationDDom,
960 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
963 ddc::Coordinate<CoordsDims...>
const,
964 BatchedInterpolationDDom,
966 memory_space>
const coords_eval,
969 batched_spline_domain_type<BatchedInterpolationDDom>,
971 memory_space>
const spline_coef)
const
976 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
977 || std::is_same_v<InterestDim, continuous_dimension_type2>);
978 if constexpr (std::is_same_v<
980 typename evaluation_discrete_dimension_type1::
981 continuous_dimension_type>) {
982 return deriv_dim_1(spline_eval, coords_eval, spline_coef);
983 }
else if constexpr (std::is_same_v<
985 typename evaluation_discrete_dimension_type2::
986 continuous_dimension_type>) {
987 return deriv_dim_2(spline_eval, coords_eval, spline_coef);
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005 template <
class InterestDim,
class Layout1,
class Layout2,
class BatchedInterpolationDDom>
1007 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
1011 batched_spline_domain_type<BatchedInterpolationDDom>,
1013 memory_space>
const spline_coef)
const
1018 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1019 || std::is_same_v<InterestDim, continuous_dimension_type2>);
1020 if constexpr (std::is_same_v<
1022 typename evaluation_discrete_dimension_type1::
1023 continuous_dimension_type>) {
1024 return deriv_dim_1(spline_eval, spline_coef);
1025 }
else if constexpr (std::is_same_v<
1027 typename evaluation_discrete_dimension_type2::
1028 continuous_dimension_type>) {
1029 return deriv_dim_2(spline_eval, spline_coef);
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1062 class BatchedInterpolationDDom,
1063 class... CoordsDims>
1065 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
1068 ddc::Coordinate<CoordsDims...>
const,
1069 BatchedInterpolationDDom,
1071 memory_space>
const coords_eval,
1074 batched_spline_domain_type<BatchedInterpolationDDom>,
1076 memory_space>
const spline_coef)
const
1081 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1082 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1085 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1086 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
1087 return deriv_1_and_2(spline_eval, coords_eval, spline_coef);
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1113 class BatchedInterpolationDDom>
1115 ddc::
ChunkSpan<
double, BatchedInterpolationDDom, Layout1, memory_space>
const
1119 batched_spline_domain_type<BatchedInterpolationDDom>,
1121 memory_space>
const spline_coef)
const
1126 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1127 && std::is_same_v<InterestDim2, continuous_dimension_type2>)
1130 typename evaluation_discrete_dimension_type1::continuous_dimension_type>
1131 && std::is_same_v<InterestDim1, continuous_dimension_type2>));
1132 return deriv_1_and_2(spline_eval, spline_coef);
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148 template <
class Layout1,
class Layout2,
class BatchedDDom,
class BatchedSplineDDom>
1150 ddc::
ChunkSpan<
double, BatchedDDom, Layout1, memory_space>
const integrals,
1151 ddc::
ChunkSpan<
double const, BatchedSplineDDom, Layout2, memory_space>
const
1155 ddc::type_seq_contains_v<
1156 ddc::detail::TypeSeq<bsplines_type1, bsplines_type2>,
1157 to_type_seq_t<BatchedSplineDDom>>,
1158 "The spline coefficients domain must contain the bsplines dimensions");
1159 using batch_domain_type
1160 =
ddc::remove_dims_of_t<BatchedSplineDDom, bsplines_type1, bsplines_type2>;
1162 std::is_same_v<batch_domain_type, BatchedDDom>,
1163 "The integrals domain must only contain the batch dimensions");
1165 batch_domain_type batch_domain(integrals.domain());
1170 ddc::integrals(exec_space(), values1);
1175 ddc::integrals(exec_space(), values2);
1177 ddc::parallel_for_each(
1178 "ddc_splines_integrate_bsplines",
1181 KOKKOS_LAMBDA(
typename batch_domain_type::discrete_element_type
const j) {
1183 for (
typename spline_domain_type1::discrete_element_type
const i1 :
1185 for (
typename spline_domain_type2::discrete_element_type
const i2 :
1187 integrals(j) += spline_coef(i1, i2, j) * values1(i1) * values2(i2);
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209 template <
class Layout,
class... CoordsDims>
1210 KOKKOS_INLINE_FUNCTION
double eval(
1211 ddc::Coordinate<CoordsDims...> coord_eval,
1212 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
1215 using Dim1 = continuous_dimension_type1;
1216 using Dim2 = continuous_dimension_type2;
1217 if constexpr (bsplines_type1::is_periodic()) {
1218 if (
ddc::get<Dim1>(coord_eval) <
ddc::discrete_space<bsplines_type1>().rmin()
1219 ||
ddc::get<Dim1>(coord_eval) >
ddc::discrete_space<bsplines_type1>().rmax()) {
1220 ddc::get<Dim1>(coord_eval)
1222 (
ddc::get<Dim1>(coord_eval)
1223 -
ddc::discrete_space<bsplines_type1>().rmin())
1224 /
ddc::discrete_space<bsplines_type1>().length())
1225 *
ddc::discrete_space<bsplines_type1>().length();
1228 if constexpr (bsplines_type2::is_periodic()) {
1229 if (
ddc::get<Dim2>(coord_eval) <
ddc::discrete_space<bsplines_type2>().rmin()
1230 ||
ddc::get<Dim2>(coord_eval) >
ddc::discrete_space<bsplines_type2>().rmax()) {
1231 ddc::get<Dim2>(coord_eval)
1233 (
ddc::get<Dim2>(coord_eval)
1234 -
ddc::discrete_space<bsplines_type2>().rmin())
1235 /
ddc::discrete_space<bsplines_type2>().length())
1236 *
ddc::discrete_space<bsplines_type2>().length();
1239 if constexpr (!bsplines_type1::is_periodic()) {
1240 if (
ddc::get<Dim1>(coord_eval) <
ddc::discrete_space<bsplines_type1>().rmin()) {
1241 return m_lower_extrap_rule_1(coord_eval, spline_coef);
1243 if (
ddc::get<Dim1>(coord_eval) >
ddc::discrete_space<bsplines_type1>().rmax()) {
1244 return m_upper_extrap_rule_1(coord_eval, spline_coef);
1247 if constexpr (!bsplines_type2::is_periodic()) {
1248 if (
ddc::get<Dim2>(coord_eval) <
ddc::discrete_space<bsplines_type2>().rmin()) {
1249 return m_lower_extrap_rule_2(coord_eval, spline_coef);
1251 if (
ddc::get<Dim2>(coord_eval) >
ddc::discrete_space<bsplines_type2>().rmax()) {
1252 return m_upper_extrap_rule_2(coord_eval, spline_coef);
1255 return eval_no_bc<eval_type, eval_type>(
1256 ddc::Coordinate<continuous_dimension_type1, continuous_dimension_type2>(
1257 ddc::get<Dim1>(coord_eval),
1258 ddc::get<Dim2>(coord_eval)),
1263
1264
1265
1266
1267
1268
1269
1270 template <
class EvalType1,
class EvalType2,
class Layout,
class... CoordsDims>
1271 KOKKOS_INLINE_FUNCTION
double eval_no_bc(
1272 ddc::Coordinate<CoordsDims...>
const& coord_eval,
1273 ddc::
ChunkSpan<
double const, spline_domain_type, Layout, memory_space>
const
1277 std::is_same_v<EvalType1, eval_type> || std::is_same_v<EvalType1, eval_deriv_type>);
1279 std::is_same_v<EvalType2, eval_type> || std::is_same_v<EvalType2, eval_deriv_type>);
1280 ddc::DiscreteElement<bsplines_type1> jmin1;
1281 ddc::DiscreteElement<bsplines_type2> jmin2;
1283 std::array<
double, bsplines_type1::degree() + 1> vals1_ptr;
1284 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type1::degree() + 1>>
const
1285 vals1(vals1_ptr.data());
1286 std::array<
double, bsplines_type2::degree() + 1> vals2_ptr;
1287 Kokkos::mdspan<
double, Kokkos::extents<std::size_t, bsplines_type2::degree() + 1>>
const
1288 vals2(vals2_ptr.data());
1289 ddc::Coordinate<continuous_dimension_type1>
const coord_eval_interest1(coord_eval);
1290 ddc::Coordinate<continuous_dimension_type2>
const coord_eval_interest2(coord_eval);
1292 if constexpr (std::is_same_v<EvalType1, eval_type>) {
1293 jmin1 =
ddc::discrete_space<bsplines_type1>().eval_basis(vals1, coord_eval_interest1);
1294 }
else if constexpr (std::is_same_v<EvalType1, eval_deriv_type>) {
1295 jmin1 =
ddc::discrete_space<bsplines_type1>().eval_deriv(vals1, coord_eval_interest1);
1297 if constexpr (std::is_same_v<EvalType2, eval_type>) {
1298 jmin2 =
ddc::discrete_space<bsplines_type2>().eval_basis(vals2, coord_eval_interest2);
1299 }
else if constexpr (std::is_same_v<EvalType2, eval_deriv_type>) {
1300 jmin2 =
ddc::discrete_space<bsplines_type2>().eval_deriv(vals2, coord_eval_interest2);
1304 for (std::size_t i = 0; i < bsplines_type1::degree() + 1; ++i) {
1305 for (std::size_t j = 0; j < bsplines_type2::degree() + 1; ++j) {
1307 ddc::DiscreteElement<
1309 bsplines_type2>(jmin1 + i, jmin2 + j))
1310 * vals1[i] * vals2[j];
friend class DiscreteDomain
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.
void deriv(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 spline function (described by its spline coefficients) on a mesh along a specified dime...
SplineEvaluator2D & operator=(SplineEvaluator2D const &x)=default
Copy-assigns.
KOKKOS_FUNCTION double deriv2(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Double-differentiate 2D spline function (described by its spline coefficients) at a given coordinate ...
KOKKOS_FUNCTION double deriv_dim_2(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 s...
void deriv_dim_1(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 first dimensi...
SplineEvaluator2D(SplineEvaluator2D const &x)=default
Copy-constructs.
~SplineEvaluator2D()=default
Destructs.
void deriv2(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
Double-differentiate 2D spline function (described by its spline coefficients) on a mesh along specif...
KOKKOS_FUNCTION double deriv(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 a...
void deriv_dim_1(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 first dimensi...
void deriv2(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
Double-differentiate 2D spline function (described by its spline coefficients) on a mesh along specif...
void deriv_1_and_2(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
Cross-differentiate 2D spline function (described by its spline coefficients) on a mesh along dimensi...
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...
void deriv(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 spline function (described by its spline coefficients) on a mesh along a specified dime...
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 deriv_1_and_2(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
Cross-differentiate 2D spline function (described by its spline coefficients) on a mesh along dimensi...
KOKKOS_FUNCTION double deriv_1_and_2(ddc::Coordinate< CoordsDims... > const &coord_eval, ddc::ChunkSpan< double const, spline_domain_type, Layout, memory_space > const spline_coef) const
Cross-differentiate 2D spline function (described by its spline coefficients) at a given coordinate.
void deriv_dim_2(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 second dimens...
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.
KOKKOS_FUNCTION double deriv_dim_1(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 f...
void deriv_dim_2(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 second dimens...
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.
The top-level namespace of DDC.