OGS
ProcessLib::LIE Namespace Reference

Namespaces

namespace  anonymous_namespace{MeshUtils.cpp}
 
namespace  anonymous_namespace{PostUtils.cpp}
 
namespace  HydroMechanics
 
namespace  SmallDeformation
 

Classes

struct  BranchProperty
 
struct  FractureProperty
 
struct  FracturePropertyHM
 
struct  JunctionProperty
 
class  PostProcessTool
 

Functions

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)
 
bool levelsetFracture (FractureProperty const &frac, Eigen::Vector3d const &x)
 
bool levelsetBranch (BranchProperty const &branch, Eigen::Vector3d const &x)
 
std::vector< double > uGlobalEnrichments (std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
 
std::vector< double > duGlobalEnrichments (std::size_t this_frac_id, std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
 
void getFractureMatrixDataInMesh (MeshLib::Mesh const &mesh, std::vector< MeshLib::Element * > &vec_matrix_elements, std::vector< int > &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_elements, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_matrix_elements, std::vector< std::vector< MeshLib::Node * > > &vec_fracture_nodes, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)
 
template<typename Derived >
Eigen::Vector3d computePhysicalCoordinates (MeshLib::Element const &e, Eigen::MatrixBase< Derived > const &shape)
 

Function Documentation

◆ computePhysicalCoordinates()

template<typename Derived >
Eigen::Vector3d ProcessLib::LIE::computePhysicalCoordinates ( MeshLib::Element const & e,
Eigen::MatrixBase< Derived > const & shape )

compute physical coordinates from the given shape vector, i.e. from the natural coordinates

Definition at line 24 of file Utils.h.

26{
27 Eigen::Vector3d pt = Eigen::Vector3d::Zero();
28 for (unsigned i = 0; i < e.getNumberOfNodes(); i++)
29 {
30 MeshLib::Node const& node = *e.getNode(i);
31 for (unsigned j = 0; j < 3; j++)
32 {
33 pt[j] += shape[i] * node[j];
34 }
35 }
36 return pt;
37}

References MeshLib::Element::getNode(), and MeshLib::Element::getNumberOfNodes().

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector().

◆ createBranchProperty()

BranchProperty ProcessLib::LIE::createBranchProperty ( MeshLib::Node const & branchNode,
FractureProperty const & master_frac,
FractureProperty const & slave_frac )
inline

Definition at line 100 of file FractureProperty.h.

103{
104 BranchProperty branch{branchNode, master_frac.fracture_id,
105 slave_frac.fracture_id};
106
107 // set a normal vector from the master to the slave fracture
108 Eigen::Vector3d branch_vector =
109 slave_frac.point_on_fracture - branch.coords;
110 double sign = (branch_vector.dot(master_frac.normal_vector) < 0) ? -1 : 1;
111 branch.normal_vector_branch = sign * master_frac.normal_vector;
112 return branch;
113}
Eigen::Vector3d const coords

References ProcessLib::LIE::FractureProperty::fracture_id, ProcessLib::LIE::FractureProperty::normal_vector, and ProcessLib::LIE::FractureProperty::point_on_fracture.

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess().

◆ duGlobalEnrichments()

std::vector< double > ProcessLib::LIE::duGlobalEnrichments ( std::size_t this_frac_id,
std::vector< FractureProperty * > const & frac_props,
std::vector< JunctionProperty * > const & junction_props,
std::unordered_map< int, int > const & fracID_to_local,
Eigen::Vector3d const & x )

calculate the enrichment function for fracture relative displacements Remarks:

  • branch/junction intersections of two fractures are supported in 2D
Parameters
this_frac_idthe fracture ID
frac_propsfracture properties
junction_propsjunction properties
fracID_to_locala mapping table from a fracture ID to a local index in frac_props
xevaluating point coordinates
Returns
a vector of enrichment values for fracture relative displacements

Definition at line 82 of file LevelSetFunction.cpp.

88{
89 auto this_frac_local_index = fracID_to_local.at(this_frac_id);
90 auto const& this_frac = *frac_props[this_frac_local_index];
91 // pre-calculate levelsets for all fractures
92 std::vector<bool> levelsets(frac_props.size());
93 for (std::size_t i = 0; i < frac_props.size(); i++)
94 {
95 levelsets[i] = levelsetFracture(*frac_props[i], x);
96 }
97
98 std::vector<double> enrichments(frac_props.size() + junction_props.size());
99 enrichments[this_frac_local_index] = 1.0;
100
101 // fractures possibly with branches
102 if (frac_props.size() > 1)
103 {
104 for (auto const& branch : this_frac.branches_master)
105 {
106 if (branch.master_fracture_id != this_frac.fracture_id)
107 {
108 continue;
109 }
110
111 if (fracID_to_local.find(branch.slave_fracture_id) ==
112 fracID_to_local.end())
113 {
114 continue;
115 }
116
117 double sign = boost::math::sign(
118 this_frac.normal_vector.dot(branch.normal_vector_branch));
119 auto slave_fid = fracID_to_local.at(branch.slave_fracture_id);
120 double const enrich = levelsets[slave_fid] ? 1. : 0.;
121 enrichments[slave_fid] = sign * enrich;
122 }
123 }
124
125 // junctions
126 for (unsigned i = 0; i < junction_props.size(); i++)
127 {
128 auto const* junction = junction_props[i];
129 if (!ranges::contains(junction->fracture_ids, this_frac.fracture_id))
130 {
131 continue;
132 }
133
134 auto another_frac_id =
135 (junction->fracture_ids[0] == this_frac.fracture_id)
136 ? junction->fracture_ids[1]
137 : junction->fracture_ids[0];
138 auto fid = fracID_to_local.at(another_frac_id);
139 double const enrich = levelsets[fid] ? 1. : 0.;
140 enrichments[i + frac_props.size()] = enrich;
141 }
142
143 return enrichments;
144}
bool levelsetFracture(FractureProperty const &frac, Eigen::Vector3d const &x)

References levelsetFracture().

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector().

◆ getFractureMatrixDataInMesh()

void ProcessLib::LIE::getFractureMatrixDataInMesh ( MeshLib::Mesh const & mesh,
std::vector< MeshLib::Element * > & vec_matrix_elements,
std::vector< int > & vec_fracture_mat_IDs,
std::vector< std::vector< MeshLib::Element * > > & vec_fracture_elements,
std::vector< std::vector< MeshLib::Element * > > & vec_fracture_matrix_elements,
std::vector< std::vector< MeshLib::Node * > > & vec_fracture_nodes,
std::vector< std::pair< std::size_t, std::vector< int > > > & vec_branch_nodeID_matIDs,
std::vector< std::pair< std::size_t, std::vector< int > > > & vec_junction_nodeID_matIDs )

get data about fracture and matrix elements/nodes from a mesh

Parameters
meshA mesh which includes fracture elements, i.e. lower-dimensional elements. It is assumed that elements forming a fracture have a distinct material ID.
vec_matrix_elementsa vector of matrix elements
vec_fracture_mat_IDsfracture material IDs found in the mesh
vec_fracture_elementsa vector of fracture elements (grouped by fracture IDs)
vec_fracture_matrix_elementsa vector of fracture elements and matrix elements connecting to the fracture (grouped by fracture IDs)
vec_fracture_nodesa vector of fracture element nodes (grouped by fracture IDs)
vec_branch_nodeID_matIDsa vector of branch node IDs found in the mesh (and corresponding fracture material IDs)
vec_junction_nodeID_matIDsa vector of junction node IDs found in the mesh (and corresponding fracture material IDs)

Definition at line 213 of file MeshUtils.cpp.

224{
225 IsCrackTip isCrackTip(mesh);
226
227 // get vectors of matrix elements and fracture elements
228 vec_matrix_elements.reserve(mesh.getNumberOfElements());
229 std::vector<MeshLib::Element*> all_fracture_elements;
230 for (MeshLib::Element* e : mesh.getElements())
231 {
232 if (e->getDimension() == mesh.getDimension())
233 {
234 vec_matrix_elements.push_back(e);
235 }
236 else
237 {
238 all_fracture_elements.push_back(e);
239 }
240 }
241 DBUG("-> found total {:d} matrix elements and {:d} fracture elements",
242 vec_matrix_elements.size(), all_fracture_elements.size());
243
244 // get fracture material IDs
245 auto const material_ids = materialIDs(mesh);
246 if (!material_ids)
247 {
248 OGS_FATAL("Could not access MaterialIDs property from mesh.");
249 }
250 transform(cbegin(all_fracture_elements), cend(all_fracture_elements),
251 back_inserter(vec_fracture_mat_IDs),
252 [&material_ids](auto const* e)
253 { return (*material_ids)[e->getID()]; });
254
255 BaseLib::makeVectorUnique(vec_fracture_mat_IDs);
256 DBUG("-> found {:d} fracture material groups", vec_fracture_mat_IDs.size());
257
258 // create a vector of fracture elements for each material
259 vec_fracture_elements.resize(vec_fracture_mat_IDs.size());
260 for (unsigned frac_id = 0; frac_id < vec_fracture_mat_IDs.size(); frac_id++)
261 {
262 const auto frac_mat_id = vec_fracture_mat_IDs[frac_id];
263 std::vector<MeshLib::Element*>& vec_elements =
264 vec_fracture_elements[frac_id];
265 std::copy_if(all_fracture_elements.begin(), all_fracture_elements.end(),
266 std::back_inserter(vec_elements),
267 [&](MeshLib::Element const* const e)
268 { return (*material_ids)[e->getID()] == frac_mat_id; });
269 DBUG("-> found {:d} elements on the fracture {:d}", vec_elements.size(),
270 frac_id);
271 }
272
273 // get a vector of fracture nodes for each material
274 vec_fracture_nodes.resize(vec_fracture_mat_IDs.size());
275 for (unsigned frac_id = 0; frac_id < vec_fracture_mat_IDs.size(); frac_id++)
276 {
277 std::vector<MeshLib::Node*>& vec_nodes = vec_fracture_nodes[frac_id];
278 for (MeshLib::Element const* const e : vec_fracture_elements[frac_id])
279 {
280 for (unsigned i = 0; i < e->getNumberOfNodes(); i++)
281 {
282 if (isCrackTip(*e->getNode(i)))
283 {
284 continue;
285 }
286 vec_nodes.push_back(const_cast<MeshLib::Node*>(e->getNode(i)));
287 }
288 }
291 DBUG("-> found {:d} nodes on the fracture {:d}", vec_nodes.size(),
292 frac_id);
293 }
294
295 // find branch/junction nodes which connect to multiple fractures
296 std::vector<std::vector<MeshLib::Element*>> intersected_fracture_elements;
297 findFracutreIntersections(mesh, vec_fracture_mat_IDs, vec_fracture_nodes,
298 intersected_fracture_elements,
299 vec_branch_nodeID_matIDs,
300 vec_junction_nodeID_matIDs);
301
302 // create a vector fracture elements and connected matrix elements,
303 // which are passed to a DoF table
304 for (unsigned fid = 0; fid < vec_fracture_elements.size(); fid++)
305 {
306 auto const& fracture_elements = vec_fracture_elements[fid];
307 std::vector<MeshLib::Element*> vec_ele;
308 // first, collect matrix elements
309 for (MeshLib::Element const* const e : fracture_elements)
310 {
311 // it is sufficient to iterate over base nodes, because they are
312 // already connected to all neighbours
313 for (unsigned i = 0; i < e->getNumberOfBaseNodes(); i++)
314 {
315 MeshLib::Node const* node = e->getNode(i);
316 if (isCrackTip(*node))
317 {
318 continue;
319 }
320 auto const& elements_connected_to_node =
321 mesh.getElementsConnectedToNode(*node);
322 for (unsigned j = 0; j < elements_connected_to_node.size(); j++)
323 {
324 // only matrix elements
325 if (elements_connected_to_node[j]->getDimension() <
326 mesh.getDimension())
327 {
328 continue;
329 }
330 vec_ele.push_back(const_cast<MeshLib::Element*>(
331 elements_connected_to_node[j]));
332 }
333 }
334 }
337
338 // second, append fracture elements
339 std::copy(fracture_elements.begin(), fracture_elements.end(),
340 std::back_inserter(vec_ele));
341 // thirdly, append intersected fracture elements
342 std::copy(intersected_fracture_elements[fid].begin(),
343 intersected_fracture_elements[fid].end(),
344 std::back_inserter(vec_ele));
345
346 vec_fracture_matrix_elements.push_back(vec_ele);
347 }
348}
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void makeVectorUnique(std::vector< T > &v)
Definition Algorithm.h:176
bool idsComparator(T const a, T const b)
Definition Mesh.h:206
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition Mesh.cpp:268
void findFracutreIntersections(MeshLib::Mesh const &mesh, std::vector< int > const &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Node * > > const &vec_fracture_nodes, std::vector< std::vector< MeshLib::Element * > > &intersected_fracture_elements, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)
Definition MeshUtils.cpp:70
unsigned getDimension(MeshLib::MeshElemType eleType)

References DBUG(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getElementsConnectedToNode(), MeshLib::Mesh::getNumberOfElements(), MeshLib::idsComparator(), BaseLib::makeVectorUnique(), and OGS_FATAL.

Referenced by ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::HydroMechanicsProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess(), and anonymous_namespace{postLIE.cpp}::postVTU().

◆ levelsetBranch()

bool ProcessLib::LIE::levelsetBranch ( BranchProperty const & branch,
Eigen::Vector3d const & x )

Definition at line 39 of file LevelSetFunction.cpp.

40{
41 return branch.normal_vector_branch.dot(x - branch.coords) > 0;
42}

References ProcessLib::LIE::BranchProperty::coords, and ProcessLib::LIE::BranchProperty::normal_vector_branch.

◆ levelsetFracture()

bool ProcessLib::LIE::levelsetFracture ( FractureProperty const & frac,
Eigen::Vector3d const & x )

Definition at line 34 of file LevelSetFunction.cpp.

35{
36 return frac.normal_vector.dot(x - frac.point_on_fracture) > 0;
37}

References ProcessLib::LIE::FractureProperty::normal_vector, and ProcessLib::LIE::FractureProperty::point_on_fracture.

Referenced by duGlobalEnrichments(), and uGlobalEnrichments().

◆ setFractureProperty()

void ProcessLib::LIE::setFractureProperty ( int const dim,
MeshLib::Element const & e,
FractureProperty & frac_prop )
inline

configure fracture property based on a fracture element assuming a fracture is a straight line/flat plane

Definition at line 76 of file FractureProperty.h.

78{
79 auto& n = frac_prop.normal_vector;
80 // 1st node is used but using other node is also possible, because
81 // a fracture is not curving
82 for (int j = 0; j < 3; j++)
83 {
84 frac_prop.point_on_fracture[j] = getCenterOfGravity(e).data()[j];
85 }
86
87 const MeshLib::ElementCoordinatesMappingLocal ele_local_coord(e, dim);
88
89 // Global to local rotation matrix:
90 Eigen::MatrixXd const global2local_rotation =
91 ele_local_coord.getRotationMatrixToGlobal().transpose();
92 n = global2local_rotation.row(dim - 1);
93
94 frac_prop.R = global2local_rotation.topLeftCorner(dim, dim);
95
96 DBUG("Normal vector of the fracture element {:d}: [{:g}, {:g}, {:g}]",
97 e.getID(), n[0], n[1], n[2]);
98}
const double * data() const
Definition Point3d.h:60
Eigen::MatrixXd R
Rotation matrix from global to local coordinates.

References MathLib::Point3d::data(), DBUG(), MeshLib::Element::getID(), MeshLib::ElementCoordinatesMappingLocal::getRotationMatrixToGlobal(), ProcessLib::LIE::FractureProperty::normal_vector, ProcessLib::LIE::FractureProperty::point_on_fracture, and ProcessLib::LIE::FractureProperty::R.

Referenced by ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::HydroMechanicsProcess(), and ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess().

◆ uGlobalEnrichments()

std::vector< double > ProcessLib::LIE::uGlobalEnrichments ( std::vector< FractureProperty * > const & frac_props,
std::vector< JunctionProperty * > const & junction_props,
std::unordered_map< int, int > const & fracID_to_local,
Eigen::Vector3d const & x )

calculate the enrichment function for displacements at the given point Remarks:

  • branch/junction intersections of two fractures are supported in 2D
Parameters
frac_propsfracture properties
junction_propsjunction properties
fracID_to_locala mapping table from a fracture ID to a local index in frac_props
xevaluating point coordinates
Returns
a vector of enrichment values for displacements

Definition at line 44 of file LevelSetFunction.cpp.

49{
50 // pre-calculate levelsets for all fractures
51 std::vector<bool> levelsets(frac_props.size());
52 for (std::size_t i = 0; i < frac_props.size(); i++)
53 {
54 levelsets[i] = levelsetFracture(*frac_props[i], x);
55 }
56
57 std::vector<double> enrichments(frac_props.size() + junction_props.size());
58 // fractures possibly with branches
59 for (std::size_t i = 0; i < frac_props.size(); i++)
60 {
61 auto const* frac = frac_props[i];
62 enrichments[i] = Heaviside(
63 std::accumulate(cbegin(frac->branches_slave),
64 cend(frac->branches_slave), levelsets[i],
65 [&](bool const enrich, auto const& branch)
66 { return enrich && levelsetBranch(branch, x); }));
67 }
68
69 // junctions
70 for (std::size_t i = 0; i < junction_props.size(); i++)
71 {
72 auto const* junction = junction_props[i];
73 auto fid1 = fracID_to_local.at(junction->fracture_ids[0]);
74 auto fid2 = fracID_to_local.at(junction->fracture_ids[1]);
75 bool const enrich = levelsets[fid1] && levelsets[fid2];
76 enrichments[i + frac_props.size()] = Heaviside(enrich);
77 }
78
79 return enrichments;
80}

References levelsetFracture().

Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerMatrixNearFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian(), ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::assembleWithJacobianConcrete(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::initializeConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), and ProcessLib::LIE::HydroMechanics::HydroMechanicsLocalAssemblerMatrixNearFracture< ShapeFunctionDisplacement, ShapeFunctionPressure, GlobalDim >::postTimestepConcreteWithVector().