https://mooseframework.inl.gov
Public Types | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
FaceCenteredMapFunctor< T, Map > Class Template Reference

A functor whose evaluation relies on querying a map where the keys are face info ids and the values correspond to the face values. More...

#include <FaceCenteredMapFunctor.h>

Inheritance diagram for FaceCenteredMapFunctor< T, Map >:
[legend]

Public Types

using ElemArg = Moose::ElemArg
 
using FaceArg = Moose::FaceArg
 
using ElemQpArg = Moose::ElemQpArg
 
using ElemSideQpArg = Moose::ElemSideQpArg
 
using ElemPointArg = Moose::ElemPointArg
 
using StateArg = Moose::StateArg
 
using NodeArg = Moose::NodeArg
 
typedef FunctorBase< T > FunctorType
 
typedef T ValueType
 
typedef typename FunctorReturnType< T, FunctorEvaluationKind::Gradient >::type GradientType
 
typedef ValueType DotType
 

Public Member Functions

 FaceCenteredMapFunctor (const MooseMesh &mesh, const std::string &name)
 
 FaceCenteredMapFunctor (const MooseMesh &mesh, const std::set< SubdomainID > &sub_ids, const std::string &name)
 
bool hasBlocks (SubdomainID sub_id) const override
 
ValueType evaluate (const FaceInfo *const fi) const
 Evaluate the face functor using a FaceInfo argument. More...
 
bool supportsFaceArg () const override final
 
bool supportsElemSideQpArg () const override final
 
FunctorReturnType< T, FET >::type genericEvaluate (const Space &r, const State &state) const
 
const MooseFunctorName & functorName () const
 
virtual void residualSetup () override
 
virtual void jacobianSetup () override
 
virtual void timestepSetup () override
 
virtual void customSetup (const ExecFlagType &exec_type) override
 
void setCacheClearanceSchedule (const std::set< ExecFlagType > &clearance_schedule)
 
virtual bool isExtrapolatedBoundaryFace (const FaceInfo &, const Elem *, const StateArg &) const
 
bool isInternalFace (const FaceInfo &) const
 
virtual bool isConstant () const
 
virtual bool hasFaceSide (const FaceInfo &fi, const bool fi_elem_side) const override
 
void checkFace (const Moose::FaceArg &face) const
 
ValueType operator() (const ElemArg &elem, const StateArg &state) const
 
ValueType operator() (const FaceArg &face, const StateArg &state) const
 
ValueType operator() (const ElemQpArg &qp, const StateArg &state) const
 
ValueType operator() (const ElemSideQpArg &qp, const StateArg &state) const
 
ValueType operator() (const ElemPointArg &elem_point, const StateArg &state) const
 
ValueType operator() (const NodeArg &node, const StateArg &state) const
 
ValueType operator() (const ElemArg &elem, const StateArg &state) const
 
ValueType operator() (const FaceArg &face, const StateArg &state) const
 
ValueType operator() (const ElemQpArg &qp, const StateArg &state) const
 
ValueType operator() (const ElemSideQpArg &qp, const StateArg &state) const
 
ValueType operator() (const ElemPointArg &elem_point, const StateArg &state) const
 
ValueType operator() (const NodeArg &node, const StateArg &state) const
 
GradientType gradient (const ElemArg &elem, const StateArg &state) const
 
GradientType gradient (const FaceArg &face, const StateArg &state) const
 
GradientType gradient (const ElemQpArg &qp, const StateArg &state) const
 
GradientType gradient (const ElemSideQpArg &qp, const StateArg &state) const
 
GradientType gradient (const ElemPointArg &elem_point, const StateArg &state) const
 
GradientType gradient (const NodeArg &node, const StateArg &state) const
 
GradientType gradient (const ElemArg &elem, const StateArg &state) const
 
GradientType gradient (const FaceArg &face, const StateArg &state) const
 
GradientType gradient (const ElemQpArg &qp, const StateArg &state) const
 
