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

Upwind numerical flux scheme for the advection equation using a cell-centered finite volume method. More...

#include <AEFVUpwindInternalSideFlux.h>

Inheritance diagram for AEFVUpwindInternalSideFlux:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 AEFVUpwindInternalSideFlux (const InputParameters &parameters)
 
virtual ~AEFVUpwindInternalSideFlux ()
 
virtual void calcFlux (unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave, std::vector< Real > &flux) const override
 Solve the Riemann problem. More...
 
virtual void calcJacobian (unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave, DenseMatrix< Real > &jac1, DenseMatrix< Real > &jac2) const override
 Compute the Jacobian matrix. More...
 
virtual void execute () override
 
virtual void initialize () override
 
virtual void finalize () override
 
virtual void threadJoin (const UserObject &) override
 
virtual const std::vector< Real > & getFlux (unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave) const
 Get the flux vector. More...
 
virtual const DenseMatrix< Real > & getJacobian (Moose::DGResidualType type, unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave) const
 Get the Jacobian matrix. More...
 
virtual void subdomainSetup () override
 
bool needThreadedCopy () const override final
 
SubProblemgetSubProblem () const
 
bool shouldDuplicateInitialExecution () const
 
virtual Real spatialValue (const Point &) const
 
virtual const std::vector< Point > spatialPoints () const
 
void gatherSum (T &value)
 
void gatherMax (T &value)
 
void gatherMin (T &value)
 
void gatherProxyValueMax (T1 &proxy, T2 &value)
 
void gatherProxyValueMin (T1 &proxy, T2 &value)
 
void setPrimaryThreadCopy (UserObject *primary)
 
UserObjectprimaryThreadCopy ()
 
std::set< UserObjectName > getDependObjects () const
 
const std::set< std::string > & getRequestedItems () override
 
const std::set< std::string > & getSuppliedItems () override
 
unsigned int systemNumber () const
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () 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
 
virtual void initialSetup ()
 
virtual void timestepSetup ()
 
virtual void jacobianSetup ()
 
virtual void residualSetup ()
 
virtual void customSetup (const ExecFlagType &)
 
const ExecFlagEnumgetExecuteOnEnum () const
 
UserObjectName getUserObjectName (const std::string &param_name) const
 
const T & getUserObject (const std::string &param_name, bool is_dependency=true) const
 
const T & getUserObjectByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBase (const std::string &param_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBaseByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const std::vector< MooseVariableScalar *> & getCoupledMooseScalarVars ()
 
const std::set< TagID > & getScalarVariableCoupleableVectorTags () const
 
const std::set< TagID > & getScalarVariableCoupleableMatrixTags () const
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, MaterialData &material_data, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, MaterialData &material_data, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, const unsigned int state=0)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name, MaterialData &material_data)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data, const unsigned int state)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name)
 
