11 #include "PanzerDiscFE_config.hpp"
15 #include "Kokkos_ViewFactory.hpp"
18 #include "Intrepid2_Utils.hpp"
19 #include "Intrepid2_FunctionSpaceTools.hpp"
20 #include "Intrepid2_Orientation.hpp"
21 #include "Intrepid2_OrientationTools.hpp"
24 #include "Phalanx_GetNonConstDynRankViewFromConstMDField.hpp"
29 template<
typename Scalar>
31 applyOrientationsImpl(
const int num_cells,
32 Kokkos::DynRankView<Scalar, PHX::Device> view,
33 Kokkos::DynRankView<Intrepid2::Orientation, PHX::Device> orientations,
34 const typename BasisValues2<Scalar>::IntrepidBasis & basis)
36 using ots=Intrepid2::OrientationTools<PHX::Device>;
38 auto sub_orientations = Kokkos::subview(orientations, std::make_pair(0,num_cells));
43 auto sub_view = Kokkos::subview(view, std::make_pair(0,num_cells), Kokkos::ALL(), Kokkos::ALL());
44 auto sub_view_clone =
Kokkos::createDynRankView(view,
"sub_view_clone", sub_view.extent(0), sub_view.extent(1), sub_view.extent(2));
45 Kokkos::deep_copy(sub_view_clone, sub_view);
48 ots::modifyBasisByOrientation(sub_view, sub_view_clone, sub_orientations, &basis);
49 }
else if (view.rank() == 4){
51 auto sub_view = Kokkos::subview(view, std::make_pair(0,num_cells), Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
52 auto sub_view_clone =
Kokkos::createDynRankView(view,
"sub_view_clone", sub_view.extent(0), sub_view.extent(1), sub_view.extent(2), sub_view.extent(3));
53 Kokkos::deep_copy(sub_view_clone, sub_view);
56 ots::modifyBasisByOrientation(sub_view, sub_view_clone, sub_orientations, &basis);
58 throw std::logic_error(
"applyOrientationsImpl : Unknown view of rank " + std::to_string(view.rank()));
61 template<
typename Scalar>
63 applyOrientationsImpl(
const int num_cells,
64 Kokkos::DynRankView<Scalar, PHX::Device> view,
65 const std::vector<Intrepid2::Orientation> & orientations,
66 const typename BasisValues2<Scalar>::IntrepidBasis & basis)
70 Kokkos::DynRankView<Intrepid2::Orientation,PHX::Device> device_orientations(
"drv_orts", num_cells);
71 auto host_orientations = Kokkos::create_mirror_view(device_orientations);
72 for(
int i=0; i < num_cells; ++i)
73 host_orientations(i) = orientations[i];
74 Kokkos::deep_copy(device_orientations,host_orientations);
77 applyOrientationsImpl(num_cells, view, device_orientations, basis);
83 template <
typename Scalar>
86 const bool allocArrays,
87 const bool buildWeighted)
88 : compute_derivatives(true)
89 , build_weighted(buildWeighted)
90 , alloc_arrays(allocArrays)
93 , references_evaluated(false)
94 , orientations_applied_(false)
96 , num_evaluate_cells_(0)
98 , num_orientations_cells_(0)
124 template <
typename Scalar>
130 const int in_num_cells)
132 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::evaluateValues(5 arg)",bv_ev_5);
134 build_weighted =
false;
136 setupUniform(basis_layout, cub_points, jac, jac_det, jac_inv, in_num_cells);
138 const auto elmtspace = getElementSpace();
139 const int num_dims = jac.extent(2);
143 getBasisValuesRef(
true,
true);
146 getVectorBasisValuesRef(
true,
true);
149 getGradBasisValuesRef(
true,
true);
153 getCurl2DVectorBasisRef(
true,
true);
154 else if(num_dims == 3)
155 getCurlVectorBasisRef(
true,
true);
159 getDivVectorBasisRef(
true,
true);
161 references_evaluated =
true;
165 getBasisValues(
false,
true,
true);
168 getVectorBasisValues(
false,
true,
true);
171 getGradBasisValues(
false,
true,
true);
175 getCurl2DVectorBasis(
false,
true,
true);
176 else if(num_dims == 3)
177 getCurlVectorBasis(
false,
true,
true);
181 getDivVectorBasis(
false,
true,
true);
184 orientations_applied_ = (orientations_.size()>0);
187 template <
typename Scalar>
195 bool use_node_coordinates,
196 const int in_num_cells)
198 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::evaluateValues(8 arg, uniform cub pts)",bv_ev_8u);
201 evaluateValues(cub_points, jac, jac_det, jac_inv, in_num_cells);
202 if(weighted_measure.size() > 0)
203 setWeightedMeasure(weighted_measure);
205 cell_node_coordinates_ = node_coordinates;
209 const auto elmtspace = getElementSpace();
210 const int num_dims = jac.extent(2);
213 getBasisValues(
true,
true,
true);
216 getVectorBasisValues(
true,
true,
true);
219 getGradBasisValues(
true,
true,
true);
223 getCurl2DVectorBasis(
true,
true,
true);
224 else if(num_dims == 3)
225 getCurlVectorBasis(
true,
true,
true);
229 getDivVectorBasis(
true,
true,
true);
232 if(use_node_coordinates){
233 getBasisCoordinatesRef(
true,
true);
234 getBasisCoordinates(
true,
true);
238 orientations_applied_ = (orientations_.size()>0);
241 template <
typename Scalar>
249 bool use_node_coordinates,
250 const int in_num_cells)
252 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::evaluateValues(8 arg,nonuniform cub pts)",bv_ev_8nu);
254 cell_node_coordinates_ = node_coordinates;
256 setup(basis_layout, cub_points, jac, jac_det, jac_inv, in_num_cells);
257 if(weighted_measure.size() > 0)
258 setWeightedMeasure(weighted_measure);
260 const auto elmtspace = getElementSpace();
261 const int num_dims = jac.extent(2);
264 getBasisValues(
false,
true,
true);
265 if(build_weighted) getBasisValues(
true,
true,
true);
269 getVectorBasisValues(
false,
true,
true);
270 if(build_weighted) getVectorBasisValues(
true,
true,
true);
274 getGradBasisValues(
false,
true,
true);
275 if(build_weighted) getGradBasisValues(
true,
true,
true);
280 getCurl2DVectorBasis(
false,
true,
true);
281 if(build_weighted) getCurl2DVectorBasis(
true,
true,
true);
282 }
else if(num_dims == 3) {
283 getCurlVectorBasis(
false,
true,
true);
284 if(build_weighted) getCurlVectorBasis(
true,
true,
true);
289 getDivVectorBasis(
false,
true,
true);
290 if(build_weighted) getDivVectorBasis(
true,
true,
true);
294 if(use_node_coordinates){
295 getBasisCoordinatesRef(
true,
true);
296 getBasisCoordinates(
true,
true);
300 orientations_applied_ = (orientations_.size()>0);
303 template <
typename Scalar>
310 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::evaluateValuesCV(5 arg)",bv_ev_cv_5);
314 evaluateValues(cub_points,jac,jac_det,jac_inv,weighted_measure,node_coordinates,
false,jac.extent(0));
317 template <
typename Scalar>
324 bool use_node_coordinates,
325 const int in_num_cells)
327 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::evaluateValuesCV(7 arg)",bv_ev_cv_7);
329 evaluateValues(cell_cub_points,jac,jac_det,jac_inv,weighted_measure,node_coordinates,use_node_coordinates,in_num_cells);
332 template <
typename Scalar>
335 const int in_num_cells)
337 num_evaluate_cells_ = in_num_cells < 0 ? node_coordinates.extent(0) : in_num_cells;
338 cell_node_coordinates_ = node_coordinates;
340 getBasisCoordinates(
true,
true);
344 template <
typename Scalar>
347 const int in_num_cells)
349 if (!intrepid_basis->requireOrientation())
352 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::applyOrientations()",bv_ev_app_orts);
359 const int num_cell_basis_layout = in_num_cells < 0 ? basis_layout->numCells() : in_num_cells;
360 const int num_cell_orientation = orientations.size();
361 const int num_cells = num_cell_basis_layout < num_cell_orientation ? num_cell_basis_layout : num_cell_orientation;
362 const int num_dim = basis_layout->dimension();
365 Kokkos::DynRankView<Intrepid2::Orientation,PHX::Device> device_orientations(
"device_orientations", num_cells);
366 auto host_orientations = Kokkos::create_mirror_view(device_orientations);
367 for(
int i=0; i < num_cells; ++i)
368 host_orientations(i) = orientations[i];
369 Kokkos::deep_copy(device_orientations,host_orientations);
372 applyOrientationsImpl<Scalar>(num_cells, basis_scalar.get_view(), device_orientations, *intrepid_basis);
373 if(build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_basis_scalar.get_view(), device_orientations, *intrepid_basis);
374 if(compute_derivatives) applyOrientationsImpl<Scalar>(num_cells, grad_basis.get_view(), device_orientations, *intrepid_basis);
375 if(compute_derivatives and build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_grad_basis.get_view(), device_orientations, *intrepid_basis);
379 applyOrientationsImpl<Scalar>(num_cells, basis_vector.get_view(), device_orientations, *intrepid_basis);
380 if(build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_basis_vector.get_view(), device_orientations, *intrepid_basis);
382 if(compute_derivatives) applyOrientationsImpl<Scalar>(num_cells, curl_basis_scalar.get_view(), device_orientations, *intrepid_basis);
383 if(compute_derivatives and build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_curl_basis_scalar.get_view(), device_orientations, *intrepid_basis);
386 if(compute_derivatives) applyOrientationsImpl<Scalar>(num_cells, curl_basis_vector.get_view(), device_orientations, *intrepid_basis);
387 if(compute_derivatives and build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_curl_basis_vector.get_view(), device_orientations, *intrepid_basis);
392 applyOrientationsImpl<Scalar>(num_cells, basis_vector.get_view(), device_orientations, *intrepid_basis);
393 if(build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_basis_vector.get_view(), device_orientations, *intrepid_basis);
394 if(compute_derivatives) applyOrientationsImpl<Scalar>(num_cells, div_basis.get_view(), device_orientations, *intrepid_basis);
395 if(compute_derivatives and build_weighted) applyOrientationsImpl<Scalar>(num_cells, weighted_div_basis.get_view(), device_orientations, *intrepid_basis);
398 orientations_applied_ =
true;
402 template <
typename Scalar>
409 template <
typename Scalar>
411 {
return basis_layout->getBasis()->getElementSpace(); }
413 template <
typename Scalar>
416 bool computeDerivatives)
420 compute_derivatives = computeDerivatives;
421 basis_layout = layout;
422 num_cells_ = basis_layout->
numCells();
429 int numcells = basisDesc->
numCells();
433 intrepid_basis = basisDesc->
getIntrepid2Basis<PHX::Device::execution_space,Scalar,Scalar>();
448 weighted_basis_scalar = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"weighted_basis",numcells,card,num_quad);
453 if(compute_derivatives) {
454 grad_basis_ref = af.
buildStaticArray<Scalar,BASIS,IP,
Dim>(
"grad_basis_ref",card,num_quad,dim);
455 grad_basis = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"grad_basis",numcells,card,num_quad,dim);
458 weighted_grad_basis = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"weighted_grad_basis",numcells,card,num_quad,dim);
473 basis_vector = af.
buildStaticArray<Scalar,
Cell,BASIS,IP,Dim>(
"basis",numcells,card,num_quad,dim);
476 weighted_basis_vector = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"weighted_basis",numcells,card,num_quad,dim);
486 if(compute_derivatives) {
489 curl_basis_ref_scalar = af.
buildStaticArray<Scalar,BASIS,IP>(
"curl_basis_ref",card,num_quad);
490 curl_basis_scalar = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"curl_basis",numcells,card,num_quad);
493 weighted_curl_basis_scalar = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"weighted_curl_basis",numcells,card,num_quad);
496 curl_basis_ref_vector = af.
buildStaticArray<Scalar,BASIS,IP,Dim>(
"curl_basis_ref",card,num_quad,dim);
497 curl_basis_vector = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"curl_basis",numcells,card,num_quad,dim);
500 weighted_curl_basis_vector = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"weighted_curl_basis",numcells,card,num_quad,dim);
512 basis_vector = af.
buildStaticArray<Scalar,
Cell,BASIS,IP,Dim>(
"basis",numcells,card,num_quad,dim);
515 weighted_basis_vector = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"weighted_basis",numcells,card,num_quad,dim);
530 if(compute_derivatives) {
531 div_basis_ref = af.
buildStaticArray<Scalar,BASIS,IP>(
"div_basis_ref",card,num_quad);
532 div_basis = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"div_basis",numcells,card,num_quad);
535 weighted_div_basis = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"weighted_div_basis",numcells,card,num_quad);
547 weighted_basis_scalar = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"weighted_basis",numcells,card,num_quad);
567 basis_coordinates = af.
buildStaticArray<Scalar,
Cell,BASIS,Dim>(
"basis_coordinates",numcells,card,dim);
575 template <
typename Scalar>
583 const int num_evaluated_cells)
585 basis_layout = basis;
586 intrepid_basis = basis->
getBasis()->getIntrepid2Basis<PHX::Device::execution_space,Scalar,Scalar>();
588 num_cells_ = basis_layout->numCells();
589 num_evaluate_cells_ = num_evaluated_cells >= 0 ? num_evaluated_cells : num_cells_;
590 build_weighted =
false;
593 cubature_points_ref_ = reference_points;
594 cubature_jacobian_ = point_jacobian;
595 cubature_jacobian_determinant_ = point_jacobian_determinant;
596 cubature_jacobian_inverse_ = point_jacobian_inverse;
602 template <
typename Scalar>
610 const int num_evaluated_cells)
612 basis_layout = basis;
613 intrepid_basis = basis->
getBasis()->getIntrepid2Basis<PHX::Device::execution_space,Scalar,Scalar>();
615 num_cells_ = basis_layout->numCells();
616 num_evaluate_cells_ = num_evaluated_cells >= 0 ? num_evaluated_cells : num_cells_;
617 cubature_points_uniform_ref_ = reference_points;
618 build_weighted =
false;
621 cubature_jacobian_ = point_jacobian;
622 cubature_jacobian_determinant_ = point_jacobian_determinant;
623 cubature_jacobian_inverse_ = point_jacobian_inverse;
629 template <
typename Scalar>
633 const int num_orientations_cells)
635 if(num_orientations_cells < 0)
636 num_orientations_cells_ = num_evaluate_cells_;
638 num_orientations_cells_ = num_orientations_cells;
639 if(orientations.size() == 0){
640 orientations_applied_ =
false;
643 orientations_ = orientations;
644 orientations_applied_ =
true;
650 template <
typename Scalar>
655 cell_node_coordinates_ = node_coordinates;
658 template <
typename Scalar>
661 auto pure_basis = basis_layout->getBasis();
662 return {pure_basis->order(),pure_basis->type()};
665 template <
typename Scalar>
671 basis_ref_scalar_evaluated_ =
false;
672 basis_scalar_evaluated_ =
false;
673 basis_ref_vector_evaluated_ =
false;
674 basis_vector_evaluated_ =
false;
675 grad_basis_ref_evaluated_ =
false;
676 grad_basis_evaluated_ =
false;
677 curl_basis_ref_scalar_evaluated_ =
false;
678 curl_basis_scalar_evaluated_ =
false;
679 curl_basis_ref_vector_evaluated_ =
false;
680 curl_basis_vector_evaluated_ =
false;
681 div_basis_ref_evaluated_ =
false;
682 div_basis_evaluated_ =
false;
683 weighted_basis_scalar_evaluated_ =
false;
684 weighted_basis_vector_evaluated_ =
false;
685 weighted_grad_basis_evaluated_ =
false;
686 weighted_curl_basis_scalar_evaluated_ =
false;
687 weighted_curl_basis_vector_evaluated_ =
false;
688 weighted_div_basis_evaluated_ =
false;
689 basis_coordinates_ref_evaluated_ =
false;
690 basis_coordinates_evaluated_ =
false;
718 template <
typename Scalar>
724 "BasisValues2::setWeightedMeasure : Weighted measure already set. Can only set weighted measure once after setup or setupUniform have beens called.");
725 cubature_weights_ = weighted_measure;
726 build_weighted =
true;
734 #define PANZER_CACHE_DATA(name) \
736 if(name.size()==tmp_##name.size()){ \
737 Kokkos::deep_copy(name.get_view(), tmp_##name.get_view()); \
741 name##_evaluated_ = true; \
744 template <
typename Scalar>
745 typename BasisValues2<Scalar>::ConstArray_BasisDim
748 const bool force)
const
751 if(basis_coordinates_ref_evaluated_ and not force)
752 return basis_coordinates_ref;
754 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getBasisCoordinatesRef()",bv_get_bc_ref);
758 const int num_card = basis_layout->cardinality();
759 const int num_dim = basis_layout->dimension();
761 using coordsScalarType =
typename Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>::scalarType;
762 auto tmp_basis_coordinates_ref = af.
buildStaticArray<coordsScalarType,
BASIS,
Dim>(
"basis_coordinates_ref", num_card, num_dim);
763 intrepid_basis->getDofCoords(tmp_basis_coordinates_ref.get_view());
764 PHX::Device().fence();
769 return tmp_basis_coordinates_ref;
772 template <
typename Scalar>
776 const bool force)
const
779 if(basis_ref_scalar_evaluated_ and not force)
780 return basis_ref_scalar;
782 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getBasisValuesRef()",bv_get_bV_ref);
793 const int num_quad = basis_layout->numPoints();
794 const int num_card = basis_layout->cardinality();
797 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
799 intrepid_basis->getValues(tmp_basis_ref_scalar.get_view(), cubature_points_uniform_ref, Intrepid2::OPERATOR_VALUE);
800 PHX::Device().fence();
805 return tmp_basis_ref_scalar;
808 template <
typename Scalar>
812 const bool force)
const
815 if(basis_ref_vector_evaluated_ and not force)
816 return basis_ref_vector;
818 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getVectorBasisValuesRef()",bv_get_vec_bv_ref);
829 const int num_quad = basis_layout->numPoints();
830 const int num_card = basis_layout->cardinality();
831 const int num_dim = basis_layout->dimension();
834 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
836 intrepid_basis->getValues(tmp_basis_ref_vector.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_VALUE);
837 PHX::Device().fence();
842 return tmp_basis_ref_vector;
845 template <
typename Scalar>
849 const bool force)
const
852 if(grad_basis_ref_evaluated_ and not force)
853 return grad_basis_ref;
855 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getGradBasisValuesRef()",bv_get_grad_bv_ref);
866 const int num_quad = basis_layout->numPoints();
867 const int num_card = basis_layout->cardinality();
868 const int num_dim = basis_layout->dimension();
871 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
873 intrepid_basis->getValues(tmp_grad_basis_ref.get_view(), cubature_points_uniform_ref, Intrepid2::OPERATOR_GRAD);
874 PHX::Device().fence();
879 return tmp_grad_basis_ref;
882 template <
typename Scalar>
886 const bool force)
const
889 if(curl_basis_ref_scalar_evaluated_ and not force)
890 return curl_basis_ref_scalar;
892 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getCurl2DVectorBasisRef()",bv_get_curl2_bv_ref);
904 const int num_quad = basis_layout->numPoints();
905 const int num_card = basis_layout->cardinality();
907 auto tmp_curl_basis_ref_scalar = af.
buildStaticArray<Scalar,
BASIS,
IP>(
"dyn_curl_basis_ref_scalar",num_card,num_quad);
908 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
910 intrepid_basis->getValues(tmp_curl_basis_ref_scalar.get_view(), cubature_points_uniform_ref, Intrepid2::OPERATOR_CURL);
911 PHX::Device().fence();
916 return tmp_curl_basis_ref_scalar;
919 template <
typename Scalar>
923 const bool force)
const
926 if(curl_basis_ref_vector_evaluated_ and not force)
927 return curl_basis_ref_vector;
929 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getCurlVectorBasisRef()",bv_get_curl_vec_bv_ref);
941 const int num_quad = basis_layout->numPoints();
942 const int num_card = basis_layout->cardinality();
943 const int num_dim = basis_layout->dimension();
945 auto tmp_curl_basis_ref_vector = af.
buildStaticArray<Scalar,
BASIS,
IP,
Dim>(
"dyn_curl_basis_ref_vector",num_card,num_quad,num_dim);
946 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
948 intrepid_basis->getValues(tmp_curl_basis_ref_vector.get_view(), cubature_points_uniform_ref, Intrepid2::OPERATOR_CURL);
949 PHX::Device().fence();
954 return tmp_curl_basis_ref_vector;
957 template <
typename Scalar>
961 const bool force)
const
964 if(div_basis_ref_evaluated_ and not force)
965 return div_basis_ref;
967 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getDivVectorBasisRef()",bv_get_div_vec_bv_ref);
978 const int num_quad = basis_layout->numPoints();
979 const int num_card = basis_layout->cardinality();
982 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
984 intrepid_basis->getValues(tmp_div_basis_ref.get_view(), cubature_points_uniform_ref, Intrepid2::OPERATOR_DIV);
985 PHX::Device().fence();
990 return tmp_div_basis_ref;
993 template <
typename Scalar>
997 const bool force)
const
999 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getBasisCoordinates()",bv_get_bc);
1002 if(basis_coordinates_evaluated_ and not force)
1003 return basis_coordinates;
1007 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1008 const auto s_node_coordinates = Kokkos::subview(cell_node_coordinates_.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
1012 const int num_card = basis_layout->cardinality();
1013 const int num_dim = basis_layout->dimension();
1015 auto tmp_basis_coordinates = af.buildStaticArray<Scalar,
Cell,
BASIS,
IP>(
"basis_coordinates", num_cells_, num_card, num_dim);
1016 auto s_aux = Kokkos::subview(tmp_basis_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
1019 auto const_bcr = getBasisCoordinatesRef(
false);
1020 auto bcr = PHX::getNonConstDynRankViewFromConstMDField(const_bcr);
1022 Intrepid2::CellTools<PHX::Device::execution_space> cell_tools;
1023 cell_tools.mapToPhysicalFrame(s_aux, bcr, s_node_coordinates, *cell_topology_);
1024 PHX::Device().fence();
1029 return tmp_basis_coordinates;
1032 template <
typename Scalar>
1037 const bool force)
const
1040 if(weighted_basis_scalar_evaluated_ and not force)
1041 return weighted_basis_scalar;
1043 if(basis_scalar_evaluated_ and not force)
1044 return basis_scalar;
1046 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getBasisValues()",bv_get_bv);
1050 const int num_cells = num_cells_;
1051 const int num_points = basis_layout->numPoints();
1052 const int num_card = basis_layout->cardinality();
1058 const auto bv = getBasisValues(
false, force);
1061 auto tmp_weighted_basis_scalar = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP>(
"weighted_basis_scalar", num_cells, num_card, num_points);
1063 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1064 auto s_aux = Kokkos::subview(tmp_weighted_basis_scalar.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
1065 auto s_cw = Kokkos::subview(cubature_weights_.get_view(),cell_range,Kokkos::ALL());
1066 auto s_bv = Kokkos::subview(bv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
1068 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1069 fst::multiplyMeasure(s_aux,s_cw,s_bv);
1077 return tmp_weighted_basis_scalar;
1081 const auto element_space = getElementSpace();
1091 if(hasUniformReferenceSpace()){
1093 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1096 auto cell_basis_ref_scalar = af.
buildStaticArray<Scalar,BASIS,IP>(
"cell_basis_ref_scalar",num_card,num_points);
1097 intrepid_basis->getValues(cell_basis_ref_scalar.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_VALUE);
1099 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1100 auto s_aux = Kokkos::subview(tmp_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1103 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1105 auto s_cjd = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1106 fst::HVOLtransformVALUE(s_aux,s_cjd,cell_basis_ref_scalar.get_view());
1108 fst::HGRADtransformVALUE(s_aux,cell_basis_ref_scalar.get_view());
1110 PHX::Device().fence();
1123 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1124 #ifdef KOKKOS_ENABLE_CUDA
1125 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_basis_scalar.get_view())::memory_space>::value) {
1127 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_basis_scalar.get_view())::memory_space>::value) {
1129 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1130 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1131 auto tmp_basis_scalar_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_basis_scalar.get_view());
1132 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1134 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1135 auto my_cell_basis_host = Kokkos::subview(tmp_basis_scalar_host,cell,Kokkos::ALL(),Kokkos::ALL());
1136 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1137 intrepid_basis_host->getValues(my_cell_basis_host,my_cell_cub_points_ref_host);
1139 auto tmp_basis_scalar_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_basis_scalar_ref",num_cells,num_card,num_points);
1140 Kokkos::deep_copy(tmp_basis_scalar_ref.get_view(),tmp_basis_scalar_host);
1141 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1142 auto s_aux = Kokkos::subview(tmp_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1143 auto s_ref = Kokkos::subview(tmp_basis_scalar_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1145 using fst=Intrepid2::FunctionSpaceTools<PHX::Device>;
1147 auto s_cjd = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1148 fst::HVOLtransformVALUE(s_aux,s_cjd,s_ref);
1150 fst::HGRADtransformVALUE(s_aux,s_ref);
1154 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1155 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1156 auto tmp_basis_scalar_host = Kokkos::create_mirror_view(tmp_basis_scalar.get_view());
1157 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1159 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1160 auto my_cell_basis_host = Kokkos::subview(tmp_basis_scalar_host,cell,Kokkos::ALL(),Kokkos::ALL());
1161 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1162 intrepid_basis_host->getValues(my_cell_basis_host,my_cell_cub_points_ref_host);
1164 auto tmp_basis_scalar_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_basis_scalar_ref",num_cells,num_card,num_points);
1165 Kokkos::deep_copy(tmp_basis_scalar_ref.get_view(),tmp_basis_scalar_host);
1166 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1167 auto s_aux = Kokkos::subview(tmp_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1168 auto s_ref = Kokkos::subview(tmp_basis_scalar_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1170 using fst=Intrepid2::FunctionSpaceTools<PHX::Device>;
1172 auto s_cjd = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1173 fst::HVOLtransformVALUE(s_aux,s_cjd,s_ref);
1175 fst::HGRADtransformVALUE(s_aux,s_ref);
1177 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1180 PHX::Device().fence();
1187 if(orientations_.size() > 0)
1188 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_basis_scalar.get_view(), orientations_, *intrepid_basis);
1193 return tmp_basis_scalar;
1199 template <
typename Scalar>
1204 const bool force)
const
1207 if(weighted_basis_vector_evaluated_ and not force)
1208 return weighted_basis_vector;
1210 if(basis_vector_evaluated_ and not force)
1211 return basis_vector;
1213 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getVectorBasisValues()",bv_get_vec_bv);
1217 const int num_cells = num_cells_;
1218 const int num_points = basis_layout->numPoints();
1219 const int num_card = basis_layout->cardinality();
1220 const int num_dim = basis_layout->dimension();
1227 const auto bv = getVectorBasisValues(
false, force);
1230 auto tmp_weighted_basis_vector = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP,
Dim>(
"weighted_basis_vector", num_cells, num_card, num_points, num_dim);
1232 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1233 auto s_aux = Kokkos::subview(tmp_weighted_basis_vector.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1234 auto s_cw = Kokkos::subview(cubature_weights_.get_view(),cell_range,Kokkos::ALL());
1235 auto s_bv = Kokkos::subview(bv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1237 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1238 fst::multiplyMeasure(s_aux, s_cw, s_bv);
1243 return tmp_weighted_basis_vector;
1247 const auto element_space = getElementSpace();
1255 TEUCHOS_ASSERT(cubature_jacobian_.size() > 0 && cubature_jacobian_determinant_.size() > 0);
1260 if(hasUniformReferenceSpace()){
1262 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1265 auto cell_basis_ref_vector = af.
buildStaticArray<Scalar,BASIS,IP,Dim>(
"cell_basis_ref_scalar",num_card,num_points,num_dim);
1266 intrepid_basis->getValues(cell_basis_ref_vector.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_VALUE);
1268 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1269 auto s_aux = Kokkos::subview(tmp_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1272 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1274 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1275 fst::HCURLtransformVALUE(s_aux,s_jac_inv,cell_basis_ref_vector.get_view());
1277 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1278 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1279 fst::HDIVtransformVALUE(s_aux,s_jac, s_jac_det, cell_basis_ref_vector.get_view());
1281 PHX::Device().fence();
1295 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1296 #ifdef KOKKOS_ENABLE_CUDA
1297 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_basis_vector.get_view())::memory_space>::value) {
1299 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_basis_vector.get_view())::memory_space>::value) {
1301 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1302 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1303 auto tmp_basis_vector_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_basis_vector.get_view());
1305 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1306 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1307 auto my_cell_basis_host = Kokkos::subview(tmp_basis_vector_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1308 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1309 intrepid_basis_host->getValues(my_cell_basis_host,my_cell_cub_points_ref_host);
1311 auto tmp_basis_vector_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_basis_vector_ref",num_cells,num_card,num_points,num_dim);
1312 Kokkos::deep_copy(tmp_basis_vector_ref.get_view(),tmp_basis_vector_host);
1314 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1315 auto s_aux = Kokkos::subview(tmp_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1316 auto s_ref = Kokkos::subview(tmp_basis_vector_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1318 using fst=Intrepid2::FunctionSpaceTools<PHX::Device>;
1320 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1321 fst::HCURLtransformVALUE(s_aux,s_jac_inv,s_ref);
1323 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1324 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1325 fst::HDIVtransformVALUE(s_aux,s_jac, s_jac_det, s_ref);
1329 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1330 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1331 auto tmp_basis_vector_host = Kokkos::create_mirror_view(tmp_basis_vector.get_view());
1333 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1334 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1335 auto my_cell_basis_host = Kokkos::subview(tmp_basis_vector_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1336 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1337 intrepid_basis_host->getValues(my_cell_basis_host,my_cell_cub_points_ref_host);
1339 auto tmp_basis_vector_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_basis_vector_ref",num_cells,num_card,num_points,num_dim);
1340 Kokkos::deep_copy(tmp_basis_vector_ref.get_view(),tmp_basis_vector_host);
1342 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1343 auto s_aux = Kokkos::subview(tmp_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1344 auto s_ref = Kokkos::subview(tmp_basis_vector_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1346 using fst=Intrepid2::FunctionSpaceTools<PHX::Device>;
1348 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1349 fst::HCURLtransformVALUE(s_aux,s_jac_inv,s_ref);
1351 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1352 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1353 fst::HDIVtransformVALUE(s_aux,s_jac, s_jac_det, s_ref);
1355 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1358 PHX::Device().fence();
1362 if(orientations_.size() > 0)
1363 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_basis_vector.get_view(), orientations_, *intrepid_basis);
1368 return tmp_basis_vector;
1374 template <
typename Scalar>
1379 const bool force)
const
1382 if(weighted_grad_basis_evaluated_ and not force)
1383 return weighted_grad_basis;
1385 if(grad_basis_evaluated_ and not force)
1388 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getGradBasisValues()",bv_get_grad_bv);
1392 const int num_cells = num_cells_;
1393 const int num_points = basis_layout->numPoints();
1394 const int num_card = basis_layout->cardinality();
1395 const int num_dim = basis_layout->dimension();
1402 const auto bv = getGradBasisValues(
false, force);
1405 auto tmp_weighted_grad_basis = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP,
Dim>(
"weighted_grad_basis", num_cells, num_card, num_points, num_dim);
1407 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1408 auto s_aux = Kokkos::subview(tmp_weighted_grad_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1409 auto s_cw = Kokkos::subview(cubature_weights_.get_view(), cell_range, Kokkos::ALL());
1410 auto s_bv = Kokkos::subview(bv.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1412 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1413 fst::multiplyMeasure(s_aux,s_cw,s_bv);
1418 return tmp_weighted_grad_basis;
1424 const auto element_space = getElementSpace();
1428 auto tmp_grad_basis = af.
buildStaticArray<Scalar,
Cell,BASIS,IP,Dim>(
"basis_scalar",num_cells,num_card,num_points,num_dim);
1430 if(hasUniformReferenceSpace()){
1432 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1435 intrepid_basis->getValues(cell_grad_basis_ref.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_GRAD);
1437 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1438 auto s_aux = Kokkos::subview(tmp_grad_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1439 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1442 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1443 fst::HGRADtransformGRAD(s_aux, s_jac_inv,cell_grad_basis_ref.get_view());
1445 PHX::Device().fence();
1459 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1460 #ifdef KOKKOS_ENABLE_CUDA
1461 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_grad_basis.get_view())::memory_space>::value) {
1463 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_grad_basis.get_view())::memory_space>::value) {
1465 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1466 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1467 auto tmp_grad_basis_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_grad_basis.get_view());
1469 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1470 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1471 auto my_cell_grad_basis_host = Kokkos::subview(tmp_grad_basis_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1472 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1473 intrepid_basis_host->getValues(my_cell_grad_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_GRAD);
1475 auto tmp_grad_basis_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_grad_basis_ref",num_cells,num_card,num_points,num_dim);
1476 Kokkos::deep_copy(tmp_grad_basis_ref.get_view(),tmp_grad_basis_host);
1478 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1479 auto s_aux = Kokkos::subview(tmp_grad_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1480 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1481 auto s_ref = Kokkos::subview(tmp_grad_basis_ref.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1484 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1485 fst::HGRADtransformGRAD(s_aux, s_jac_inv, s_ref);
1488 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1489 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1490 auto tmp_grad_basis_host = Kokkos::create_mirror_view(tmp_grad_basis.get_view());
1492 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1493 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1494 auto my_cell_grad_basis_host = Kokkos::subview(tmp_grad_basis_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1495 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1496 intrepid_basis_host->getValues(my_cell_grad_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_GRAD);
1498 auto tmp_grad_basis_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_grad_basis_ref",num_cells,num_card,num_points,num_dim);
1499 Kokkos::deep_copy(tmp_grad_basis_ref.get_view(),tmp_grad_basis_host);
1501 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1502 auto s_aux = Kokkos::subview(tmp_grad_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1503 auto s_jac_inv = Kokkos::subview(cubature_jacobian_inverse_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1504 auto s_ref = Kokkos::subview(tmp_grad_basis_ref.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1507 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1508 fst::HGRADtransformGRAD(s_aux, s_jac_inv, s_ref);
1509 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1512 PHX::Device().fence();
1516 if(orientations_.size() > 0)
1517 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_grad_basis.get_view(), orientations_, *intrepid_basis);
1522 return tmp_grad_basis;
1528 template <
typename Scalar>
1533 const bool force)
const
1536 if(weighted_curl_basis_scalar_evaluated_ and not force)
1537 return weighted_curl_basis_scalar;
1539 if(curl_basis_scalar_evaluated_ and not force)
1540 return curl_basis_scalar;
1542 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getCurl2DVectorBasis()",bv_get_curl2d_vec_bv);
1546 const int num_cells = num_cells_;
1547 const int num_points = basis_layout->numPoints();
1548 const int num_card = basis_layout->cardinality();
1549 const int num_dim = basis_layout->dimension();
1556 const auto bv = getCurl2DVectorBasis(
false, force);
1559 auto tmp_weighted_curl_basis_scalar = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP>(
"weighted_curl_basis_scalar", num_cells, num_card, num_points);
1561 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1562 auto s_aux = Kokkos::subview(tmp_weighted_curl_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1563 auto s_cw = Kokkos::subview(cubature_weights_.get_view(), cell_range, Kokkos::ALL());
1564 auto s_bv = Kokkos::subview(bv.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1566 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1567 fst::multiplyMeasure(s_aux,s_cw,s_bv);
1572 return tmp_weighted_curl_basis_scalar;
1579 const auto element_space = getElementSpace();
1584 if(hasUniformReferenceSpace()){
1586 auto cell_curl_basis_ref_scalar = af.
buildStaticArray<Scalar,BASIS,IP>(
"cell_curl_basis_ref_scalar",num_card,num_points);
1587 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1589 intrepid_basis->getValues(cell_curl_basis_ref_scalar.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_CURL);
1591 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1592 auto s_aux = Kokkos::subview(tmp_curl_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1593 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1598 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1599 fst::HDIVtransformDIV(s_aux,s_jac_det,cell_curl_basis_ref_scalar.get_view());
1600 PHX::Device().fence();
1614 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1615 #ifdef KOKKOS_ENABLE_CUDA
1616 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_curl_basis_scalar.get_view())::memory_space>::value) {
1618 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_curl_basis_scalar.get_view())::memory_space>::value) {
1620 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1621 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1622 auto tmp_curl_basis_scalar_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_curl_basis_scalar.get_view());
1624 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1625 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1626 auto my_cell_curl_basis_host = Kokkos::subview(tmp_curl_basis_scalar_host,cell,Kokkos::ALL(),Kokkos::ALL());
1627 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1628 intrepid_basis_host->getValues(my_cell_curl_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_CURL);
1630 auto tmp_curl_basis_scalar_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_curl_basis_scalar_ref",num_cells,num_card,num_points);
1631 Kokkos::deep_copy(tmp_curl_basis_scalar_ref.get_view(),tmp_curl_basis_scalar_host);
1633 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1634 auto s_aux = Kokkos::subview(tmp_curl_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1635 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1636 auto s_ref = Kokkos::subview(tmp_curl_basis_scalar_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1641 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1642 fst::HDIVtransformDIV(s_aux,s_jac_det,s_ref);
1645 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1646 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1647 auto tmp_curl_basis_scalar_host = Kokkos::create_mirror_view(tmp_curl_basis_scalar.get_view());
1649 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1650 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1651 auto my_cell_curl_basis_host = Kokkos::subview(tmp_curl_basis_scalar_host,cell,Kokkos::ALL(),Kokkos::ALL());
1652 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1653 intrepid_basis_host->getValues(my_cell_curl_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_CURL);
1655 auto tmp_curl_basis_scalar_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_curl_basis_scalar_ref",num_cells,num_card,num_points);
1656 Kokkos::deep_copy(tmp_curl_basis_scalar_ref.get_view(),tmp_curl_basis_scalar_host);
1658 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1659 auto s_aux = Kokkos::subview(tmp_curl_basis_scalar.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1660 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1661 auto s_ref = Kokkos::subview(tmp_curl_basis_scalar_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1666 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1667 fst::HDIVtransformDIV(s_aux,s_jac_det,s_ref);
1668 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1671 PHX::Device().fence();
1674 if(orientations_.size() > 0)
1675 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_curl_basis_scalar.get_view(), orientations_, *intrepid_basis);
1680 return tmp_curl_basis_scalar;
1686 template <
typename Scalar>
1691 const bool force)
const
1694 if(weighted_curl_basis_vector_evaluated_ and not force)
1695 return weighted_curl_basis_vector;
1697 if(curl_basis_vector_evaluated_ and not force)
1698 return curl_basis_vector;
1700 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getCurlVectorBasis()",bv_get_curl_vec_bv);
1704 const int num_cells = num_cells_;
1705 const int num_points = basis_layout->numPoints();
1706 const int num_card = basis_layout->cardinality();
1707 const int num_dim = basis_layout->dimension();
1714 const auto bv = getCurlVectorBasis(
false, force);
1717 auto tmp_weighted_curl_basis_vector = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP,
Dim>(
"weighted_curl_basis_vector", num_cells, num_card, num_points, num_dim);
1719 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1720 auto s_aux = Kokkos::subview(tmp_weighted_curl_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1721 auto s_cw = Kokkos::subview(cubature_weights_.get_view(), cell_range, Kokkos::ALL());
1722 auto s_bv = Kokkos::subview(bv.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1724 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1725 fst::multiplyMeasure(s_aux, s_cw, s_bv);
1730 return tmp_weighted_curl_basis_vector;
1737 const auto element_space = getElementSpace();
1742 if(hasUniformReferenceSpace()){
1744 auto cell_curl_basis_ref_vector = af.
buildStaticArray<Scalar,BASIS,IP,Dim>(
"cell_curl_basis_ref_vector",num_card,num_points,num_dim);
1745 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1747 intrepid_basis->getValues(cell_curl_basis_ref_vector.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_CURL);
1749 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1750 auto s_aux = Kokkos::subview(tmp_curl_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1751 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1752 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1754 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1755 fst::HCURLtransformCURL(s_aux, s_jac, s_jac_det,cell_curl_basis_ref_vector.get_view());
1756 PHX::Device().fence();
1770 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1771 #ifdef KOKKOS_ENABLE_CUDA
1772 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_curl_basis_vector.get_view())::memory_space>::value) {
1774 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_curl_basis_vector.get_view())::memory_space>::value) {
1776 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1777 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1778 auto tmp_curl_basis_vector_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_curl_basis_vector.get_view());
1780 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1781 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1782 auto my_cell_curl_basis_host = Kokkos::subview(tmp_curl_basis_vector_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1783 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1784 intrepid_basis_host->getValues(my_cell_curl_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_CURL);
1786 auto tmp_curl_basis_vector_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_curl_basis_scalar_ref",num_cells,num_card,num_points,num_dim);
1787 Kokkos::deep_copy(tmp_curl_basis_vector_ref.get_view(),tmp_curl_basis_vector_host);
1789 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1790 auto s_aux = Kokkos::subview(tmp_curl_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1791 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1792 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1793 auto s_ref = Kokkos::subview(tmp_curl_basis_vector_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1795 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1796 fst::HCURLtransformCURL(s_aux, s_jac, s_jac_det, s_ref);
1799 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1800 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1801 auto tmp_curl_basis_vector_host = Kokkos::create_mirror_view(tmp_curl_basis_vector.get_view());
1803 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1804 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1805 auto my_cell_curl_basis_host = Kokkos::subview(tmp_curl_basis_vector_host,cell,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
1806 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1807 intrepid_basis_host->getValues(my_cell_curl_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_CURL);
1809 auto tmp_curl_basis_vector_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP,Dim>(
"tmp_curl_basis_scalar_ref",num_cells,num_card,num_points,num_dim);
1810 Kokkos::deep_copy(tmp_curl_basis_vector_ref.get_view(),tmp_curl_basis_vector_host);
1812 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1813 auto s_aux = Kokkos::subview(tmp_curl_basis_vector.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1814 auto s_jac = Kokkos::subview(cubature_jacobian_.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1815 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1816 auto s_ref = Kokkos::subview(tmp_curl_basis_vector_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
1818 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1819 fst::HCURLtransformCURL(s_aux, s_jac, s_jac_det, s_ref);
1820 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1823 PHX::Device().fence();
1827 if(orientations_.size() > 0)
1828 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_curl_basis_vector.get_view(), orientations_, *intrepid_basis);
1833 return tmp_curl_basis_vector;
1839 template <
typename Scalar>
1844 const bool force)
const
1847 if(weighted_div_basis_evaluated_ and not force)
1848 return weighted_div_basis;
1850 if(div_basis_evaluated_ and not force)
1853 PANZER_FUNC_TIME_MONITOR_DIFF(
"panzer::basisValues2::getDevVectorBasis()",bv_get_div_vec_bv);
1857 const int num_cells = num_cells_;
1858 const int num_points = basis_layout->numPoints();
1859 const int num_card = basis_layout->cardinality();
1866 const auto bv = getDivVectorBasis(
false, force);
1869 auto tmp_weighted_div_basis = af.
buildStaticArray<Scalar,
Cell,
BASIS,
IP>(
"weighted_div_basis", num_cells, num_card, num_points);
1871 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1872 auto s_aux = Kokkos::subview(tmp_weighted_div_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1873 auto s_cw = Kokkos::subview(cubature_weights_.get_view(), cell_range, Kokkos::ALL());
1874 auto s_bv = Kokkos::subview(bv.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1876 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1877 fst::multiplyMeasure(s_aux, s_cw, s_bv);
1882 return tmp_weighted_div_basis;
1888 const auto element_space = getElementSpace();
1893 if(hasUniformReferenceSpace()){
1895 auto cell_div_basis_ref = af.
buildStaticArray<Scalar,BASIS,IP>(
"cell_div_basis_ref",num_card,num_points);
1896 auto cubature_points_uniform_ref = PHX::getNonConstDynRankViewFromConstMDField(cubature_points_uniform_ref_);
1898 intrepid_basis->getValues(cell_div_basis_ref.get_view(),cubature_points_uniform_ref,Intrepid2::OPERATOR_DIV);
1900 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1901 auto s_aux = Kokkos::subview(tmp_div_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1902 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1904 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1905 fst::HDIVtransformDIV(s_aux,s_jac_det,cell_div_basis_ref.get_view());
1906 PHX::Device().fence();
1920 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1921 #ifdef KOKKOS_ENABLE_CUDA
1922 if constexpr (std::is_same<Kokkos::CudaUVMSpace,
typename decltype(tmp_div_basis.get_view())::memory_space>::value) {
1924 if constexpr (std::is_same<Kokkos::HIPSpace,
typename decltype(tmp_div_basis.get_view())::memory_space>::value) {
1926 auto cubature_points_ref_host = Kokkos::create_mirror(Kokkos::HostSpace{},cubature_points_ref_.get_view());
1927 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1928 auto tmp_div_basis_host = Kokkos::create_mirror(Kokkos::HostSpace{},tmp_div_basis.get_view());
1930 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1931 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1932 auto my_cell_div_basis_host = Kokkos::subview(tmp_div_basis_host,cell,Kokkos::ALL(),Kokkos::ALL());
1933 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1934 intrepid_basis_host->getValues(my_cell_div_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_DIV);
1936 auto tmp_div_basis_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_div_basis_ref",num_cells,num_card,num_points);
1937 Kokkos::deep_copy(tmp_div_basis_ref.get_view(),tmp_div_basis_host);
1939 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1940 auto s_aux = Kokkos::subview(tmp_div_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1941 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1942 auto s_ref = Kokkos::subview(tmp_div_basis_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1944 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1945 fst::HDIVtransformDIV(s_aux,s_jac_det,s_ref);
1948 auto cubature_points_ref_host = Kokkos::create_mirror_view(cubature_points_ref_.get_view());
1949 Kokkos::deep_copy(cubature_points_ref_host,cubature_points_ref_.get_view());
1950 auto tmp_div_basis_host = Kokkos::create_mirror_view(tmp_div_basis.get_view());
1952 auto intrepid_basis_host = intrepid_basis->getHostBasis();
1953 for(
int cell=0; cell<num_evaluate_cells_; ++cell) {
1954 auto my_cell_div_basis_host = Kokkos::subview(tmp_div_basis_host,cell,Kokkos::ALL(),Kokkos::ALL());
1955 auto my_cell_cub_points_ref_host = Kokkos::subview(cubature_points_ref_host,cell,Kokkos::ALL(),Kokkos::ALL());
1956 intrepid_basis_host->getValues(my_cell_div_basis_host,my_cell_cub_points_ref_host,Intrepid2::OPERATOR_DIV);
1958 auto tmp_div_basis_ref = af.
buildStaticArray<Scalar,Cell,BASIS,IP>(
"tmp_div_basis_ref",num_cells,num_card,num_points);
1959 Kokkos::deep_copy(tmp_div_basis_ref.get_view(),tmp_div_basis_host);
1961 const std::pair<int,int> cell_range(0,num_evaluate_cells_);
1962 auto s_aux = Kokkos::subview(tmp_div_basis.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1963 auto s_jac_det = Kokkos::subview(cubature_jacobian_determinant_.get_view(), cell_range, Kokkos::ALL());
1964 auto s_ref = Kokkos::subview(tmp_div_basis_ref.get_view(), cell_range, Kokkos::ALL(), Kokkos::ALL());
1966 using fst=Intrepid2::FunctionSpaceTools<PHX::Device::execution_space>;
1967 fst::HDIVtransformDIV(s_aux,s_jac_det,s_ref);
1968 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_IMPL_HIP_UNIFIED_MEMORY)
1971 PHX::Device().fence();
1975 if(orientations_.size() > 0)
1976 applyOrientationsImpl<Scalar>(num_orientations_cells_, tmp_div_basis.get_view(), orientations_, *intrepid_basis);
1981 return tmp_div_basis;
void setOrientations(const std::vector< Intrepid2::Orientation > &orientations, const int num_orientations_cells=-1)
Set the orientations object for applying orientations using the lazy evaluation path - required for c...
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
bool grad_basis_evaluated_
ConstArray_BasisIP getCurl2DVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
bool basis_ref_vector_evaluated_
void evaluateValuesCV(const PHX::MDField< Scalar, Cell, IP, Dim > &cell_cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac, const PHX::MDField< Scalar, Cell, IP > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac_inv)
bool basis_ref_scalar_evaluated_
Used to check if arrays have been cached.
bool curl_basis_ref_scalar_evaluated_
EElementSpace getElementSpace() const
ConstArray_BasisIPDim getVectorBasisValuesRef(const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at reference points.
Teuchos::RCP< const CellTopologyInfo > getCellTopologyInfo() const
int cardinality() const
Returns the number of basis coefficients.
bool weighted_curl_basis_scalar_evaluated_
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
ConstArray_BasisIP getBasisValuesRef(const bool cache=true, const bool force=false) const
Get the basis values evaluated at reference points.
panzer::BasisDescriptor getBasisDescriptor() const
Return the basis descriptor.
ConstArray_CellBasisIPDim getCurlVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 3D vector basis evaluated at mesh points.
bool curl_basis_vector_evaluated_
Teuchos::RCP< const PureBasis > getBasis() const
bool weighted_basis_vector_evaluated_
int dimension() const
Returns the dimension of the basis from the topology.
bool weighted_grad_basis_evaluated_
void applyOrientations(const PHX::MDField< const Scalar, Cell, BASIS > &orientations)
Method to apply orientations to a basis values container.
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
bool basis_coordinates_evaluated_
bool weighted_div_basis_evaluated_
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
void setup(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, PHX::MDField< const Scalar, Cell, IP, Dim > reference_points, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian, PHX::MDField< const Scalar, Cell, IP > point_jacobian_determinant, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian_inverse, const int num_evaluated_cells=-1)
Setup for lazy evaluation for non-uniform point layout.
ConstArray_CellBasisDim getBasisCoordinates(const bool cache=true, const bool force=false) const
Carterisan coordinates for basis coefficients in mesh space.
PureBasis::EElementSpace getElementSpace() const
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, double, double > > getIntrepid2Basis() const
ConstArray_BasisIP getDivVectorBasisRef(const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at reference points.
ConstArray_CellBasisIP getDivVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at mesh points.
#define PANZER_CACHE_DATA(name)
void setupArrays(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true)
Sizes/allocates memory for arrays.
int numCells() const
Returns the number of cells in the data layouts.
bool basis_scalar_evaluated_
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
ConstArray_BasisIPDim getCurlVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
BasisValues2(const std::string &prefix="", const bool allocArrays=false, const bool buildWeighted=false)
Main constructor.
ConstArray_BasisDim getBasisCoordinatesRef(const bool cache=true, const bool force=false) const
Get the reference coordinates for basis.
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, const int in_num_cells=-1)
bool basis_coordinates_ref_evaluated_
bool curl_basis_ref_vector_evaluated_
#define TEUCHOS_ASSERT(assertion_test)
void setupUniform(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, PHX::MDField< const Scalar, IP, Dim > reference_points, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian, PHX::MDField< const Scalar, Cell, IP > point_jacobian_determinant, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian_inverse, const int num_evaluated_cells=-1)
Setup for lazy evaluation for uniform point layout.
bool grad_basis_ref_evaluated_
void setCellNodeCoordinates(PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates)
Set the cell node coordinates (required for getBasisCoordinates())
bool weighted_curl_basis_vector_evaluated_
void evaluateValues(const PHX::MDField< Scalar, IP, Dim > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac, const PHX::MDField< Scalar, Cell, IP > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac_inv, const int in_num_cells=-1)
ConstArray_CellBasisIPDim getGradBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at mesh points.
void setWeightedMeasure(PHX::MDField< const Scalar, Cell, IP > weighted_measure)
Set the cubature weights (weighted measure) for the basis values object - required to get weighted ba...
bool curl_basis_scalar_evaluated_
ConstArray_CellBasisIP getCurl2DVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 2D vector basis evaluated at mesh points.
bool basis_vector_evaluated_
ConstArray_BasisIPDim getGradBasisValuesRef(const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at reference points.
bool div_basis_ref_evaluated_
bool div_basis_evaluated_
bool weighted_basis_scalar_evaluated_