GradientType gradient (const ElemSideQpArg &qp, const StateArg &state) const
 
GradientType gradient (const ElemPointArg &elem_point, const StateArg &state) const
 
GradientType gradient (const NodeArg &node, const StateArg &state) const
 
DotType dot (const ElemArg &elem, const StateArg &state) const
 
DotType dot (const FaceArg &face, const StateArg &state) const
 
DotType dot (const ElemQpArg &qp, const StateArg &state) const
 
DotType dot (const ElemSideQpArg &qp, const StateArg &state) const
 
DotType dot (const ElemPointArg &elem_point, const StateArg &state) const
 
DotType dot (const NodeArg &node, const StateArg &state) const
 
DotType dot (const ElemArg &elem, const StateArg &state) const
 
DotType dot (const FaceArg &face, const StateArg &state) const
 
DotType dot (const ElemQpArg &qp, const StateArg &state) const
 
DotType dot (const ElemSideQpArg &qp, const StateArg &state) const
 
DotType dot (const ElemPointArg &elem_point, const StateArg &state) const
 
DotType dot (const NodeArg &node, const StateArg &state) const
 
GradientType gradDot (const ElemArg &elem, const StateArg &state) const
 
GradientType gradDot (const FaceArg &face, const StateArg &state) const
 
GradientType gradDot (const ElemQpArg &qp, const StateArg &state) const
 
GradientType gradDot (const ElemSideQpArg &qp, const StateArg &state) const
 
GradientType gradDot (const ElemPointArg &elem_point, const StateArg &state) const
 
GradientType gradDot (const NodeArg &node, const StateArg &state) const
 
GradientType gradDot (const ElemArg &elem, const StateArg &state) const
 
GradientType gradDot (const FaceArg &face, const StateArg &state) const
 
GradientType gradDot (const ElemQpArg &qp, const StateArg &state) const
 
GradientType gradDot (const ElemSideQpArg &qp, const StateArg &state) const
 
GradientType gradDot (const ElemPointArg &elem_point, const StateArg &state) const
 
GradientType gradDot (const NodeArg &node, const StateArg &state) const
 

Protected Member Functions

virtual GradientType evaluateGradient (const ElemArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const FaceArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemSideQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemPointArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const NodeArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const FaceArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemSideQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const ElemPointArg &, const StateArg &) const
 
virtual GradientType evaluateGradient (const NodeArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemArg &, const StateArg &) const
 
virtual DotType evaluateDot (const FaceArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemQpArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemSideQpArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemPointArg &, const StateArg &) const
 
virtual DotType evaluateDot (const NodeArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemArg &, const StateArg &) const
 
virtual DotType evaluateDot (const FaceArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemQpArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemSideQpArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemPointArg &, const StateArg &) const
 
virtual DotType evaluateDot (const NodeArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const FaceArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemSideQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemPointArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const NodeArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const FaceArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemSideQpArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const ElemPointArg &, const StateArg &) const
 
virtual GradientType evaluateGradDot (const NodeArg &, const StateArg &) const
 

Private Member Functions

ValueType evaluate (const ElemArg &elem_arg, const StateArg &state) const override final
 
ValueType evaluate (const FaceArg &face, const StateArg &state) const override final
 
ValueType evaluate (const ElemPointArg &, const StateArg &) const override
 
ValueType evaluate (const ElemQpArg &, const StateArg &) const override
 
ValueType evaluate (const ElemSideQpArg &, const StateArg &) const override
 
ValueType evaluate (const NodeArg &node_arg, const StateArg &state) const override final
 

Private Attributes

const MooseMesh_mesh
 The mesh that this functor lives on. More...
 
const std::set< SubdomainID_sub_ids
 The subdomain IDs that this functor lives on. More...
 

Detailed Description

template<typename T, typename Map>
class FaceCenteredMapFunctor< T, Map >

