43#ifndef PANZER_BASIS_VALUES2_HPP
44#define PANZER_BASIS_VALUES2_HPP
46#include "Teuchos_RCP.hpp"
48#include "Intrepid2_Basis.hpp"
49#include "Intrepid2_Orientation.hpp"
57 class OrientationsInterface;
99 template <
typename Scalar>
103 using IntrepidBasis = Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>;
134 const bool allocArrays=
false,
135 const bool buildWeighted=
false);
138 void setupArrays(
const Teuchos::RCP<const panzer::BasisIRLayout>& basis,
139 bool computeDerivatives=
true);
141 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim> & cub_points,
142 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
143 const PHX::MDField<Scalar,Cell,IP> & jac_det,
144 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
145 const int in_num_cells = -1);
147 void evaluateValues(
const PHX::MDField<Scalar,IP,Dim> & cub_points,
148 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
149 const PHX::MDField<Scalar,Cell,IP> & jac_det,
150 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
151 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
152 const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
153 bool use_vertex_coordinates=
true,
154 const int in_num_cells = -1);
156 void evaluateValuesCV(
const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
157 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
158 const PHX::MDField<Scalar,Cell,IP> & jac_det,
159 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv);
162 void evaluateValuesCV(
const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
163 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
164 const PHX::MDField<Scalar,Cell,IP> & jac_det,
165 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
166 const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
167 bool use_vertex_coordinates=
true,
168 const int in_num_cells = -1);
171 void evaluateValues(
const PHX::MDField<Scalar,Cell,IP,Dim> & cub_points,
172 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> &
jac,
173 const PHX::MDField<Scalar,Cell,IP> & jac_det,
174 const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
175 const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
176 const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
177 bool use_vertex_coordinates=
true,
178 const int in_num_cells = -1);
183 void applyOrientations(
const PHX::MDField<const Scalar,Cell,BASIS> & orientations);
187 const int in_num_cells = -1);
241 std::vector<PHX::index_size_type>
ddims_;
247 const int in_num_cells = -1);
331 setup(
const Teuchos::RCP<const panzer::BasisIRLayout> & basis,
332 PHX::MDField<const Scalar, Cell, IP, Dim> reference_points,
333 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian,
334 PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
335 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
336 const int num_evaluated_cells = -1);
352 setupUniform(
const Teuchos::RCP<const panzer::BasisIRLayout> & basis,
353 PHX::MDField<const Scalar, IP, Dim> reference_points,
354 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian,
355 PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
356 PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
357 const int num_evaluated_cells = -1);
361 setOrientations(
const std::vector<Intrepid2::Orientation> & orientations,
362 const int num_orientations_cells = -1);
378 const std::vector<PHX::index_size_type> &
398 const bool force =
false)
const;
413 const bool force =
false)
const;
428 const bool force =
false)
const;
443 const bool force =
false)
const;
459 const bool force =
false)
const;
475 const bool force =
false)
const;
491 const bool force =
false)
const;
505 const bool force =
false)
const;
520 const bool cache =
true,
521 const bool force =
false)
const;
536 const bool cache =
true,
537 const bool force =
false)
const;
552 const bool cache =
true,
553 const bool force =
false)
const;
569 const bool cache =
true,
570 const bool force =
false)
const;
586 const bool cache =
true,
587 const bool force =
false)
const;
603 const bool cache =
true,
604 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
BasisValues2(const std::string &prefix="", const bool allocArrays=false, const bool buildWeighted=false)
Main constructor.
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
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
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_
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int in_num_cells=-1)
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_
PHX::MDField< const Scalar, Cell, NODE, Dim > cell_vertex_coordinates_
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
void setCellVertexCoordinates(PHX::MDField< Scalar, Cell, NODE, Dim > vertex_coordinates)
Set the cell vertex coordinates (required for getBasisCoordinates())
Array_BasisIPDim basis_ref_vector
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
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