43 #ifndef PANZER_BASIS_VALUES2_HPP
44 #define PANZER_BASIS_VALUES2_HPP
48 #include "Intrepid2_Basis.hpp"
49 #include "Intrepid2_Orientation.hpp"
62 template <
typename Scalar>
68 typedef PHX::MDField<Scalar,BASIS,IP,void,void,void,void,void,void>
Array_BasisIP;
72 typedef PHX::MDField<Scalar,BASIS,Dim,void,void,void,void,void,void>
Array_BasisDim;
75 BasisValues2(
const std::string & pre=
"",
bool allocArrays=
false,
bool buildWeighted=
false)
81 bool computeDerivatives=
true);
83 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
84 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
85 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
86 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
88 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
89 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
90 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
91 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
92 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
93 const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
94 bool use_vertex_coordinates=
true);
96 void evaluateValuesCV(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cell_cub_points,
97 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
98 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
99 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
101 void evaluateValues(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
102 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
103 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
104 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
105 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
106 const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
107 bool use_vertex_coordinates=
true);
112 void applyOrientations(
const PHX::MDField<const Scalar,Cell,BASIS> & orientations);
115 void applyOrientations(
const std::vector<Intrepid2::Orientation> & orientations);
169 std::vector<PHX::index_size_type>
ddims_;
174 void evaluateValues_Const(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
175 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
176 const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
178 void evaluateValues_HGrad(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
179 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
180 const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
182 void evaluateValues_HCurl(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
183 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
184 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
185 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
186 const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
188 void evaluateValues_HDiv(
const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
189 const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> &
jac,
190 const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
191 const PHX::MDField<Scalar,Cell,IP> & weighted_measure);