A functor whose evaluation relies on querying a map where the keys are face info ids and the values correspond to the face values.

The primary purpose of this functor is to store face based fields (face flux, face velocity, normal gradient) often encountered in fluid dynamics problems using the finite volume method

Definition at line 25 of file FaceCenteredMapFunctor.h.

Member Typedef Documentation

◆ ElemArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::ElemArg = Moose::ElemArg

Definition at line 31 of file FaceCenteredMapFunctor.h.

◆ ElemPointArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::ElemPointArg = Moose::ElemPointArg

Definition at line 35 of file FaceCenteredMapFunctor.h.

◆ ElemQpArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::ElemQpArg = Moose::ElemQpArg

Definition at line 33 of file FaceCenteredMapFunctor.h.

◆ ElemSideQpArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::ElemSideQpArg = Moose::ElemSideQpArg

Definition at line 34 of file FaceCenteredMapFunctor.h.

◆ FaceArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::FaceArg = Moose::FaceArg

Definition at line 32 of file FaceCenteredMapFunctor.h.

◆ NodeArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::NodeArg = Moose::NodeArg

Definition at line 37 of file FaceCenteredMapFunctor.h.

◆ StateArg

template<typename T, typename Map>
using FaceCenteredMapFunctor< T, Map >::StateArg = Moose::StateArg

Definition at line 36 of file FaceCenteredMapFunctor.h.

Constructor & Destructor Documentation

◆ FaceCenteredMapFunctor() [1/2]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::FaceCenteredMapFunctor ( const MooseMesh mesh,
const std::string &  name 
)

Definition at line 73 of file FaceCenteredMapFunctor.h.

76 {
77 }
const std::string name
Definition: Setup.h:20
const MooseMesh & _mesh
The mesh that this functor lives on.

◆ FaceCenteredMapFunctor() [2/2]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::FaceCenteredMapFunctor ( const MooseMesh mesh,
const std::set< SubdomainID > &  sub_ids,
const std::string &  name 
)

Definition at line 80 of file FaceCenteredMapFunctor.h.

84  _mesh(mesh),
85  _sub_ids(sub_ids == mesh.meshSubdomains() ? std::set<SubdomainID>() : sub_ids)
86 {
87 }
MeshBase & mesh
const std::string name
Definition: Setup.h:20
const std::set< SubdomainID > _sub_ids
The subdomain IDs that this functor lives on.
const MooseMesh & _mesh
The mesh that this functor lives on.

Member Function Documentation

◆ evaluate() [1/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const FaceInfo *const  fi) const

Evaluate the face functor using a FaceInfo argument.

Parameters
fiThe object containing the face information

Definition at line 132 of file FaceCenteredMapFunctor.C.

Referenced by RhieChowMassFlux::getMassFlux().

133 {
134  try
135  {
136  return libmesh_map_find(*this, fi->id());
137  }
138  catch (libMesh::LogicError &)
139  {
140  if (!_sub_ids.empty() && !_sub_ids.count(fi->elem().subdomain_id()))
141  {
142  if (fi->neighborPtr() && !_sub_ids.count(fi->neighborPtr()->subdomain_id()))
143  mooseError("Attempted to evaluate FaceCenteredMapFunctor '",
144  this->functorName(),
145  "' with an element subdomain id of '",
146  fi->elem().subdomain_id(),
147  fi->neighborPtr() ? " or neighbor subdomain id of '" +
148  std::to_string(fi->neighborPtr()->subdomain_id()) + "'"
149  : "",
150  "' but that subdomain id is not one of the subdomain ids the functor is "
151  "restricted to.");
152  }
153  else
154  mooseError("Attempted access into FaceCenteredMapFunctor '",
155  this->functorName(),
156  "' with a key that does not yet exist in the map. Make sure to fill your "
157  "FaceCenteredMapFunctor for all elements you will attempt to access later.");
158 
160  }
161 }
void mooseError(Args &&... args)
const Elem & elem() const
const Elem * neighborPtr() const
const MooseFunctorName & functorName() const
subdomain_id_type subdomain_id() const
const std::set< SubdomainID > _sub_ids
The subdomain IDs that this functor lives on.
dof_id_type id() const

