Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
panzer::IntegrationValues2< Scalar > Class Template Reference

#include <Panzer_IntegrationValues2.hpp>

Public Types

typedef ArrayTraits< Scalar,
PHX::MDField< Scalar >
>::size_type 
size_type
 
typedef PHX::MDField< Scalar > ArrayDynamic
 
typedef PHX::MDField< double > DblArrayDynamic
 
typedef PHX::MDField< Scalar, IPArray_IP
 
typedef PHX::MDField< Scalar,
IP, Dim
Array_IPDim
 
typedef PHX::MDField< Scalar,
Point
Array_Point
 
typedef PHX::MDField< Scalar,
Cell, IP
Array_CellIP
 
typedef PHX::MDField< Scalar,
Cell, IP, Dim
Array_CellIPDim
 
typedef PHX::MDField< Scalar,
Cell, IP, Dim, Dim
Array_CellIPDimDim
 
typedef PHX::MDField< Scalar,
Cell, BASIS, Dim
Array_CellBASISDim
 
typedef PHX::MDField< const
Scalar, IP
ConstArray_IP
 
typedef PHX::MDField< const
Scalar, IP, Dim
ConstArray_IPDim
 
typedef PHX::MDField< const
Scalar, Point
ConstArray_Point
 
typedef PHX::MDField< const
Scalar, Cell, IP
ConstArray_CellIP
 
typedef PHX::MDField< const
Scalar, Cell, IP, Dim
ConstArray_CellIPDim
 
typedef PHX::MDField< const
Scalar, Cell, IP, Dim, Dim
ConstArray_CellIPDimDim
 
typedef PHX::MDField< const
Scalar, Cell, BASIS, Dim
ConstArray_CellBASISDim
 

Public Member Functions

 IntegrationValues2 (const std::string &pre="", const bool allocArrays=false)
 Base constructor. More...
 
void setupArrays (const Teuchos::RCP< const panzer::IntegrationRule > &ir)
 Sizes/allocates memory for arrays. More...
 
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. More...
 
void evaluateValues (const PHX::MDField< Scalar, Cell, NODE, Dim > &cell_node_coordinates, const PHX::MDField< Scalar, Cell, IP, Dim > &other_ip_coordinates, const int num_cells=-1)
 Match IP. More...
 
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. More...
 
void setupPermutations (const Teuchos::RCP< const SubcellConnectivity > &face_connectivity, const int num_virtual_cells)
 Initialize the permutation arrays given a face connectivity. More...
 
void setupPermutations (const PHX::MDField< Scalar, Cell, IP, Dim > &other_ip_coordinates)
 Initialize the permutation arrays given another IntegrationValues2<Scalar>::getCubaturePoints() array. More...
 
ConstArray_IPDim getUniformCubaturePointsRef (const bool cache=true, const bool force=false, const bool apply_permutation=true) const
 Get the uniform cubature points. More...
 
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. More...
 
ConstArray_IP getUniformCubatureWeightsRef (const bool cache=true, const bool force=false, const bool apply_permutation=true) const
 Get the uniform cubature weights. More...
 
ConstArray_CellBASISDim getNodeCoordinates () const
 Get the node coordinates describing the geometry of the mesh. More...
 
ConstArray_CellIPDimDim getJacobian (const bool cache=true, const bool force=false) const
 Get the Jacobian matrix evaluated at the cubature points. More...
 
ConstArray_CellIPDimDim getJacobianInverse (const bool cache=true, const bool force=false) const
 Get the inverse of the Jacobian matrix evaluated at the cubature points. More...
 
ConstArray_CellIP getJacobianDeterminant (const bool cache=true, const bool force=false) const
 Get the determinant of the Jacobian matrix evaluated at the cubature points. More...
 
ConstArray_CellIP getWeightedMeasure (const bool cache=true, const bool force=false) const
 Get the weighted measure (integration weights) More...
 
