https://mooseframework.inl.gov
FaceInfo.C
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 #include "FaceInfo.h"
11 #include "MooseTypes.h"
12 
13 #include "libmesh/elem.h"
14 #include "libmesh/node.h"
15 #include "libmesh/fe_base.h"
16 #include "libmesh/quadrature_gauss.h"
17 #include "libmesh/remote_elem.h"
18 
19 FaceInfo::FaceInfo(const ElemInfo * elem_info, unsigned int side, const dof_id_type id)
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),
25  _neighbor_side_id(libMesh::invalid_uint),
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<libMesh::FEBase> fe(
34  fe->attach_quadrature_rule(&qface);
35  const std::vector<Point> & normals = fe->get_normals();
36  fe->reinit(_elem_info->elem(), _elem_side_id);
37  mooseAssert(normals.size() == 1, "FaceInfo construction broken w.r.t. computing face normals");
38  _normal = normals[0];
39 
40  _face_area = face->volume();
41  _face_centroid = face->vertex_average();
42 }
43 
44 void
45 FaceInfo::computeInternalCoefficients(const ElemInfo * const neighbor_info)
46 {
47  mooseAssert(neighbor_info,
48  "We need a neighbor if we want to compute interpolation coefficients!");
49  _neighbor_info = neighbor_info;
50  _neighbor_side_id = _neighbor_info->elem()->which_neighbor_am_i(_elem_info->elem());
51 
52  // Setup quantities used for the approximation of the spatial derivatives
54  _d_cn_mag = _d_cn.norm();
55  _e_cn = _d_cn / _d_cn_mag;
56 
57  Point r_intersection =
58  _elem_info->centroid() +
60 
61  // For interpolation coefficients
62  _gc = (_neighbor_info->centroid() - r_intersection).norm() / _d_cn_mag;
63 }
64 
65 void
67 {
68  mooseAssert(!_neighbor_info, "This functions shall only be called on a boundary!");
69 
70  // Setup quantities used for the approximation of the spatial derivatives
72  _d_cn_mag = _d_cn.norm();
73  _e_cn = _d_cn / _d_cn_mag;
74 
75  // For interpolation coefficients
76  _gc = 1.0;
77 }
78 
79 Point
81 {
82  const Point r_intersection =
83  _elem_info->centroid() +
85 
86  return _face_centroid - r_intersection;
87 }
const ElemInfo * _neighbor_info
Definition: FaceInfo.h:177
Point _e_cn
Definition: FaceInfo.h:195
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
std::unique_ptr< FEGenericBase< Real > > build(const unsigned int dim, const FEType &fet)
const unsigned int invalid_uint
Point skewnessCorrectionVector() const
Returns the skewness-correction vector (vector between the approximate and real face centroids)...
Definition: FaceInfo.C:80
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:153
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
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
const Point & centroid() const
Definition: ElemInfo.h:36
Real _gc
Geometric weighting factor for face value interpolation.
Definition: FaceInfo.h:201
processor_id_type _processor_id
auto norm(const T &a) -> decltype(std::abs(a))
void computeBoundaryCoefficients()
Computes the interpolation weights and similar quantities with the assumption that the face is on a b...
Definition: FaceInfo.C:66
const unsigned int _elem_side_id
the elem local side id
Definition: FaceInfo.h:187
unsigned int _neighbor_side_id
Definition: FaceInfo.h:188
const ElemInfo *const _elem_info
the elem and neighbor elems
Definition: FaceInfo.h:176
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
uint8_t dof_id_type
Real _face_area
Definition: FaceInfo.h:190