◆ evaluate() [2/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const ElemArg elem_arg,
const StateArg state 
) const
finaloverrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 45 of file FaceCenteredMapFunctor.C.

46 {
47  // The following reconstruction is based on Weller's method. For more information on this,
48  // we recommend:
49  //
50  // Weller, Hilary. "Non-orthogonal version of the arbitrary polygonal C-grid and a new diamond
51  // grid." Geoscientific Model Development 7.3 (2014): 779-797.
52  //
53  // and
54  //
55  // Aguerre, Horacio J., et al. "An oscillation-free flow solver based on flux reconstruction."
56  // Journal of Computational Physics 365 (2018): 135-148.
57  //
58  // This basically reconstructs the cell value based on flux values as follows:
59  //
60  // $\left( \sum_f n_f \outer S_f \right)^{-1} \sum_f (\phi_f \cdot S_f)n_f$
61  //
62  // where $S_f$ is the surface normal vector, $n_f$ is the unit surface vector and $\phi_f$ is a
63  // vector value on the field. Hence the restriction to vector values.
64 
66 
67  // Primitive type one rank below the type of the stored data (T). If we are storing a rank two
68  // tensor, this is a vector, if we are storing a vector this is just a number
69  using PrimitiveType = typename MetaPhysicL::ReplaceAlgebraicType<
70  T,
72 
74  {
75  const auto dim = _mesh.dimension();
76  // The reason why these are DenseVector/Matrix is that when the mesh dimension is lower than 3,
77  // we get singular matrixes if we try to invert TensorValues.
78  DenseMatrix<PrimitiveType> n_x_Sf(dim, dim);
79  DenseVector<PrimitiveType> sum_normal_flux(dim);
80 
81  const Elem * const elem = elem_arg.elem;
82 
83  for (const auto side : make_range(elem->n_sides()))
84  {
85  const Elem * const neighbor = elem->neighbor_ptr(side);
86 
87  // We need to check if the faceinfo belongs to the element or the neighbor. Based on that we
88  // query the faceinfo and adjust the normal to point outward of the current cell
89  const bool elem_has_fi = Moose::FV::elemHasFaceInfo(*elem, neighbor);
90  const FaceInfo * const fi = _mesh.faceInfo(
91  elem_has_fi ? elem : neighbor, elem_has_fi ? side : neighbor->which_neighbor_am_i(elem));
92  const Point & normal = elem_has_fi ? fi->normal() : Point(-fi->normal());
93 
94  const Point area_vector = normal * fi->faceArea();
95  const ValueType face_value = this->evaluate(fi);
96 
97  const auto product = Moose::outer_product(normal, area_vector);
98 
99  const auto flux_contrib = normal * (face_value * area_vector);
100  for (const auto i : make_range(dim))
101  {
102  sum_normal_flux(i) += flux_contrib(i);
103  for (const auto j : make_range(dim))
104  n_x_Sf(i, j) += product(i, j);
105  }
106  }
107 
108  // We do the inversion of the surface vector matrix here. It is symmetric
109  // and small so we can do it using a Cholesky decomposition.
110  DenseVector<PrimitiveType> dense_result(dim);
111  n_x_Sf.cholesky_solve(sum_normal_flux, dense_result);
112 
113  ValueType result;
114  for (const auto i : make_range(dim))
115  result(i) = dense_result(i);
116 
117  return result;
118  }
119  else
120  mooseError("Cell center reconstruction is not implemented!");
121 }
bool elemHasFaceInfo(const Elem &elem, const Elem *const neighbor)
void mooseError(Args &&... args)
unsigned int dim
TypeVector< typename CompareTypes< T, T2 >::supertype > outer_product(const T &a, const TypeVector< T2 > &b)
Real faceArea() const
const std::vector< const FaceInfo *> & faceInfo() const
ValueType evaluate(const FaceInfo *const fi) const
Evaluate the face functor using a FaceInfo argument.
virtual unsigned int dimension() const
const Point & normal() const
virtual unsigned int n_sides() const=0
const Elem * neighbor_ptr(unsigned int i) const
IntRange< T > make_range(T beg, T end)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MooseMesh & _mesh
The mesh that this functor lives on.