ConstArray_CellIPDim getWeightedNormals (const bool cache=true, const bool force=false) const
 Get the weighted normals. More...
 
ConstArray_CellIPDim getSurfaceNormals (const bool cache=true, const bool force=false) const
 Get the surface normals. More...
 
ConstArray_CellIPDimDim getSurfaceRotationMatrices (const bool cache=true, const bool force=false) const
 Get the surface rotation matrices. More...
 
ConstArray_CellIPDimDim getCovarientMatrix (const bool cache=true, const bool force=false) const
 Get the covarient matrix. More...
 
ConstArray_CellIPDimDim getContravarientMatrix (const bool cache=true, const bool force=false) const
 Get the contravarient matrix. More...
 
ConstArray_CellIP getNormContravarientMatrix (const bool cache=true, const bool force=false) const
 Get the contravarient matrix. More...
 
ConstArray_CellIPDim getCubaturePoints (const bool cache=true, const bool force=false) const
 Get the cubature points in physical space. More...
 
ConstArray_CellIPDim getCubaturePointsRef (const bool cache=true, const bool force=false) const
 Get the cubature points in the reference space. More...
 
Teuchos::RCP< const
panzer::IntegrationRule
getIntegrationRule () const
 Returns the IntegrationRule. More...
 

Public Attributes

Array_IPDim cub_points
 
Array_IPDim side_cub_points
 
Array_IP cub_weights
 
Array_CellBASISDim node_coordinates
 
Array_CellIPDimDim jac
 
Array_CellIPDimDim jac_inv
 
Array_CellIP jac_det
 
Array_CellIP weighted_measure
 
Array_CellIPDim weighted_normals
 
Array_CellIPDim surface_normals
 
Array_CellIPDimDim surface_rotation_matrices
 
Array_CellIPDimDim covarient
 
Array_CellIPDimDim contravarient
 
Array_CellIP norm_contravarient
 
Array_CellIPDim ip_coordinates
 
Array_CellIPDim ref_ip_coordinates
 
Teuchos::RCP< const
panzer::IntegrationRule
int_rule
 
Teuchos::RCP
< Intrepid2::Cubature
< PHX::Device::execution_space,
double, double > > 
intrepid_cubature
 

Protected Member Functions

void resetArrays ()
 
void evaluateEverything ()
 

Protected Attributes

int num_cells_
 
int num_evaluate_cells_
 
int num_virtual_cells_
 
bool requires_permutation_
 
PHX::MDField< const int, Cell, IPpermutations_
 
Teuchos::RCP< const
SubcellConnectivity
side_connectivity_
 
bool cub_points_evaluated_
 
bool side_cub_points_evaluated_
 
bool cub_weights_evaluated_
 
bool node_coordinates_evaluated_
 
bool jac_evaluated_
 
bool jac_inv_evaluated_
 
bool jac_det_evaluated_
 
bool weighted_measure_evaluated_
 
bool weighted_normals_evaluated_
 
bool surface_normals_evaluated_
 
bool surface_rotation_matrices_evaluated_
 
bool covarient_evaluated_
 
bool contravarient_evaluated_
 
bool norm_contravarient_evaluated_
 
bool ip_coordinates_evaluated_
 
bool ref_ip_coordinates_evaluated_
 

Private Attributes

bool alloc_arrays_
 
std::string prefix_
 
std::vector< PHX::index_size_type > ddims_
 

Detailed Description

template<typename Scalar>
class panzer::IntegrationValues2< Scalar >

Definition at line 61 of file Panzer_IntegrationValues2.hpp.

Member Typedef Documentation

template<typename Scalar>
typedef ArrayTraits<Scalar,PHX::MDField<Scalar> >::size_type panzer::IntegrationValues2< Scalar >::size_type

Definition at line 63 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar> panzer::IntegrationValues2< Scalar >::ArrayDynamic

Definition at line 65 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<double> panzer::IntegrationValues2< Scalar >::DblArrayDynamic

