17 #include <unordered_map> 19 #include <unordered_set> 21 #include "libmesh/petsc_vector.h" 52 bool subtract_mesh_velocity)
const override;
78 const std::vector<unsigned int> & momentum_system_numbers);
84 void computeHbyA(
const bool with_updated_pressure,
const bool verbose);
88 std::vector<std::unique_ptr<NumericVector<Number>>> &
102 template <
typename VarType>
120 std::vector<const MooseLinearVariableFVReal *>
_vel;
137 std::vector<std::unique_ptr<NumericVector<Number>>>
_HbyA_raw;
149 std::vector<std::unique_ptr<NumericVector<Number>>>
_Ainv_raw;
197 template <
typename VarType>
201 const auto & var_blocks = var.blockIDs();
202 const auto & uo_blocks =
blockIDs();
205 std::set<SubdomainID> uo_blocks_minus_var_blocks;
206 std::set_difference(uo_blocks.begin(),
210 std::inserter(uo_blocks_minus_var_blocks, uo_blocks_minus_var_blocks.end()));
211 if (uo_blocks_minus_var_blocks.size() > 0)
212 mooseError(
"Block restriction of interpolator user object '",
216 ") includes blocks not in the block restriction of variable '",
223 std::set<SubdomainID> var_blocks_minus_uo_blocks;
224 std::set_difference(var_blocks.begin(),
228 std::inserter(var_blocks_minus_uo_blocks, var_blocks_minus_uo_blocks.end()));
232 for (
auto & block_id : var_blocks_minus_uo_blocks)
235 std::set<SubdomainID> connected_blocks_on_uo;
236 std::set_intersection(connected_blocks.begin(),
237 connected_blocks.end(),
240 std::inserter(connected_blocks_on_uo, connected_blocks_on_uo.end()));
241 if (connected_blocks_on_uo.size() > 0)
242 mooseError(
"Block restriction of interpolator user object '",
246 ") doesn't match the block restriction of variable '",
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
virtual const std::string & name() const
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...
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< 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 ¶ms)
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)
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.