35 bool const is_equation_type_volume)
39 std::array
const required_medium_properties = {
45 std::array
const required_liquid_properties = {
50 if (is_equation_type_volume)
52 for (
auto const& medium : media_map.
media())
55 auto const& fluid_phase_density =
57 if (
typeid(fluid_phase_density) !=
61 "Since `equation_balance_type` is set to `volume`,the "
62 "phase density type must be `Constant`. Note: by "
63 "default, `equation_balance_type` is set to `volume`.");
68 for (
auto const& medium : media_map.
media())
70 checkRequiredProperties(*medium, required_medium_properties);
71 checkRequiredProperties(fluidPhase(*medium),
72 required_liquid_properties);
74 DBUG(
"Media properties verified.");
78 std::string_view
const type_in_string)
80 if (type_in_string ==
"volume")
84 if (type_in_string ==
"mass")
89 "The value of `equation_balance_type` must be either `volume` or "
94 std::string
const& name,
96 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
97 std::vector<ProcessVariable>
const& variables,
98 std::vector<std::unique_ptr<ParameterLib::ParameterBase>>
const& parameters,
99 unsigned const integration_order,
101 std::vector<std::unique_ptr<MeshLib::Mesh>>
const& meshes,
102 std::map<
int, std::shared_ptr<MaterialPropertyLib::Medium>>
const& media)
107 DBUG(
"Create LiquidFlowProcess.");
116 variables, pv_config,
118 "process_variable"});
119 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
121 process_variables.push_back(std::move(per_process_variables));
128 std::vector<double>
const b =
131 int const mesh_space_dimension =
133 if (
static_cast<int>(b.size()) != mesh_space_dimension)
136 "specific body force (gravity vector) has {:d} components, but the "
137 "space dimension is {:d}.",
138 b.size(), mesh_space_dimension);
141 Eigen::VectorXd specific_body_force(b.size());
145 std::copy_n(b.data(), b.size(), specific_body_force.data());
148 std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
149 auto calculatesurfaceflux_config =
152 if (calculatesurfaceflux_config)
155 *calculatesurfaceflux_config, meshes);
161 auto const is_linear =
166 INFO(
"LiquidFlow process is set to be linear.");
169 auto const equation_balance_type_str =
174 DBUG(
"Check the media properties of LiquidFlow process ...");
176 DBUG(
"Media properties verified.");
181 auto const aperture_config =
188 *aperture_config,
"parameter", parameters, 1);
194 (*ranges::begin(media_map.media()))->hasPhase(
"Gas")
200 std::move(media_map),
203 mesh_space_dimension,
204 std::move(specific_body_force),
206 *aperture_size_parameter,
210 return std::make_unique<LiquidFlowProcess>(
211 std::move(name), mesh, std::move(jacobian_assembler), parameters,
212 integration_order, std::move(process_variables),
213 std::move(process_data), std::move(secondary_variables),
214 std::move(surfaceflux), is_linear);
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
std::optional< ConfigTree > getConfigSubtreeOptional(std::string const &root) const
T getConfigParameter(std::string const ¶m) const
ConfigTree getConfigSubtree(std::string const &root) const
void checkConfigParameter(std::string const ¶m, std::string_view const value) const
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name
Handles configuration of several secondary variables from the project file.
void checkMPLPhasesForSinglePhaseFlow(MeshLib::Mesh const &mesh, MaterialPropertyLib::MaterialSpatialDistributionMap const &media_map)
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 ¶meter_name, std::vector< std::unique_ptr< ParameterBase > > const ¶meters, 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
Governing equation balance type.
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)
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)