https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
BernoulliPressureVariable Class Referenceabstract

A special variable class for pressure which flags faces at which porosity jumps occur as extrapolated or Dirichlet boundary faces. More...

#include <BernoulliPressureVariable.h>

Inheritance diagram for BernoulliPressureVariable:
[legend]

Public Types

typedef typename MooseVariableField< Real >::OutputGradient OutputGradient
 
typedef typename MooseVariableField< Real >::OutputSecond OutputSecond
 
typedef typename MooseVariableField< Real >::OutputDivergence OutputDivergence
 
typedef typename MooseVariableField< Real >::FieldVariableValue FieldVariableValue
 
typedef typename MooseVariableField< Real >::FieldVariableGradient FieldVariableGradient
 
typedef typename MooseVariableField< Real >::FieldVariableSecond FieldVariableSecond
 
typedef typename MooseVariableField< Real >::FieldVariableCurl FieldVariableCurl
 
typedef typename MooseVariableField< Real >::FieldVariableDivergence FieldVariableDivergence
 
typedef typename MooseVariableField< Real >::OutputShape OutputShape
 
typedef typename MooseVariableField< Real >::OutputShapeGradient OutputShapeGradient
 
typedef typename MooseVariableField< Real >::OutputShapeSecond OutputShapeSecond
 
typedef typename MooseVariableField< Real >::OutputShapeDivergence OutputShapeDivergence
 
typedef typename MooseVariableField< Real >::OutputData OutputData
 
typedef typename MooseVariableField< Real >::DoFValue DoFValue
 
typedef typename MooseVariableField< Real >::FieldVariablePhiValue FieldVariablePhiValue
 
typedef typename MooseVariableField< Real >::FieldVariablePhiDivergence FieldVariablePhiDivergence
 
typedef typename MooseVariableField< Real >::FieldVariablePhiGradient FieldVariablePhiGradient
 
typedef typename MooseVariableField< Real >::FieldVariablePhiSecond FieldVariablePhiSecond
 
typedef Moose::ElemQpArg ElemQpArg
 
typedef Moose::ElemSideQpArg ElemSideQpArg
 
typedef Moose::ElemArg ElemArg
 
typedef Moose::FaceArg FaceArg
 
typedef Moose::StateArg StateArg
 
typedef Moose::NodeArg NodeArg
 
typedef Moose::ElemPointArg ElemPointArg
 
typedef MooseArray< std::vector< OutputShape > > FieldVariablePhiCurl
 
typedef MooseArray< std::vector< OutputShape > > FieldVariableTestValue
 
typedef MooseArray< std::vector< OutputShapeGradient > > FieldVariableTestGradient
 
typedef MooseArray< std::vector< OutputShapeSecond > > FieldVariableTestSecond
 
typedef MooseArray< std::vector< OutputShape > > FieldVariableTestCurl
 
typedef MooseArray< std::vector< OutputShapeDivergence > > FieldVariableTestDivergence
 
typedef DataFileName DataFileParameterType
 
typedef FunctorBase< T > FunctorType
 
typedef T ValueType
 
typedef typename FunctorReturnType< T, FunctorEvaluationKind::Gradient >::type GradientType
 
typedef ValueType DotType
 

Public Member Functions

 BernoulliPressureVariable (const InputParameters &params)
 
