OGS
CreateOrthotropicEmbeddedFracturePermeability.cpp
Go to the documentation of this file.
1
11#include <Eigen/Geometry>
12
13#include "BaseLib/ConfigTree.h"
15#include "ParameterLib/Utils.h"
16
17namespace MaterialPropertyLib
18{
20 int const geometry_dimension, BaseLib::ConfigTree const& config,
21 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters)
22{
23 if ((geometry_dimension != 2) && (geometry_dimension != 3))
24 {
26 "The OrthotropicEmbeddedFracturePermeability is implemented only "
27 "for 2D or 3D problems");
28 }
29
31 config.checkConfigParameter("type",
32 "OrthotropicEmbeddedFracturePermeability");
33
34 // Second access for storage.
36 auto property_name = config.peekConfigParameter<std::string>("name");
37
38 DBUG("Create OrthotropicEmbeddedFracturePermeability medium property");
39
40 auto const a_i =
42 config.getConfigParameter<std::vector<double>>("mean_frac_distances");
43 if (a_i.size() != 3)
44 {
46 "The size of the mean fracture distances vector must be 3, but is "
47 "{}.",
48 a_i.size());
49 }
50
51 auto const e_i0 =
53 config.getConfigParameter<std::vector<double>>("threshold_strains");
54 if (e_i0.size() != 3)
55 {
57 "The size of the mean threshold strains vector must be 3, but is "
58 "{}.",
59 e_i0.size());
60 }
61
62 auto const n =
64 config.getConfigParameter<std::vector<double>>("fracture_normals");
65 if (n.size() != 6)
66 {
68 "The size of the fracture normals vector must be 6, but is {}.",
69 n.size());
70 }
71 Eigen::Vector3d const n1 = Eigen::Vector3d({n[0], n[1], n[2]}).normalized();
72 Eigen::Vector3d const n2 = Eigen::Vector3d({n[3], n[4], n[5]}).normalized();
73
74 if (n1.dot(n2) > std::numeric_limits<double>::epsilon())
75 {
77 "The given fracture normals are not orthogonal. Please provide two "
78 "orthogonal fracture normals");
79 }
80
81 Eigen::Matrix3d const n_i =
82 (Eigen::Matrix3d() << n1, n2, n1.cross(n2)).finished();
83
84 std::string const intrinsic_permeability_param_name =
86 config.getConfigParameter<std::string>("intrinsic_permeability");
87
89 intrinsic_permeability_param_name, parameters, 0, nullptr);
90
91 std::string const fracture_rotation_xy_param_name =
93 config.getConfigParameter<std::string>("fracture_rotation_xy");
94
95 auto const& phi_xy = ParameterLib::findParameter<double>(
96 fracture_rotation_xy_param_name, parameters, 0, nullptr);
97
98 std::string const fracture_rotation_yz_param_name =
100 config.getConfigParameter<std::string>("fracture_rotation_yz");
101
102 auto const& phi_yz = ParameterLib::findParameter<double>(
103 fracture_rotation_yz_param_name, parameters, 0, nullptr);
104
105 auto const jf =
107 config.getConfigParameter<double>("jacobian_factor", 0.);
108
109 if (geometry_dimension == 2)
110 {
111 return std::make_unique<OrthotropicEmbeddedFracturePermeability<2>>(
112 std::move(property_name), a_i, e_i0, n_i, k, phi_xy, phi_yz, jf);
113 }
114 return std::make_unique<OrthotropicEmbeddedFracturePermeability<3>>(
115 std::move(property_name), a_i, e_i0, n_i, k, phi_xy, phi_yz, jf);
116}
117} // namespace MaterialPropertyLib
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
T peekConfigParameter(std::string const &param) const
T getConfigParameter(std::string const &param) const
void checkConfigParameter(std::string const &param, std::string_view const value) const
std::unique_ptr< Property > createOrthotropicEmbeddedFracturePermeability(int const geometry_dimension, BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters)
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)
Definition Utils.h:102