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

Functions

void checkMPLProperties (MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map, bool const is_equation_type_volume)
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

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 25 of file CreateLiquidFlowProcess.cpp.

29{
31
32 std::array const required_medium_properties = {
37
38 std::array const required_liquid_properties = {
41
42 // Check Constant-type density.
43 if (is_equation_type_volume)
44 {
46 media_map);
47 }
48
49 for (auto const& medium : media_map.media())
50 {
51 checkRequiredProperties(*medium, required_medium_properties);
53 required_liquid_properties);
54 }
55 DBUG("Media properties verified.");
56}
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:95
void checkVolumeBalanceEquationSetting(MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)

References MaterialPropertyLib::checkMPLPhasesForSinglePhaseFlow(), ProcessLib::Common::HydraulicProcess::checkVolumeBalanceEquationSetting(), DBUG(), MaterialPropertyLib::density, MaterialPropertyLib::MaterialSpatialDistributionMap::media(), MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::storage, and MaterialPropertyLib::viscosity.

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 58 of file CreateLiquidFlowProcess.cpp.

68{
70 config.checkConfigParameter("type", "LIQUID_FLOW");
71
72 DBUG("Create LiquidFlowProcess.");
73
75
77 auto const pv_config = config.getConfigSubtree("process_variables");
78
80 auto per_process_variables = findProcessVariables(
81 variables, pv_config,
82 {
83 "process_variable"});
84 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
85 process_variables;
86 process_variables.push_back(std::move(per_process_variables));
87
88 SecondaryVariableCollection secondary_variables;
89
90 ProcessLib::createSecondaryVariables(config, secondary_variables);
91
93 std::vector<double> const b =
95 config.getConfigParameter<std::vector<double>>("specific_body_force");
96 int const mesh_space_dimension =
98 if (static_cast<int>(b.size()) != mesh_space_dimension)
99 {
100 OGS_FATAL(
101 "specific body force (gravity vector) has {:d} components, but the "
102 "space dimension is {:d}.",
103 b.size(), mesh_space_dimension);
104 }
105
106 Eigen::VectorXd specific_body_force(b.size());
107 bool const has_gravity = MathLib::toVector(b).norm() > 0;
108 if (has_gravity)
109 {
110 std::copy_n(b.data(), b.size(), specific_body_force.data());
111 }
112
113 std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
114 auto calculatesurfaceflux_config =
116 config.getConfigSubtreeOptional("calculatesurfaceflux");
117 if (calculatesurfaceflux_config)
118 {
120 *calculatesurfaceflux_config, meshes);
121 }
122
123 auto media_map =
125
126 auto const is_linear =
128 config.getConfigParameter<bool>("linear", false);
129 if (is_linear)
130 {
131 INFO("LiquidFlow process is set to be linear.");
132 }
133
134 auto const equation_balance_type_str =
136 config.getConfigParameter<std::string>("equation_balance_type",
137 "volume");
138
139 DBUG("Check the media properties of LiquidFlow process ...");
140 checkMPLProperties(mesh, media_map, equation_balance_type_str == "volume");
141 DBUG("Media properties verified.");
142
143 auto const* aperture_size_parameter = &ParameterLib::findParameter<double>(
145
146 auto const aperture_config =
148 config.getConfigSubtreeOptional("aperture_size");
149 if (aperture_config)
150 {
151 aperture_size_parameter = &ParameterLib::findParameter<double>(
153 *aperture_config, "parameter", parameters, 1);
154 }
155
156 // The uniqueness of phase has already been checked in
157 // `checkMPLProperties`.
158 MaterialPropertyLib::Variable const phase_variable =
159 (*ranges::begin(media_map.media()))
163
164 LiquidFlowData process_data{
165 equation_balance_type_str ==
166 "volume", // is_volume_balance_equation_type
167 std::move(media_map),
169 mesh_space_dimension, mesh.getDimension(), mesh.getElements()),
170 mesh_space_dimension,
171 std::move(specific_body_force),
172 has_gravity,
173 *aperture_size_parameter,
174 NumLib::ShapeMatrixCache{integration_order},
175 phase_variable};
176
177 return std::make_unique<LiquidFlowProcess>(
178 std::move(name), mesh, std::move(jacobian_assembler), parameters,
179 integration_order, std::move(process_variables),
180 std::move(process_data), std::move(secondary_variables),
181 std::move(surfaceflux), is_linear);
182}
#define OGS_FATAL(...)
Definition Error.h:19
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:98
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:101
unsigned getDimension() const
Definition Mesh.h:80
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)
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, MaterialPropertyLib::createMaterialSpatialDistributionMap(), ProcessLib::createSecondaryVariables(), ProcessLib::SurfaceFluxData::createSurfaceFluxData(), DBUG(), ParameterLib::findParameter(), ProcessLib::findProcessVariables(), MaterialPropertyLib::Gas, 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.