https://mooseframework.inl.gov
FVUtils.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 "MooseMesh.h"
13 #include "MooseError.h"
14 #include "MooseMeshUtils.h"
15 #include "FaceInfo.h"
16 #include "libmesh/elem.h"
17 
18 #include <tuple>
19 
20 template <typename>
21 class MooseVariableFV;
22 
23 template <typename>
25 
26 namespace Moose
27 {
28 namespace FV
29 {
30 
38 bool elemHasFaceInfo(const Elem & elem, const Elem * const neighbor);
39 
40 template <typename ActionFunctor>
41 void
42 loopOverElemFaceInfo(const Elem & elem,
43  const MooseMesh & mesh,
44  ActionFunctor & act,
45  const Moose::CoordinateSystemType coord_type,
46  const unsigned int rz_radial_coord = libMesh::invalid_uint)
47 {
48  mooseAssert(elem.active(), "We should never call this method with an inactive element");
49 
50  for (const auto side : elem.side_index_range())
51  {
52  const Elem * const candidate_neighbor = elem.neighbor_ptr(side);
53 
54  bool elem_has_info = elemHasFaceInfo(elem, candidate_neighbor);
55 
56  std::set<const Elem *> neighbors;
57 
58  const bool inactive_neighbor_detected =
59  candidate_neighbor ? !candidate_neighbor->active() : false;
60 
61  // See MooseMesh::buildFaceInfo for corresponding checks/additions of FaceInfo
62  if (inactive_neighbor_detected)
63  {
64  // We must be next to an element that has been refined
65  mooseAssert(candidate_neighbor->has_children(), "We should have children");
66 
67  const auto candidate_neighbor_side = candidate_neighbor->which_neighbor_am_i(&elem);
68 
69  for (const auto child_num : make_range(candidate_neighbor->n_children()))
70  if (candidate_neighbor->is_child_on_side(child_num, candidate_neighbor_side))
71  {
72  const Elem * const child = candidate_neighbor->child_ptr(child_num);
73  mooseAssert(child->level() - elem.level() == 1, "The math doesn't work out here.");
74  mooseAssert(child->has_neighbor(&elem), "Elem should be a neighbor of this child.");
75  mooseAssert(child->active(),
76  "We shouldn't have greater than a face mismatch level of one");
77  neighbors.insert(child);
78  }
79  }
80  else
81  neighbors.insert(candidate_neighbor);
82 
83  for (const Elem * const neighbor : neighbors)
84  {
85  const FaceInfo * const fi =
86  elem_has_info ? mesh.faceInfo(&elem, side)
87  : mesh.faceInfo(neighbor, neighbor->which_neighbor_am_i(&elem));
88 
89  mooseAssert(fi, "We should have found a FaceInfo");
90  mooseAssert(elem_has_info ? &elem == &fi->elem() : &elem == fi->neighborPtr(),
91  "Doesn't seem like we understand how this FaceInfo thing is working");
92  if (neighbor)
93  {
94  mooseAssert(neighbor != libMesh::remote_elem,
95  "Remote element detected. This indicates that you have insufficient geometric "
96  "ghosting. Please contact your application developers.");
97  mooseAssert(elem_has_info ? neighbor == fi->neighborPtr() : neighbor == &fi->elem(),
98  "Doesn't seem like we understand how this FaceInfo thing is working");
99  }
100 
101  const Point elem_normal = elem_has_info ? fi->normal() : Point(-fi->normal());
102 
103  Real coord;
104  MooseMeshUtils::coordTransformFactor(fi->faceCentroid(), coord, coord_type, rz_radial_coord);
105 
106  const Point surface_vector = elem_normal * fi->faceArea() * coord;
107 
108  act(elem, neighbor, fi, surface_vector, coord, elem_has_info);
109  }
110  }
111 }
112 
128 template <typename FVVar>
129 std::tuple<const Elem *, const Elem *, bool>
130 determineElemOneAndTwo(const FaceInfo & fi, const FVVar & var)
131 {
132  auto ft = fi.faceType(std::make_pair(var.number(), var.sys().number()));
133  mooseAssert(ft == FaceInfo::VarFaceNeighbors::BOTH
134  ? var.hasBlocks(fi.elem().subdomain_id()) && fi.neighborPtr() &&
135  var.hasBlocks(fi.neighborPtr()->subdomain_id())
136  : true,
137  "Finite volume variable " << var.name()
138  << " does not exist on both sides of the face despite "
139  "what the FaceInfo is telling us.");
140  mooseAssert(ft == FaceInfo::VarFaceNeighbors::ELEM
141  ? var.hasBlocks(fi.elem().subdomain_id()) &&
142  (!fi.neighborPtr() || !var.hasBlocks(fi.neighborPtr()->subdomain_id()))
143  : true,
144  "Finite volume variable " << var.name()
145  << " does not exist on or only on the elem side of the "
146  "face despite what the FaceInfo is telling us.");
147  mooseAssert(ft == FaceInfo::VarFaceNeighbors::NEIGHBOR
148  ? fi.neighborPtr() && var.hasBlocks(fi.neighborPtr()->subdomain_id()) &&
149  !var.hasBlocks(fi.elem().subdomain_id())
150  : true,
151  "Finite volume variable " << var.name()
152  << " does not exist on or only on the neighbor side of the "
153  "face despite what the FaceInfo is telling us.");
154 
155  const bool one_is_elem =
157  const Elem * const elem_one = one_is_elem ? &fi.elem() : fi.neighborPtr();
158  mooseAssert(elem_one, "This elem should be non-null!");
159  const Elem * const elem_two = one_is_elem ? fi.neighborPtr() : &fi.elem();
160 
161  return std::make_tuple(elem_one, elem_two, one_is_elem);
162 }
163 }
164 }
void loopOverElemFaceInfo(const Elem &elem, const MooseMesh &mesh, ActionFunctor &act, const Moose::CoordinateSystemType coord_type, const unsigned int rz_radial_coord=libMesh::invalid_uint)
Definition: FVUtils.h:42
bool elemHasFaceInfo(const Elem &elem, const Elem *const neighbor)
This function infers based on elements if the faceinfo between them belongs to the element or not...
Definition: FVUtils.C:21
const unsigned int invalid_uint
void coordTransformFactor(const P &point, C &factor, const Moose::CoordinateSystemType coord_type, const unsigned int rz_radial_coord=libMesh::invalid_uint)
compute a coordinate transformation factor
std::tuple< const Elem *, const Elem *, bool > determineElemOneAndTwo(const FaceInfo &fi, const FVVar &var)
This utility determines element one and element two given a FaceInfo fi and variable var...
Definition: FVUtils.h:130
const Elem & elem() const
Definition: FaceInfo.h:81
const Point & faceCentroid() const
Returns the coordinates of the face centroid.
Definition: FaceInfo.h:71
MeshBase & mesh
Real faceArea() const
Returns the face area of face id.
Definition: FaceInfo.h:60
This data structure is used to store geometric and variable related metadata about each cell face in ...
Definition: FaceInfo.h:36
const Elem * neighborPtr() const
Definition: FaceInfo.h:84
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
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
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
CoordinateSystemType
Definition: MooseTypes.h:809
IntRange< T > make_range(T beg, T end)
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
This class provides variable solution values for other classes/objects to bind to when looping over f...
This class provides variable solution interface for linear finite volume problems.
Definition: FVUtils.h:24
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...
Definition: FaceInfo.h:225
const RemoteElem * remote_elem