OGS
CreateLiquidFlowProcess.cpp
Go to the documentation of this file.
1 
13 
14 #include <algorithm>
15 
16 #include "LiquidFlowProcess.h"
22 #include "ParameterLib/Utils.h"
25 
26 namespace ProcessLib
27 {
28 namespace LiquidFlow
29 {
31  MeshLib::Mesh const& mesh,
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 
48  mesh, media_map, required_medium_properties, required_solid_properties,
49  required_liquid_properties);
50 }
51 
52 std::unique_ptr<Process> createLiquidFlowProcess(
53  std::string name,
54  MeshLib::Mesh& mesh,
55  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
56  std::vector<ProcessVariable> const& variables,
57  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
58  unsigned const integration_order,
59  BaseLib::ConfigTree const& config,
60  std::vector<std::unique_ptr<MeshLib::Mesh>> const& meshes,
61  std::map<int, std::shared_ptr<MaterialPropertyLib::Medium>> const& media)
62 {
64  config.checkConfigParameter("type", "LIQUID_FLOW");
65 
66  DBUG("Create LiquidFlowProcess.");
67 
68  // Process variable.
69 
71  auto const pv_config = config.getConfigSubtree("process_variables");
72 
73  auto per_process_variables = findProcessVariables(
74  variables, pv_config,
75  {
76  "process_variable"});
77  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
78  process_variables;
79  process_variables.push_back(std::move(per_process_variables));
80 
81  SecondaryVariableCollection secondary_variables;
82 
83  ProcessLib::createSecondaryVariables(config, secondary_variables);
84 
85  std::vector<double> const b =
87  config.getConfigParameter<std::vector<double>>("specific_body_force");
88  int const mesh_space_dimension =
90  if (static_cast<int>(b.size()) != mesh_space_dimension)
91  {
92  OGS_FATAL(
93  "specific body force (gravity vector) has {:d} components, but the "
94  "space dimension is {:d}.",
95  b.size(), mesh_space_dimension);
96  }
97 
98  Eigen::VectorXd specific_body_force(b.size());
99  bool const has_gravity = MathLib::toVector(b).norm() > 0;
100  if (has_gravity)
101  {
102  std::copy_n(b.data(), b.size(), specific_body_force.data());
103  }
104 
105  std::unique_ptr<ProcessLib::SurfaceFluxData> surfaceflux;
106  auto calculatesurfaceflux_config =
108  config.getConfigSubtreeOptional("calculatesurfaceflux");
109  if (calculatesurfaceflux_config)
110  {
112  *calculatesurfaceflux_config, meshes);
113  }
114 
115  auto media_map =
117 
118  DBUG("Check the media properties of LiquidFlow process ...");
119  checkMPLProperties(mesh, *media_map);
120  DBUG("Media properties verified.");
121 
122  LiquidFlowData process_data{
123  std::move(media_map),
125  mesh_space_dimension, mesh.getDimension(), mesh.getElements()),
126  mesh_space_dimension, std::move(specific_body_force), has_gravity};
127 
128  return std::make_unique<LiquidFlowProcess>(
129  std::move(name), mesh, std::move(jacobian_assembler), parameters,
130  integration_order, std::move(process_variables),
131  std::move(process_data), std::move(secondary_variables),
132  std::move(surfaceflux));
133 }
134 
135 } // namespace LiquidFlow
136 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
void checkConfigParameter(std::string const &param, T const &value) const
std::optional< ConfigTree > getConfigSubtreeOptional(std::string const &root) const
Definition: ConfigTree.cpp:155
T getConfigParameter(std::string const &param) const
ConfigTree getConfigSubtree(std::string const &root) const
Definition: ConfigTree.cpp:146
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:95
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
Handles configuration of several secondary variables from the project file.
std::unique_ptr< MaterialSpatialDistributionMap > createMaterialSpatialDistributionMap(std::map< int, std::shared_ptr< Medium >> const &media, MeshLib::Mesh const &mesh)
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)
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::unique_ptr< Process > createLiquidFlowProcess(std::string 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)
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)