Panzer
Version of the Day
|
#include <Panzer_BasisValues2.hpp>
Public Member Functions | |
BasisValues2 (const std::string &prefix="", const bool allocArrays=false, const bool buildWeighted=false) | |
Main constructor. More... | |
void | setupArrays (const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true) |
Sizes/allocates memory for arrays. More... | |
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) |
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 PHX::MDField< Scalar, Cell, IP > &weighted_measure, const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, bool use_node_coordinates=true, const int in_num_cells=-1) |
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) |
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, const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, bool use_node_coordinates=true, const int in_num_cells=-1) |
void | evaluateValues (const PHX::MDField< Scalar, Cell, 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 PHX::MDField< Scalar, Cell, IP > &weighted_measure, const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, bool use_node_coordinates=true, const int in_num_cells=-1) |
void | applyOrientations (const PHX::MDField< const Scalar, Cell, BASIS > &orientations) |
Method to apply orientations to a basis values container. More... | |
void | applyOrientations (const std::vector< Intrepid2::Orientation > &orientations, const int in_num_cells=-1) |
void | setExtendedDimensions (const std::vector< PHX::index_size_type > &ddims) |
PureBasis::EElementSpace | getElementSpace () const |
bool | orientationsApplied () const |
void | evaluateBasisCoordinates (const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, const int in_num_cells=-1) |
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. More... | |
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. More... | |
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 certain bases. More... | |
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 basis objects. More... | |
void | setCellVertexCoordinates (PHX::MDField< Scalar, Cell, NODE, Dim > vertex_coordinates) |
void | setCellNodeCoordinates (PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates) |
Set the cell node coordinates (required for getBasisCoordinates()) More... | |
bool | hasUniformReferenceSpace () const |
Check if reference point space is uniform across all cells (faster evaluation) More... | |
panzer::BasisDescriptor | getBasisDescriptor () const |
Return the basis descriptor. More... | |
const std::vector < PHX::index_size_type > & | getExtendedDimensions () const |
Get the extended dimensions used by sacado AD allocations. More... | |
ConstArray_BasisDim | getBasisCoordinatesRef (const bool cache=true, const bool force=false) const |
Get the reference coordinates for basis. More... | |
ConstArray_BasisIP | getBasisValuesRef (const bool cache=true, const bool force=false) const |
Get the basis values evaluated at reference points. More... | |
ConstArray_BasisIPDim | getVectorBasisValuesRef (const bool cache=true, const bool force=false) const |
Get the vector basis values evaluated at reference points. More... | |
ConstArray_BasisIPDim | getGradBasisValuesRef (const bool cache=true, const bool force=false) const |
Get the gradient of the basis evaluated at reference points. More... | |
ConstArray_BasisIP | getCurl2DVectorBasisRef (const bool cache=true, const bool force=false) const |
Get the curl of a vector basis evaluated at reference points. More... | |
ConstArray_BasisIPDim | getCurlVectorBasisRef (const bool cache=true, const bool force=false) const |
Get the curl of a vector basis evaluated at reference points. More... | |
ConstArray_BasisIP | getDivVectorBasisRef (const bool cache=true, const bool force=false) const |
Get the divergence of a vector basis evaluated at reference points. More... | |
ConstArray_CellBasisDim | getBasisCoordinates (const bool cache=true, const bool force=false) const |
Carterisan coordinates for basis coefficients in mesh space. More... | |
ConstArray_CellBasisIP | getBasisValues (const bool weighted, const bool cache=true, const bool force=false) const |
Get the basis values evaluated at mesh points. More... | |
ConstArray_CellBasisIPDim | getVectorBasisValues (const bool weighted, const bool cache=true, const bool force=false) const |
Get the vector basis values evaluated at mesh points. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
Protected Member Functions | |
void | resetArrays () |
Private Attributes | |
bool | references_evaluated |
bool | orientations_applied_ |
Data structure that holds all evaluated fields associated with a basis function and integration rule. This class will allocate the memory and evaluate the basis functions. The orientations must be applied using the applyOrientations
method.
Lazy evaluation path:
BasisValues can be defined in standard or uniform forms. Uniform means that the reference space is the same for all cells, while standard means that there is a different reference space per cell. If the standard method is used, all get*Ref calls will throw errors.
The construction path for the lazy evaluation form of the BasisValues2 object is as follows:
auto basis_values = Teuchos::rcp(new BasisValues<Scalar>("prefix"));
basis_values->setup[Uniform](basis, reference_points, jacobian, jacobian_determinant, jacobian_inverse);
basis_values->setOrientations(orientations);
basis_values->setExtendedDimensions(ddims);
basis_values->setWeightedMeasure(weighted_measure);
Note that the above optional/required calls should happen in this order, but all must happen before any 'get' calls can be called.
Definition at line 68 of file Panzer_BasisValues2.hpp.
typedef ArrayTraits<Scalar,PHX::MDField<Scalar> >::size_type panzer::BasisValues2< Scalar >::size_type |
Definition at line 70 of file Panzer_BasisValues2.hpp.
using panzer::BasisValues2< Scalar >::IntrepidBasis = Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar> |
Definition at line 71 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar> panzer::BasisValues2< Scalar >::ArrayDynamic |
Definition at line 73 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,BASIS,IP> panzer::BasisValues2< Scalar >::Array_BasisIP |
Definition at line 74 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,BASIS,IP,Dim> panzer::BasisValues2< Scalar >::Array_BasisIPDim |
Definition at line 75 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,BASIS,Dim> panzer::BasisValues2< Scalar >::Array_BasisDim |
Definition at line 76 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,Cell,BASIS,IP> panzer::BasisValues2< Scalar >::Array_CellBasisIP |
Definition at line 77 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,Cell,BASIS,IP,Dim> panzer::BasisValues2< Scalar >::Array_CellBasisIPDim |
Definition at line 78 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<Scalar,Cell,BASIS,Dim> panzer::BasisValues2< Scalar >::Array_CellBasisDim |
Definition at line 79 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar> panzer::BasisValues2< Scalar >::ConstArrayDynamic |
Definition at line 81 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,BASIS,IP> panzer::BasisValues2< Scalar >::ConstArray_BasisIP |
Definition at line 82 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,BASIS,IP,Dim> panzer::BasisValues2< Scalar >::ConstArray_BasisIPDim |
Definition at line 83 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,BASIS,Dim> panzer::BasisValues2< Scalar >::ConstArray_BasisDim |
Definition at line 84 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,Cell,BASIS,IP> panzer::BasisValues2< Scalar >::ConstArray_CellBasisIP |
Definition at line 85 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,Cell,BASIS,IP,Dim> panzer::BasisValues2< Scalar >::ConstArray_CellBasisIPDim |
Definition at line 86 of file Panzer_BasisValues2.hpp.
typedef PHX::MDField<const Scalar,Cell,BASIS,Dim> panzer::BasisValues2< Scalar >::ConstArray_CellBasisDim |
Definition at line 87 of file Panzer_BasisValues2.hpp.
panzer::BasisValues2< Scalar >::BasisValues2 | ( | const std::string & | prefix = "" , |
const bool | allocArrays = false , |
||
const bool | buildWeighted = false |
||
) |
Main constructor.
[in] | prefix | Name to prefix all arrays with |
[in] | allocArrays | If true we will allocate all arrays in the setupArrays call |
[in] | buildWeighted | Builds the weighted components for the basis (i.e. multiplied by weighted measure for integration purposes) |
Definition at line 85 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::setupArrays | ( | const Teuchos::RCP< const panzer::BasisIRLayout > & | basis, |
bool | computeDerivatives = true |
||
) |
Sizes/allocates memory for arrays.
Definition at line 415 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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 |
||
) |
Definition at line 126 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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 PHX::MDField< Scalar, Cell, IP > & | weighted_measure, | ||
const PHX::MDField< Scalar, Cell, NODE, Dim > & | node_coordinates, | ||
bool | use_node_coordinates = true , |
||
const int | in_num_cells = -1 |
||
) |
Definition at line 189 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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 | ||
) |
Definition at line 305 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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, | ||
const PHX::MDField< Scalar, Cell, NODE, Dim > & | node_coordinates, | ||
bool | use_node_coordinates = true , |
||
const int | in_num_cells = -1 |
||
) |
Definition at line 319 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::evaluateValues | ( | const PHX::MDField< Scalar, Cell, 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 PHX::MDField< Scalar, Cell, IP > & | weighted_measure, | ||
const PHX::MDField< Scalar, Cell, NODE, Dim > & | node_coordinates, | ||
bool | use_node_coordinates = true , |
||
const int | in_num_cells = -1 |
||
) |
Definition at line 243 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::applyOrientations | ( | const PHX::MDField< const Scalar, Cell, BASIS > & | orientations | ) |
Method to apply orientations to a basis values container.
Definition at line 404 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::applyOrientations | ( | const std::vector< Intrepid2::Orientation > & | orientations, |
const int | in_num_cells = -1 |
||
) |
Definition at line 346 of file Panzer_BasisValues2_impl.hpp.
|
inline |
Definition at line 157 of file Panzer_BasisValues2.hpp.
PureBasis::EElementSpace panzer::BasisValues2< Scalar >::getElementSpace | ( | ) | const |
Definition at line 410 of file Panzer_BasisValues2_impl.hpp.
|
inline |
Definition at line 211 of file Panzer_BasisValues2.hpp.
void panzer::BasisValues2< Scalar >::evaluateBasisCoordinates | ( | const PHX::MDField< Scalar, Cell, NODE, Dim > & | node_coordinates, |
const int | in_num_cells = -1 |
||
) |
Definition at line 334 of file Panzer_BasisValues2_impl.hpp.
|
protected |
Definition at line 668 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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.
[in] | basis | Basis layout - contains information for intrepid |
[in] | reference_points | Points (e.g. cubature points) in the reference space of the cell |
[in] | point_jacobian | Cell jacobian evaluated at the reference points |
[in] | point_jacobian_determinant | Determinant of point_jacobian array |
[in] | point_jacobian_inverse | Inverse of point_jacobian array |
[in] | num_evaluated_cells | Used to force evaluation of arrays over subset of cells (default: all cells) |
Definition at line 578 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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.
[in] | basis | Basis layout - contains information for intrepid |
[in] | reference_points | Points (e.g. cubature points) in the reference space of the cell |
[in] | point_jacobian | Cell jacobian evaluated at the reference points |
[in] | point_jacobian_determinant | Determinant of point_jacobian array |
[in] | point_jacobian_inverse | Inverse of point_jacobian array |
[in] | num_evaluated_cells | Used to force evaluation of arrays over subset of cells (default: all cells) |
Definition at line 605 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::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 certain bases.
Definition at line 632 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::setWeightedMeasure | ( | PHX::MDField< const Scalar, Cell, IP > | weighted_measure | ) |
Set the cubature weights (weighted measure) for the basis values object - required to get weighted basis objects.
Definition at line 721 of file Panzer_BasisValues2_impl.hpp.
void panzer::BasisValues2< Scalar >::setCellVertexCoordinates | ( | PHX::MDField< Scalar, Cell, NODE, Dim > | vertex_coordinates | ) |
TO BE DEPRECATED..... Set the cell vertex coordinates (required for getBasisCoordinates())
void panzer::BasisValues2< Scalar >::setCellNodeCoordinates | ( | PHX::MDField< Scalar, Cell, NODE, Dim > | node_coordinates | ) |
Set the cell node coordinates (required for getBasisCoordinates())
Definition at line 653 of file Panzer_BasisValues2_impl.hpp.
|
inline |
Check if reference point space is uniform across all cells (faster evaluation)
Definition at line 352 of file Panzer_BasisValues2.hpp.
panzer::BasisDescriptor panzer::BasisValues2< Scalar >::getBasisDescriptor | ( | ) | const |
Return the basis descriptor.
Definition at line 659 of file Panzer_BasisValues2_impl.hpp.
|
inline |
Get the extended dimensions used by sacado AD allocations.
Definition at line 360 of file Panzer_BasisValues2.hpp.
BasisValues2< Scalar >::ConstArray_BasisDim panzer::BasisValues2< Scalar >::getBasisCoordinatesRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the reference coordinates for basis.
If | reference points are not uniform |
If | basis does not support a coordinate space |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
return Array <basis, dim>
Definition at line 747 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIP panzer::BasisValues2< Scalar >::getBasisValuesRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the basis values evaluated at reference points.
If | not a scalar basis |
If | reference points are not uniform |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 775 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIPDim panzer::BasisValues2< Scalar >::getVectorBasisValuesRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the vector basis values evaluated at reference points.
If | not a vector basis |
If | reference points are not uniform |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 811 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIPDim panzer::BasisValues2< Scalar >::getGradBasisValuesRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the gradient of the basis evaluated at reference points.
If | not a scalar basis |
If | reference points are not uniform |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 848 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIP panzer::BasisValues2< Scalar >::getCurl2DVectorBasisRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the curl of a vector basis evaluated at reference points.
If | reference points are not uniform |
If | not a 2D space |
If | not a HCurl basis |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 885 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIPDim panzer::BasisValues2< Scalar >::getCurlVectorBasisRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the curl of a vector basis evaluated at reference points.
If | reference points are not uniform |
If | not a 3D space |
If | not a HCurl basis |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 922 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_BasisIP panzer::BasisValues2< Scalar >::getDivVectorBasisRef | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Get the divergence of a vector basis evaluated at reference points.
If | reference points are not uniform |
If | not a 3D space |
If | not a HDiv basis |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 960 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisDim panzer::BasisValues2< Scalar >::getBasisCoordinates | ( | const bool | cache = true , |
const bool | force = false |
||
) | const |
Carterisan coordinates for basis coefficients in mesh space.
If | basis does not support a coordinate space |
return Array <cell, basis, dim>
Definition at line 996 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIP panzer::BasisValues2< Scalar >::getBasisValues | ( | const bool | weighted, |
const bool | cache = true , |
||
const bool | force = false |
||
) | const |
Get the basis values evaluated at mesh points.
If | not a scalar basis |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1035 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIPDim panzer::BasisValues2< Scalar >::getVectorBasisValues | ( | const bool | weighted, |
const bool | cache = true , |
||
const bool | force = false |
||
) | const |
Get the vector basis values evaluated at mesh points.
If | not a vector basis |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1202 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIPDim panzer::BasisValues2< Scalar >::getGradBasisValues | ( | const bool | weighted, |
const bool | cache = true , |
||
const bool | force = false |
||
) | const |
Get the gradient of the basis evaluated at mesh points.
If | not a scalar basis |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1377 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIP panzer::BasisValues2< Scalar >::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.
If | not a HCurl basis |
If | not a 2D space |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1531 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIPDim panzer::BasisValues2< Scalar >::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.
If | not a HCurl basis |
If | not a 3D space |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1689 of file Panzer_BasisValues2_impl.hpp.
BasisValues2< Scalar >::ConstArray_CellBasisIP panzer::BasisValues2< Scalar >::getDivVectorBasis | ( | const bool | weighted, |
const bool | cache = true , |
||
const bool | force = false |
||
) | const |
Get the divergence of a vector basis evaluated at mesh points.
If | not a HDiv basis |
If | not a 3D space |
[in] | weighted | Add cubature weighting for integration purposes |
[in] | cache | If true, the returned object will be cached for later use |
[in] | force | Force re-evaluation of cached array |
Definition at line 1842 of file Panzer_BasisValues2_impl.hpp.
|
mutable |
Definition at line 162 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 163 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 165 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 166 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 168 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 169 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 171 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 172 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 174 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 175 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 177 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 178 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 180 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 181 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 182 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 183 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 184 of file Panzer_BasisValues2.hpp.
|
mutable |
Definition at line 185 of file Panzer_BasisValues2.hpp.
|
mutable |
Carterisan coordinates for basis coefficients
NOTE: This quantity is not always available. Certain bases may not have a corresponding coordiante value
Definition at line 192 of file Panzer_BasisValues2.hpp.
|
mutable |
Carterisan coordinates for basis coefficients
NOTE: This quantity is not always available. Certain bases may not have a corresponding coordiante value
Definition at line 199 of file Panzer_BasisValues2.hpp.
Teuchos::RCP<const panzer::BasisIRLayout> panzer::BasisValues2< Scalar >::basis_layout |
Definition at line 201 of file Panzer_BasisValues2.hpp.
Teuchos::RCP<IntrepidBasis> panzer::BasisValues2< Scalar >::intrepid_basis |
Definition at line 203 of file Panzer_BasisValues2.hpp.
bool panzer::BasisValues2< Scalar >::compute_derivatives |
Definition at line 205 of file Panzer_BasisValues2.hpp.
bool panzer::BasisValues2< Scalar >::build_weighted |
Definition at line 206 of file Panzer_BasisValues2.hpp.
bool panzer::BasisValues2< Scalar >::alloc_arrays |
Definition at line 207 of file Panzer_BasisValues2.hpp.
std::string panzer::BasisValues2< Scalar >::prefix |
Definition at line 208 of file Panzer_BasisValues2.hpp.
std::vector<PHX::index_size_type> panzer::BasisValues2< Scalar >::ddims_ |
Definition at line 209 of file Panzer_BasisValues2.hpp.
|
private |
Definition at line 219 of file Panzer_BasisValues2.hpp.
|
private |
Definition at line 222 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 233 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 236 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 239 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 242 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 245 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 248 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 249 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 250 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 251 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 254 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 257 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 260 of file Panzer_BasisValues2.hpp.
|
protected |
Definition at line 263 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Used to check if arrays have been cached.
Definition at line 266 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 267 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 268 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 269 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 270 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 271 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 272 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 273 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 274 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 275 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 276 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 277 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 278 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 279 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 280 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 281 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 282 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 283 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 284 of file Panzer_BasisValues2.hpp.
|
mutableprotected |
Definition at line 285 of file Panzer_BasisValues2.hpp.