OGS
ProcessLib::LiquidFlow Namespace Reference

Classes

struct  IntegrationPointData
 
struct  LiquidFlowData
 
class  LiquidFlowLocalAssembler
 
class  LiquidFlowLocalAssemblerInterface
 
class  LiquidFlowProcess
 A class to simulate the liquid flow process in porous media described by. More...
 

Functions

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

Definition at line 30 of file CreateLiquidFlowProcess.cpp.

33{
34 std::array const required_medium_properties = {
39
40 std::array const required_liquid_properties = {
43
44 std::array<MaterialPropertyLib::PropertyType, 0> const
45 required_solid_properties{};
46
47 std::array<MaterialPropertyLib::PropertyType, 0> const
48 required_gas_properties{};
49
51 mesh, media_map, required_medium_properties, required_solid_properties,
52 required_liquid_properties, required_gas_properties);
53}
void checkMaterialSpatialDistributionMap(MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map, ContainerMedium const &required_properties_medium, ContainerSolid const &required_properties_solid_phase, ContainerLiquid const &required_properties_liquid_phase, ContainerGas const &required_properties_gas_phase)

References MaterialPropertyLib::checkMaterialSpatialDistributionMap(), MaterialPropertyLib::density, 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__LIQUID_FLOW__linear
Input File Parameter
prj__processes__process__LIQUID_FLOW__aperture_size
Input File Parameter
prj__processes__process__LIQUID_FLOW__aperture_size__parameter

Definition at line 55 of file CreateLiquidFlowProcess.cpp.

65{
67 config.checkConfigParameter("type", "LIQUID_FLOW");
68
69 DBUG("Create LiquidFlowProcess.");
70
72
74 auto const pv_config = config.getConfigSubtree("process_variables");
75
77 auto per_process_variables = findProcessVariables(
78 variables, pv_config,
79 {
80 "process_variable"});
81 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
82 process_variables;
83 process_variables.push_back(std::move(per_process_variables));
84
85 SecondaryVariableCollection secondary_variables;
86
87 ProcessLib::createSecondaryVariables(config, secondary_variables);
88
90 std::vector<double> const b =
92 config.getConfigParameter<std::vector<double>>("specific_body_force");
93 int const mesh_space_dimension =
95 if (static_cast<int>(b.size()) != mesh_space_dimension)
96 {
98 "specific body force (gravity vector) has {:d} components, but the "
99 "space dimension is {:d}.",
100 b.size(), mesh_space_dimension);
101 }
102
103 Eigen::VectorXd specific_body_force(b.size());
104 bool const has_gravity = MathLib::toVector(b).norm() > 0;
105 if (has_gravity)
106 {
107 std::copy_n(b.data(), b.size(), specific_body_force.data());
108 }
109
110 std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
111 auto calculatesurfaceflux_config =
113 config.getConfigSubtreeOptional("calculatesurfaceflux");
114 if (calculatesurfaceflux_config)
115 {
117 *calculatesurfaceflux_config, meshes);
118 }
119
120 auto media_map =
122
123 auto const is_linear =
125 config.getConfigParameter<bool>("linear", false);
126 if (is_linear)
127 {
128 INFO("LiquidFlow process is set to be linear.");
129 }
130
131 DBUG("Check the media properties of LiquidFlow process ...");
132 checkMPLProperties(mesh, media_map);
133 DBUG("Media properties verified.");
134
135 auto const* aperture_size_parameter = &ParameterLib::findParameter<double>(
137
138 auto const aperture_config =
140 config.getConfigSubtreeOptional("aperture_size");
141 if (aperture_config)
142 {
143 aperture_size_parameter = &ParameterLib::findParameter<double>(
145 *aperture_config, "parameter", parameters, 1);
146 }
147
148 LiquidFlowData process_data{
149 std::move(media_map),
151 mesh_space_dimension, mesh.getDimension(), mesh.getElements()),
152 mesh_space_dimension,
153 std::move(specific_body_force),
154 has_gravity,
155 *aperture_size_parameter,
156 NumLib::ShapeMatrixCache{integration_order}};
157
158 return std::make_unique<LiquidFlowProcess>(
159 std::move(name), mesh, std::move(jacobian_assembler), parameters,
160 integration_order, std::move(process_variables),
161 std::move(process_data), std::move(secondary_variables),
162 std::move(surfaceflux), is_linear);
163}
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Definition Process.h:46
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.
void checkMPLProperties(MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)
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(), ProcessLib::findProcessVariables(), BaseLib::ConfigTree::getConfigParameter(), BaseLib::ConfigTree::getConfigSubtree(), BaseLib::ConfigTree::getConfigSubtreeOptional(), MeshLib::Mesh::getDimension(), MeshLib::getElementRotationMatrices(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getNodes(), MeshLib::getSpaceDimension(), INFO(), 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 52 of file LiquidFlowLocalAssembler.h.