std::pair< const MaterialProperty< T > *, std::set< SubdomainID > > getBlockMaterialProperty (const MaterialPropertyName &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialProperty ()
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialPropertyByName (const std::string &prop_name)
 
const MaterialProperty< T > & getZeroMaterialProperty (Ts... args)
 
std::set< SubdomainIDgetMaterialPropertyBlocks (const std::string &name)
 
std::vector< SubdomainName > getMaterialPropertyBlockNames (const std::string &name)
 
std::set< BoundaryIDgetMaterialPropertyBoundaryIDs (const std::string &name)
 
std::vector< BoundaryName > getMaterialPropertyBoundaryNames (const std::string &name)
 
void checkBlockAndBoundaryCompatibility (std::shared_ptr< MaterialBase > discrete)
 
std::unordered_map< SubdomainID, std::vector< MaterialBase *> > buildRequiredMaterials (bool allow_stateful=true)
 
void statefulPropertiesAllowed (bool)
 
bool getMaterialPropertyCalled () const
 
virtual const std::unordered_set< unsigned int > & getMatPropDependencies () const
 
virtual void resolveOptionalProperties ()
 
const GenericMaterialProperty< T, is_ad > & getPossiblyConstantGenericMaterialPropertyByName (const MaterialPropertyName &prop_name, MaterialData &material_data, const unsigned int state)
 
bool isImplicit ()
 
Moose::StateArg determineState () const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
const GenericOptionalMaterialProperty< T, is_ad > & getGenericOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const GenericOptionalMaterialProperty< T, is_ad > & getGenericOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalMaterialProperty< T > & getOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalMaterialProperty< T > & getOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalADMaterialProperty< T > & getOptionalADMaterialProperty (const std::string &name)
 
const OptionalADMaterialProperty< T > & getOptionalADMaterialProperty (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOld (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOld (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOlder (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOlder (const std::string &name)
 
MaterialBasegetMaterial (const std::string &name)
 
MaterialBasegetMaterial (const std::string &name)
 
MaterialBasegetMaterialByName (const std::string &name, bool no_warn=false)
 
MaterialBasegetMaterialByName (const std::string &name, bool no_warn=false)
 
bool hasMaterialProperty (const std::string &name)
 
bool hasMaterialProperty (const std::string &name)
 
bool hasMaterialPropertyByName (const std::string &name)
 
bool hasMaterialPropertyByName (const std::string &name)
 
bool hasADMaterialProperty (const std::string &name)
 
bool hasADMaterialProperty (const std::string &name)
 
bool hasADMaterialPropertyByName (const std::string &name)
 
bool hasADMaterialPropertyByName (const std::string &name)
 
bool hasGenericMaterialProperty (const std::string &name)
 
bool hasGenericMaterialProperty (const std::string &name)
 
bool hasGenericMaterialPropertyByName (const std::string &name)
 
bool hasGenericMaterialPropertyByName (const std::string &name)
 
const FunctiongetFunction (const std::string &name) const
 
const FunctiongetFunctionByName (const FunctionName &name) const
 
bool hasFunction (const std::string &param_name) const
 
bool hasFunctionByName (const FunctionName &name) const
 
bool isDefaultPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessor (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessorByName (const PostprocessorName &name) const
 
std::size_t coupledPostprocessors (const std::string &param_name) const
 
const PostprocessorName & getPostprocessorName (const std::string &param_name, const unsigned int index=0) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name) const
 
const VectorPostprocessorName & getVectorPostprocessorName (const std::string &param_name) const
 
T & getSampler (const std::string &name)
 
SamplergetSampler (const std::string &name)
 
T & getSamplerByName (const SamplerName &name)
 
SamplergetSamplerByName (const SamplerName &name)
 
virtual void meshChanged ()
 
virtual void meshDisplaced ()
 
PerfGraphperfGraph ()
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) 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 void sort (typename std::vector< T > &vector)
 
static void sortDFS (typename std::vector< T > &vector)
 
static void cyclicDependencyError (CyclicDependencyException< T2 > &e, const std::string &header)
 

Public Attributes

const ConsoleStream _console
 

Static Public Attributes

static constexpr PropertyValue::id_type default_property_id
 
static constexpr PropertyValue::id_type zero_property_id
 
static constexpr auto SYSTEM
 
static constexpr auto NAME
 

Protected Member Functions

virtual void addPostprocessorDependencyHelper (const PostprocessorName &name) const override
 
virtual void addVectorPostprocessorDependencyHelper (const VectorPostprocessorName &name) const override
 
virtual void addUserObjectDependencyHelper (const UserObject &uo) const override
 
void addReporterDependencyHelper (const ReporterName &reporter_name) override
 
const ReporterNamegetReporterName (const std::string &param_name) const
 
T & declareRestartableData (const std::string &data_name, Args &&... args)
 
ManagedValue< T > declareManagedRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
const T & getRestartableData (const std::string &data_name) const
 
T & declareRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
T & declareRecoverableData (const std::string &data_name, Args &&... args)
 
T & declareRestartableDataWithObjectName (const std::string &data_name, const std::string &object_name, Args &&... args)
 
T & declareRestartableDataWithObjectNameWithContext (const std::string &data_name, const std::string &object_name, void *context, Args &&... args)
 
std::string restartableName (const std::string &data_name) const
 
const T & getMeshProperty (const std::string &data_name, const std::string &prefix)
 
const T & getMeshProperty (const std::string &data_name)
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name) const
 
bool hasMeshProperty (const std::string &data_name) const
 
std::string meshPropertyName (const std::string &data_name) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 
bool isCoupledScalar (const std::string &var_name, unsigned int i=0) const
 
unsigned int coupledScalarComponents (const std::string &var_name) const
 
unsigned int coupledScalar (const std::string &var_name, unsigned int comp=0) const
 
libMesh::Order coupledScalarOrder (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const ADVariableValueadCoupledScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const GenericVariableValue< is_ad > & coupledGenericScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const GenericVariableValue< false > & coupledGenericScalarValue (const std::string &var_name, const unsigned int comp) const
 
const GenericVariableValue< true > & coupledGenericScalarValue (const std::string &var_name, const unsigned int comp) const
 
const VariableValuecoupledVectorTagScalarValue (const std::string &var_name, TagID tag, unsigned int comp=0) const
 
const VariableValuecoupledMatrixTagScalarValue (const std::string &var_name, TagID tag, unsigned int comp=0) const
 
const VariableValuecoupledScalarValueOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarValueOlder (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDot (const std::string &var_name, unsigned int comp=0) const
 
const ADVariableValueadCoupledScalarDot (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDot (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDotOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDu (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDotDu (const std::string &var_name, unsigned int comp=0) const
 
const MooseVariableScalargetScalarVar (const std::string &var_name, unsigned int comp) const
 
virtual void checkMaterialProperty (const std::string &name, const unsigned int state)
 
void markMatPropRequested (const std::string &)
 
MaterialPropertyName getMaterialPropertyName (const std::string &name) const
 
void checkExecutionStage ()
 
const T & getReporterValue (const std::string &param_name, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, ReporterMode mode, const std::size_t time_index=0)
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
const GenericMaterialProperty< T, is_ad > * defaultGenericMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > * defaultGenericMaterialProperty (const std::string &name)
 
const MaterialProperty< T > * defaultMaterialProperty (const std::string &name)
 
const MaterialProperty< T > * defaultMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > * defaultADMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > * defaultADMaterialProperty (const std::string &name)
 

Static Protected Member Functions

static std::string meshPropertyName (const std::string &data_name, const std::string &prefix)
 

Protected Attributes

const Real _velocity
 advective velocity More...
 
unsigned int _cached_flux_elem_id
 element ID of the cached flux values More...
 
unsigned int _cached_flux_neig_id
 neighbor element ID of the cached flux values More...
 
unsigned int _cached_jacobian_elem_id
 element ID of the cached Jacobian values More...
 
unsigned int _cached_jacobian_neig_id
 neighbor element ID of the cached Jacobian values More...
 
std::vector< Real_flux
 flux vector of this side More...
 
DenseMatrix< Real_jac1
 Jacobian matrix contribution to the "left" cell. More...
 
DenseMatrix< Real_jac2
 Jacobian matrix contribution to the "right" cell. More...
 
SubProblem_subproblem
 
FEProblemBase_fe_problem
 
SystemBase_sys
 
const THREAD_ID _tid
 
Assembly_assembly
 
const Moose::CoordinateSystemType_coord_sys
 
const bool _duplicate_initial_execution
 
std::set< std::string > _depend_uo
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const ExecFlagEnum_execute_enum
 
const ExecFlagType_current_execute_flag
 
MooseApp_restartable_app
 
const std::string _restartable_system_name
 
const THREAD_ID _restartable_tid
 
const bool _restartable_read_only
 
FEProblemBase_mci_feproblem
 
FEProblemBase_mdi_feproblem
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
FEProblemBase_sc_fe_problem
 
const THREAD_ID _sc_tid
 
const Real_real_zero
 
const VariableValue_scalar_zero
 
const Point & _point_zero
 
const InputParameters_mi_params
 
const std::string _mi_name
 
const MooseObjectName _mi_moose_object_name
 
FEProblemBase_mi_feproblem
 
SubProblem_mi_subproblem
 
const THREAD_ID _mi_tid
 
const Moose::MaterialDataType _material_data_type
 
MaterialData_material_data
 
bool _stateful_allowed
 
bool _get_material_property_called
 
std::vector< std::unique_ptr< PropertyValue > > _default_properties
 
std::unordered_set< unsigned int_material_property_dependencies
 
const MaterialPropertyName _get_suffix
 
const bool _use_interpolated_state
 
const InputParameters_ti_params
 
FEProblemBase_ti_feproblem
 
bool _is_implicit
 
Real_t
 
const Real_t_old
 
int_t_step
 
Real_dt
 
Real_dt_old
 
bool _is_transient
 
const Parallel::Communicator & _communicator
 

Static Protected Attributes

static const std::string _interpolated_old
 
static const std::string _interpolated_older
 

Detailed Description

Upwind numerical flux scheme for the advection equation using a cell-centered finite volume method.

Definition at line 21 of file AEFVUpwindInternalSideFlux.h.

Constructor & Destructor Documentation

◆ AEFVUpwindInternalSideFlux()

AEFVUpwindInternalSideFlux::AEFVUpwindInternalSideFlux ( const InputParameters parameters)

Definition at line 24 of file AEFVUpwindInternalSideFlux.C.

25  : InternalSideFluxBase(parameters), _velocity(getParam<Real>("velocity"))
26 {
27 }
InternalSideFluxBase(const InputParameters &parameters)
const Real _velocity
advective velocity
const InputParameters & parameters() const

◆ ~AEFVUpwindInternalSideFlux()

AEFVUpwindInternalSideFlux::~AEFVUpwindInternalSideFlux ( )
virtual

Definition at line 29 of file AEFVUpwindInternalSideFlux.C.

29 {}

Member Function Documentation

◆ calcFlux()

void AEFVUpwindInternalSideFlux::calcFlux ( unsigned int  iside,
dof_id_type  ielem,
dof_id_type  ineig,
const std::vector< Real > &  uvec1,
const std::vector< Real > &  uvec2,
const RealVectorValue dwave,
std::vector< Real > &  flux 
) const
overridevirtual

Solve the Riemann problem.

Parameters
[in]isidelocal index of current side
[in]ielemglobal index of the current element
[in]ineigglobal index of the neighbor element
[in]uvec1vector of variables on the "left"
[in]uvec2vector of variables on the "right"
[in]dwavevector of unit normal
[out]fluxflux vector across the side

Implements InternalSideFluxBase.

Definition at line 32 of file AEFVUpwindInternalSideFlux.C.

39 {
40  mooseAssert(uvec1.size() == 1, "Invalid size for uvec1. Must be single variable coupling.");
41  mooseAssert(uvec2.size() == 1, "Invalid size for uvec1. Must be single variable coupling.");
42 
43  // assign the size of flux vector, e.g. = 1 for the advection equation
44  flux.resize(1);
45 
46  // assume a constant velocity on the left
47  RealVectorValue uadv1(_velocity, 0.0, 0.0);
48 
49  // assume a constant velocity on the right
50  RealVectorValue uadv2(_velocity, 0.0, 0.0);
51 
52  // normal velocity on the left and right
53  Real vdon1 = uadv1 * dwave;
54  Real vdon2 = uadv2 * dwave;
55 
56  // calculate the so-called a^plus and a^minus
57  Real aplus = 0.5 * (vdon1 + std::abs(vdon1));
58  Real amins = 0.5 * (vdon2 - std::abs(vdon2));
59 
60  // finally calculate the flux
61  flux[0] = aplus * uvec1[0] + amins * uvec2[0];
62 }
const Real _velocity
advective velocity
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ calcJacobian()

void AEFVUpwindInternalSideFlux::calcJacobian ( unsigned int  iside,
dof_id_type  ielem,
dof_id_type  ineig,
const std::vector< Real > &  uvec1,
const std::vector< Real > &  uvec2,
const RealVectorValue dwave,
DenseMatrix< Real > &  jac1,
DenseMatrix< Real > &  jac2 
) const
overridevirtual

Compute the Jacobian matrix.

Parameters
[in]isidelocal index of current side
[in]ielemglobal index of the current element
[in]ineigglobal index of the neighbor element
[in]uvec1vector of variables on the "left"
[in]uvec2vector of variables on the "right"
[in]dwavevector of unit normal
[out]jac1Jacobian matrix contribution to the "left" cell
[out]jac2Jacobian matrix contribution to the "right" cell

Implements InternalSideFluxBase.

Definition at line 65 of file AEFVUpwindInternalSideFlux.C.

73 {
74  mooseAssert(uvec1.size() == 1, "Invalid size for uvec1. Must be single variable coupling.");
75  mooseAssert(uvec2.size() == 1, "Invalid size for uvec1. Must be single variable coupling.");
76 
77  // assign the size of Jacobian matrix, e.g. = (1, 1) for the advection equation
78  jac1.resize(1, 1);
79  jac2.resize(1, 1);
80 
81  // assume a constant velocity on the left
82  RealVectorValue uadv1(_velocity, 0.0, 0.0);
83 
84  // assume a constant velocity on the right
85  RealVectorValue uadv2(_velocity, 0.0, 0.0);
86 
87  // normal velocity on the left and right
88  Real vdon1 = uadv1 * dwave;
89  Real vdon2 = uadv2 * dwave;
90 
91  // calculate the so-called a^plus and a^minus
92  Real aplus = 0.5 * (vdon1 + std::abs(vdon1));
93  Real amins = 0.5 * (vdon2 - std::abs(vdon2));
94 
95  // finally calculate the Jacobian matrix
96  jac1(0, 0) = aplus;
97  jac2(0, 0) = amins;
98 }
const Real _velocity
advective velocity
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void resize(const unsigned int new_m, const unsigned int new_n)

◆ execute()

void InternalSideFluxBase::execute ( )
overridevirtualinherited

Implements ThreadedGeneralUserObject.

Definition at line 39 of file InternalSideFluxBase.C.

40 {
41 }

◆ finalize()

void InternalSideFluxBase::finalize ( )
overridevirtualinherited

Implements ThreadedGeneralUserObject.

Definition at line 44 of file InternalSideFluxBase.C.

45 {
46 }

◆ getFlux()

const std::vector< Real > & InternalSideFluxBase::getFlux ( unsigned int  iside,
dof_id_type  ielem,
dof_id_type  ineig,
const std::vector< Real > &  uvec1,
const std::vector< Real > &  uvec2,
const RealVectorValue dwave 
) const
virtualinherited

Get the flux vector.

Parameters
[in]isidelocal index of current side
[in]ielemglobal index of the current element
[in]ineigglobal index of the neighbor element
[in]uvec1vector of variables on the "left"
[in]uvec2vector of variables on the "right"
[in]dwavevector of unit normal

Definition at line 54 of file InternalSideFluxBase.C.

Referenced by AEFVKernel::computeQpResidual().

60 {
61  if (_cached_flux_elem_id != ielem || _cached_flux_neig_id != ineig)
62  {
63  _cached_flux_elem_id = ielem;
64  _cached_flux_neig_id = ineig;
65 
66  calcFlux(iside, ielem, ineig, uvec1, uvec2, dwave, _flux);
67  }
68  return _flux;
69 }
unsigned int _cached_flux_elem_id
element ID of the cached flux values
virtual void calcFlux(unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave, std::vector< Real > &flux) const =0
Solve the Riemann problem.
std::vector< Real > _flux
flux vector of this side
unsigned int _cached_flux_neig_id
neighbor element ID of the cached flux values

◆ getJacobian()

const DenseMatrix< Real > & InternalSideFluxBase::getJacobian ( Moose::DGResidualType  type,
unsigned int  iside,
dof_id_type  ielem,
dof_id_type  ineig,
const std::vector< Real > &  uvec1,
const std::vector< Real > &  uvec2,
const RealVectorValue dwave 
) const
virtualinherited

Get the Jacobian matrix.

Parameters
[in]isidelocal index of current side
[in]ielemglobal index of the current element
[in]ineigglobal index of the neighbor element
[in]uvec1vector of variables on the "left"
[in]uvec2vector of variables on the "right"
[in]dwavevector of unit normal

Definition at line 72 of file InternalSideFluxBase.C.

Referenced by AEFVKernel::computeQpJacobian().

79 {
80  if (_cached_jacobian_elem_id != ielem || _cached_jacobian_neig_id != ineig)
81  {
84 
85  calcJacobian(iside, ielem, ineig, uvec1, uvec2, dwave, _jac1, _jac2);
86  }
87 
88  if (type == Moose::Element)
89  return _jac1;
90  else
91  return _jac2;
92 }
virtual void calcJacobian(unsigned int iside, dof_id_type ielem, dof_id_type ineig, const std::vector< Real > &uvec1, const std::vector< Real > &uvec2, const RealVectorValue &dwave, DenseMatrix< Real > &jac1, DenseMatrix< Real > &jac2) const =0
Compute the Jacobian matrix.
unsigned int _cached_jacobian_neig_id
neighbor element ID of the cached Jacobian values
DenseMatrix< Real > _jac2
Jacobian matrix contribution to the "right" cell.
const std::string & type() const
unsigned int _cached_jacobian_elem_id
element ID of the cached Jacobian values
DenseMatrix< Real > _jac1
Jacobian matrix contribution to the "left" cell.

◆ initialize()

void InternalSideFluxBase::initialize ( )
overridevirtualinherited

Implements ThreadedGeneralUserObject.

Definition at line 30 of file InternalSideFluxBase.C.

31 {
36 }
unsigned int _cached_flux_elem_id
element ID of the cached flux values
const unsigned int invalid_uint
unsigned int _cached_jacobian_neig_id
neighbor element ID of the cached Jacobian values
unsigned int _cached_jacobian_elem_id
element ID of the cached Jacobian values
unsigned int _cached_flux_neig_id
neighbor element ID of the cached flux values

◆ threadJoin()

void InternalSideFluxBase::threadJoin ( const UserObject )
overridevirtualinherited

Reimplemented from ThreadedGeneralUserObject.

Definition at line 49 of file InternalSideFluxBase.C.

50 {
51 }

◆ validParams()

InputParameters AEFVUpwindInternalSideFlux::validParams ( )
static

Definition at line 15 of file AEFVUpwindInternalSideFlux.C.

16 {
18  params.addClassDescription("Upwind numerical flux scheme for the advection equation using a "
19  "cell-centered finite volume method.");
20  params.addParam<Real>("velocity", 1.0, "Advective velocity");
21  return params;
22 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)
static InputParameters validParams()

Member Data Documentation

◆ _cached_flux_elem_id

unsigned int InternalSideFluxBase::_cached_flux_elem_id
mutableprotectedinherited

element ID of the cached flux values

Definition at line 112 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getFlux(), and InternalSideFluxBase::initialize().

◆ _cached_flux_neig_id

unsigned int InternalSideFluxBase::_cached_flux_neig_id
mutableprotectedinherited

neighbor element ID of the cached flux values

Definition at line 114 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getFlux(), and InternalSideFluxBase::initialize().

◆ _cached_jacobian_elem_id

unsigned int InternalSideFluxBase::_cached_jacobian_elem_id
mutableprotectedinherited

element ID of the cached Jacobian values

Definition at line 117 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getJacobian(), and InternalSideFluxBase::initialize().

◆ _cached_jacobian_neig_id

unsigned int InternalSideFluxBase::_cached_jacobian_neig_id
mutableprotectedinherited

neighbor element ID of the cached Jacobian values

Definition at line 119 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getJacobian(), and InternalSideFluxBase::initialize().

◆ _flux

std::vector<Real> InternalSideFluxBase::_flux
mutableprotectedinherited

flux vector of this side

Definition at line 122 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getFlux().

◆ _jac1

DenseMatrix<Real> InternalSideFluxBase::_jac1
mutableprotectedinherited

Jacobian matrix contribution to the "left" cell.

Definition at line 124 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getJacobian().

◆ _jac2

DenseMatrix<Real> InternalSideFluxBase::_jac2
mutableprotectedinherited

Jacobian matrix contribution to the "right" cell.

Definition at line 126 of file InternalSideFluxBase.h.

Referenced by InternalSideFluxBase::getJacobian().

◆ _velocity

const Real AEFVUpwindInternalSideFlux::_velocity
protected

advective velocity

Definition at line 48 of file AEFVUpwindInternalSideFlux.h.

Referenced by calcFlux(), and calcJacobian().


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