bool isExtrapolatedBoundaryFace (const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
 
void initialSetup () override
 
void computeFaceValues (const FaceInfo &) override
 
void computeElemValues () override
 
void computeElemValuesFace () override
 
void computeNeighborValuesFace () override
 
void computeNeighborValues () override
 
void requireQpComputations () const override
 
virtual void timestepSetup () override
 
virtual void meshChanged () override
 
bool isSeparatorBoundary (const FaceInfo &fi) const
 
virtual bool isExtrapolatedBoundaryFace (const FaceInfo &, const Elem *, const StateArg &) const
 
virtual bool isFV () const override
 
virtual void prepare () override final
 
virtual void prepareNeighbor () override final
 
virtual void prepareAux () override final
 
virtual void reinitNode () override final
 
virtual void reinitNodes (const std::vector< dof_id_type > &) override final
 
virtual void reinitNodesNeighbor (const std::vector< dof_id_type > &) override final
 
virtual void reinitAux () override final
 
virtual void reinitAuxNeighbor () override final
 
virtual void prepareLowerD () override final
 
virtual const dof_id_typenodalDofIndex () const override final
 
virtual const dof_id_typenodalDofIndexNeighbor () const override final
 
virtual std::size_t phiSize () const override final
 
virtual std::size_t phiFaceSize () const override final
 
virtual std::size_t phiNeighborSize () const override final
 
virtual std::size_t phiFaceNeighborSize () const override final
 
virtual std::size_t phiLowerSize () const override final
 
virtual void computeLowerDValues () override final
 
virtual void computeNodalNeighborValues () override final
 
virtual void computeNodalValues () override final
 
virtual const std::vector< dof_id_type > & dofIndicesLower () const override final
 
unsigned int numberOfDofs () const override final
 
virtual unsigned int numberOfDofsNeighbor () override final
 
virtual bool isNodal () const override final
 
bool hasDoFsOnNodes () const override final
 
libMesh::FEContinuity getContinuity () const override final
 
virtual bool isNodalDefined () const override final
 
virtual void setNodalValue (const Real &value, unsigned int idx=0) override
 
virtual void setDofValue (const OutputData &value, unsigned int index) override
 
virtual void setDofValue (const OutputData &value, unsigned int index)=0
 
void clearDofIndices () override
 
virtual void prepareIC () override
 
virtual const Elem *const & currentElem () const override
 
virtual void getDofIndices (const Elem *elem, std::vector< dof_id_type > &dof_indices) const override
 
virtual const std::vector< dof_id_type > & dofIndices () const final
 
virtual const std::vector< dof_id_type > & dofIndicesNeighbor () const final
 
Moose::FV::InterpMethod faceInterpolationMethod () const
 
void clearAllDofIndices () final
 
const DoFValuenodalVectorTagValue (TagID) const override
 
const DoFValuenodalMatrixTagValue (TagID) const override
 
const FieldVariableValuevectorTagValue (TagID tag) const override
 
const DoFValuevectorTagDofValue (TagID tag) const override
 
const FieldVariableValuematrixTagValue (TagID tag) const override
 
const FieldVariableValuevectorTagValueNeighbor (TagID tag)
 
const FieldVariableValuematrixTagValueNeighbor (TagID tag)
 
const FieldVariableValueuDot () const
 
const FieldVariableValuesln () const override
 
const FieldVariableValueslnOld () const override
 
const FieldVariableValueslnOlder () const override
 
const FieldVariableGradientgradSln () const override
 
const FieldVariableGradientgradSlnOld () const override
 
const FieldVariableValueuDotNeighbor () const
 
const FieldVariableValueslnNeighbor () const override
 
const FieldVariableValueslnOldNeighbor () const override
 
const FieldVariableGradientgradSlnNeighbor () const override
 
const FieldVariableGradientgradSlnOldNeighbor () const override
 
const VariableValueduDotDu () const
 
const VariableValueduDotDotDu () const
 
const VariableValueduDotDuNeighbor () const
 
const VariableValueduDotDotDuNeighbor () const
 
const ADTemplateVariableValue< Real > & adSln () const override
 
const ADTemplateVariableGradient< Real > & adGradSln () const override
 
virtual const VectorValue< ADReal > & adGradSln (const Elem *const elem, const StateArg &state, const bool correct_skewness=false) const
 
virtual VectorValue< ADRealadGradSln (const FaceInfo &fi, const StateArg &state, const bool correct_skewness=false) const
 
virtual VectorValue< ADRealuncorrectedAdGradSln (const FaceInfo &fi, const StateArg &state, const bool correct_skewness=false) const
 
ADReal getBoundaryFaceValue (const FaceInfo &fi, const StateArg &state, bool correct_skewness=false) const
 
const ADTemplateVariableSecond< Real > & adSecondSln () const override
 
const ADTemplateVariableValue< Real > & adUDot () const override
 
const ADTemplateVariableValue< Real > & adUDotDot () const override
 
const ADTemplateVariableGradient< Real > & adGradSlnDot () const override
 
const ADTemplateVariableCurl< Real > & adCurlSln () const override
 
const ADTemplateVariableValue< Real > & adSlnNeighbor () const override
 
const ADTemplateVariableGradient< Real > & adGradSlnNeighbor () const override
 
const ADTemplateVariableSecond< Real > & adSecondSlnNeighbor () const override
 
const ADTemplateVariableValue< Real > & adUDotNeighbor () const override
 
const ADTemplateVariableValue< Real > & adUDotDotNeighbor () const override
 
const ADTemplateVariableGradient< Real > & adGradSlnNeighborDot () const override
 
const ADTemplateVariableCurl< Real > & adCurlSlnNeighbor () const override
 
virtual void setDofValues (const DenseVector< OutputData > &values) override
 
virtual void setDofValues (const DenseVector< OutputData > &values)=0
 
virtual void setLowerDofValues (const DenseVector< OutputData > &values) override
 
virtual void setLowerDofValues (const DenseVector< OutputData > &values)=0
 
OutputData getElementalValue (const Elem *elem, unsigned int idx=0) const
 
OutputData getElementalValueOld (const Elem *elem, unsigned int idx=0) const
 
OutputData getElementalValueOlder (const Elem *elem, unsigned int idx=0) const
 
virtual void insert (libMesh::NumericVector< libMesh::Number > &vector) override
 
virtual void insertLower (libMesh::NumericVector< libMesh::Number > &vector) override
 
virtual void add (libMesh::NumericVector< libMesh::Number > &vector) override
 
const DoFValuedofValues () const override
 
const DoFValuedofValuesOld () const override
 
const DoFValuedofValuesOlder () const override
 
const DoFValuedofValuesPreviousNL () const override
 
const DoFValuedofValuesNeighbor () const override
 
const DoFValuedofValuesOldNeighbor () const override
 
const DoFValuedofValuesOlderNeighbor () const override
 
const DoFValuedofValuesPreviousNLNeighbor () const override
 
const DoFValuedofValuesDot () const override
 
const DoFValuedofValuesDotNeighbor () const override
 
const DoFValuedofValuesDotOld () const override
 
const DoFValuedofValuesDotOldNeighbor () const override
 
const DoFValuedofValuesDotDot () const override
 
const DoFValuedofValuesDotDotNeighbor () const override
 
const DoFValuedofValuesDotDotOld () const override
 
const DoFValuedofValuesDotDotOldNeighbor () const override
 
const MooseArray< libMesh::Number > & dofValuesDuDotDu () const override
 
const MooseArray< libMesh::Number > & dofValuesDuDotDuNeighbor () const override
 
const MooseArray< libMesh::Number > & dofValuesDuDotDotDu () const override
 
const MooseArray< libMesh::Number > & dofValuesDuDotDotDuNeighbor () const override
 
const MooseArray< ADReal > & adDofValues () const override
 
const MooseArray< ADReal > & adDofValuesNeighbor () const override
 
const MooseArray< ADReal > & adDofValuesDot () const override
 
Real getValue (const Elem *elem) const
 
OutputTools< Real >::OutputGradient getGradient (const Elem *elem) const
 
bool hasDirichletBC () const
 
std::pair< bool, const FVDirichletBCBase * > getDirichletBC (const FaceInfo &fi) const
 
std::pair< bool, std::vector< const FVFluxBC * > > getFluxBCs (const FaceInfo &fi) const
 
virtual void residualSetup () override
 
virtual void jacobianSetup () override
 
ADReal getElemValue (const Elem *elem, const StateArg &state) const
 
void setActiveTags (const std::set< TagID > &vtags) override
 
bool supportsFaceArg () const override final
 
bool supportsElemSideQpArg () const override final
 
const MooseArray< Real > & nodalValueArray () const override
 
const MooseArray< Real > & nodalValueOldArray () const override
 
const MooseArray< Real > & nodalValueOlderArray () const override
 
bool computingSecond () const override final
 
bool computingCurl () const override final
 
bool computingDiv () const override final
 
bool usesSecondPhiNeighbor () const override final
 
const FieldVariablePhiValuephi () const override final
 
const FieldVariablePhiGradientgradPhi () const override final
 
const FieldVariablePhiSecondsecondPhi () const override final
 
const FieldVariablePhiValuecurlPhi () const override final
 
const FieldVariablePhiDivergencedivPhi () const override final
 
const FieldVariablePhiValuephiFace () const override final
 
const FieldVariablePhiGradientgradPhiFace () const override final
 
const FieldVariablePhiSecondsecondPhiFace () const override final
 
const FieldVariablePhiValuephiFaceNeighbor () const override final
 
const FieldVariablePhiGradientgradPhiFaceNeighbor () const override final
 
const FieldVariablePhiSecondsecondPhiFaceNeighbor () const override final
 
const FieldVariablePhiValuephiNeighbor () const override final
 
const FieldVariablePhiGradientgradPhiNeighbor () const override final
 
const FieldVariablePhiSecondsecondPhiNeighbor () const override final
 
virtual const FieldVariablePhiValuephiLower () const override
 
unsigned int oldestSolutionStateRequested () const override final
 
virtual ADReal getExtrapolatedBoundaryFaceValue (const FaceInfo &fi, bool two_term_expansion, bool correct_skewness, const Elem *elem_side_to_extrapolate_from, const StateArg &state) const
 
const bool & getTwoTermBoundaryExpansion () const
 
virtual Moose::VarFieldType fieldType () const override
 
virtual bool isArray () const override
 
virtual bool isVector () const override
 
bool usesPhiNeighbor () const
 
bool usesGradPhiNeighbor () const
 
bool hasBlocks (const SubdomainName &name) const
 
bool hasBlocks (const std::vector< SubdomainName > &names) const
 
bool hasBlocks (const std::set< SubdomainName > &names) const
 
bool hasBlocks (const std::vector< SubdomainID > &ids) const
 
bool hasBlocks (const std::set< SubdomainID > &ids) const
 
bool hasBlocks (const SubdomainID id) const override
 
const std::string & componentName (const unsigned int comp) const
 
const std::set< SubdomainID > & activeSubdomains () const
 
bool activeOnSubdomain (SubdomainID subdomain) const
 
bool activeOnSubdomains (const std::set< SubdomainID > &subdomains) const
 
virtual bool needsGradientVectorStorage () const
 
const std::string & arrayVariableComponent (const unsigned int i) const
 
unsigned int number () const
 
const libMesh::FETypefeType () const
 
SystemBasesys ()
 
const SystemBasesys () const
 
const std::string & name () const override
 
bool useDual () const
 
const std::vector< dof_id_type > & allDofIndices () const
 
unsigned int totalVarDofs ()
 
Moose::VarKindType kind () const
 
void scalingFactor (const std::vector< Real > &factor)
 
Real scalingFactor () const
 
const std::vector< Real > & arrayScalingFactor () const
 
libMesh::Order order () const
 
unsigned int count () const
 
const libMesh::DofMapdofMap () const
 
std::vector< dof_id_typecomponentDofIndices (const std::vector< dof_id_type > &dof_indices, unsigned int component) const
 
bool eigen () const
 
void eigen (bool eigen)
 
bool isLowerD () const
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
const std::vector< SubdomainName > & blocks () const
 
unsigned int numBlocks () const
 
virtual const std::set< SubdomainID > & blockIDs () const
 
unsigned int blocksMaxDimension () const
 
bool isBlockSubset (const std::set< SubdomainID > &ids) const
 
bool isBlockSubset (const std::vector< SubdomainID > &ids) const
 
bool hasBlockMaterialProperty (const std::string &prop_name)
 
const std::set< SubdomainID > & meshBlockIDs () const
 
virtual bool blockRestricted () const
 
virtual void checkVariable (const MooseVariableFieldBase &variable) const
 
void buildOutputHideVariableList (std::set< std::string > variable_names)
 
const std::set< OutputName > & getOutputs ()
 
virtual void subdomainSetup ()
 
virtual void customSetup (const ExecFlagType &)
 
virtual void customSetup (const ExecFlagType &exec_type) override
 
const ExecFlagEnumgetExecuteOnEnum () const
 
FunctorReturnType< T, FET >::type genericEvaluate (const Space &r, const State &state) const
 
const MooseFunctorName & functorName () const
 
void setCacheClearanceSchedule (const std::set< ExecFlagType > &clearance_schedule)
 
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
 
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
 
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
 
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
 
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
 
const MooseArray< GenericReal< is_ad > > & genericDofValues () const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static std::string deduceFunctorName (const std::string &name, const InputParameters &params)
 

Public Attributes

const ConsoleStream _console
 

Protected Member Functions

bool isDirichletBoundaryFace (const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
 
ADReal getDirichletBoundaryFaceValue (const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
 
std::pair< bool, ADRealVectorValueelemIsUpwind (const Elem &elem, const FaceInfo &fi, const Moose::StateArg &time) const
 Checks to see whether the provided element is upwind of the provided face. More...
 
bool isFullyDevelopedFlowFace (const FaceInfo &fi) const
 Returns whether the passed-in FaceInfo corresponds to a fully-developed flow face. More...
 
void cacheSeparatorBoundaries ()
 Caches the separator boundaries. More...
 
void clearCaches ()
 
bool doDerivatives () const
 
virtual bool hasBlockMaterialPropertyHelper (const std::string &prop_name)
 
void initializeBlockRestrictable (const MooseObject *moose_object)
 
Moose::CoordinateSystemType getBlockCoordSystem ()
 
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 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 ElemSideQpArg &, const StateArg &) const
 
virtual DotType evaluateDot (const ElemPointArg &, const StateArg &) const
 
virtual DotType evaluateDot (const NodeArg &, 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
 
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 ValueType evaluate (const ElemArg &elem, const StateArg &state) const=0
 
virtual ValueType evaluate (const FaceArg &face, const StateArg &state) const=0
 
virtual ValueType evaluate (const ElemQpArg &qp, const StateArg &state) const=0
 
virtual ValueType evaluate (const ElemSideQpArg &side_qp, const StateArg &state) const=0
 
virtual ValueType evaluate (const ElemPointArg &elem_point, const StateArg &state) const=0
 
virtual ValueType evaluate (const NodeArg &node, const StateArg &state) const=0
 
const libMesh::NumericVector< libMesh::Number > & getSolution (const Moose::StateArg &state) const
 
std::string deduceFunctorName (const std::string &name) const
 
const Moose::Functor< T > & getFunctor (const std::string &name)
 
const Moose::Functor< T > & getFunctor (const std::string &name, THREAD_ID tid)
 
const Moose::Functor< T > & getFunctor (const std::string &name, SubProblem &subproblem)
 
const Moose::Functor< T > & getFunctor (const std::string &name, SubProblem &subproblem, THREAD_ID tid)
 
bool isFunctor (const std::string &name) const
 
bool isFunctor (const std::string &name, const SubProblem &subproblem) const
 
Moose::ElemArg makeElemArg (const Elem *elem, bool correct_skewnewss=false) const
 
void checkFunctorSupportsSideIntegration (const std::string &name, bool qp_integration)
 

Protected Attributes

const Moose::Functor< ADReal > * _u
 The x-component of velocity. More...
 
const Moose::Functor< ADReal > * _v
 The y-component of velocity. More...
 
const Moose::Functor< ADReal > * _w
 The z-component of velocity. More...
 
const Moose::Functor< ADReal > * _eps
 The porosity. More...
 
const Moose::Functor< ADReal > * _rho
 The density. More...
 
std::vector< BoundaryName > _pressure_drop_sidesets
 The names of the sidesets which will have associated form loss coefficients. More...
 
const std::vector< BoundaryID_pressure_drop_sideset_ids
 The boundary IDs corresponding to the form loss sidesets. More...
 
std::vector< Real_pressure_drop_form_factors
 The form loss coefficients corresponding to the sidesets. More...
 
std::unordered_map< BoundaryID, const INSFVHydraulicSeparatorInterface * > _boundary_id_to_separator
 A container for quick access of hydraulic separator BCs associated with this variable. More...
 
 usingMooseVariableFieldMembers
 
std::unique_ptr< MooseVariableDataFV< Real > > _element_data
 
std::unique_ptr< MooseVariableDataFV< Real > > _neighbor_data
 
std::unordered_map< const Elem *, VectorValue< ADReal > > _elem_to_grad
 
bool _two_term_boundary_expansion
 
VectorValue< ADReal_temp_cell_gradient
 
const bool _cache_cell_gradients
 
Moose::FV::InterpMethod _face_interp_method
 
SystemBase_sys
 
libMesh::FEType _fe_type
 
unsigned int _var_num
 
unsigned int _index
 
bool _is_eigen
 
Moose::VarKindType _var_kind
 
SubProblem_subproblem
 
const libMesh::Variable_variable
 
Assembly_assembly
 
const libMesh::DofMap_dof_map
 
std::vector< dof_id_type_dof_indices
 
MooseMesh_mesh
 
THREAD_ID _tid
 
const unsigned int _count
 
std::vector< Real_scaling_factor
 
std::string _var_name
 
bool _use_dual
 
bool _is_lower_d
 
std::vector< std::string > _array_var_component_names
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const MaterialData_blk_material_data
 
const ExecFlagEnum_execute_enum
 
const ExecFlagType_current_execute_flag
 
FEProblemBase_mci_feproblem
 
const TimeIntegrator *const _time_integrator
 
ADReal _ad_real_dummy
 
const Parallel::Communicator & _communicator
 

Private Attributes

const bool _allow_two_term_expansion_on_bernoulli_faces
 Switch to enable the two-term extrapolation on porosity jump faces. More...
 

Detailed Description

A special variable class for pressure which flags faces at which porosity jumps occur as extrapolated or Dirichlet boundary faces.

The downwind is flagged as an extrapolated boundary face while the upwind face is flagged as a Dirichlet face. The upwind Dirichlet value is computed using the downwind extrapolated pressure value and the Bernoulli equation

Definition at line 23 of file BernoulliPressureVariable.h.

Constructor & Destructor Documentation

◆ BernoulliPressureVariable()

BernoulliPressureVariable::BernoulliPressureVariable ( const InputParameters params)

Definition at line 52 of file BernoulliPressureVariable.C.

53  : INSFVPressureVariable(params),
54  ADFunctorInterface(this),
55  _u(nullptr),
56  _v(nullptr),
57  _w(nullptr),
58  _eps(nullptr),
59  _rho(nullptr),
60  _pressure_drop_sidesets(getParam<std::vector<BoundaryName>>("pressure_drop_sidesets")),
62  _pressure_drop_form_factors(getParam<std::vector<Real>>("pressure_drop_form_factors")),
64  getParam<bool>("allow_two_term_expansion_on_bernoulli_faces"))
65 {
68  paramError(
69  "allow_two_term_expansion_on_bernoulli_faces",
70  "Skewness correction with two term extrapolation on Bernoulli faces is not supported!");
71 
73  paramError("pressure_drop_sidesets",
74  "The number of sidesets and the number of supplied form losses are not equal!");
75 }
const Moose::Functor< ADReal > * _eps
The porosity.
const std::vector< BoundaryID > _pressure_drop_sideset_ids
The boundary IDs corresponding to the form loss sidesets.
Moose::FV::InterpMethod _face_interp_method
const Moose::Functor< ADReal > * _rho
The density.
INSFVPressureVariable(const InputParameters &params)
const Moose::Functor< ADReal > * _v
The y-component of velocity.
std::vector< Real > _pressure_drop_form_factors
The form loss coefficients corresponding to the sidesets.
const Moose::Functor< ADReal > * _u
The x-component of velocity.
const Moose::Functor< ADReal > * _w
The z-component of velocity.
const T & getParam(const std::string &name) const
const bool _allow_two_term_expansion_on_bernoulli_faces
Switch to enable the two-term extrapolation on porosity jump faces.
void paramError(const std::string &param, Args... args) const
ADFunctorInterface(const MooseObject *moose_object)
std::vector< BoundaryName > _pressure_drop_sidesets
The names of the sidesets which will have associated form loss coefficients.
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const

Member Function Documentation

◆ cacheSeparatorBoundaries()

void INSFVVariable::cacheSeparatorBoundaries ( )
protectedinherited

Caches the separator boundaries.

Definition at line 35 of file INSFVVariable.C.

Referenced by INSFVVariable::initialSetup(), INSFVVariable::meshChanged(), and INSFVVariable::timestepSetup().

36 {
38  std::vector<FVFluxBC *> bcs;
39 
40  const auto base_query = this->_subproblem.getMooseApp()
41  .theWarehouse()
42  .query()
43  .template condition<AttribSystem>("FVFluxBC")
44  .template condition<AttribThread>(_tid);
45 
46  for (const auto bnd_id : this->_mesh.getBoundaryIDs())
47  {
48  auto base_query_copy = base_query;
49  base_query_copy.template condition<AttribBoundaries>(std::set<BoundaryID>({bnd_id}))
50  .queryInto(bcs);
51  for (const auto bc : bcs)
52  {
53  const auto separator = dynamic_cast<const INSFVHydraulicSeparatorInterface *>(bc);
54  if (separator)
55  _boundary_id_to_separator.emplace(bnd_id, separator);
56  }
57  }
58 }
MooseApp & getMooseApp() const
std::unordered_map< BoundaryID, const INSFVHydraulicSeparatorInterface * > _boundary_id_to_separator
A container for quick access of hydraulic separator BCs associated with this variable.
Definition: INSFVVariable.h:56
SubProblem & _subproblem
Query query()
A base class which serves as a tag for hydraulic separators.
TheWarehouse & theWarehouse()
std::vector< BoundaryID > getBoundaryIDs(const Elem *const elem, const unsigned short int side) const

◆ computeElemValues()

void INSFVVariable::computeElemValues ( )
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 73 of file INSFVVariable.h.

74 {
75  if (_qp_calculations)
77  else
78  _element_data->setGeometry(Moose::Volume);
79 }
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62
std::unique_ptr< MooseVariableDataFV< Real > > _element_data
virtual void computeElemValues() override

◆ computeElemValuesFace()

void INSFVVariable::computeElemValuesFace ( )
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 82 of file INSFVVariable.h.

83 {
84  if (_qp_calculations)
86 }
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62
virtual void computeElemValuesFace() override

◆ computeFaceValues()

void INSFVVariable::computeFaceValues ( const FaceInfo fi)
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 66 of file INSFVVariable.h.

67 {
68  if (_qp_calculations)
70 }
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62
virtual void computeFaceValues(const FaceInfo &fi) override

◆ computeNeighborValues()

void INSFVVariable::computeNeighborValues ( )
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 96 of file INSFVVariable.h.

97 {
98  if (_qp_calculations)
100 }
virtual void computeNeighborValues() override
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62

◆ computeNeighborValuesFace()

void INSFVVariable::computeNeighborValuesFace ( )
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 89 of file INSFVVariable.h.

90 {
91  if (_qp_calculations)
93 }
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62
virtual void computeNeighborValuesFace() override

◆ elemIsUpwind()

std::pair< bool, ADRealVectorValue > BernoulliPressureVariable::elemIsUpwind ( const Elem &  elem,
const FaceInfo fi,
const Moose::StateArg time 
) const
protected

Checks to see whether the provided element is upwind of the provided face.

Parameters
elemthe element to check whether it is upwind of the face
fithe face
timeThe time at which to evaluate the velocity
Returns
a pair in which the first member is whether the element is upwind of the face and the second member is the evaluated face superficial velocity

Definition at line 90 of file BernoulliPressureVariable.C.

Referenced by getDirichletBoundaryFaceValue(), isDirichletBoundaryFace(), and isExtrapolatedBoundaryFace().

93 {
94  const Moose::FaceArg face{
95  &fi, Moose::FV::LimiterType::CentralDifference, true, false, nullptr, nullptr};
96 
97  const VectorValue<ADReal> vel_face{(*_u)(face, time), (*_v)(face, time), (*_w)(face, time)};
98  const bool fi_elem_is_upwind = vel_face * fi.normal() > 0;
99  const bool elem_is_upwind = &elem == &fi.elem() ? fi_elem_is_upwind : !fi_elem_is_upwind;
100 
101  return {elem_is_upwind, vel_face};
102 }
const Elem & elem() const
const Moose::Functor< ADReal > * _v
The y-component of velocity.
const Moose::Functor< ADReal > * _w
The z-component of velocity.
const Point & normal() const

◆ getDirichletBoundaryFaceValue()

ADReal BernoulliPressureVariable::getDirichletBoundaryFaceValue ( const FaceInfo fi,
const Elem *  elem,
const Moose::StateArg time 
) const
overrideprotectedvirtual

Reimplemented from MooseVariableFV< Real >.

Definition at line 147 of file BernoulliPressureVariable.C.

150 {
151  mooseAssert(isDirichletBoundaryFace(fi, elem, time), "This better be a Dirichlet face");
152 
155 
156  const auto [is_jump_face, eps_elem, eps_neighbor] = NS::isPorosityJumpFace(*_eps, fi, time);
157 #ifndef NDEBUG
158  mooseAssert(is_jump_face,
159  "If we are not a traditional Dirichlet face, then we must be a jump face");
160 #else
161  libmesh_ignore(is_jump_face);
162 #endif
163 
164  const Moose::FaceArg face_elem{
165  &fi, Moose::FV::LimiterType::CentralDifference, true, false, &fi.elem(), nullptr};
166  const Moose::FaceArg face_neighbor{
167  &fi, Moose::FV::LimiterType::CentralDifference, true, false, fi.neighborPtr(), nullptr};
168 
169  const auto [elem_is_upwind, vel_face] = elemIsUpwind(*elem, fi, time);
170 
171  const bool fi_elem_is_upwind = elem == &fi.elem() ? elem_is_upwind : !elem_is_upwind;
172  const auto & downwind_face = fi_elem_is_upwind ? face_neighbor : face_elem;
173 
174  const auto & upwind_elem = makeElemArg(fi_elem_is_upwind ? fi.elemPtr() : fi.neighborPtr());
175  const auto & downwind_elem = makeElemArg(fi_elem_is_upwind ? fi.neighborPtr() : fi.elemPtr());
176 
177  const auto & vel_upwind = vel_face;
178  const auto & vel_downwind = vel_face;
179 
180  const auto & eps_upwind = fi_elem_is_upwind ? eps_elem : eps_neighbor;
181  const auto & eps_downwind = fi_elem_is_upwind ? eps_neighbor : eps_elem;
182 
183  /*
184  If a weakly compressible material is used where the density slightly
185  depends on pressure. Given that the two-term expansion on Bernoulli faces is
186  enabled, we take use the downwind face assuming the continuity (minimal jump) of the
187  density. This protects against an infinite recursion between pressure and
188  density.
189  */
190  const auto & rho_upwind = (*_rho)(upwind_elem, time);
191  const auto & rho_downwind = (*_rho)(downwind_elem, time);
192 
193  const VectorValue<ADReal> interstitial_vel_upwind = vel_upwind * (1 / eps_upwind);
194  const VectorValue<ADReal> interstitial_vel_downwind = vel_downwind * (1 / eps_downwind);
195 
196  const auto v_dot_n_upwind = interstitial_vel_upwind * fi.normal();
197  const auto v_dot_n_downwind = interstitial_vel_downwind * fi.normal();
198 
199  // Iterate through sidesets to see if we have associated irreversible pressure losses
200  // or not.
201  ADReal factor = 0.0;
202  for (const auto & bd_id : fi.boundaryIDs())
203  for (const auto i : index_range(_pressure_drop_sideset_ids))
204  if (_pressure_drop_sideset_ids[i] == bd_id)
205  factor += _pressure_drop_form_factors[i];
206 
207  auto upwind_bernoulli_vel_chunk = 0.5 * rho_upwind * v_dot_n_upwind * v_dot_n_upwind;
208  auto downwind_bernoulli_vel_chunk = 0.5 * rho_downwind * v_dot_n_downwind * v_dot_n_downwind;
209 
210  // We add the additional, irreversible pressure loss here.
211  // If it is a contraction we have to use the downwind values, otherwise the upwind values.
212  if (eps_downwind < eps_upwind)
213  downwind_bernoulli_vel_chunk +=
214  0.5 * factor * rho_downwind * v_dot_n_downwind * v_dot_n_downwind;
215  else
216  upwind_bernoulli_vel_chunk += -0.5 * factor * rho_upwind * v_dot_n_upwind * v_dot_n_upwind;
217 
218  ADReal p_downwind;
220  p_downwind = (*this)(downwind_face, time);
221  else
222  p_downwind = (*this)(downwind_elem, time);
223 
224  return p_downwind + downwind_bernoulli_vel_chunk - upwind_bernoulli_vel_chunk;
225 }
const Moose::Functor< ADReal > * _eps
The porosity.
const std::vector< BoundaryID > _pressure_drop_sideset_ids
The boundary IDs corresponding to the form loss sidesets.
const std::set< BoundaryID > & boundaryIDs() const
const Elem & elem() const
std::tuple< bool, ADReal, ADReal > isPorosityJumpFace(const Moose::Functor< ADReal > &porosity, const FaceInfo &fi, const Moose::StateArg &time)
Checks to see whether the porosity value jumps from one side to the other of the provided face...
Definition: NSFVUtils.C:58
bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
std::vector< Real > _pressure_drop_form_factors
The form loss coefficients corresponding to the sidesets.
DualNumber< Real, DNDerivativeType, true > ADReal
virtual bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
Moose::ElemArg makeElemArg(const Elem *elem, bool correct_skewnewss=false) const
void libmesh_ignore(const Args &...)
const Elem * neighborPtr() const
const bool _allow_two_term_expansion_on_bernoulli_faces
Switch to enable the two-term extrapolation on porosity jump faces.
const Point & normal() const
std::pair< bool, ADRealVectorValue > elemIsUpwind(const Elem &elem, const FaceInfo &fi, const Moose::StateArg &time) const
Checks to see whether the provided element is upwind of the provided face.
const Elem * elemPtr() const
virtual ADReal getDirichletBoundaryFaceValue(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
auto index_range(const T &sizable)

◆ initialSetup()

void BernoulliPressureVariable::initialSetup ( )
overridevirtual

Reimplemented from INSFVVariable.

Definition at line 78 of file BernoulliPressureVariable.C.

79 {
80  _u = &getFunctor<ADReal>("u", _subproblem);
81  _v = &getFunctor<ADReal>("v", _subproblem);
82  _w = &getFunctor<ADReal>("w", _subproblem);
83  _eps = &getFunctor<ADReal>(NS::porosity, _subproblem);
84  _rho = &getFunctor<ADReal>(NS::density, _subproblem);
85 
87 }
const Moose::Functor< ADReal > * _eps
The porosity.
const Moose::Functor< ADReal > * _rho
The density.
static const std::string density
Definition: NS.h:33
const Moose::Functor< ADReal > * _v
The y-component of velocity.
static const std::string porosity
Definition: NS.h:104
const Moose::Functor< ADReal > * _u
The x-component of velocity.
virtual void initialSetup() override
Definition: INSFVVariable.C:75
SubProblem & _subproblem
const Moose::Functor< ADReal > * _w
The z-component of velocity.

◆ isDirichletBoundaryFace()

bool BernoulliPressureVariable::isDirichletBoundaryFace ( const FaceInfo fi,
const Elem *  elem,
const Moose::StateArg time 
) const
overrideprotectedvirtual

Reimplemented from MooseVariableFV< Real >.

Definition at line 126 of file BernoulliPressureVariable.C.

Referenced by getDirichletBoundaryFaceValue(), and isExtrapolatedBoundaryFace().

129 {
131  return true;
132 
133  if (isInternalFace(fi))
134  {
135  if (isSeparatorBoundary(fi))
136  return false;
137 
138  if (std::get<0>(NS::isPorosityJumpFace(*_eps, fi, time)))
139  // We choose to apply the Dirichlet condition for the upwind element
140  return std::get<0>(elemIsUpwind(*elem, fi, time));
141  }
142 
143  return false;
144 }
const Moose::Functor< ADReal > * _eps
The porosity.
std::tuple< bool, ADReal, ADReal > isPorosityJumpFace(const Moose::Functor< ADReal > &porosity, const FaceInfo &fi, const Moose::StateArg &time)
Checks to see whether the porosity value jumps from one side to the other of the provided face...
Definition: NSFVUtils.C:58
bool isSeparatorBoundary(const FaceInfo &fi) const
Definition: INSFVVariable.C:82
virtual bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &state) const
bool isInternalFace(const FaceInfo &) const
std::pair< bool, ADRealVectorValue > elemIsUpwind(const Elem &elem, const FaceInfo &fi, const Moose::StateArg &time) const
Checks to see whether the provided element is upwind of the provided face.

◆ isExtrapolatedBoundaryFace()

bool BernoulliPressureVariable::isExtrapolatedBoundaryFace ( const FaceInfo fi,
const Elem *  elem,
const Moose::StateArg time 
) const
overridevirtual

Reimplemented from MooseVariableFV< Real >.

Definition at line 105 of file BernoulliPressureVariable.C.

108 {
109  if (isDirichletBoundaryFace(fi, elem, time))
110  return false;
111  if (!isInternalFace(fi))
112  // We are neither a Dirichlet nor an internal face
113  return true;
114 
115  if (isSeparatorBoundary(fi))
116  return true;
117 
118  if (std::get<0>(NS::isPorosityJumpFace(*_eps, fi, time)))
119  // We choose to extrapolate for the element that is downwind
120  return !std::get<0>(elemIsUpwind(*elem, fi, time));
121 
122  return false;
123 }
const Moose::Functor< ADReal > * _eps
The porosity.
std::tuple< bool, ADReal, ADReal > isPorosityJumpFace(const Moose::Functor< ADReal > &porosity, const FaceInfo &fi, const Moose::StateArg &time)
Checks to see whether the porosity value jumps from one side to the other of the provided face...
Definition: NSFVUtils.C:58
bool isDirichletBoundaryFace(const FaceInfo &fi, const Elem *elem, const Moose::StateArg &time) const override
bool isSeparatorBoundary(const FaceInfo &fi) const
Definition: INSFVVariable.C:82
bool isInternalFace(const FaceInfo &) const
std::pair< bool, ADRealVectorValue > elemIsUpwind(const Elem &elem, const FaceInfo &fi, const Moose::StateArg &time) const
Checks to see whether the provided element is upwind of the provided face.

◆ isFullyDevelopedFlowFace()

bool INSFVVariable::isFullyDevelopedFlowFace ( const FaceInfo fi) const
protectedinherited

Returns whether the passed-in FaceInfo corresponds to a fully-developed flow face.

Definition at line 110 of file INSFVVariable.C.

Referenced by INSFVVelocityVariable::adGradSln(), INSFVVelocityVariable::getExtrapolatedBoundaryFaceValue(), and INSFVVelocityVariable::uncorrectedAdGradSln().

111 {
112  const auto & face_type = fi.faceType(std::make_pair(this->number(), this->sys().number()));
113 
114  mooseAssert(face_type != FaceInfo::VarFaceNeighbors::NEITHER,
115  "I'm concerned that if you're calling this method with a FaceInfo that doesn't have "
116  "this variable defined on either side, that you are doing something dangerous.");
117 
118  // If we're defined on both sides of the face, then we're not a boundary face
119  if (face_type == FaceInfo::VarFaceNeighbors::BOTH)
120  return false;
121 
122  std::vector<INSFVFullyDevelopedFlowBC *> bcs;
123 
124  this->_subproblem.getMooseApp()
125  .theWarehouse()
126  .query()
127  .template condition<AttribINSFVBCs>(INSFVBCs::INSFVFullyDevelopedFlowBC)
128  .template condition<AttribBoundaries>(fi.boundaryIDs())
129  .queryInto(bcs);
130 
131  return !bcs.empty();
132 }
const std::set< BoundaryID > & boundaryIDs() const
unsigned int number() const
MooseApp & getMooseApp() const
SubProblem & _subproblem
Query query()
TheWarehouse & theWarehouse()
VarFaceNeighbors faceType(const std::pair< unsigned int, unsigned int > &var_sys) const

◆ isSeparatorBoundary()

bool INSFVVariable::isSeparatorBoundary ( const FaceInfo fi) const
inherited

Definition at line 82 of file INSFVVariable.C.

Referenced by isDirichletBoundaryFace(), isExtrapolatedBoundaryFace(), and INSFVVariable::isExtrapolatedBoundaryFace().

83 {
84  for (const auto bid : fi.boundaryIDs())
85  if (_boundary_id_to_separator.count(bid))
86  return true;
87 
88  return false;
89 }
const std::set< BoundaryID > & boundaryIDs() const
std::unordered_map< BoundaryID, const INSFVHydraulicSeparatorInterface * > _boundary_id_to_separator
A container for quick access of hydraulic separator BCs associated with this variable.
Definition: INSFVVariable.h:56

◆ meshChanged()

void INSFVVariable::meshChanged ( )
overridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 68 of file INSFVVariable.C.

69 {
72 }
virtual void meshChanged() override
void cacheSeparatorBoundaries()
Caches the separator boundaries.
Definition: INSFVVariable.C:35

◆ requireQpComputations()

void INSFVVariable::requireQpComputations ( ) const
inlineoverridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Reimplemented in INSFVPressureNoQpComputation, and INSFVVelocityNoQpComputation.

Definition at line 30 of file INSFVVariable.h.

30 { _qp_calculations = true; }
bool _qp_calculations
Whether to pre-initialize variable data for use in traditional MOOSE quadrature point based objects...
Definition: INSFVVariable.h:62

◆ timestepSetup()

void INSFVVariable::timestepSetup ( )
overridevirtualinherited

Reimplemented from MooseVariableFV< Real >.

Definition at line 61 of file INSFVVariable.C.

62 {
65 }
virtual void timestepSetup() override
void cacheSeparatorBoundaries()
Caches the separator boundaries.
Definition: INSFVVariable.C:35

◆ validParams()

InputParameters BernoulliPressureVariable::validParams ( )
static

Definition at line 17 of file BernoulliPressureVariable.C.

18 {
19  auto params = INSFVPressureVariable::validParams();
21  params.addRequiredParam<MooseFunctorName>("u", "The x-component of velocity");
22  params.addParam<MooseFunctorName>("v", 0, "The y-component of velocity");
23  params.addParam<MooseFunctorName>("w", 0, "The z-component of velocity");
24  params.addRequiredParam<MooseFunctorName>(NS::porosity, "The porosity");
25  params.addRequiredParam<MooseFunctorName>(NS::density, "The density");
26  params.addParam<std::vector<BoundaryName>>(
27  "pressure_drop_sidesets", {}, "Sidesets over which form loss coefficients are to be applied");
28  params.addParam<std::vector<Real>>(
29  "pressure_drop_form_factors",
30  {},
31  "User-supplied form loss coefficients to be applied over the sidesets listed above");
32  params.addParam<bool>(
33  "allow_two_term_expansion_on_bernoulli_faces",
34  false,
35  "Switch to enable the two-term extrapolation on porosity jump faces. "
36  "WARNING: This might lead to crushes in parallel runs if porosity jump faces are connected "
37  "with one cell (usually corners) due to the insufficient number of ghosted "
38  "layers.");
39 
40  // Assuming that this variable is used for advection problems, due to the
41  // utilization of the pressure gradient in the advecting velocity
42  // through the Rhie-Chow interpolation, we have to extend the ghosted layers.
43  params.addRelationshipManager(
44  "ElementSideNeighborLayers",
47  [](const InputParameters & /*obj_params*/, InputParameters & rm_params)
48  { rm_params.set<unsigned short>("layers") = 3; });
49  return params;
50 }
T & set(const std::string &name, bool quiet_mode=false)
static InputParameters validParams()
static const std::string density
Definition: NS.h:33
static InputParameters validParams()
static const std::string porosity
Definition: NS.h:104

Member Data Documentation

◆ _allow_two_term_expansion_on_bernoulli_faces

const bool BernoulliPressureVariable::_allow_two_term_expansion_on_bernoulli_faces
private

Switch to enable the two-term extrapolation on porosity jump faces.

Definition at line 78 of file BernoulliPressureVariable.h.

Referenced by BernoulliPressureVariable(), and getDirichletBoundaryFaceValue().

◆ _boundary_id_to_separator

std::unordered_map<BoundaryID, const INSFVHydraulicSeparatorInterface *> INSFVVariable::_boundary_id_to_separator
protectedinherited

A container for quick access of hydraulic separator BCs associated with this variable.

Definition at line 56 of file INSFVVariable.h.

Referenced by INSFVVariable::cacheSeparatorBoundaries(), and INSFVVariable::isSeparatorBoundary().

◆ _eps

const Moose::Functor<ADReal>* BernoulliPressureVariable::_eps
protected

◆ _pressure_drop_form_factors

std::vector<Real> BernoulliPressureVariable::_pressure_drop_form_factors
protected

The form loss coefficients corresponding to the sidesets.

Definition at line 74 of file BernoulliPressureVariable.h.

Referenced by BernoulliPressureVariable(), and getDirichletBoundaryFaceValue().

◆ _pressure_drop_sideset_ids

const std::vector<BoundaryID> BernoulliPressureVariable::_pressure_drop_sideset_ids
protected

The boundary IDs corresponding to the form loss sidesets.

Definition at line 71 of file BernoulliPressureVariable.h.

Referenced by getDirichletBoundaryFaceValue().

◆ _pressure_drop_sidesets

std::vector<BoundaryName> BernoulliPressureVariable::_pressure_drop_sidesets
protected

The names of the sidesets which will have associated form loss coefficients.

Definition at line 68 of file BernoulliPressureVariable.h.

Referenced by BernoulliPressureVariable().

◆ _rho

const Moose::Functor<ADReal>* BernoulliPressureVariable::_rho
protected

The density.

Definition at line 65 of file BernoulliPressureVariable.h.

Referenced by initialSetup().

◆ _u

const Moose::Functor<ADReal>* BernoulliPressureVariable::_u
protected

The x-component of velocity.

Definition at line 57 of file BernoulliPressureVariable.h.

Referenced by initialSetup().

◆ _v

const Moose::Functor<ADReal>* BernoulliPressureVariable::_v
protected

The y-component of velocity.

Definition at line 59 of file BernoulliPressureVariable.h.

Referenced by elemIsUpwind(), and initialSetup().

◆ _w

const Moose::Functor<ADReal>* BernoulliPressureVariable::_w
protected

The z-component of velocity.

Definition at line 61 of file BernoulliPressureVariable.h.

Referenced by elemIsUpwind(), and initialSetup().


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