https://mooseframework.inl.gov
RhieChowInterpolatorBase.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 
11 #include "INSFVAttributes.h"
13 #include "GatherRCDataFaceThread.h"
14 #include "SubProblem.h"
15 #include "MooseMesh.h"
16 #include "SystemBase.h"
17 #include "NS.h"
18 #include "Assembly.h"
19 #include "INSFVVelocityVariable.h"
20 #include "INSFVPressureVariable.h"
22 #include "VectorCompositeFunctor.h"
23 #include "FVElementalKernel.h"
24 
25 using namespace libMesh;
26 
29 {
33 
34  params.addClassDescription(
35  "Computes the Rhie-Chow velocity based on gathered 'a' coefficient data.");
36 
37  // Avoid uninitialized residual objects
38  params.suppressParameter<bool>("force_preic");
39 
40  params.addRequiredParam<VariableName>(NS::pressure, "The pressure variable.");
41  params.addRequiredParam<VariableName>("u", "The x-component of velocity");
42  params.addParam<VariableName>("v", "The y-component of velocity");
43  params.addParam<VariableName>("w", "The z-component of velocity");
44 
45  MooseEnum velocity_interp_method("average rc", "rc");
46  params.addParam<MooseEnum>(
47  "velocity_interp_method",
48  velocity_interp_method,
49  "The interpolation to use for the velocity. Options are "
50  "'average' and 'rc' which stands for Rhie-Chow. The default is Rhie-Chow.");
51 
52  return params;
53 }
54 
56  : RhieChowFaceFluxProvider(params),
57  TaggingInterface(this),
58  ADFunctorInterface(this),
59  _moose_mesh(UserObject::_subproblem.mesh()),
60  _mesh(_moose_mesh.getMesh()),
61  _dim(blocksMaxDimension()),
62  _p(dynamic_cast<INSFVPressureVariable *>(
63  &UserObject::_subproblem.getVariable(0, getParam<VariableName>(NS::pressure)))),
64  _u(dynamic_cast<INSFVVelocityVariable *>(
65  &UserObject::_subproblem.getVariable(0, getParam<VariableName>("u")))),
66  _v(isParamValid("v") ? dynamic_cast<INSFVVelocityVariable *>(
67  &UserObject::_subproblem.getVariable(0, getParam<VariableName>("v")))
68  : nullptr),
69  _w(isParamValid("w") ? dynamic_cast<INSFVVelocityVariable *>(
70  &UserObject::_subproblem.getVariable(0, getParam<VariableName>("w")))
71  : nullptr),
72  _ps(libMesh::n_threads(), nullptr),
73  _us(libMesh::n_threads(), nullptr),
74  _vs(libMesh::n_threads(), nullptr),
75  _ws(libMesh::n_threads(), nullptr),
76  _displaced(dynamic_cast<DisplacedProblem *>(&(UserObject::_subproblem)))
77 {
78  if (!_p)
79  paramError(NS::pressure, "the pressure must be a INSFVPressureVariable.");
80  checkBlocks(*_p);
81 
82  if (!_u)
83  paramError("u", "the u velocity must be an INSFVVelocityVariable.");
84  checkBlocks(*_u);
85  _var_numbers.push_back(_u->number());
86 
87  if (_dim >= 2)
88  {
89  if (!_v)
90  mooseError("In two or more dimensions, the v velocity must be supplied and it must be an "
91  "INSFVVelocityVariable.");
92  checkBlocks(*_v);
93  _var_numbers.push_back(_v->number());
94  }
95 
96  if (_dim >= 3)
97  {
98  if (!_w)
99  mooseError("In three-dimensions, the w velocity must be supplied and it must be an "
100  "INSFVVelocityVariable.");
101  checkBlocks(*_w);
102  _var_numbers.push_back(_w->number());
103  }
104 
106  fillContainer("u", _us);
107 
108  if (_dim >= 2)
109  {
110  fillContainer("v", _vs);
112  mooseError("x and y velocity component face interpolation methods do not match");
113 
114  if (_dim >= 3)
115  {
116  fillContainer("w", _ws);
118  mooseError("x and z velocity component face interpolation methods do not match");
119  }
120  }
121 
123  mooseError("Different subproblems in RhieChowInterpolatorBase!");
124 
125  const auto & velocity_interp_method = params.get<MooseEnum>("velocity_interp_method");
126  if (velocity_interp_method == "average")
128  else if (velocity_interp_method == "rc")
130 }
131 
132 Real
134  const FaceInfo & fi,
135  const Moose::StateArg & time,
136  const THREAD_ID tid,
137  bool subtract_mesh_velocity) const
138 {
139  return raw_value(this->getVelocity(m, fi, time, tid, subtract_mesh_velocity)) * fi.normal();
140 }
INSFVPressureVariable *const _p
The thread 0 copy of the pressure variable.
std::vector< MooseVariableFVReal * > _ps
All the thread copies of the pressure variable.
SubProblem & _subproblem
unsigned int n_threads()
void paramError(const std::string &param, Args... args) const
T & getMesh(MooseMesh &mesh)
function to cast mesh
Definition: SCM.h:35
INSFVVelocityVariable *const _u
The thread 0 copy of the x-velocity variable.
unsigned int number() const
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
INSFVVelocityVariable *const _v
The thread 0 copy of the y-velocity variable (null if the problem is 1D)
static InputParameters validParams()
MeshBase & mesh
auto raw_value(const Eigen::Map< T > &in)
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
static InputParameters validParams()
SubProblem & _subproblem
virtual VectorValue< ADReal > getVelocity(const Moose::FV::InterpMethod m, const FaceInfo &fi, const Moose::StateArg &time, const THREAD_ID tid, bool subtract_mesh_velocity) const =0
Retrieve a face velocity.
std::vector< MooseVariableFVReal * > _us
All the thread copies of the x-velocity variable.
RhieChowInterpolatorBase(const InputParameters &params)
const unsigned int _dim
The dimension of the mesh, e.g. 3 for hexes and tets, 2 for quads and tris.
const Point & normal() const
Moose::FV::InterpMethod _velocity_interp_method
The interpolation method to use for the velocity.
void fillContainer(const std::string &var_name, Container &container)
Fill the passed-in variable container with the thread copies of var_name.
virtual Real getVolumetricFaceFlux(const Moose::FV::InterpMethod m, const FaceInfo &fi, const Moose::StateArg &time, const THREAD_ID tid, bool subtract_mesh_velocity) const override
Retrieve the volumetric face flux, will not include derivatives.
void checkBlocks(const VarType &var) const
Check the block consistency between the passed in var and us.
Moose::FV::InterpMethod faceInterpolationMethod() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
INSFVVelocityVariable *const _w
The thread 0 copy of the z-velocity variable (null if the problem is not 3D)
std::vector< unsigned int > _var_numbers
The velocity variable numbers.
static const std::string pressure
Definition: NS.h:56
void mooseError(Args &&... args) const
std::vector< MooseVariableFVReal * > _ws
All the thread copies of the z-velocity variable.
static InputParameters validParams()
std::vector< MooseVariableFVReal * > _vs
All the thread copies of the y-velocity variable.
static InputParameters validParams()
unsigned int THREAD_ID