This data structure is used to store geometric and variable related metadata about each cell face in the mesh. More...
#include <FaceInfo.h>
Public Types | |
enum | VarFaceNeighbors { VarFaceNeighbors::BOTH, VarFaceNeighbors::NEITHER, VarFaceNeighbors::ELEM, VarFaceNeighbors::NEIGHBOR } |
This enum is used to indicate which side(s) of a face a particular variable is defined on. More... | |
Public Member Functions | |
FaceInfo (const ElemInfo *const elem_info, const unsigned int side, const dof_id_type id) | |
dof_id_type | id () const |
Return the ID of the face. More... | |
Real | faceArea () const |
Returns the face area of face id. More... | |
Real & | faceCoord () |
Sets/gets the coordinate transformation factor (for e.g. More... | |
Real | faceCoord () const |
const Point & | normal () const |
Returns the unit normal vector for the face oriented outward from the face's elem element. More... | |
const Point & | faceCentroid () const |
Returns the coordinates of the face centroid. More... | |
Point | skewnessCorrectionVector () const |
Returns the skewness-correction vector (vector between the approximate and real face centroids). More... | |
const Point & | elemCentroid () const |
Returns the element centroids of the elements on the elem and neighbor sides of the face. More... | |
const Point & | neighborCentroid () const |
VarFaceNeighbors | faceType (const std::pair< unsigned int, unsigned int > &var_sys) const |
Returns which side(s) the given variable-system number pair is defined on for this face. More... | |
VarFaceNeighbors & | faceType (const std::pair< unsigned int, unsigned int > &var_sys) |
Mutably returns which side(s) the given variable-system number pair is defined on for this face. More... | |
const std::set< BoundaryID > & | boundaryIDs () const |
Const getter for every associated boundary ID. More... | |
std::set< BoundaryID > & | boundaryIDs () |
Returns the set of boundary ids for all boundaries that include this face. More... | |
Real | elemVolume () const |
Return the element volume. More... | |
Real | neighborVolume () const |
Return the neighbor volume. More... | |
Real | gC () const |
Return the geometric weighting factor. More... | |
const Point & | dCN () const |
Real | dCNMag () const |
const Point & | eCN () const |
processor_id_type | processor_id () const |
void | computeInternalCoefficients (const ElemInfo *const neighbor_info) |
Takes the ElemInfo of the neighbor cell and computes interpolation weights together with other quantities used to generate spatial operators. More... | |
void | computeBoundaryCoefficients () |
Computes the interpolation weights and similar quantities with the assumption that the face is on a boundary. More... | |
const Elem & | elem () const |
const Elem * | elemPtr () const |
const Elem & | neighbor () const |
const Elem * | neighborPtr () const |
const ElemInfo * | elemInfo () const |
const ElemInfo * | neighborInfo () const |
SubdomainID | elemSubdomainID () const |
SubdomainID | neighborSubdomainID () const |
unsigned int | elemSideID () const |
unsigned int | neighborSideID () const |
Private Member Functions | |
std::vector< std::vector< VarFaceNeighbors > > & | faceType () |
Getter for the face type for every stored variable. More... | |
Private Attributes | |
const ElemInfo *const | _elem_info |
the elem and neighbor elems More... | |
const ElemInfo * | _neighbor_info |
const dof_id_type | _id |
Real | _face_coord = 0 |
Point | _normal |
const processor_id_type | _processor_id |
const unsigned int | _elem_side_id |
the elem local side id More... | |
unsigned int | _neighbor_side_id |
Real | _face_area |
Point | _face_centroid |
Point | _d_cn |
the distance vector between neighbor and element centroids More... | |
Point | _e_cn |
Real | _d_cn_mag |
the distance norm between neighbor and element centroids More... | |
Real | _gc |
Geometric weighting factor for face value interpolation. More... | |
std::vector< std::vector< VarFaceNeighbors > > | _face_types_by_var |
A vector that provides the information about what face type this is for each variable. More... | |
std::set< BoundaryID > | _boundary_ids |
the set of boundary ids that this face is associated with More... | |
friend | MooseMesh |
Allows access to private members from moose mesh only. More... | |
This data structure is used to store geometric and variable related metadata about each cell face in the mesh.
This info is used by face loops (e.g. for finite volumes method numerical flux loop). These objects can be created and cached up front. Since it only stores information that changes when the mesh is modified it only needs an update whenever the mesh changes.
Definition at line 36 of file FaceInfo.h.
|
strong |
This enum is used to indicate which side(s) of a face a particular variable is defined on.
This is important for certain BC-related finite volume calculations. Because of the way side-sets and variable block-restriction work in MOOSE, there may be boundary conditions applied to internal faces on the mesh where a variable is only active on one or even zero sides of the face. For such faces, FV needs to know which sides (if any) to add BC residual contributions to.
Enumerator | |
---|---|
BOTH | |
NEITHER | |
ELEM | |
NEIGHBOR |
Definition at line 48 of file FaceInfo.h.
FaceInfo::FaceInfo | ( | const ElemInfo *const | elem_info, |
const unsigned int | side, | ||
const dof_id_type | id | ||
) |
Definition at line 19 of file FaceInfo.C.
|
inline |
Const getter for every associated boundary ID.
Definition at line 120 of file FaceInfo.h.
Referenced by LinearFVFluxKernel::addMatrixContribution(), LinearFVFluxKernel::addRightHandSideContribution(), FVFluxKernel::avoidBoundary(), MooseLinearVariableFV< Real >::evaluate(), MooseVariableFV< Real >::getDirichletBC(), MooseVariableFV< Real >::getFluxBCs(), MooseLinearVariableFV< Real >::isDirichletBoundaryFace(), and FVFluxKernel::skipForBoundary().
|
inline |
Returns the set of boundary ids for all boundaries that include this face.
Definition at line 123 of file FaceInfo.h.
void FaceInfo::computeBoundaryCoefficients | ( | ) |
Computes the interpolation weights and similar quantities with the assumption that the face is on a boundary.
Definition at line 66 of file FaceInfo.C.
Takes the ElemInfo of the neighbor cell and computes interpolation weights together with other quantities used to generate spatial operators.
Definition at line 45 of file FaceInfo.C.
|
inline |
Definition at line 138 of file FaceInfo.h.
Referenced by LinearFVBoundaryCondition::computeCellToFaceVector(), LinearFVDiffusion::computeFluxMatrixContribution(), LinearFVAnisotropicDiffusion::computeFluxMatrixContribution(), Moose::FV::fullLimitedInterpolation(), and Moose::FV::interpCoeffs().
|
inline |
Definition at line 144 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVAnisotropicDiffusion::computeFluxMatrixContribution(), LinearFVDiffusion::computeFluxMatrixContribution(), and FVOrthogonalDiffusion::computeQpResidual().
|
inline |
Definition at line 151 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVAnisotropicDiffusion::computeFluxRHSContribution(), and LinearFVDiffusion::computeFluxRHSContribution().
|
inline |
Returns the elem and neighbor elements adjacent to the face. If a face is on a mesh boundary, the neighborPtr will return nullptr - the elem will never be null.
Definition at line 81 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), FVDiffusionInterface::computeQpResidual(), Moose::FV::determineElemOneAndTwo(), FVFluxBC::elemArg(), Function::evaluate(), MooseVariableFE< Real >::faceEvaluate(), MooseLinearVariableFV< Real >::gradSln(), SideIntegralFunctorPostprocessorTempl< false >::hasFaceSide(), FVFluxKernel::hasFaceSide(), Moose::FunctorBase< libMesh::VectorValue >::hasFaceSide(), PiecewiseByBlockLambdaFunctor< T >::isExtrapolatedBoundaryFace(), Moose::FunctorBase< libMesh::VectorValue >::isInternalFace(), Moose::FV::loopOverElemFaceInfo(), Moose::FaceArg::makeElem(), FaceArgProducerInterface::makeFace(), Moose::FV::onBoundary(), Assembly::reinitFVFace(), FVInterfaceKernel::setupData(), FVInterfaceKernel::singleSidedFaceArg(), and MooseVariableFV< Real >::uncorrectedAdGradSln().
|
inline |
Returns the element centroids of the elements on the elem and neighbor sides of the face.
If no neighbor face is defined, a "ghost" neighbor centroid is calculated by reflecting/extrapolating from the elem centroid through the face centroid
Definition at line 95 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVBoundaryCondition::computeCellToFaceVector(), FVDiffusionInterface::computeQpResidual(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), Moose::FV::fullLimitedInterpolation(), MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue(), and Moose::FV::Limiter< T >::rf_minmax().
|
inline |
Definition at line 85 of file FaceInfo.h.
Referenced by LinearFVFluxKernel::addMatrixContribution(), LinearFVFluxKernel::addRightHandSideContribution(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryGradientRHSContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryNormalGradient(), LinearFVAnisotropicDiffusion::computeBoundaryRHSContribution(), LinearFVDiffusion::computeBoundaryRHSContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValue(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryValueRHSContribution(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValueRHSContribution(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValueRHSContribution(), LinearFVAnisotropicDiffusion::computeFluxRHSContribution(), LinearFVDiffusion::computeFluxRHSContribution(), and MooseLinearVariableFV< Real >::gradSln().
|
inline |
Definition at line 82 of file FaceInfo.h.
Referenced by FVFunctorDirichletBCTempl< is_ad >::boundaryValue(), Moose::FunctorBase< libMesh::VectorValue >::checkFace(), LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryNormalGradient(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValue(), SideAdvectiveFluxIntegralTempl< is_ad >::computeFaceInfoIntegral(), FVFluxKernel::elemArg(), FVInterfaceKernel::elemArg(), Function::evaluate(), MooseLinearVariableFV< Real >::evaluate(), MooseVariableFV< Real >::evaluate(), MooseVariableFE< Real >::faceEvaluate(), MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue(), and FVInterfaceKernel::singleSidedFaceArg().
|
inline |
Returns the elem and neighbor centroids. If no neighbor element exists, then the maximum unsigned int is returned for the neighbor side ID.
Definition at line 109 of file FaceInfo.h.
Referenced by MooseVariableFE< Real >::faceEvaluate(), and Assembly::reinitFVFace().
|
inline |
Returns the elem and neighbor subdomain IDs. If no neighbor element exists, then an invalid ID is returned for the neighbor subdomain ID.
Definition at line 102 of file FaceInfo.h.
Referenced by Moose::FV::onBoundary().
|
inline |
Return the element volume.
Definition at line 126 of file FaceInfo.h.
|
inline |
Returns the face area of face id.
Definition at line 60 of file FaceInfo.h.
Referenced by FVScalarLagrangeMultiplierInterface::computeJacobian(), FVBoundaryScalarLagrangeMultiplierConstraint::computeJacobian(), FVFluxBC::computeJacobian(), FVFluxKernel::computeJacobian(), FVInterfaceKernel::computeJacobian(), FVScalarLagrangeMultiplierInterface::computeResidual(), FVBoundaryScalarLagrangeMultiplierConstraint::computeResidual(), FVFluxBC::computeResidual(), FVFluxKernel::computeResidual(), FVInterfaceKernel::computeResidual(), and Moose::FV::loopOverElemFaceInfo().
|
inline |
Returns the coordinates of the face centroid.
Definition at line 71 of file FaceInfo.h.
Referenced by LinearFVFluxKernel::addMatrixContribution(), MooseVariableFV< Real >::adGradSln(), FVFunctionDirichletBC::boundaryValue(), Moose::FunctorBase< libMesh::VectorValue >::checkFace(), LinearFVDiffusion::computeBoundaryRHSContribution(), LinearFVAnisotropicDiffusion::computeBoundaryRHSContribution(), LinearFVBoundaryCondition::computeCellToFaceVector(), MooseMeshUtils::computeFiniteVolumeCoords(), FVFunctionNeumannBC::computeQpResidual(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), FVFluxBC::computeResidual(), Moose::FV::fullLimitedInterpolation(), MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue(), Moose::FaceArg::getPoint(), Moose::FV::loopOverElemFaceInfo(), and Moose::FV::Limiter< T >::rf_minmax().
|
inline |
Sets/gets the coordinate transformation factor (for e.g.
rz, spherical coords) to be used for integration over faces.
Definition at line 64 of file FaceInfo.h.
Referenced by MooseMeshUtils::computeFiniteVolumeCoords(), FVScalarLagrangeMultiplierInterface::computeJacobian(), FVBoundaryScalarLagrangeMultiplierConstraint::computeJacobian(), FVFluxBC::computeJacobian(), FVFluxKernel::computeJacobian(), FVInterfaceKernel::computeJacobian(), FVScalarLagrangeMultiplierInterface::computeResidual(), FVBoundaryScalarLagrangeMultiplierConstraint::computeResidual(), FVFluxBC::computeResidual(), FVFluxKernel::computeResidual(), and FVInterfaceKernel::computeResidual().
|
inline |
Definition at line 65 of file FaceInfo.h.
|
inline |
Returns which side(s) the given variable-system number pair is defined on for this face.
Definition at line 225 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::computeFaceValues(), FVBoundaryScalarLagrangeMultiplierConstraint::computeJacobian(), FVFluxBC::computeJacobian(), FVFluxKernel::computeJacobian(), FVBoundaryScalarLagrangeMultiplierConstraint::computeResidual(), FVFluxKernel::computeResidual(), Moose::FV::determineElemOneAndTwo(), MooseLinearVariableFV< Real >::evaluate(), LinearFVFluxKernel::hasFaceSide(), LinearFVBoundaryCondition::hasFaceSide(), FVBoundaryCondition::hasFaceSide(), FVInterfaceKernel::setupData(), LinearFVFluxKernel::setupFaceData(), FVFluxBC::uOnGhost(), FVFluxBC::uOnUSub(), and FVFluxBC::updateCurrentFace().
|
inline |
Mutably returns which side(s) the given variable-system number pair is defined on for this face.
Definition at line 234 of file FaceInfo.h.
|
inlineprivate |
Getter for the face type for every stored variable.
This will be a friend of MooseMesh to make sure we can only access it from there.
Definition at line 173 of file FaceInfo.h.
|
inline |
Return the geometric weighting factor.
Definition at line 132 of file FaceInfo.h.
Referenced by Moose::FV::interpCoeffs().
|
inline |
|
inline |
Definition at line 216 of file FaceInfo.h.
Referenced by Function::evaluate(), FVFluxKernel::hasFaceSide(), Moose::FunctorBase< libMesh::VectorValue >::hasFaceSide(), PiecewiseByBlockLambdaFunctor< T >::isExtrapolatedBoundaryFace(), Moose::FV::onBoundary(), and FVInterfaceKernel::singleSidedFaceArg().
|
inline |
Definition at line 243 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVBoundaryCondition::computeCellToFaceVector(), FVDiffusionInterface::computeQpResidual(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), Moose::FV::fullLimitedInterpolation(), MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue(), and Moose::FV::Limiter< T >::rf_minmax().
|
inline |
Definition at line 86 of file FaceInfo.h.
Referenced by LinearFVFluxKernel::addMatrixContribution(), LinearFVFluxKernel::addRightHandSideContribution(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryGradientRHSContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryNormalGradient(), LinearFVAnisotropicDiffusion::computeBoundaryRHSContribution(), LinearFVDiffusion::computeBoundaryRHSContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValueRHSContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryValueRHSContribution(), LinearFVDiffusion::computeFluxRHSContribution(), LinearFVAnisotropicDiffusion::computeFluxRHSContribution(), and MooseLinearVariableFV< Real >::gradSln().
|
inline |
Definition at line 84 of file FaceInfo.h.
Referenced by MooseVariableFV< Real >::adGradSln(), FVFunctorDirichletBCTempl< is_ad >::boundaryValue(), Moose::FunctorBase< libMesh::VectorValue >::checkFace(), LinearFVAdvectionDiffusionFunctorDirichletBC::computeBoundaryNormalGradient(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValue(), LinearFVBoundaryCondition::computeCellToFaceVector(), FVDiffusionInterface::computeQpResidual(), Moose::FV::determineElemOneAndTwo(), Function::evaluate(), MooseLinearVariableFV< Real >::evaluate(), MooseVariableFV< Real >::evaluate(), MooseVariableFE< Real >::faceEvaluate(), SideIntegralFunctorPostprocessorTempl< false >::hasFaceSide(), FVFluxKernel::hasFaceSide(), Moose::FunctorBase< libMesh::VectorValue >::hasFaceSide(), PiecewiseByBlockLambdaFunctor< T >::isExtrapolatedBoundaryFace(), Moose::FunctorBase< libMesh::VectorValue >::isInternalFace(), Moose::FV::loopOverElemFaceInfo(), FaceArgProducerInterface::makeFace(), Moose::FaceArg::makeNeighbor(), FVFluxBC::neighborArg(), FVFluxKernel::neighborArg(), FVInterfaceKernel::neighborArg(), Moose::FV::onBoundary(), Assembly::reinitFVFace(), ComputeFVFluxThread< RangeType, AttribMatrixTags >::reinitVariables(), FVInterfaceKernel::singleSidedFaceArg(), and MooseVariableFV< Real >::uncorrectedAdGradSln().
|
inline |
Definition at line 110 of file FaceInfo.h.
Referenced by MooseVariableFE< Real >::faceEvaluate(), and Assembly::reinitFVFace().
|
inline |
Definition at line 252 of file FaceInfo.h.
Referenced by Moose::FV::onBoundary().
|
inline |
|
inline |
Returns the unit normal vector for the face oriented outward from the face's elem element.
Definition at line 68 of file FaceInfo.h.
Referenced by LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryGradientMatrixContribution(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryGradientRHSContribution(), LinearFVAdvection::computeBoundaryMatrixContribution(), LinearFVAnisotropicDiffusion::computeBoundaryMatrixContribution(), LinearFVAdvectionDiffusionExtrapolatedBC::computeBoundaryNormalGradient(), LinearFVAdvection::computeBoundaryRHSContribution(), LinearFVAnisotropicDiffusion::computeBoundaryRHSContribution(), LinearFVDiffusion::computeBoundaryRHSContribution(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValue(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValueMatrixContribution(), LinearFVAdvectionDiffusionFunctorRobinBC::computeBoundaryValueRHSContribution(), LinearFVAdvectionDiffusionFunctorNeumannBC::computeBoundaryValueRHSContribution(), LinearFVBoundaryCondition::computeCellToFaceDistance(), LinearFVAdvection::computeElemMatrixContribution(), SideAdvectiveFluxIntegralTempl< is_ad >::computeFaceInfoIntegral(), LinearFVAnisotropicDiffusion::computeFluxMatrixContribution(), LinearFVDiffusion::computeFluxMatrixContribution(), LinearFVDiffusion::computeFluxRHSContribution(), LinearFVAnisotropicDiffusion::computeFluxRHSContribution(), FVBoundaryScalarLagrangeMultiplierConstraint::computeJacobian(), FVFluxBC::computeJacobian(), FVFluxKernel::computeJacobian(), LinearFVAdvection::computeNeighborMatrixContribution(), FVBoundaryScalarLagrangeMultiplierConstraint::computeResidual(), FVFluxKernel::computeResidual(), Function::evaluate(), Moose::FV::gradUDotNormal(), Moose::FV::interpolate(), Moose::FV::loopOverElemFaceInfo(), FVInterfaceKernel::setupData(), and FVFluxBC::updateCurrentFace().
|
inline |
Point FaceInfo::skewnessCorrectionVector | ( | ) | const |
Returns the skewness-correction vector (vector between the approximate and real face centroids).
Definition at line 80 of file FaceInfo.C.
Referenced by Moose::FV::skewCorrectedLinearInterpolation().
|
private |
the set of boundary ids that this face is associated with
Definition at line 209 of file FaceInfo.h.
Referenced by boundaryIDs().
|
private |
the distance vector between neighbor and element centroids
Definition at line 194 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and dCN().
|
private |
the distance norm between neighbor and element centroids
Definition at line 198 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and dCNMag().
|
private |
Definition at line 195 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), eCN(), and skewnessCorrectionVector().
|
private |
the elem and neighbor elems
Definition at line 176 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), elem(), elemCentroid(), elemInfo(), elemPtr(), elemSubdomainID(), elemVolume(), FaceInfo(), and skewnessCorrectionVector().
|
private |
the elem local side id
Definition at line 187 of file FaceInfo.h.
Referenced by elemSideID(), and FaceInfo().
|
private |
Definition at line 190 of file FaceInfo.h.
Referenced by faceArea(), and FaceInfo().
|
private |
Definition at line 191 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), faceCentroid(), FaceInfo(), and skewnessCorrectionVector().
|
private |
Definition at line 181 of file FaceInfo.h.
Referenced by faceCoord().
|
private |
A vector that provides the information about what face type this is for each variable.
The first index is the system number; the second index of the key is the variable number within the system.
Definition at line 206 of file FaceInfo.h.
Referenced by faceType().
|
private |
Geometric weighting factor for face value interpolation.
Definition at line 201 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and gC().
|
private |
Definition at line 179 of file FaceInfo.h.
Referenced by id().
|
private |
Definition at line 177 of file FaceInfo.h.
Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), neighbor(), neighborCentroid(), neighborInfo(), neighborPtr(), neighborSubdomainID(), and neighborVolume().
|
private |
Definition at line 188 of file FaceInfo.h.
Referenced by computeInternalCoefficients(), and neighborSideID().
|
private |
Definition at line 182 of file FaceInfo.h.
Referenced by computeInternalCoefficients(), FaceInfo(), normal(), and skewnessCorrectionVector().
|
private |
Definition at line 184 of file FaceInfo.h.
Referenced by processor_id().
|
private |
Allows access to private members from moose mesh only.
Definition at line 212 of file FaceInfo.h.