11 #ifndef PANZER_INTEGRATION_VALUES2_HPP
12 #define PANZER_INTEGRATION_VALUES2_HPP
16 #include "PanzerDiscFE_config.hpp"
20 #include "Phalanx_MDField.hpp"
21 #include "Intrepid2_Cubature.hpp"
25 class SubcellConnectivity;
27 template <
typename Scalar>
58 const bool allocArrays=
false);
79 const int num_cells = -1,
81 const int num_virtual_cells = -1);
99 const int num_cells = -1);
160 const int num_cells = -1);
174 const int num_virtual_cells);
201 const bool force =
false,
202 const bool apply_permutation =
true)
const;
218 const bool force =
false,
219 const bool apply_permutation =
true)
const;
235 const bool force =
false,
236 const bool apply_permutation =
true)
const;
259 const bool force =
false)
const;
274 const bool force =
false)
const;
288 const bool force =
false)
const;
302 const bool force =
false)
const;
316 const bool force =
false)
const;
330 const bool force =
false)
const;
344 const bool force =
false)
const;
360 const bool force =
false)
const;
376 const bool force =
false)
const;
392 const bool force =
false)
const;
406 const bool force =
false)
const;
420 const bool force =
false)
const;
484 std::vector<PHX::index_size_type>
ddims_;
bool requires_permutation_
ConstArray_IPDim getUniformSideCubaturePointsRef(const bool cache=true, const bool force=false, const bool apply_permutation=true) const
Get the uniform cubature points for a side.
void setupPermutations(const Teuchos::RCP< const SubcellConnectivity > &face_connectivity, const int num_virtual_cells)
Initialize the permutation arrays given a face connectivity.
Array_CellIPDimDim covarient
bool weighted_measure_evaluated_
bool ref_ip_coordinates_evaluated_
ConstArray_CellIPDimDim getContravarientMatrix(const bool cache=true, const bool force=false) const
Get the contravarient matrix.
PHX::MDField< Scalar, Cell, IP > Array_CellIP
PHX::MDField< Scalar, Point > Array_Point
void setup(const Teuchos::RCP< const panzer::IntegrationRule > &ir, const PHX::MDField< Scalar, Cell, NODE, Dim > &cell_node_coordinates, const int num_cells=-1)
Main setup call for the lazy evaluation interface.
Array_CellIPDimDim jac_inv
ConstArray_CellIP getWeightedMeasure(const bool cache=true, const bool force=false) const
Get the weighted measure (integration weights)
PHX::MDField< const int, Cell, IP > permutations_
bool covarient_evaluated_
PHX::MDField< const Scalar, Point > ConstArray_Point
Array_CellIPDim ref_ip_coordinates
PHX::MDField< Scalar, IP > Array_IP
PHX::MDField< const Scalar, IP > ConstArray_IP
Array_CellBASISDim node_coordinates
PHX::MDField< Scalar > ArrayDynamic
PHX::MDField< const Scalar, Cell, IP > ConstArray_CellIP
Array_CellIPDimDim surface_rotation_matrices
void evaluateEverything()
bool surface_normals_evaluated_
bool contravarient_evaluated_
bool cub_points_evaluated_
ConstArray_CellIPDimDim getJacobianInverse(const bool cache=true, const bool force=false) const
Get the inverse of the Jacobian matrix evaluated at the cubature points.
PHX::MDField< Scalar, Cell, IP, Dim, Dim > Array_CellIPDimDim
Array_CellIPDim weighted_normals
Teuchos::RCP< Intrepid2::Cubature< PHX::Device::execution_space, double, double > > intrepid_cubature
ConstArray_CellIPDimDim getCovarientMatrix(const bool cache=true, const bool force=false) const
Get the covarient matrix.
PHX::MDField< const Scalar, IP, Dim > ConstArray_IPDim
bool node_coordinates_evaluated_
ConstArray_CellIPDim getWeightedNormals(const bool cache=true, const bool force=false) const
Get the weighted normals.
bool cub_weights_evaluated_
ConstArray_CellIP getNormContravarientMatrix(const bool cache=true, const bool force=false) const
Get the contravarient matrix.
Array_CellIPDimDim contravarient
Teuchos::RCP< const panzer::IntegrationRule > int_rule
ConstArray_CellIPDimDim getSurfaceRotationMatrices(const bool cache=true, const bool force=false) const
Get the surface rotation matrices.
Array_CellIPDim ip_coordinates
PHX::MDField< Scalar, Cell, BASIS, Dim > Array_CellBASISDim
Array_IPDim side_cub_points
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
bool surface_rotation_matrices_evaluated_
ConstArray_CellIPDim getCubaturePoints(const bool cache=true, const bool force=false) const
Get the cubature points in physical space.
void setupArrays(const Teuchos::RCP< const panzer::IntegrationRule > &ir)
Sizes/allocates memory for arrays.
ConstArray_CellIPDim getCubaturePointsRef(const bool cache=true, const bool force=false) const
Get the cubature points in the reference space.
Array_CellIP weighted_measure
bool side_cub_points_evaluated_
bool weighted_normals_evaluated_
PHX::MDField< Scalar, Cell, IP, Dim > Array_CellIPDim
Array_CellIP norm_contravarient
ConstArray_IPDim getUniformCubaturePointsRef(const bool cache=true, const bool force=false, const bool apply_permutation=true) const
Get the uniform cubature points.
ConstArray_CellIPDimDim getJacobian(const bool cache=true, const bool force=false) const
Get the Jacobian matrix evaluated at the cubature points.
PHX::MDField< double > DblArrayDynamic
Teuchos::RCP< const panzer::IntegrationRule > getIntegrationRule() const
Returns the IntegrationRule.
ConstArray_CellIP getJacobianDeterminant(const bool cache=true, const bool force=false) const
Get the determinant of the Jacobian matrix evaluated at the cubature points.
PHX::MDField< Scalar, IP, Dim > Array_IPDim
Array_CellIPDim surface_normals
IntegrationValues2(const std::string &pre="", const bool allocArrays=false)
Base constructor.
PHX::MDField< const Scalar, Cell, IP, Dim > ConstArray_CellIPDim
ConstArray_CellBASISDim getNodeCoordinates() const
Get the node coordinates describing the geometry of the mesh.
std::vector< PHX::index_size_type > ddims_
bool ip_coordinates_evaluated_
void evaluateValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &cell_node_coordinates, const int num_cells=-1, const Teuchos::RCP< const SubcellConnectivity > &face_connectivity=Teuchos::null, const int num_virtual_cells=-1)
Evaluate basis values.
ConstArray_IP getUniformCubatureWeightsRef(const bool cache=true, const bool force=false, const bool apply_permutation=true) const
Get the uniform cubature weights.
ConstArray_CellIPDim getSurfaceNormals(const bool cache=true, const bool force=false) const
Get the surface normals.
PHX::MDField< const Scalar, Cell, IP, Dim, Dim > ConstArray_CellIPDimDim
Teuchos::RCP< const SubcellConnectivity > side_connectivity_
bool norm_contravarient_evaluated_
PHX::MDField< const Scalar, Cell, BASIS, Dim > ConstArray_CellBASISDim