Definition at line 66 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,IP> panzer::IntegrationValues2< Scalar >::Array_IP

Definition at line 68 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,IP,Dim> panzer::IntegrationValues2< Scalar >::Array_IPDim

Definition at line 69 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,Point> panzer::IntegrationValues2< Scalar >::Array_Point

Definition at line 70 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,Cell,IP> panzer::IntegrationValues2< Scalar >::Array_CellIP

Definition at line 71 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,Cell,IP,Dim> panzer::IntegrationValues2< Scalar >::Array_CellIPDim

Definition at line 72 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,Cell,IP,Dim,Dim> panzer::IntegrationValues2< Scalar >::Array_CellIPDimDim

Definition at line 73 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<Scalar,Cell,BASIS,Dim> panzer::IntegrationValues2< Scalar >::Array_CellBASISDim

Definition at line 74 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,IP> panzer::IntegrationValues2< Scalar >::ConstArray_IP

Definition at line 76 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,IP,Dim> panzer::IntegrationValues2< Scalar >::ConstArray_IPDim

Definition at line 77 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,Point> panzer::IntegrationValues2< Scalar >::ConstArray_Point

Definition at line 78 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,Cell,IP> panzer::IntegrationValues2< Scalar >::ConstArray_CellIP

Definition at line 79 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,Cell,IP,Dim> panzer::IntegrationValues2< Scalar >::ConstArray_CellIPDim

Definition at line 80 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,Cell,IP,Dim,Dim> panzer::IntegrationValues2< Scalar >::ConstArray_CellIPDimDim

Definition at line 81 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
typedef PHX::MDField<const Scalar,Cell,BASIS,Dim> panzer::IntegrationValues2< Scalar >::ConstArray_CellBASISDim

Definition at line 82 of file Panzer_IntegrationValues2.hpp.

Constructor & Destructor Documentation

template<typename Scalar >
panzer::IntegrationValues2< Scalar >::IntegrationValues2 ( const std::string &  pre = "",
const bool  allocArrays = false 
)

Base constructor.

Parameters
[in]prePrefix to apply to all internal field names
[in]allocArrays(Classic Interface Only) Allocate array data in 'setupArrays' call

Definition at line 587 of file Panzer_IntegrationValues2.cpp.

Member Function Documentation

template<typename Scalar >
void panzer::IntegrationValues2< Scalar >::setupArrays ( const Teuchos::RCP< const panzer::IntegrationRule > &  ir)

Sizes/allocates memory for arrays.

Definition at line 627 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar>
void panzer::IntegrationValues2< Scalar >::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.

Parameters
cell_node_coordinates[in] Cell node coordinates, not basis coordinates.
num_cells[in] (optional) number of cells in the workset. This can be less than the workset size. If set to zero, extent(0) of the evaluated array is used which equates to the workset size.
face_connectivity[in] (optional) connectivity used to enforce quadrature alignment for surface integration.

Definition at line 690 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar>
void panzer::IntegrationValues2< Scalar >::evaluateValues ( const PHX::MDField< Scalar, Cell, NODE, Dim > &  cell_node_coordinates,
const PHX::MDField< Scalar, Cell, IP, Dim > &  other_ip_coordinates,
const int  num_cells = -1 
)

Match IP.

Optionally provide IP coordinates for an element 'other' that shares the same side. If provided, a permutation of the cubature points is calculated so that the integration values are ordered according to the other element's. This permutation is then applied so that all fields are ordered accordingly in their IP dimension.

Parameters
num_cells[in] (optional) number of cells in the workset. This can be less than the workset size. If set to zero, extent(0) of the evaluated array is used which equates to the workset size.

Definition at line 711 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar>
void panzer::IntegrationValues2< Scalar >::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.

The lazy evaluation construction path is designed such that you only allocate and fill arrays on demand, with an option of caching those generated fields. This is useful for when we are worried about a code's memory footprint.

The process for setting up one of these objects is to initialize the IntegrationValues2 class as follows:

