https://mooseframework.inl.gov
RhieChowMassFlux.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 
13 #include "CellCenteredMapFunctor.h"
14 #include "FaceCenteredMapFunctor.h"
15 #include "VectorComponentFunctor.h"
18 #include <unordered_map>
19 #include <set>
20 #include <unordered_set>
21 
22 #include "libmesh/petsc_vector.h"
23 
24 class MooseMesh;
27 namespace libMesh
28 {
29 class Elem;
30 class MeshBase;
31 }
32 
38 {
39 public:
41  RhieChowMassFlux(const InputParameters & params);
42 
44  Real getMassFlux(const FaceInfo & fi) const;
45 
47  Real getVolumetricFaceFlux(const FaceInfo & fi) const;
48 
50  const FaceInfo & fi,
51  const Moose::StateArg & time,
52  const THREAD_ID tid,
53  bool subtract_mesh_velocity) const override;
54 
56  void initFaceMassFlux();
58  void initCouplingField();
60  void computeFaceMassFlux();
62  void computeCellVelocity();
63 
64  virtual void meshChanged() override;
65  virtual void initialize() override;
66  virtual void execute() override {}
67  virtual void finalize() override {}
68  virtual void initialSetup() override;
69 
77  void linkMomentumPressureSystems(const std::vector<LinearSystem *> & momentum_systems,
78  const LinearSystem & pressure_system,
79  const std::vector<unsigned int> & momentum_system_numbers);
80 
85  void computeHbyA(const bool with_updated_pressure, const bool verbose);
86 
87 protected:
89  std::vector<std::unique_ptr<NumericVector<Number>>> &
90  selectPressureGradient(const bool updated_pressure);
91 
93  void setupMeshInformation();
94 
96  void
97  populateCouplingFunctors(const std::vector<std::unique_ptr<NumericVector<Number>>> & raw_hbya,
98  const std::vector<std::unique_ptr<NumericVector<Number>>> & raw_Ainv);
99 
103  template <typename VarType>
104  void checkBlocks(const VarType & var) const;
105 
106  virtual bool supportMeshVelocity() const override { return false; }
107 
110 
113 
115  const unsigned int _dim;
116 
119 
121  std::vector<const MooseLinearVariableFVReal *> _vel;
122 
125 
133 
138  std::vector<std::unique_ptr<NumericVector<Number>>> _HbyA_raw;
139 
145 
150  std::vector<std::unique_ptr<NumericVector<Number>>> _Ainv_raw;
151 
156 
158  std::vector<std::vector<LinearFVElementalKernel *>> _body_force_kernels;
160  std::vector<std::vector<std::string>> _body_force_kernel_names;
161 
166  std::vector<std::unique_ptr<NumericVector<Number>>> _grad_p_current;
167 
172 
174  std::vector<LinearSystem *> _momentum_systems;
175 
177  std::vector<unsigned int> _momentum_system_numbers;
178 
180  std::vector<unsigned int> _global_momentum_system_numbers;
181 
183  std::vector<libMesh::LinearImplicitSystem *> _momentum_implicit_systems;
184 
187 
190 
192  std::unique_ptr<NumericVector<Number>> _cell_volumes;
193 
196 
197 private:
200  std::vector<const FaceInfo *> _flow_face_info;
201 };
202 
203 template <typename VarType>
204 void
205 RhieChowMassFlux::checkBlocks(const VarType & var) const
206 {
207  const auto & var_blocks = var.blockIDs();
208  const auto & uo_blocks = blockIDs();
209 
210  // Error if this UO has any blocks that the variable does not
211  std::set<SubdomainID> uo_blocks_minus_var_blocks;
212  std::set_difference(uo_blocks.begin(),
213  uo_blocks.end(),
214  var_blocks.begin(),
215  var_blocks.end(),
216  std::inserter(uo_blocks_minus_var_blocks, uo_blocks_minus_var_blocks.end()));
217  if (uo_blocks_minus_var_blocks.size() > 0)
218  mooseError("Block restriction of interpolator user object '",
219  this->name(),
220  "' (",
222  ") includes blocks not in the block restriction of variable '",
223  var.name(),
224  "' (",
225  Moose::stringify(var.blocks()),
226  ")");
227 
228  // Get the blocks in the variable but not this UO
229  std::set<SubdomainID> var_blocks_minus_uo_blocks;
230  std::set_difference(var_blocks.begin(),
231  var_blocks.end(),
232  uo_blocks.begin(),
233  uo_blocks.end(),
234  std::inserter(var_blocks_minus_uo_blocks, var_blocks_minus_uo_blocks.end()));
235 
236  // For each block in the variable but not this UO, error if there is connection
237  // to any blocks on the UO.
238  for (auto & block_id : var_blocks_minus_uo_blocks)
239  {
240  const auto connected_blocks = _moose_mesh.getBlockConnectedBlocks(block_id);
241  std::set<SubdomainID> connected_blocks_on_uo;
242  std::set_intersection(connected_blocks.begin(),
243  connected_blocks.end(),
244  uo_blocks.begin(),
245  uo_blocks.end(),
246  std::inserter(connected_blocks_on_uo, connected_blocks_on_uo.end()));
247  if (connected_blocks_on_uo.size() > 0)
248  mooseError("Block restriction of interpolator user object '",
249  this->name(),
250  "' (",
251  Moose::stringify(uo_blocks),
252  ") doesn't match the block restriction of variable '",
253  var.name(),
254  "' (",
255  Moose::stringify(var_blocks),
256  ")");
257  }
258 }
std::vector< LinearSystem * > _momentum_systems
Pointers to the linear system(s) in moose corresponding to the momentum equation(s) ...
const MooseEnum _pressure_projection_method
Enumerator for the method used for pressure projection.
virtual void initialize() override
void setupMeshInformation()
Compute the cell volumes on the mesh.
User object responsible for determining the face fluxes using the Rhie-Chow interpolation in a segreg...
A functor whose evaluation relies on querying a map where the keys are face info ids and the values c...
void checkBlocks(const VarType &var) const
Check the block consistency between the passed in var and us.
std::vector< libMesh::LinearImplicitSystem * > _momentum_implicit_systems
Pointers to the momentum equation implicit system(s) from libmesh.
Real getMassFlux(const FaceInfo &fi) const
Get the face velocity times density (used in advection terms)
void computeFaceMassFlux()
Update the values of the face velocities in the containers.
const std::set< SubdomainID > & getBlockConnectedBlocks(const SubdomainID subdomain_id) const
std::vector< std::unique_ptr< NumericVector< Number > > > _Ainv_raw
We hold on to the cell-based 1/A vectors so that we can easily reconstruct the cell velocities as wel...
const LinearSystem * _pressure_system
Pointer to the pressure system.
FaceCenteredMapFunctor< RealVectorValue, std::unordered_map< dof_id_type, RealVectorValue > > _Ainv
A map functor from faces to $(1/A)_f$.
virtual const std::set< SubdomainID > & blockIDs() const
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
std::vector< unsigned int > _momentum_system_numbers
Numbers of the momentum system(s)
virtual void initialSetup() override
virtual void meshChanged() override
std::vector< unsigned int > _global_momentum_system_numbers
Global numbers of the momentum system(s)
void computeHbyA(const bool with_updated_pressure, const bool verbose)
Computes the inverse of the diagonal (1/A) of the system matrix plus the H/A components for the press...
const std::string & name() const
std::vector< std::vector< std::string > > _body_force_kernel_names
Vector of body force term names.
FaceCenteredMapFunctor< Real, std::unordered_map< dof_id_type, Real > > _HbyA_flux
A map functor from faces to $HbyA_{ij} = (A_{offdiag}*{(predicted~velocity)} - {Source})_{ij}/A_{ij}$...
virtual bool supportMeshVelocity() const override
Returns whether the UO can support mesh velocity advection.
std::vector< std::vector< LinearFVElementalKernel * > > _body_force_kernels
Pointer to the body force terms.
std::vector< const MooseLinearVariableFVReal * > _vel
The thread 0 copy of the x-velocity variable.
const MooseLinearVariableFVReal *const _p
The thread 0 copy of the pressure variable.
unsigned int _global_pressure_system_number
Global number of the pressure system.
RhieChowMassFlux(const InputParameters &params)
void linkMomentumPressureSystems(const std::vector< LinearSystem *> &momentum_systems, const LinearSystem &pressure_system, const std::vector< unsigned int > &momentum_system_numbers)
Update the momentum system-related information.
const Moose::Functor< Real > & _rho
Functor describing the density of the fluid.
static InputParameters validParams()
void computeCellVelocity()
Update the cell values of the velocity variables.
const MooseMesh & _moose_mesh
The MooseMesh that this user object operates on.
Real getVolumetricFaceFlux(const FaceInfo &fi) const
Get the volumetric face flux (used in advection terms)
std::string stringify(const T &t)
std::vector< std::unique_ptr< NumericVector< Number > > > & selectPressureGradient(const bool updated_pressure)
Select the right pressure gradient field and return a reference to the container. ...
virtual void execute() override
FaceCenteredMapFunctor< Real, std::unordered_map< dof_id_type, Real > > & _face_mass_flux
A map functor from faces to mass fluxes which are used in the advection terms.
std::vector< std::unique_ptr< NumericVector< Number > > > _grad_p_current
for a PISO iteration we need to hold on to the original pressure gradient field.
const std::vector< SubdomainName > & blocks() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void mooseError(Args &&... args) const
std::vector< const FaceInfo * > _flow_face_info
The subset of the FaceInfo objects that actually cover the subdomains which the flow field is defined...
void populateCouplingFunctors(const std::vector< std::unique_ptr< NumericVector< Number >>> &raw_hbya, const std::vector< std::unique_ptr< NumericVector< Number >>> &raw_Ainv)
Populate the face values of the H/A and 1/A fields.
LinearFVAnisotropicDiffusion * _p_diffusion_kernel
Pointer to the pressure diffusion term in the pressure Poisson equation.
std::vector< std::unique_ptr< NumericVector< Number > > > _HbyA_raw
We hold on to the cell-based HbyA vectors so that we can easily reconstruct the cell velocities as we...
void initCouplingField()
Initialize the coupling fields (HbyA and Ainv)
unsigned int THREAD_ID
const libMesh::MeshBase & _mesh
The libMesh mesh that this object acts on.
void initFaceMassFlux()
Initialize the container for face velocities.
std::unique_ptr< NumericVector< Number > > _cell_volumes
We will hold a vector of cell volumes to make sure we can do volume corrections rapidly.
virtual void finalize() override
const unsigned int _dim
The dimension of the mesh, e.g. 3 for hexes and tets, 2 for quads and tris.