OGS
ProcessLib::LiquidFlow Namespace Reference

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< ProcesscreateLiquidFlowProcess (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)

Variables

const unsigned NUM_NODAL_DOF = 1

Enumeration Type Documentation

◆ EquationBalanceType

Governing equation balance type.

Enumerator
volume 
mass 

Definition at line 19 of file LiquidFlowData.h.

Function Documentation

◆ checkMPLProperties()

void ProcessLib::LiquidFlow::checkMPLProperties ( MeshLib::Mesh const & mesh,
MaterialPropertyLib::MaterialSpatialDistributionMap const & media_map,
bool const is_equation_type_volume )

Definition at line 24 of file CreateLiquidFlowProcess.cpp.

28{
30
31 std::array const required_medium_properties = {
36
37 std::array const required_liquid_properties = {
40
41 // Check Constant-type density.
42 if (is_equation_type_volume)
43 {
44 for (auto const& medium : media_map.media())
45 {
46 // auto const& medium = *media_map.getMedium(element_id);
47 auto const& fluid_phase_density =
49 if (typeid(fluid_phase_density) !=
51 {
53 "Since `equation_balance_type` is set to `volume`,the "
54 "phase density type must be `Constant`. Note: by "
55 "default, `equation_balance_type` is set to `volume`.");
56 }
57 }
58 }
59
60 for (auto const& medium : media_map.media())
61 {
62 checkRequiredProperties(*medium, required_medium_properties);
64 required_liquid_properties);
65 }
66 DBUG("Media properties verified.");
67}
#define OGS_FATAL(...)
Definition Error.h:19
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void checkMPLPhasesForSinglePhaseFlow(MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)
void checkRequiredProperties(Component const &c, std::span< PropertyType const > const required_properties)
Definition Component.cpp:51
Phase const & fluidPhase(Medium const &medium)
Returns a gas or aqueous liquid phase of the given medium.
Definition Medium.cpp:91

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().

◆ covertEquationBalanceTypeFromString()

EquationBalanceType ProcessLib::LiquidFlow::covertEquationBalanceTypeFromString ( std::string_view const type_in_string)

Definition at line 69 of file CreateLiquidFlowProcess.cpp.

71{
72 if (type_in_string == "volume")
73 {
75 }
76 if (type_in_string == "mass")
77 {
79 }
81 "The value of `equation_balance_type` must be either `volume` or "
82 "`mass`");
83};

References mass, OGS_FATAL, and volume.

Referenced by createLiquidFlowProcess().

◆ 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 )
Input File Parameter
prj__processes__process__type

Process Variables

Input File Parameter
prj__processes__process__LIQUID_FLOW__process_variables

Primary process variables as they appear in the global component vector:

Input File Parameter
prj__processes__process__LIQUID_FLOW__process_variables__process_variable

Process Parameters

Input File Parameter
prj__processes__process__LIQUID_FLOW__specific_body_force
Input File Parameter
prj__processes__process__calculatesurfaceflux
Input File Parameter
prj__processes__process__linear
Input File Parameter
prj__processes__process__LIQUID_FLOW__equation_balance_type
Input File Parameter
prj__processes__process__LIQUID_FLOW__aperture_size
Input File Parameter
prj__processes__process__LIQUID_FLOW__aperture_size__parameter

Definition at line 85 of file CreateLiquidFlowProcess.cpp.

