www.mooseframework.org
Public Types | Public Member Functions | Private Attributes | List of all members
FaceInfo Class Reference

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...
 
RealfaceCoord ()
 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...
 
VarFaceNeighborsfaceType (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
 
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 ElemInfoelemInfo () const
 
const ElemInfoneighborInfo () const
 
SubdomainID elemSubdomainID () const
 
SubdomainID neighborSubdomainID () const
 
unsigned int elemSideID () const
 
unsigned int neighborSideID () const
 

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::map< std::pair< unsigned int, unsigned int >, VarFaceNeighbors_face_types_by_var
 a map 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...
 

Detailed Description

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 35 of file FaceInfo.h.

Member Enumeration Documentation

◆ VarFaceNeighbors

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 47 of file FaceInfo.h.

48  {
49  BOTH,
50  NEITHER,
51  ELEM,
52  NEIGHBOR
53  };
BOTH

Constructor & Destructor Documentation

◆ FaceInfo()

FaceInfo::FaceInfo ( const ElemInfo *const  elem_info,
const unsigned int  side,
const dof_id_type  id 
)

Definition at line 19 of file FaceInfo.C.

20  : _elem_info(elem_info),
21  _neighbor_info(nullptr),
22  _id(id),
23  _processor_id(_elem_info->elem()->processor_id()),
24  _elem_side_id(side),
26  _gc(0.5)
27 {
28  // Compute face-related quantities
29  unsigned int dim = _elem_info->elem()->dim();
30  const std::unique_ptr<const Elem> face = _elem_info->elem()->build_side_ptr(_elem_side_id);
31  std::unique_ptr<FEBase> fe(FEBase::build(dim, FEType(_elem_info->elem()->default_order())));
32  QGauss qface(dim - 1, CONSTANT);
33  fe->attach_quadrature_rule(&qface);
34  const std::vector<Point> & normals = fe->get_normals();
35  fe->reinit(_elem_info->elem(), _elem_side_id);
36  mooseAssert(normals.size() == 1, "FaceInfo construction broken w.r.t. computing face normals");
37  _normal = normals[0];
38 
39  _face_area = face->volume();
40  _face_centroid = face->vertex_average();
41 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
const unsigned int invalid_uint
const Elem * elem() const
Definition: ElemInfo.h:34
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:148
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:195
const dof_id_type _id
Definition: FaceInfo.h:173
const processor_id_type _processor_id
Definition: FaceInfo.h:178
const unsigned int _elem_side_id
the elem local side id
Definition: FaceInfo.h:181
unsigned int _neighbor_side_id
Definition: FaceInfo.h:182
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
Point _face_centroid
Definition: FaceInfo.h:185
Point _normal
Definition: FaceInfo.h:176
Real _face_area
Definition: FaceInfo.h:184

Member Function Documentation

◆ boundaryIDs() [1/2]

const std::set<BoundaryID>& FaceInfo::boundaryIDs ( ) const
inline

◆ boundaryIDs() [2/2]

std::set<BoundaryID>& FaceInfo::boundaryIDs ( )
inline

Returns the set of boundary ids for all boundaries that include this face.

Definition at line 121 of file FaceInfo.h.

121 { return _boundary_ids; }
std::set< BoundaryID > _boundary_ids
the set of boundary ids that this face is associated with
Definition: FaceInfo.h:203

◆ computeBoundaryCoefficients()

void FaceInfo::computeBoundaryCoefficients ( )

Computes the interpolation weights and similar quantities with the assumption that the face is on a boundary.

Definition at line 65 of file FaceInfo.C.

66 {
67  mooseAssert(!_neighbor_info, "This functions shall only be called on a boundary!");
68 
69  // Setup quantities used for the approximation of the spatial derivatives
71  _d_cn_mag = _d_cn.norm();
72  _e_cn = _d_cn / _d_cn_mag;
73 
74  // For interpolation coefficients
75  _gc = 1.0;
76 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
Point _e_cn
Definition: FaceInfo.h:189
Real _d_cn_mag
the distance norm between neighbor and element centroids
Definition: FaceInfo.h:192
const Point & centroid() const
Definition: ElemInfo.h:36
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:195
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
Point _d_cn
the distance vector between neighbor and element centroids
Definition: FaceInfo.h:188
Point _face_centroid
Definition: FaceInfo.h:185

◆ computeInternalCoefficients()

void FaceInfo::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.

Definition at line 44 of file FaceInfo.C.

45 {
46  mooseAssert(neighbor_info,
47  "We need a neighbor if we want to compute interpolation coefficients!");
48  _neighbor_info = neighbor_info;
49  _neighbor_side_id = _neighbor_info->elem()->which_neighbor_am_i(_elem_info->elem());
50 
51  // Setup quantities used for the approximation of the spatial derivatives
53  _d_cn_mag = _d_cn.norm();
54  _e_cn = _d_cn / _d_cn_mag;
55 
56  Point r_intersection =
57  _elem_info->centroid() +
59 
60  // For interpolation coefficients
61  _gc = (_neighbor_info->centroid() - r_intersection).norm() / _d_cn_mag;
62 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
Point _e_cn
Definition: FaceInfo.h:189
const Elem * elem() const
Definition: ElemInfo.h:34
Real _d_cn_mag
the distance norm between neighbor and element centroids
Definition: FaceInfo.h:192
const Point & centroid() const
Definition: ElemInfo.h:36
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:195
auto norm(const T &a) -> decltype(std::abs(a))
unsigned int _neighbor_side_id
Definition: FaceInfo.h:182
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
Point _d_cn
the distance vector between neighbor and element centroids
Definition: FaceInfo.h:188
Point _face_centroid
Definition: FaceInfo.h:185
Point _normal
Definition: FaceInfo.h:176

◆ dCN()

const Point& FaceInfo::dCN ( ) const
inline
Returns
the distance vector drawn from centroid C to N, or in terms of MOOSE implementation, the distance vector obtained from subtracting the element centroid from the neighbor centroid

Definition at line 136 of file FaceInfo.h.

Referenced by LinearFVBoundaryCondition::computeCellToFaceVector(), LinearFVDiffusion::computeFluxMatrixContribution(), and Moose::FV::interpCoeffs().

136 { return _d_cn; }
Point _d_cn
the distance vector between neighbor and element centroids
Definition: FaceInfo.h:188

◆ dCNMag()

Real FaceInfo::dCNMag ( ) const
inline
Returns
the magnitude of the distance vector between centroids C and N, or in terms of MOOSE implementation, the magnitude of the distance vector between neighbor and element centroids

Definition at line 142 of file FaceInfo.h.

Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVDiffusion::computeFluxMatrixContribution(), and FVOrthogonalDiffusion::computeQpResidual().

142 { return _d_cn_mag; }
Real _d_cn_mag
the distance norm between neighbor and element centroids
Definition: FaceInfo.h:192

◆ eCN()

const Point& FaceInfo::eCN ( ) const
inline
Returns
the normalized (e.g. unit) distance vector drawn from centroid C to N, or in terms of MOOSE implementation, the normalized (e.g. unit) distance vector obtained from subtracting the element centroid from the neighbor centroid

Definition at line 149 of file FaceInfo.h.

Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVDiffusion::computeBoundaryRHSContribution(), and LinearFVDiffusion::computeFluxRHSContribution().

149 { return _e_cn; }
Point _e_cn
Definition: FaceInfo.h:189

◆ elem()

const Elem& FaceInfo::elem ( ) const
inline

◆ elemCentroid()

const Point& FaceInfo::elemCentroid ( ) const
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

  • i.e. the vector from the elem element centroid to the face centroid is doubled in length. The tip of this new vector is the neighbor centroid. This is important for FV dirichlet BCs.

Definition at line 94 of file FaceInfo.h.

Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVBoundaryCondition::computeCellToFaceVector(), FVDiffusionInterface::computeQpResidual(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), and MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue().

94 { return _elem_info->centroid(); }
const Point & centroid() const
Definition: ElemInfo.h:36
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170

◆ elemInfo()

const ElemInfo* FaceInfo::elemInfo ( ) const
inline

◆ elemPtr()

const Elem* FaceInfo::elemPtr ( ) const
inline

◆ elemSideID()

unsigned int FaceInfo::elemSideID ( ) const
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 108 of file FaceInfo.h.

Referenced by MooseVariableFE< Real >::faceEvaluate(), and Assembly::reinitFVFace().

108 { return _elem_side_id; }
const unsigned int _elem_side_id
the elem local side id
Definition: FaceInfo.h:181

◆ elemSubdomainID()

SubdomainID FaceInfo::elemSubdomainID ( ) const
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 101 of file FaceInfo.h.

Referenced by Moose::FV::onBoundary().

101 { return _elem_info->subdomain_id(); }
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
SubdomainID subdomain_id() const
We return the subdomain ID of the corresponding libmesh element.
Definition: ElemInfo.h:43

◆ elemVolume()

Real FaceInfo::elemVolume ( ) const
inline

Return the element volume.

Definition at line 124 of file FaceInfo.h.

124 { return _elem_info->volume(); }
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
Real volume() const
Definition: ElemInfo.h:35

◆ faceArea()

Real FaceInfo::faceArea ( ) const
inline

◆ faceCentroid()

const Point& FaceInfo::faceCentroid ( ) const
inline

◆ faceCoord() [1/2]

Real& FaceInfo::faceCoord ( )
inline

◆ faceCoord() [2/2]

Real FaceInfo::faceCoord ( ) const
inline

Definition at line 64 of file FaceInfo.h.

64 { return _face_coord; }
Real _face_coord
Definition: FaceInfo.h:175

◆ faceType() [1/2]

FaceInfo::VarFaceNeighbors FaceInfo::faceType ( const std::pair< unsigned int, unsigned int > &  var_sys) const
inline

Returns which side(s) the given variable-system number pair is defined on for this face.

Definition at line 216 of file FaceInfo.h.

Referenced by MooseVariableFV< Real >::computeFaceValues(), FVBoundaryScalarLagrangeMultiplierConstraint::computeJacobian(), FVFluxBC::computeJacobian(), FVFluxKernel::computeJacobian(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), FVBoundaryScalarLagrangeMultiplierConstraint::computeResidual(), FVFluxBC::computeResidual(), FVFluxKernel::computeResidual(), Moose::FV::determineElemOneAndTwo(), MooseLinearVariableFV< Real >::evaluate(), LinearFVFluxKernel::hasFaceSide(), LinearFVBoundaryCondition::hasFaceSide(), FVBoundaryCondition::hasFaceSide(), ComputeLinearFVGreenGaussGradientFaceThread::onBoundaryFace(), LinearFVFluxKernel::setCurrentFaceInfo(), FVInterfaceKernel::setupData(), FVFluxBC::uOnGhost(), and FVFluxBC::uOnUSub().

217 {
218  auto it = _face_types_by_var.find(var_sys);
219  if (it == _face_types_by_var.end())
220  mooseError("Variable number ",
221  var_sys.first,
222  " in system number ",
223  var_sys.second,
224  " not found in variable to VarFaceNeighbors map");
225  return it->second;
226 }
std::map< std::pair< unsigned int, unsigned int >, VarFaceNeighbors > _face_types_by_var
a map that provides the information about what face type this is for each variable.
Definition: FaceInfo.h:200
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299

◆ faceType() [2/2]

FaceInfo::VarFaceNeighbors & FaceInfo::faceType ( const std::pair< unsigned int, unsigned int > &  var_sys)
inline

Mutably returns which side(s) the given variable-system number pair is defined on for this face.

Definition at line 229 of file FaceInfo.h.

230 {
231  return _face_types_by_var[var_sys];
232 }
std::map< std::pair< unsigned int, unsigned int >, VarFaceNeighbors > _face_types_by_var
a map that provides the information about what face type this is for each variable.
Definition: FaceInfo.h:200

◆ gC()

Real FaceInfo::gC ( ) const
inline

Return the geometric weighting factor.

Definition at line 130 of file FaceInfo.h.

Referenced by Moose::FV::interpCoeffs().

130 { return _gc; }
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:195

◆ id()

dof_id_type FaceInfo::id ( ) const
inline

Return the ID of the face.

Definition at line 56 of file FaceInfo.h.

56 { return _id; }
const dof_id_type _id
Definition: FaceInfo.h:173

◆ neighbor()

const Elem & FaceInfo::neighbor ( ) const
inline

Definition at line 207 of file FaceInfo.h.

Referenced by Function::evaluate(), FVFluxKernel::hasFaceSide(), Moose::FunctorBase< libMesh::VectorValue >::hasFaceSide(), PiecewiseByBlockLambdaFunctor< T >::isExtrapolatedBoundaryFace(), and Moose::FV::onBoundary().

208 {
209  mooseAssert(_neighbor_info,
210  "FaceInfo object 'const Elem & neighbor()' is called but neighbor element pointer "
211  "is null. This occurs for faces at the domain boundary");
212  return *(_neighbor_info->elem());
213 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
const Elem * elem() const
Definition: ElemInfo.h:34

◆ neighborCentroid()

const Point & FaceInfo::neighborCentroid ( ) const
inline

Definition at line 235 of file FaceInfo.h.

Referenced by MooseVariableFV< Real >::adGradSln(), LinearFVBoundaryCondition::computeCellToFaceVector(), FVDiffusionInterface::computeQpResidual(), FVOrthogonalBoundaryDiffusion::computeQpResidual(), and MooseVariableFV< Real >::getExtrapolatedBoundaryFaceValue().

236 {
237  mooseAssert(_neighbor_info,
238  "The neighbor is not defined on this faceInfo! A possible explanation is that the "
239  "face is a (physical/processor) boundary face.");
240  return _neighbor_info->centroid();
241 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
const Point & centroid() const
Definition: ElemInfo.h:36

◆ neighborInfo()

const ElemInfo* FaceInfo::neighborInfo ( ) const
inline

◆ neighborPtr()

const Elem* FaceInfo::neighborPtr ( ) const
inline

◆ neighborSideID()

unsigned int FaceInfo::neighborSideID ( ) const
inline

Definition at line 109 of file FaceInfo.h.

Referenced by MooseVariableFE< Real >::faceEvaluate(), and Assembly::reinitFVFace().

109 { return _neighbor_side_id; }
unsigned int _neighbor_side_id
Definition: FaceInfo.h:182

◆ neighborSubdomainID()

SubdomainID FaceInfo::neighborSubdomainID ( ) const
inline

Definition at line 244 of file FaceInfo.h.

Referenced by Moose::FV::onBoundary().

245 {
247 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
const SubdomainID INVALID_BLOCK_ID
Definition: MooseTypes.C:22
SubdomainID subdomain_id() const
We return the subdomain ID of the corresponding libmesh element.
Definition: ElemInfo.h:43

◆ neighborVolume()

Real FaceInfo::neighborVolume ( ) const
inline

Return the neighbor volume.

Definition at line 250 of file FaceInfo.h.

251 {
252  mooseAssert(_neighbor_info,
253  "The neighbor is not defined on this faceInfo! A possible explanation is that the "
254  "face is a (physical/processor) boundary face.");
255  return _neighbor_info->volume();
256 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:171
Real volume() const
Definition: ElemInfo.h:35

◆ normal()

const Point& FaceInfo::normal ( ) const
inline

◆ processor_id()

processor_id_type FaceInfo::processor_id ( ) const
inline
Returns
the ID of the processor that owns this object

Definition at line 154 of file FaceInfo.h.

154 { return _processor_id; }
const processor_id_type _processor_id
Definition: FaceInfo.h:178

◆ skewnessCorrectionVector()

Point FaceInfo::skewnessCorrectionVector ( ) const

Returns the skewness-correction vector (vector between the approximate and real face centroids).

Definition at line 79 of file FaceInfo.C.

Referenced by Moose::FV::skewCorrectedLinearInterpolation().

80 {
81  const Point r_intersection =
82  _elem_info->centroid() +
84 
85  return _face_centroid - r_intersection;
86 }
Point _e_cn
Definition: FaceInfo.h:189
const Point & centroid() const
Definition: ElemInfo.h:36
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:170
Point _face_centroid
Definition: FaceInfo.h:185
Point _normal
Definition: FaceInfo.h:176

Member Data Documentation

◆ _boundary_ids

std::set<BoundaryID> FaceInfo::_boundary_ids
private

the set of boundary ids that this face is associated with

Definition at line 203 of file FaceInfo.h.

Referenced by boundaryIDs().

◆ _d_cn

Point FaceInfo::_d_cn
private

the distance vector between neighbor and element centroids

Definition at line 188 of file FaceInfo.h.

Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and dCN().

◆ _d_cn_mag

Real FaceInfo::_d_cn_mag
private

the distance norm between neighbor and element centroids

Definition at line 192 of file FaceInfo.h.

Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and dCNMag().

◆ _e_cn

Point FaceInfo::_e_cn
private

◆ _elem_info

const ElemInfo* const FaceInfo::_elem_info
private

◆ _elem_side_id

const unsigned int FaceInfo::_elem_side_id
private

the elem local side id

Definition at line 181 of file FaceInfo.h.

Referenced by elemSideID(), and FaceInfo().

◆ _face_area

Real FaceInfo::_face_area
private

Definition at line 184 of file FaceInfo.h.

Referenced by faceArea(), and FaceInfo().

◆ _face_centroid

Point FaceInfo::_face_centroid
private

◆ _face_coord

Real FaceInfo::_face_coord = 0
private

Definition at line 175 of file FaceInfo.h.

Referenced by faceCoord().

◆ _face_types_by_var

std::map<std::pair<unsigned int, unsigned int>, VarFaceNeighbors> FaceInfo::_face_types_by_var
private

a map that provides the information about what face type this is for each variable.

The first member of the key is the variable number; the second member of the key is the number of the system that the variable lives in

Definition at line 200 of file FaceInfo.h.

Referenced by faceType().

◆ _gc

Real FaceInfo::_gc
private

Geometric weighting factor for face value interpolation.

Definition at line 195 of file FaceInfo.h.

Referenced by computeBoundaryCoefficients(), computeInternalCoefficients(), and gC().

◆ _id

const dof_id_type FaceInfo::_id
private

Definition at line 173 of file FaceInfo.h.

Referenced by id().

◆ _neighbor_info

const ElemInfo* FaceInfo::_neighbor_info
private

◆ _neighbor_side_id

unsigned int FaceInfo::_neighbor_side_id
private

Definition at line 182 of file FaceInfo.h.

Referenced by computeInternalCoefficients(), and neighborSideID().

◆ _normal

Point FaceInfo::_normal
private

◆ _processor_id

const processor_id_type FaceInfo::_processor_id
private

Definition at line 178 of file FaceInfo.h.

Referenced by processor_id().


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