OGS
|
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) |
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.
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().
|
inline |
Definition at line 100 of file FractureProperty.h.
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().
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:
this_frac_id | the fracture ID |
frac_props | fracture properties |
junction_props | junction properties |
fracID_to_local | a mapping table from a fracture ID to a local index in frac_props |
x | evaluating point coordinates |
Definition at line 82 of file LevelSetFunction.cpp.
References levelsetFracture().
Referenced by ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian(), and ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector().
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
mesh | A 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_elements | a vector of matrix elements |
vec_fracture_mat_IDs | fracture material IDs found in the mesh |
vec_fracture_elements | a vector of fracture elements (grouped by fracture IDs) |
vec_fracture_matrix_elements | a vector of fracture elements and matrix elements connecting to the fracture (grouped by fracture IDs) |
vec_fracture_nodes | a vector of fracture element nodes (grouped by fracture IDs) |
vec_branch_nodeID_matIDs | a vector of branch node IDs found in the mesh (and corresponding fracture material IDs) |
vec_junction_nodeID_matIDs | a vector of junction node IDs found in the mesh (and corresponding fracture material IDs) |
Definition at line 213 of file MeshUtils.cpp.
References DBUG(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getElementsConnectedToNode(), MeshLib::Mesh::getNumberOfElements(), 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().
bool ProcessLib::LIE::levelsetBranch | ( | BranchProperty const & | branch, |
Eigen::Vector3d const & | x ) |
Definition at line 39 of file LevelSetFunction.cpp.
References ProcessLib::LIE::BranchProperty::coords, and ProcessLib::LIE::BranchProperty::normal_vector_branch.
bool ProcessLib::LIE::levelsetFracture | ( | FractureProperty const & | frac, |
Eigen::Vector3d const & | x ) |
Definition at line 34 of file LevelSetFunction.cpp.
References ProcessLib::LIE::FractureProperty::normal_vector, and ProcessLib::LIE::FractureProperty::point_on_fracture.
Referenced by duGlobalEnrichments(), and uGlobalEnrichments().
|
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.
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().
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:
frac_props | fracture properties |
junction_props | junction properties |
fracID_to_local | a mapping table from a fracture ID to a local index in frac_props |
x | evaluating point coordinates |
Definition at line 44 of file LevelSetFunction.cpp.
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().