OGS
FractureProperty.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <Eigen/Core>
7#include <memory>
8
9#include "BranchProperty.h"
10#include "JunctionProperty.h"
13
14namespace ParameterLib
15{
16template <typename T>
17struct Parameter;
18}
19
20namespace ProcessLib
21{
22namespace LIE
23{
25{
26 int fracture_id = 0;
27 int mat_id = 0;
28 Eigen::Vector3d point_on_fracture;
29 Eigen::Vector3d normal_vector;
31 Eigen::MatrixXd R;
34 std::vector<BranchProperty> branches_master;
35 std::vector<BranchProperty> branches_slave;
36
37 FractureProperty(int const fracture_id_, int const material_id,
38 ParameterLib::Parameter<double> const& initial_aperture)
39 : fracture_id(fracture_id_),
40 mat_id(material_id),
41 aperture0(initial_aperture)
42 {
43 }
44
45 virtual ~FractureProperty() = default;
46};
47
50inline void setFractureProperty(int const dim, MeshLib::Element const& e,
51 FractureProperty& frac_prop)
52{
53 auto& n = frac_prop.normal_vector;
54 // 1st node is used but using other node is also possible, because
55 // a fracture is not curving
56 for (int j = 0; j < 3; j++)
57 {
58 frac_prop.point_on_fracture[j] = getCenterOfGravity(e).data()[j];
59 }
60
61 const MeshLib::ElementCoordinatesMappingLocal ele_local_coord(e, dim);
62
63 // Global to local rotation matrix:
64 Eigen::MatrixXd const global2local_rotation =
65 ele_local_coord.getRotationMatrixToGlobal().transpose();
66 n = global2local_rotation.row(dim - 1);
67
68 frac_prop.R = global2local_rotation.topLeftCorner(dim, dim);
69
70 DBUG("Normal vector of the fracture element {:d}: [{:g}, {:g}, {:g}]",
71 e.getID(), n[0], n[1], n[2]);
72}
73
75 FractureProperty const& master_frac,
76 FractureProperty const& slave_frac)
77{
78 BranchProperty branch{branchNode, master_frac.fracture_id,
79 slave_frac.fracture_id};
80
81 // set a normal vector from the master to the slave fracture
82 Eigen::Vector3d branch_vector =
83 slave_frac.point_on_fracture - branch.coords;
84 double sign = (branch_vector.dot(master_frac.normal_vector) < 0) ? -1 : 1;
85 branch.normal_vector_branch = sign * master_frac.normal_vector;
86 return branch;
87}
88} // namespace LIE
89} // namespace ProcessLib
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
const double * data() const
Definition Point3d.h:51
const RotationMatrix & getRotationMatrixToGlobal() const
return a rotation matrix converting to global coordinates
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:80
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
BranchProperty createBranchProperty(MeshLib::Node const &branchNode, FractureProperty const &master_frac, FractureProperty const &slave_frac)
Eigen::Vector3d const coords
virtual ~FractureProperty()=default
ParameterLib::Parameter< double > const & aperture0
Initial aperture.
std::vector< BranchProperty > branches_slave
FractureProperty(int const fracture_id_, int const material_id, ParameterLib::Parameter< double > const &initial_aperture)
Eigen::MatrixXd R
Rotation matrix from global to local coordinates.
std::vector< BranchProperty > branches_master