IntegrationValues2<double> values; // Constructor values.setup(ir, cell_node_coordinates); // Required - must come before all other calls values.setupPermutations(...); // Optional - only needed for surface integration and some side integration rules - must come before get*() calls auto array = values.get*(); // Lazy evaluate whatever field you need

Todo:
Instead of IntegrationRule, we just need to load the integration descriptor and the cell topology
Parameters
[in]irIntegration rule descripting integration scheme
[in]cell_node_coordinatesNode/Vertex <cell, node, dim> coordinates describing cell geometry
[in]num_cellsIn case you need to only generate integration values for the first 'num_cells' of the node_coordinates - defaults to all cells

Definition at line 764 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
void panzer::IntegrationValues2< Scalar >::setupPermutations ( const Teuchos::RCP< const SubcellConnectivity > &  face_connectivity,
const int  num_virtual_cells 
)

Initialize the permutation arrays given a face connectivity.

Note
REQUIRED FOR SURFACE INTEGRATION
Must be called AFTER setup
Virtual cells have a unique way to generate surface normals and rotation matrices, hence we need to know how many are included
Parameters
[in]face_connectivityConnectivity describing how sides are connected to cells
[in]num_virtual_cellsNumber of virtual cells included in the node coordinates (found at end of node coordinate array)

Definition at line 729 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar>
void panzer::IntegrationValues2< Scalar >::setupPermutations ( const PHX::MDField< Scalar, Cell, IP, Dim > &  other_ip_coordinates)

Initialize the permutation arrays given another IntegrationValues2<Scalar>::getCubaturePoints() array.

Note
Required if you want points to line up between two integration values (e.g. for side integration)
Must be called AFTER setup
Parameters
[in]other_ip_coordinatesCubature points to align with

Definition at line 751 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_IPDim panzer::IntegrationValues2< Scalar >::getUniformCubaturePointsRef ( const bool  cache = true,
const bool  force = false,
const bool  apply_permutation = true 
) const

Get the uniform cubature points.

Note
DEPRECATED Please use getCubaturePointsRef call instead
Option is only supported for volume integration, and even then, may not align with the getCubaturePointsRef call
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
[in]apply_permutationADVANCED Do not change this unless you know what it does (it can break things)
Returns
Array <point, dim>

Definition at line 802 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_IPDim panzer::IntegrationValues2< Scalar >::getUniformSideCubaturePointsRef ( const bool  cache = true,
const bool  force = false,
const bool  apply_permutation = true 
) const

Get the uniform cubature points for a side.

Note
DEPRECATED Please use getCubaturePointsRef call instead
Option is only supported for side integration, and even then, may not align with the getCubaturePointsRef call
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
[in]apply_permutationADVANCED Do not change this unless you know what it does (it can break things)
Returns
Array <point, dim>

Definition at line 860 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_IP panzer::IntegrationValues2< Scalar >::getUniformCubatureWeightsRef ( const bool  cache = true,
const bool  force = false,
const bool  apply_permutation = true 
) const

Get the uniform cubature weights.

Note
DEPRECATED Please do not use
Option is only supported for some integration types
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
[in]apply_permutationADVANCED Do not change this unless you know what it does (it can break things)
Returns
Array <point>

Definition at line 912 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellBASISDim panzer::IntegrationValues2< Scalar >::getNodeCoordinates ( ) const

Get the node coordinates describing the geometry of the mesh.

Returns
Array <cell, node, dim>

Definition at line 962 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDimDim panzer::IntegrationValues2< Scalar >::getJacobian ( const bool  cache = true,
const bool  force = false 
) const

Get the Jacobian matrix evaluated at the cubature points.

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim, dim>

Definition at line 970 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDimDim panzer::IntegrationValues2< Scalar >::getJacobianInverse ( const bool  cache = true,
const bool  force = false 
) const

Get the inverse of the Jacobian matrix evaluated at the cubature points.

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim, dim>

