OGS
|
Classes | |
struct | IntegrationPointData |
struct | LiquidFlowData |
class | LiquidFlowLocalAssembler |
class | LiquidFlowLocalAssemblerInterface |
class | LiquidFlowProcess |
A class to simulate the single phase flow process in porous media described by the governing equation: More... | |
Enumerations | |
enum class | EquationBalanceType { volume , mass } |
Governing equation balance type. More... | |
Functions | |
void | checkMPLProperties (MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map, bool const is_equation_type_volume) |
EquationBalanceType | covertEquationBalanceTypeFromString (std::string_view const type_in_string) |
std::unique_ptr< Process > | createLiquidFlowProcess (std::string const &name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< ProcessVariable > const &variables, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const ¶meters, unsigned const integration_order, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) |
std::tuple< double, double > | getFluidDensityAndViscosity (double const t, double const dt, ParameterLib::SpatialPosition const &pos, MaterialPropertyLib::Phase const &fluid_phase, MaterialPropertyLib::VariableArray &vars) |
Variables | |
const unsigned | NUM_NODAL_DOF = 1 |
|
strong |
Governing equation balance type.
Enumerator | |
---|---|
volume | |
mass |
Definition at line 26 of file LiquidFlowData.h.
void ProcessLib::LiquidFlow::checkMPLProperties | ( | MeshLib::Mesh const & | mesh, |
MaterialPropertyLib::MaterialSpatialDistributionMap const & | media_map, | ||
bool const | is_equation_type_volume ) |
Definition at line 32 of file CreateLiquidFlowProcess.cpp.
References MaterialPropertyLib::checkMPLPhasesForSinglePhaseFlow(), DBUG(), MaterialPropertyLib::density, MaterialPropertyLib::MaterialSpatialDistributionMap::media(), OGS_FATAL, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::storage, and MaterialPropertyLib::viscosity.
Referenced by createLiquidFlowProcess().
EquationBalanceType ProcessLib::LiquidFlow::covertEquationBalanceTypeFromString | ( | std::string_view const | type_in_string | ) |
Definition at line 77 of file CreateLiquidFlowProcess.cpp.
References mass, OGS_FATAL, and volume.
Referenced by createLiquidFlowProcess().
std::unique_ptr< Process > ProcessLib::LiquidFlow::createLiquidFlowProcess | ( | std::string const & | name, |
MeshLib::Mesh & | mesh, | ||
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && | jacobian_assembler, | ||
std::vector< ProcessVariable > const & | variables, | ||
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & | parameters, | ||
unsigned const | integration_order, | ||
BaseLib::ConfigTree const & | config, | ||
std::vector< std::unique_ptr< MeshLib::Mesh > > const & | meshes, | ||
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & | media ) |
Primary process variables as they appear in the global component vector:
Definition at line 93 of file CreateLiquidFlowProcess.cpp.
References BaseLib::ConfigTree::checkConfigParameter(), checkMPLProperties(), ProcessLib::Process::constant_one_parameter_name, covertEquationBalanceTypeFromString(), MaterialPropertyLib::createMaterialSpatialDistributionMap(), ProcessLib::createSecondaryVariables(), ProcessLib::SurfaceFluxData::createSurfaceFluxData(), DBUG(), ParameterLib::findParameter(), ProcessLib::findProcessVariables(), MaterialPropertyLib::gas_phase_pressure, BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigSubtree(), BaseLib::ConfigTree::getConfigSubtreeOptional(), MeshLib::Mesh::getDimension(), MeshLib::getElementRotationMatrices(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::getSpaceDimension(), INFO(), MaterialPropertyLib::liquid_phase_pressure, OGS_FATAL, and MathLib::toVector().
Referenced by ProjectData::parseProcesses().
std::tuple< double, double > ProcessLib::LiquidFlow::getFluidDensityAndViscosity | ( | double const | t, |
double const | dt, | ||
ParameterLib::SpatialPosition const & | pos, | ||
MaterialPropertyLib::Phase const & | fluid_phase, | ||
MaterialPropertyLib::VariableArray & | vars ) |
Definition at line 18 of file LiquidFlowLocalAssembler.cpp.
References MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::Phase::hasProperty(), MaterialPropertyLib::molar_mass, MaterialPropertyLib::VariableArray::molar_mass, MaterialPropertyLib::Phase::property(), and MaterialPropertyLib::viscosity.
Referenced by ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::assembleMatrixAndVector(), and ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::computeProjectedDarcyVelocity().
const unsigned ProcessLib::LiquidFlow::NUM_NODAL_DOF = 1 |
Definition at line 53 of file LiquidFlowLocalAssembler.h.