https://mooseframework.inl.gov
FaceInfo.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MooseTypes.h"
13 #include "ElemInfo.h"
14 #include "MooseError.h"
15 
16 #include "libmesh/vector_value.h"
17 #include "libmesh/remote_elem.h"
18 
19 #include <map>
20 #include <set>
21 #include <memory>
22 
23 class MooseMesh;
24 namespace libMesh
25 {
26 class Elem;
27 class Node;
28 }
29 
36 class FaceInfo
37 {
38 public:
39  FaceInfo(const ElemInfo * const elem_info, const unsigned int side, const dof_id_type id);
40 
48  enum class VarFaceNeighbors
49  {
50  BOTH,
51  NEITHER,
52  ELEM,
53  NEIGHBOR
54  };
55 
57  dof_id_type id() const { return _id; }
58 
60  Real faceArea() const { return _face_area; }
61 
64  Real & faceCoord() { return _face_coord; }
65  Real faceCoord() const { return _face_coord; }
66 
68  const Point & normal() const { return _normal; }
69 
71  const Point & faceCentroid() const { return _face_centroid; }
72 
75  Point skewnessCorrectionVector() const;
76 
81  const Elem & elem() const { return *(_elem_info->elem()); }
82  const Elem * elemPtr() const { return _elem_info->elem(); }
83  const Elem & neighbor() const;
84  const Elem * neighborPtr() const { return _neighbor_info ? _neighbor_info->elem() : nullptr; }
85  const ElemInfo * elemInfo() const { return _elem_info; }
86  const ElemInfo * neighborInfo() const { return _neighbor_info; }
88 
95  const Point & elemCentroid() const { return _elem_info->centroid(); }
96  const Point & neighborCentroid() const;
98 
105 
109  unsigned int elemSideID() const { return _elem_side_id; }
110  unsigned int neighborSideID() const { return _neighbor_side_id; }
112 
114  VarFaceNeighbors faceType(const std::pair<unsigned int, unsigned int> & var_sys) const;
117  VarFaceNeighbors & faceType(const std::pair<unsigned int, unsigned int> & var_sys);
118 
120  const std::set<BoundaryID> & boundaryIDs() const { return _boundary_ids; }
121 
123  std::set<BoundaryID> & boundaryIDs() { return _boundary_ids; }
124 
126  Real elemVolume() const { return _elem_info->volume(); }
127 
129  Real neighborVolume() const;
130 
132  Real gC() const { return _gc; }
133 
138  const Point & dCN() const { return _d_cn; }
139 
144  Real dCNMag() const { return _d_cn_mag; }
145 
151  const Point & eCN() const { return _e_cn; }
152 
157 
162  void computeInternalCoefficients(const ElemInfo * const neighbor_info);
163 
169 
170 private:
173  std::vector<std::vector<VarFaceNeighbors>> & faceType() { return _face_types_by_var; }
174 
176  const ElemInfo * const _elem_info;
178 
180 
182  Point _normal;
183 
185 
187  const unsigned int _elem_side_id;
188  unsigned int _neighbor_side_id;
189 
192 
194  Point _d_cn;
195  Point _e_cn;
196 
199 
202 
206  std::vector<std::vector<VarFaceNeighbors>> _face_types_by_var;
207 
209  std::set<BoundaryID> _boundary_ids;
210 
212  friend MooseMesh;
213 };
214 
215 inline const Elem &
217 {
218  mooseAssert(_neighbor_info,
219  "FaceInfo object 'const Elem & neighbor()' is called but neighbor element pointer "
220  "is null. This occurs for faces at the domain boundary");
221  return *(_neighbor_info->elem());
222 }
223 
225 FaceInfo::faceType(const std::pair<unsigned int, unsigned int> & var_sys) const
226 {
227  mooseAssert(var_sys.second < _face_types_by_var.size(), "System number out of bounds!");
228  mooseAssert(var_sys.first < _face_types_by_var[var_sys.second].size(),
229  "Variable number out of bounds!");
230  return _face_types_by_var[var_sys.second][var_sys.first];
231 }
232 
234 FaceInfo::faceType(const std::pair<unsigned int, unsigned int> & var_sys)
235 {
236  mooseAssert(var_sys.second < _face_types_by_var.size(), "System number out of bounds!");
237  mooseAssert(var_sys.first < _face_types_by_var[var_sys.second].size(),
238  "Variable number out of bounds!");
239  return _face_types_by_var[var_sys.second][var_sys.first];
240 }
241 
242 inline const Point &
244 {
245  mooseAssert(_neighbor_info,
246  "The neighbor is not defined on this faceInfo! A possible explanation is that the "
247  "face is a (physical/processor) boundary face.");
248  return _neighbor_info->centroid();
249 }
250 
251 inline SubdomainID
253 {
255 }
256 
257 inline Real
259 {
260  mooseAssert(_neighbor_info,
261  "The neighbor is not defined on this faceInfo! A possible explanation is that the "
262  "face is a (physical/processor) boundary face.");
263  return _neighbor_info->volume();
264 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:177
Point _e_cn
Definition: FaceInfo.h:195
processor_id_type processor_id() const
Definition: FaceInfo.h:156
void computeInternalCoefficients(const ElemInfo *const neighbor_info)
Takes the ElemInfo of the neighbor cell and computes interpolation weights together with other quanti...
Definition: FaceInfo.C:45
const std::set< BoundaryID > & boundaryIDs() const
Const getter for every associated boundary ID.
Definition: FaceInfo.h:120
Point skewnessCorrectionVector() const
Returns the skewness-correction vector (vector between the approximate and real face centroids)...
Definition: FaceInfo.C:80
Real elemVolume() const
Return the element volume.
Definition: FaceInfo.h:126
const Elem & elem() const
Definition: FaceInfo.h:81
const ElemInfo * neighborInfo() const
Definition: FaceInfo.h:86
const Elem * elem() const
Definition: ElemInfo.h:34
const Point & faceCentroid() const
Returns the coordinates of the face centroid.
Definition: FaceInfo.h:71
friend MooseMesh
Allows access to private members from moose mesh only.
Definition: FaceInfo.h:212
const ElemInfo * elemInfo() const
Definition: FaceInfo.h:85
unsigned int elemSideID() const
Definition: FaceInfo.h:109
std::set< BoundaryID > _boundary_ids
the set of boundary ids that this face is associated with
Definition: FaceInfo.h:209
Real & faceCoord()
Sets/gets the coordinate transformation factor (for e.g.
Definition: FaceInfo.h:64
const Point & neighborCentroid() const
Definition: FaceInfo.h:243
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
std::set< BoundaryID > & boundaryIDs()
Returns the set of boundary ids for all boundaries that include this face.
Definition: FaceInfo.h:123
Real faceArea() const
Returns the face area of face id.
Definition: FaceInfo.h:60
unsigned int neighborSideID() const
Definition: FaceInfo.h:110
FaceInfo(const ElemInfo *const elem_info, const unsigned int side, const dof_id_type id)
Definition: FaceInfo.C:19
Real _d_cn_mag
the distance norm between neighbor and element centroids
Definition: FaceInfo.h:198
std::vector< std::vector< VarFaceNeighbors > > & faceType()
Getter for the face type for every stored variable.
Definition: FaceInfo.h:173
Real neighborVolume() const
Return the neighbor volume.
Definition: FaceInfo.h:258
const SubdomainID INVALID_BLOCK_ID
Definition: MooseTypes.C:20
SubdomainID elemSubdomainID() const
Definition: FaceInfo.h:102
uint8_t processor_id_type
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
const Point & centroid() const
Definition: ElemInfo.h:36
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:201
const Elem * neighborPtr() const
Definition: FaceInfo.h:84
const Point & elemCentroid() const
Returns the element centroids of the elements on the elem and neighbor sides of the face...
Definition: FaceInfo.h:95
Real gC() const
Return the geometric weighting factor.
Definition: FaceInfo.h:132
VarFaceNeighbors
This enum is used to indicate which side(s) of a face a particular variable is defined on...
Definition: FaceInfo.h:48
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
const Elem & neighbor() const
Definition: FaceInfo.h:216
const dof_id_type _id
Definition: FaceInfo.h:179
const Point & normal() const
Returns the unit normal vector for the face oriented outward from the face&#39;s elem element...
Definition: FaceInfo.h:68
const processor_id_type _processor_id
Definition: FaceInfo.h:184
void computeBoundaryCoefficients()
Computes the interpolation weights and similar quantities with the assumption that the face is on a b...
Definition: FaceInfo.C:66
Real dCNMag() const
Definition: FaceInfo.h:144
const Elem * elemPtr() const
Definition: FaceInfo.h:82
Real faceCoord() const
Definition: FaceInfo.h:65
const unsigned int _elem_side_id
the elem local side id
Definition: FaceInfo.h:187
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Point & eCN() const
Definition: FaceInfo.h:151
unsigned int _neighbor_side_id
Definition: FaceInfo.h:188
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:176
dof_id_type id() const
Return the ID of the face.
Definition: FaceInfo.h:57
SubdomainID neighborSubdomainID() const
Definition: FaceInfo.h:252
Real _face_coord
Definition: FaceInfo.h:181
std::vector< std::vector< VarFaceNeighbors > > _face_types_by_var
A vector that provides the information about what face type this is for each variable.
Definition: FaceInfo.h:206
Point _d_cn
the distance vector between neighbor and element centroids
Definition: FaceInfo.h:194
Point _face_centroid
Definition: FaceInfo.h:191
Class used for caching additional information for elements such as the volume and centroid...
Definition: ElemInfo.h:25
Point _normal
Definition: FaceInfo.h:182
Real volume() const
Definition: ElemInfo.h:35
SubdomainID subdomain_id() const
We return the subdomain ID of the corresponding libmesh element.
Definition: ElemInfo.h:43
const Point & dCN() const
Definition: FaceInfo.h:138
uint8_t dof_id_type
Real _face_area
Definition: FaceInfo.h:190