Definition at line 1011 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIP panzer::IntegrationValues2< Scalar >::getJacobianDeterminant ( const bool  cache = true,
const bool  force = false 
) const

Get the determinant of the Jacobian matrix evaluated at the cubature points.

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point>

Definition at line 1048 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIP panzer::IntegrationValues2< Scalar >::getWeightedMeasure ( const bool  cache = true,
const bool  force = false 
) const

Get the weighted measure (integration weights)

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point>

Definition at line 1084 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDim panzer::IntegrationValues2< Scalar >::getWeightedNormals ( const bool  cache = true,
const bool  force = false 
) const

Get the weighted normals.

Note
Support: CV_SIDE
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim>

Definition at line 1252 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDim panzer::IntegrationValues2< Scalar >::getSurfaceNormals ( const bool  cache = true,
const bool  force = false 
) const

Get the surface normals.

Note
Support: SURFACE
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim>

Definition at line 1299 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDimDim panzer::IntegrationValues2< Scalar >::getSurfaceRotationMatrices ( const bool  cache = true,
const bool  force = false 
) const

Get the surface rotation matrices.

Note
Support: SURFACE
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, 3, 3>

Definition at line 1417 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDimDim panzer::IntegrationValues2< Scalar >::getCovarientMatrix ( const bool  cache = true,
const bool  force = false 
) const

Get the covarient matrix.

cov(i,j) = jacobian(i,k) * jacobian(j,k)

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim, dim>

Definition at line 1478 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDimDim panzer::IntegrationValues2< Scalar >::getContravarientMatrix ( const bool  cache = true,
const bool  force = false 
) const

Get the contravarient matrix.

contra = (getCovarientMatrix())^{-1}

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim, dim>

Definition at line 1520 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIP panzer::IntegrationValues2< Scalar >::getNormContravarientMatrix ( const bool  cache = true,
const bool  force = false 
) const

Get the contravarient matrix.

norm = sqrt({ij} cov(i,j) * cov(i,j))

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point>

Definition at line 1555 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDim panzer::IntegrationValues2< Scalar >::getCubaturePoints ( const bool  cache = true,
const bool  force = false 
) const

Get the cubature points in physical space.

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim>

Definition at line 1596 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
IntegrationValues2< Scalar >::ConstArray_CellIPDim panzer::IntegrationValues2< Scalar >::getCubaturePointsRef ( const bool  cache = true,
const bool  force = false 
) const

Get the cubature points in the reference space.

Note
Support: VOLUME, SURFACE, SIDE, CV_VOLUME, CV_SIDE, CV_BOUNDARY
Parameters
[in]cacheIf true, the result will be stored in the IntegrationValues2 class
[in]forceForce the re-evaluation of the array
Returns
Array <cell, point, dim>

Definition at line 1665 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar>
Teuchos::RCP<const panzer::IntegrationRule> panzer::IntegrationValues2< Scalar >::getIntegrationRule ( ) const
inline

Returns the IntegrationRule.

Returns
panzer::IntegrationRule

Definition at line 461 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar >
void panzer::IntegrationValues2< Scalar >::resetArrays ( )
protected

Definition at line 603 of file Panzer_IntegrationValues2.cpp.

template<typename Scalar >
void panzer::IntegrationValues2< Scalar >::evaluateEverything ( )
protected

Definition at line 1787 of file Panzer_IntegrationValues2.cpp.

Member Data Documentation

template<typename Scalar>
Array_IPDim panzer::IntegrationValues2< Scalar >::cub_points
mutable

Definition at line 135 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_IPDim panzer::IntegrationValues2< Scalar >::side_cub_points
mutable

Definition at line 136 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_IP panzer::IntegrationValues2< Scalar >::cub_weights
mutable

Definition at line 137 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellBASISDim panzer::IntegrationValues2< Scalar >::node_coordinates
mutable

Definition at line 140 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDimDim panzer::IntegrationValues2< Scalar >::jac
mutable

