https://mooseframework.inl.gov
HLLCUserObject.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 "InternalSideUserObject.h"
13 #include "MooseMesh.h"
14 
15 using side_type = std::pair<const Elem *, unsigned int>;
16 
17 // Forward Declarations
18 class HLLCUserObject;
20 
22 {
23 public:
25 
27 
28  virtual void initialSetup() override;
29  virtual void initialize() override {}
30  virtual void execute() override;
31  virtual void finalize() override {}
32  virtual void threadJoin(const UserObject & y) override;
33 
35  std::vector<ADReal> waveSpeed(const Elem * elem, unsigned int side) const;
36 
38  bool hasData(const Elem * elem, unsigned int side) const;
39 
40 protected:
42  const FaceInfo & faceInfoHelper(const Elem * elem, unsigned int side) const;
43 
45  unsigned int _qp = 0;
46 
49 
51  const std::vector<const FaceInfo *> & _face_info;
52 
54  std::map<side_type, unsigned int> _side_to_face_info;
55 
57  std::map<side_type, std::vector<ADReal>> _wave_speed;
58 
73 };
const SinglePhaseFluidProperties & _fluid
fluid properties
const ADMaterialProperty< Real > & _speed_elem
const ADMaterialProperty< Real > & _speed_neighbor
const ADMaterialProperty< Real > & _pressure_neighbor
const ADMaterialProperty< RealVectorValue > & _vel_elem
material properties computed by VarMat that Riemann solver needs
const ADMaterialProperty< Real > & _specific_internal_energy_elem
virtual void finalize() override
const std::vector< double > y
virtual void initialSetup() override
bool hasData(const Elem *elem, unsigned int side) const
Query whether this processor has data for the provided element and side.
std::map< side_type, std::vector< ADReal > > _wave_speed
data structure storing the wave speeds SL, SM, SR
virtual void threadJoin(const UserObject &y) override
virtual void execute() override
const ADMaterialProperty< Real > & _rho_elem
const MaterialProperty< Real > *const _eps_elem
const MaterialProperty< Real > *const _eps_neighbor
static InputParameters validParams()
unsigned int _qp
quadrature point dummy
Common class for single phase fluid properties.
const ADMaterialProperty< Real > & _specific_internal_energy_neighbor
std::map< side_type, unsigned int > _side_to_face_info
face info lookup allows searching for face info entry from elem/side pair
const std::vector< const FaceInfo * > & _face_info
FV face info from MooseMesh.
HLLCUserObject(const InputParameters &parameters)
std::vector< ADReal > waveSpeed(const Elem *elem, unsigned int side) const
accessor for the wave speed
std::pair< const Elem *, unsigned int > side_type
const ADMaterialProperty< Real > & _pressure_elem
const InputParameters & parameters() const
const FaceInfo & faceInfoHelper(const Elem *elem, unsigned int side) const
helper function for returning the FaceInfo object for an elem/side pair
const ADMaterialProperty< RealVectorValue > & _vel_neighbor
virtual void initialize() override
const ADMaterialProperty< Real > & _rho_neighbor