95{
97 config.checkConfigParameter("type", "LIQUID_FLOW");
98
99 DBUG("Create LiquidFlowProcess.");
100
102
104 auto const pv_config = config.getConfigSubtree("process_variables");
105
107 auto per_process_variables = findProcessVariables(
108 variables, pv_config,
109 {
110 "process_variable"});
111 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
112 process_variables;
113 process_variables.push_back(std::move(per_process_variables));
114
115 SecondaryVariableCollection secondary_variables;
116
117 ProcessLib::createSecondaryVariables(config, secondary_variables);
118
120 std::vector<double> const b =
122 config.getConfigParameter<std::vector<double>>("specific_body_force");
123 int const mesh_space_dimension =
125 if (static_cast<int>(b.size()) != mesh_space_dimension)
126 {
127 OGS_FATAL(
128 "specific body force (gravity vector) has {:d} components, but the "
129 "space dimension is {:d}.",
130 b.size(), mesh_space_dimension);
131 }
132
133 Eigen::VectorXd specific_body_force(b.size());
134 bool const has_gravity = MathLib::toVector(b).norm() > 0;
135 if (has_gravity)
136 {
137 std::copy_n(b.data(), b.size(), specific_body_force.data());
138 }
139
140 std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
141 auto calculatesurfaceflux_config =
143 config.getConfigSubtreeOptional("calculatesurfaceflux");
144 if (calculatesurfaceflux_config)
145 {
147 *calculatesurfaceflux_config, meshes);
148 }
149
150 auto media_map =
152
153 auto const is_linear =
155 config.getConfigParameter<bool>("linear", false);
156 if (is_linear)
157 {
158 INFO("LiquidFlow process is set to be linear.");
159 }
160
161 auto const equation_balance_type_str =
163 config.getConfigParameter<std::string>("equation_balance_type",
164 "volume");
165
166 DBUG("Check the media properties of LiquidFlow process ...");
167 checkMPLProperties(mesh, media_map, equation_balance_type_str == "volume");
168 DBUG("Media properties verified.");
169
170 auto const* aperture_size_parameter = &ParameterLib::findParameter<double>(
172
173 auto const aperture_config =
175 config.getConfigSubtreeOptional("aperture_size");
176 if (aperture_config)
177 {
178 aperture_size_parameter = &ParameterLib::findParameter<double>(
180 *aperture_config, "parameter", parameters, 1);
181 }
182
183 // The uniqueness of phase has already been checked in
184 // `checkMPLProperties`.
185 MaterialPropertyLib::Variable const phase_variable =
186 (*ranges::begin(media_map.media()))->hasPhase("Gas")
189
190 LiquidFlowData process_data{
191 covertEquationBalanceTypeFromString(equation_balance_type_str),
192 std::move(media_map),
194 mesh_space_dimension, mesh.getDimension(), mesh.getElements()),
195 mesh_space_dimension,
196 std::move(specific_body_force),
197 has_gravity,
198 *aperture_size_parameter,
199 NumLib::ShapeMatrixCache{integration_order},
200 phase_variable};
201
202 return std::make_unique<LiquidFlowProcess>(
203 std::move(name), mesh, std::move(jacobian_assembler), parameters,
204 integration_order, std::move(process_variables),
205 std::move(process_data), std::move(secondary_variables),
206 std::move(surfaceflux), is_linear);
207}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:97
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:79
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:40
Handles configuration of several secondary variables from the project file.
MaterialSpatialDistributionMap createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium > > const &media, MeshLib::Mesh const &mesh)
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< Eigen::MatrixXd > getElementRotationMatrices(int const space_dimension, int const mesh_dimension, std::vector< Element * > const &elements)
Element rotation matrix computation.
int getSpaceDimension(std::vector< Node * > const &nodes)
Computes dimension of the embedding space containing the set of given points.
OGS_NO_DANGLING Parameter< ParameterDataType > & findParameter(std::string const &parameter_name, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
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::vector< std::reference_wrapper< ProcessVariable > > findProcessVariables(std::vector< ProcessVariable > const &variables, BaseLib::ConfigTree const &pv_config, std::initializer_list< std::string > tags)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)
static std::unique_ptr< ProcessLib::SurfaceFluxData > createSurfaceFluxData(BaseLib::ConfigTree const &calculatesurfaceflux_config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)

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().

Variable Documentation

◆ NUM_NODAL_DOF

const unsigned ProcessLib::LiquidFlow::NUM_NODAL_DOF = 1

Definition at line 43 of file LiquidFlowLocalAssembler.h.