Definition at line 141 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDimDim panzer::IntegrationValues2< Scalar >::jac_inv
mutable

Definition at line 142 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIP panzer::IntegrationValues2< Scalar >::jac_det
mutable

Definition at line 143 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIP panzer::IntegrationValues2< Scalar >::weighted_measure
mutable

Definition at line 144 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDim panzer::IntegrationValues2< Scalar >::weighted_normals
mutable

Definition at line 145 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDim panzer::IntegrationValues2< Scalar >::surface_normals
mutable

Definition at line 146 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDimDim panzer::IntegrationValues2< Scalar >::surface_rotation_matrices
mutable

Definition at line 147 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDimDim panzer::IntegrationValues2< Scalar >::covarient
mutable

Definition at line 152 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDimDim panzer::IntegrationValues2< Scalar >::contravarient
mutable

Definition at line 153 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIP panzer::IntegrationValues2< Scalar >::norm_contravarient
mutable

Definition at line 154 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDim panzer::IntegrationValues2< Scalar >::ip_coordinates
mutable

Definition at line 157 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Array_CellIPDim panzer::IntegrationValues2< Scalar >::ref_ip_coordinates
mutable

Definition at line 158 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Teuchos::RCP<const panzer::IntegrationRule> panzer::IntegrationValues2< Scalar >::int_rule

Definition at line 161 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Teuchos::RCP<Intrepid2::Cubature<PHX::Device::execution_space,double,double> > panzer::IntegrationValues2< Scalar >::intrepid_cubature

Definition at line 163 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
int panzer::IntegrationValues2< Scalar >::num_cells_
protected

Definition at line 473 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
int panzer::IntegrationValues2< Scalar >::num_evaluate_cells_
protected

Definition at line 476 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
int panzer::IntegrationValues2< Scalar >::num_virtual_cells_
protected

Definition at line 479 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::requires_permutation_
protected

Definition at line 482 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
PHX::MDField<const int,Cell,IP> panzer::IntegrationValues2< Scalar >::permutations_
protected

Definition at line 485 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
Teuchos::RCP<const SubcellConnectivity> panzer::IntegrationValues2< Scalar >::side_connectivity_
protected

Definition at line 489 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::cub_points_evaluated_
mutableprotected

Definition at line 492 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::side_cub_points_evaluated_
mutableprotected

Definition at line 493 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::cub_weights_evaluated_
mutableprotected

Definition at line 494 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::node_coordinates_evaluated_
mutableprotected

Definition at line 495 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::jac_evaluated_
mutableprotected

Definition at line 496 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::jac_inv_evaluated_
mutableprotected

Definition at line 497 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::jac_det_evaluated_
mutableprotected

Definition at line 498 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::weighted_measure_evaluated_
mutableprotected

Definition at line 499 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::weighted_normals_evaluated_
mutableprotected

Definition at line 500 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::surface_normals_evaluated_
mutableprotected

Definition at line 501 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::surface_rotation_matrices_evaluated_
mutableprotected

Definition at line 502 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::covarient_evaluated_
mutableprotected

Definition at line 503 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::contravarient_evaluated_
mutableprotected

Definition at line 504 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::norm_contravarient_evaluated_
mutableprotected

Definition at line 505 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::ip_coordinates_evaluated_
mutableprotected

Definition at line 506 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::ref_ip_coordinates_evaluated_
mutableprotected

Definition at line 507 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
bool panzer::IntegrationValues2< Scalar >::alloc_arrays_
private

Definition at line 515 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
std::string panzer::IntegrationValues2< Scalar >::prefix_
private

Definition at line 516 of file Panzer_IntegrationValues2.hpp.

template<typename Scalar>
std::vector<PHX::index_size_type> panzer::IntegrationValues2< Scalar >::ddims_
private

Definition at line 517 of file Panzer_IntegrationValues2.hpp.


The documentation for this class was generated from the following files: