Loading [MathJax]/extensions/tex2jax.js
OGS
CreateHydroMechanicsProcess.cpp
Go to the documentation of this file.
1 
12 
13 #include <cassert>
14 
15 #include "HydroMechanicsProcess.h"
22 #include "ParameterLib/Utils.h"
24 
25 namespace ProcessLib
26 {
27 namespace LIE
28 {
29 namespace HydroMechanics
30 {
31 template <int GlobalDim>
32 std::unique_ptr<Process> createHydroMechanicsProcess(
33  std::string name,
34  MeshLib::Mesh& mesh,
35  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
36  std::vector<ProcessVariable> const& variables,
37  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
38  std::optional<ParameterLib::CoordinateSystem> const&
39  local_coordinate_system,
40  unsigned const integration_order,
41  BaseLib::ConfigTree const& config)
42 {
44  config.checkConfigParameter("type", "HYDRO_MECHANICS_WITH_LIE");
45  DBUG("Create HydroMechanicsProcess with LIE.");
46  auto const coupling_scheme =
48  config.getConfigParameterOptional<std::string>("coupling_scheme");
49  const bool use_monolithic_scheme =
50  !(coupling_scheme && (*coupling_scheme == "staggered"));
51 
52  // Process variables
54  auto const pv_conf = config.getConfigSubtree("process_variables");
55  auto range =
57  pv_conf.getConfigParameterList<std::string>("process_variable");
58  std::vector<std::reference_wrapper<ProcessVariable>> p_u_process_variables;
59  std::vector<std::reference_wrapper<ProcessVariable>> p_process_variables;
60  std::vector<std::reference_wrapper<ProcessVariable>> u_process_variables;
61  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>
62  process_variables;
63  for (std::string const& pv_name : range)
64  {
65  if (pv_name != "pressure" && pv_name != "displacement" &&
66  pv_name.find("displacement_jump") != 0)
67  {
68  OGS_FATAL(
69  "Found a process variable name '{:s}'. It should be "
70  "'displacement' or 'displacement_jumpN' or 'pressure'");
71  }
72  auto variable = std::find_if(variables.cbegin(), variables.cend(),
73  [&pv_name](ProcessVariable const& v)
74  { return v.getName() == pv_name; });
75 
76  if (variable == variables.end())
77  {
78  OGS_FATAL(
79  "Could not find process variable '{:s}' in the provided "
80  "variables list for config tag <{:s}>.",
81  pv_name, "process_variable");
82  }
83  DBUG("Found process variable '{:s}' for config tag <{:s}>.",
84  variable->getName(), "process_variable");
85 
86  if (pv_name.find("displacement") != std::string::npos &&
87  variable->getNumberOfGlobalComponents() != GlobalDim)
88  {
89  OGS_FATAL(
90  "Number of components of the process variable '{:s}' is "
91  "different from the displacement dimension: got {:d}, expected "
92  "{:d}",
93  variable->getName(),
94  variable->getNumberOfGlobalComponents(),
95  GlobalDim);
96  }
97 
98  if (!use_monolithic_scheme)
99  {
100  if (pv_name == "pressure")
101  {
102  p_process_variables.emplace_back(
103  const_cast<ProcessVariable&>(*variable));
104  }
105  else
106  {
107  u_process_variables.emplace_back(
108  const_cast<ProcessVariable&>(*variable));
109  }
110  }
111  else
112  {
113  p_u_process_variables.emplace_back(
114  const_cast<ProcessVariable&>(*variable));
115  }
116  }
117 
118  if (p_u_process_variables.size() > 3 || u_process_variables.size() > 2)
119  {
120  OGS_FATAL("Currently only one displacement jump is supported");
121  }
122 
123  if (!use_monolithic_scheme)
124  {
125  process_variables.push_back(std::move(p_process_variables));
126  process_variables.push_back(std::move(u_process_variables));
127  }
128  else
129  {
130  process_variables.push_back(std::move(p_u_process_variables));
131  }
132 
133  auto solid_constitutive_relations =
134  MaterialLib::Solids::createConstitutiveRelations<GlobalDim>(
135  parameters, local_coordinate_system, config);
136 
137  // Intrinsic permeability
138  auto& intrinsic_permeability = ParameterLib::findParameter<double>(
139  config,
141  "intrinsic_permeability", parameters, 1, &mesh);
142 
143  DBUG("Use '{:s}' as intrinsic permeability parameter.",
144  intrinsic_permeability.name);
145 
146  // Storage coefficient
147  auto& specific_storage = ParameterLib::findParameter<double>(
148  config,
150  "specific_storage", parameters, 1, &mesh);
151 
152  DBUG("Use '{:s}' as specific storage parameter.", specific_storage.name);
153 
154  // Fluid viscosity
155  auto& fluid_viscosity = ParameterLib::findParameter<double>(
156  config,
158  "fluid_viscosity", parameters, 1, &mesh);
159  DBUG("Use '{:s}' as fluid viscosity parameter.", fluid_viscosity.name);
160 
161  // Fluid density
162  auto& fluid_density = ParameterLib::findParameter<double>(
163  config,
165  "fluid_density", parameters, 1, &mesh);
166  DBUG("Use '{:s}' as fluid density parameter.", fluid_density.name);
167 
168  // Biot coefficient
169  auto& biot_coefficient = ParameterLib::findParameter<double>(
170  config,
172  "biot_coefficient", parameters, 1, &mesh);
173  DBUG("Use '{:s}' as Biot coefficient parameter.", biot_coefficient.name);
174 
175  // Porosity
176  auto& porosity = ParameterLib::findParameter<double>(
177  config,
179  "porosity", parameters, 1, &mesh);
180  DBUG("Use '{:s}' as porosity parameter.", porosity.name);
181 
182  // Solid density
183  auto& solid_density = ParameterLib::findParameter<double>(
184  config,
186  "solid_density", parameters, 1, &mesh);
187  DBUG("Use '{:s}' as solid density parameter.", solid_density.name);
188 
189  // Specific body force
190  Eigen::Matrix<double, GlobalDim, 1> specific_body_force;
191  {
192  std::vector<double> const b =
194  config.getConfigParameter<std::vector<double>>(
195  "specific_body_force");
196  if (b.size() != GlobalDim)
197  {
198  OGS_FATAL(
199  "The size of the specific body force vector does not match the "
200  "displacement dimension. Vector size is {:d}, displacement "
201  "dimension is {:d}",
202  b.size(), GlobalDim);
203  }
204 
205  std::copy_n(b.data(), b.size(), specific_body_force.data());
206  }
207 
208  // Fracture constitutive relation.
209  std::unique_ptr<MaterialLib::Fracture::FractureModelBase<GlobalDim>>
210  fracture_model = nullptr;
211  auto const opt_fracture_model_config =
213  config.getConfigSubtreeOptional("fracture_model");
214  if (opt_fracture_model_config)
215  {
216  auto& fracture_model_config = *opt_fracture_model_config;
217 
218  auto const frac_type =
220  fracture_model_config.peekConfigParameter<std::string>("type");
221 
222  if (frac_type == "LinearElasticIsotropic")
223  {
224  fracture_model =
225  MaterialLib::Fracture::createLinearElasticIsotropic<GlobalDim>(
226  parameters, fracture_model_config);
227  }
228  else if (frac_type == "Coulomb")
229  {
230  fracture_model = MaterialLib::Fracture::createCoulomb<GlobalDim>(
231  parameters, fracture_model_config);
232  }
233  else if (frac_type == "CohesiveZoneModeI")
234  {
235  fracture_model = MaterialLib::Fracture::CohesiveZoneModeI::
236  createCohesiveZoneModeI<GlobalDim>(parameters,
237  fracture_model_config);
238  }
239  else
240  {
241  OGS_FATAL(
242  "Cannot construct fracture constitutive relation of given type "
243  "'{:s}'.",
244  frac_type);
245  }
246  }
247 
248  // Fracture properties
249  std::unique_ptr<FracturePropertyHM> frac_prop = nullptr;
250  auto opt_fracture_properties_config =
252  config.getConfigSubtreeOptional("fracture_properties");
253  if (opt_fracture_properties_config)
254  {
255  auto& fracture_properties_config = *opt_fracture_properties_config;
256 
257  frac_prop = std::make_unique<ProcessLib::LIE::FracturePropertyHM>(
258  0 /*fracture_id*/,
260  fracture_properties_config.getConfigParameter<int>("material_id"),
261  ParameterLib::findParameter<double>(
263  fracture_properties_config, "initial_aperture", parameters, 1,
264  &mesh),
265  ParameterLib::findParameter<double>(
267  fracture_properties_config, "specific_storage", parameters, 1,
268  &mesh),
269  ParameterLib::findParameter<double>(
271  fracture_properties_config, "biot_coefficient", parameters, 1,
272  &mesh));
273  if (frac_prop->aperture0.isTimeDependent())
274  {
275  OGS_FATAL(
276  "The initial aperture parameter '{:s}' must not be "
277  "time-dependent.",
278  frac_prop->aperture0.name);
279  }
280 
281  auto permeability_model_config =
283  fracture_properties_config.getConfigSubtree("permeability_model");
284  frac_prop->permeability_model =
286  permeability_model_config);
287  }
288 
289  // initial effective stress in matrix
290  auto& initial_effective_stress = ParameterLib::findParameter<double>(
291  config,
293  "initial_effective_stress", parameters,
295  DBUG("Use '{:s}' as initial effective stress parameter.",
296  initial_effective_stress.name);
297 
298  // initial effective stress in fracture
299  auto& initial_fracture_effective_stress = ParameterLib::findParameter<
300  double>(
301  config,
303  "initial_fracture_effective_stress", parameters, GlobalDim, &mesh);
304  DBUG("Use '{:s}' as initial fracture effective stress parameter.",
305  initial_fracture_effective_stress.name);
306 
307  // deactivation of matrix elements in flow
308  auto opt_deactivate_matrix_in_flow =
310  config.getConfigParameterOptional<bool>("deactivate_matrix_in_flow");
311  bool const deactivate_matrix_in_flow =
312  opt_deactivate_matrix_in_flow && *opt_deactivate_matrix_in_flow;
313  ;
314  if (deactivate_matrix_in_flow)
315  INFO("Deactivate matrix elements in flow calculation.");
316 
317  // Reference temperature
318  const auto& reference_temperature =
320  config.getConfigParameter<double>(
321  "reference_temperature", std::numeric_limits<double>::quiet_NaN());
322 
324  materialIDs(mesh),
325  std::move(solid_constitutive_relations),
326  intrinsic_permeability,
327  specific_storage,
331  porosity,
332  solid_density,
333  specific_body_force,
334  std::move(fracture_model),
335  std::move(frac_prop),
336  initial_effective_stress,
337  initial_fracture_effective_stress,
338  deactivate_matrix_in_flow,
340 
341  SecondaryVariableCollection secondary_variables;
342 
343  ProcessLib::createSecondaryVariables(config, secondary_variables);
344 
345  return std::make_unique<HydroMechanicsProcess<GlobalDim>>(
346  std::move(name), mesh, std::move(jacobian_assembler), parameters,
347  integration_order, std::move(process_variables),
348  std::move(process_data), std::move(secondary_variables),
349  use_monolithic_scheme);
350 }
351 
352 template std::unique_ptr<Process> createHydroMechanicsProcess<2>(
353  std::string name,
354  MeshLib::Mesh& mesh,
355  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
356  std::vector<ProcessVariable> const& variables,
357  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
358  std::optional<ParameterLib::CoordinateSystem> const&
359  local_coordinate_system,
360  unsigned const integration_order,
361  BaseLib::ConfigTree const& config);
362 template std::unique_ptr<Process> createHydroMechanicsProcess<3>(
363  std::string name,
364  MeshLib::Mesh& mesh,
365  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
366  std::vector<ProcessVariable> const& variables,
367  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
368  std::optional<ParameterLib::CoordinateSystem> const&
369  local_coordinate_system,
370  unsigned const integration_order,
371  BaseLib::ConfigTree const& config);
372 
373 } // namespace HydroMechanics
374 } // namespace LIE
375 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
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
std::optional< T > getConfigParameterOptional(std::string const &param) const
T getConfigParameter(std::string const &param) const
ConfigTree getConfigSubtree(std::string const &root) const
Definition: ConfigTree.cpp:146
Range< ValueIterator< T > > getConfigParameterList(std::string const &param) const
Handles configuration of several secondary variables from the project file.
std::unique_ptr< Permeability > createPermeabilityModel(BaseLib::ConfigTree const &config)
constexpr int kelvin_vector_dimensions(int const displacement_dim)
Kelvin vector dimensions for given displacement dimension.
Definition: KelvinVector.h:23
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition: Mesh.cpp:258
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)
Definition: Utils.h:102
template std::unique_ptr< Process > createHydroMechanicsProcess< 2 >(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, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
template std::unique_ptr< Process > createHydroMechanicsProcess< 3 >(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, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
std::unique_ptr< Process > createHydroMechanicsProcess(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, std::optional< ParameterLib::CoordinateSystem > const &local_coordinate_system, unsigned const integration_order, BaseLib::ConfigTree const &config)
double fluid_viscosity(const double p, const double T, const double x)
double fluid_density(const double p, const double T, const double x)
void createSecondaryVariables(BaseLib::ConfigTree const &config, SecondaryVariableCollection &secondary_variables)