11#ifndef PANZER_BASIS_VALUES2_HPP
12#define PANZER_BASIS_VALUES2_HPP
14#include "Teuchos_RCP.hpp"
16#include "Intrepid2_Basis.hpp"
17#include "Intrepid2_Orientation.hpp"
25 class OrientationsInterface;
67 template <
typename Scalar>
71 using IntrepidBasis = Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>;
102 const bool allocArrays=
false,
103 const bool buildWeighted=
false);
106 void setupArrays(
const Teuchos::RCP<const panzer::BasisIRLayout>& basis,
107 bool computeDerivatives=
true);
109 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim> & cub_points,
110 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
111 const PHX::MDField<Scalar,Cell,IP> & jac_det,
112 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
113 const int in_num_cells = -1);
115 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim> & cub_points,
116 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
117 const PHX::MDField<Scalar,Cell,IP> & jac_det,
118 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
119 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
120 const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
121 bool use_node_coordinates=
true,
122 const int in_num_cells = -1);
124 void evaluateValuesCV(
const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
125 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
126 const PHX::MDField<Scalar,Cell,IP> & jac_det,
127 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv);
130 void evaluateValuesCV(
const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
131 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
132 const PHX::MDField<Scalar,Cell,IP> & jac_det,
133 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
134 const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
135 bool use_node_coordinates=
true,
136 const int in_num_cells = -1);
139 void evaluateValues(
const PHX::MDField<Scalar,Cell,IP,Dim> & cub_points,
140 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
141 const PHX::MDField<Scalar,Cell,IP> & jac_det,
142 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
143 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
144 const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
145 bool use_node_coordinates=
true,
146 const int in_num_cells = -1);
151 void applyOrientations(
const PHX::MDField<const Scalar,Cell,BASIS> & orientations);
155 const int in_num_cells = -1);
209 std::vector<PHX::index_size_type>
ddims_;
215 const int in_num_cells = -1);
303 setup(
const Teuchos::RCP<const panzer::BasisIRLayout> & basis,
304 PHX::MDField<const Scalar, Cell, IP, Dim> reference_points,
305 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian,
306 PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
307 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
308 const int num_evaluated_cells = -1);
324 setupUniform(
const Teuchos::RCP<const panzer::BasisIRLayout> & basis,
325 PHX::MDField<const Scalar, IP, Dim> reference_points,
326 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian,
327 PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
328 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
329 const int num_evaluated_cells = -1);
333 setOrientations(
const std::vector<Intrepid2::Orientation> & orientations,
334 const int num_orientations_cells = -1);
353 const std::vector<PHX::index_size_type> &
373 const bool force =
false)
const;
388 const bool force =
false)
const;
403 const bool force =
false)
const;
418 const bool force =
false)
const;
434 const bool force =
false)
const;
450 const bool force =
false)
const;
466 const bool force =
false)
const;
480 const bool force =
false)
const;
495 const bool cache =
true,
496 const bool force =
false)
const;
511 const bool cache =
true,
512 const bool force =
false)
const;
527 const bool cache =
true,
528 const bool force =
false)
const;
544 const bool cache =
true,
545 const bool force =
false)
const;
561 const bool cache =
true,
562 const bool force =
false)
const;
578 const bool cache =
true,
579 const bool force =
false)
const;
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
PHX::MDField< Scalar, Cell, BASIS, IP, Dim > Array_CellBasisIPDim
Array_CellBasisIP weighted_div_basis
Teuchos::RCP< const panzer::BasisIRLayout > basis_layout
ConstArray_BasisIP getCurl2DVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
Array_CellBasisIP basis_scalar
bool weighted_div_basis_evaluated_
Array_BasisIPDim grad_basis_ref
Intrepid2::Basis< PHX::Device::execution_space, Scalar, Scalar > IntrepidBasis
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
PHX::MDField< const Scalar, Cell, IP > cubature_jacobian_determinant_
Array_CellBasisIPDim weighted_curl_basis_vector
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.
bool grad_basis_evaluated_
PHX::MDField< Scalar, BASIS, Dim > Array_BasisDim
Array_BasisIP curl_basis_ref_scalar
Array_BasisDim basis_coordinates_ref
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)
ConstArray_BasisIPDim getGradBasisValuesRef(const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at reference points.
PHX::MDField< Scalar, BASIS, IP, Dim > Array_BasisIPDim
bool curl_basis_vector_evaluated_
Array_CellBasisIP weighted_curl_basis_scalar
void setExtendedDimensions(const std::vector< PHX::index_size_type > &ddims)
bool basis_coordinates_evaluated_
ConstArray_BasisIPDim getVectorBasisValuesRef(const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at reference points.
Array_CellBasisIPDim weighted_grad_basis
void setCellNodeCoordinates(PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates)
Set the cell node coordinates (required for getBasisCoordinates())
bool orientations_applied_
bool basis_vector_evaluated_
PureBasis::EElementSpace getElementSpace() const
bool weighted_curl_basis_scalar_evaluated_
bool curl_basis_ref_scalar_evaluated_
Array_CellBasisIPDim curl_basis_vector
bool div_basis_ref_evaluated_
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.
int num_orientations_cells_
PHX::MDField< Scalar > ArrayDynamic
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, const int in_num_cells=-1)
Array_CellBasisIP curl_basis_scalar
bool weighted_basis_vector_evaluated_
bool weighted_grad_basis_evaluated_
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.
const std::vector< PHX::index_size_type > & getExtendedDimensions() const
Get the extended dimensions used by sacado AD allocations.
bool basis_ref_vector_evaluated_
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
bool div_basis_evaluated_
ConstArray_BasisIP getBasisValuesRef(const bool cache=true, const bool force=false) const
Get the basis values evaluated at reference points.
PHX::MDField< const Scalar, BASIS, IP, Dim > ConstArray_BasisIPDim
bool grad_basis_ref_evaluated_
bool basis_coordinates_ref_evaluated_
Teuchos::RCP< const shards::CellTopology > cell_topology_
PHX::MDField< Scalar, Cell, BASIS, Dim > Array_CellBasisDim
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...
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)
Array_BasisIPDim curl_basis_ref_vector
std::vector< PHX::index_size_type > ddims_
PHX::MDField< const Scalar, IP, Dim > cubature_points_uniform_ref_
Array_BasisIP div_basis_ref
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.
bool hasUniformReferenceSpace() const
Check if reference point space is uniform across all cells (faster evaluation)
bool orientationsApplied() const
PHX::MDField< const Scalar, Cell, IP, Dim > cubature_points_ref_
bool curl_basis_ref_vector_evaluated_
Array_CellBasisIP div_basis
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
bool weighted_basis_scalar_evaluated_
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...
Array_CellBasisIPDim basis_vector
ConstArray_BasisDim getBasisCoordinatesRef(const bool cache=true, const bool force=false) const
Get the reference coordinates for basis.
Array_CellBasisDim basis_coordinates
bool weighted_curl_basis_vector_evaluated_
PHX::MDField< const Scalar, Cell, IP, Dim, Dim > cubature_jacobian_inverse_
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.
PHX::MDField< const Scalar, BASIS, IP > ConstArray_BasisIP
bool references_evaluated
Array_BasisIPDim basis_ref_vector
panzer::BasisDescriptor getBasisDescriptor() const
Return the basis descriptor.
ConstArray_BasisIPDim getCurlVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
Array_BasisIP basis_ref_scalar
bool basis_ref_scalar_evaluated_
Used to check if arrays have been cached.
Teuchos::RCP< IntrepidBasis > intrepid_basis
ConstArray_CellBasisDim getBasisCoordinates(const bool cache=true, const bool force=false) const
Carterisan coordinates for basis coefficients in mesh space.
PHX::MDField< const Scalar, Cell, IP > cubature_weights_
PHX::MDField< const Scalar > ConstArrayDynamic
bool curl_basis_scalar_evaluated_
Array_CellBasisIPDim grad_basis
void setupArrays(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true)
Sizes/allocates memory for arrays.
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.
Array_CellBasisIP weighted_basis_scalar
PHX::MDField< const Scalar, Cell, BASIS, Dim > ConstArray_CellBasisDim
PHX::MDField< const Scalar, Cell, NODE, Dim > cell_node_coordinates_
Array_CellBasisIPDim weighted_basis_vector
PHX::MDField< Scalar, Cell, BASIS, IP > Array_CellBasisIP
PHX::MDField< const Scalar, Cell, IP, Dim, Dim > cubature_jacobian_
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
std::vector< Intrepid2::Orientation > orientations_
ConstArray_BasisIP getDivVectorBasisRef(const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at reference points.
void applyOrientations(const PHX::MDField< const Scalar, Cell, BASIS > &orientations)
Method to apply orientations to a basis values container.
PHX::MDField< Scalar, BASIS, IP > Array_BasisIP
bool basis_scalar_evaluated_
PHX::MDField< const Scalar, BASIS, Dim > ConstArray_BasisDim