◆ evaluate() [3/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const FaceArg face,
const StateArg state 
) const
finaloverrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 125 of file FaceCenteredMapFunctor.C.

126 {
127  return this->evaluate(face.fi);
128 }
ValueType evaluate(const FaceInfo *const fi) const
Evaluate the face functor using a FaceInfo argument.

◆ evaluate() [4/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const ElemPointArg ,
const StateArg  
) const
overrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 98 of file FaceCenteredMapFunctor.h.

99 {
100  mooseError("not implemented");
101 }
void mooseError(Args &&... args)

◆ evaluate() [5/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const ElemQpArg ,
const StateArg  
) const
overrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 105 of file FaceCenteredMapFunctor.h.

106 {
107  mooseError("not implemented");
108 }
void mooseError(Args &&... args)

◆ evaluate() [6/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const ElemSideQpArg ,
const StateArg  
) const
overrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 112 of file FaceCenteredMapFunctor.h.

113 {
114  mooseError("not implemented");
115 }
void mooseError(Args &&... args)

◆ evaluate() [7/7]

template<typename T , typename Map >
FaceCenteredMapFunctor< T, Map >::ValueType FaceCenteredMapFunctor< T, Map >::evaluate ( const NodeArg node_arg,
const StateArg state 
) const
finaloverrideprivatevirtual

Implements Moose::FunctorBase< T >.

Definition at line 119 of file FaceCenteredMapFunctor.h.

120 {
121  mooseError("not implemented");
122 }
void mooseError(Args &&... args)

◆ hasBlocks()

template<typename T , typename Map >
bool FaceCenteredMapFunctor< T, Map >::hasBlocks ( SubdomainID  sub_id) const
overridevirtual

Reimplemented from Moose::FunctorBase< T >.

Reimplemented in FictionalFaceCenteredMapFunctor< T, Map >.

Definition at line 91 of file FaceCenteredMapFunctor.h.

92 {
93  return _sub_ids.empty() || _sub_ids.count(sub_id);
94 }
const std::set< SubdomainID > _sub_ids
The subdomain IDs that this functor lives on.

◆ supportsElemSideQpArg()

template<typename T, typename Map>
bool FaceCenteredMapFunctor< T, Map >::supportsElemSideQpArg ( ) const
inlinefinaloverridevirtual

Implements Moose::FunctorBase< T >.

Definition at line 54 of file FaceCenteredMapFunctor.h.

54 { return false; }

◆ supportsFaceArg()

template<typename T, typename Map>
bool FaceCenteredMapFunctor< T, Map >::supportsFaceArg ( ) const
inlinefinaloverridevirtual

Implements Moose::FunctorBase< T >.

Definition at line 53 of file FaceCenteredMapFunctor.h.

53 { return true; }

Member Data Documentation

◆ _mesh

template<typename T, typename Map>
const MooseMesh& FaceCenteredMapFunctor< T, Map >::_mesh
private

The mesh that this functor lives on.

Definition at line 58 of file FaceCenteredMapFunctor.h.

◆ _sub_ids

template<typename T, typename Map>
const std::set<SubdomainID> FaceCenteredMapFunctor< T, Map >::_sub_ids
private

The subdomain IDs that this functor lives on.

If empty, then we consider the functor to live on all subdomains

Definition at line 62 of file FaceCenteredMapFunctor.h.


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