Line data Source code
1 : //* This file is part of the MOOSE framework 2 : //* https://mooseframework.inl.gov 3 : //* 4 : //* All rights reserved, see COPYRIGHT for full restrictions 5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT 6 : //* 7 : //* Licensed under LGPL 2.1, please see LICENSE for details 8 : //* https://www.gnu.org/licenses/lgpl-2.1.html 9 : 10 : #pragma once 11 : 12 : // MOOSE includes 13 : #include "CrackFrontPointsProvider.h" 14 : #include "XFEMAppTypes.h" 15 : #include "SolidMechanicsAppTypes.h" 16 : 17 : #include "libmesh/libmesh_common.h" 18 : #include "libmesh/libmesh.h" // libMesh::invalid_uint 19 : #include "libmesh/elem.h" 20 : 21 : using namespace libMesh; 22 : 23 : class XFEM; 24 : 25 : namespace Xfem 26 : { 27 : /// Data structure defining a cut on an element edge 28 : struct CutEdge 29 : { 30 : /// ID of the first node on the edge 31 : unsigned int _id1; 32 : /// ID of the second node on the edge 33 : unsigned int _id2; 34 : /// Fractional distance along the edge (from node 1 to 2) where the cut is located 35 : Real _distance; 36 : /// Local ID of this side in the host element 37 : unsigned int _host_side_id; 38 : }; 39 : 40 : /** 41 : * Operator < for two CutEdge Objects 42 : * Needed to allow the use of std::set<CutEdge> 43 : * @param lhs CutEdge object on the left side of the comparison 44 : * @param rhs CutEdge object on the right side of the comparison 45 : * @return bool true if lhs < rhs 46 : */ 47 : inline bool 48 : operator<(const CutEdge & lhs, const CutEdge & rhs) 49 : { 50 884 : if (lhs._id1 != rhs._id1) 51 764 : return lhs._id1 < rhs._id1; 52 : else 53 120 : return lhs._id2 < rhs._id2; 54 : } 55 : 56 : /// Data structure defining a cut through a node 57 : struct CutNode 58 : { 59 : /// ID of the cut node 60 : unsigned int _id; 61 : /// Local ID of this node in the host element 62 : unsigned int _host_id; 63 : }; 64 : 65 : /// Data structure defining a cut through a face 66 191666 : struct CutFace 67 : { 68 : /// ID of the cut face 69 : unsigned int _face_id; 70 : /// IDs of all cut faces 71 : std::vector<unsigned int> _face_edge; 72 : /// Fractional distance along the cut edges where the cut is located 73 : std::vector<Real> _position; 74 : }; 75 : 76 : /// Data structure describing geometrically described cut through 2D element 77 : struct GeomMarkedElemInfo2D 78 : { 79 : /// Container for data about all cut edges in this element 80 : std::vector<CutEdge> _elem_cut_edges; 81 : /// Container for data about all cut nodes in this element 82 : std::vector<CutNode> _elem_cut_nodes; 83 : /// Container for data about all cut fragments in this element 84 : std::vector<CutEdge> _frag_cut_edges; 85 : /// Container for data about all cut edges in cut fragments in this element 86 : std::vector<std::vector<Point>> _frag_edges; 87 : }; 88 : 89 : /// Data structure describing geometrically described cut through 3D element 90 : struct GeomMarkedElemInfo3D 91 : { 92 : /// Container for data about all cut faces in this element 93 : std::vector<CutFace> _elem_cut_faces; 94 : /// Container for data about all faces this element's fragment 95 : std::vector<CutFace> _frag_cut_faces; 96 : /// Container for data about all cut faces in cut fragments in this element 97 : std::vector<std::vector<Point>> _frag_faces; 98 : }; 99 : 100 : } // namespace Xfem 101 : 102 : class GeometricCutUserObject : public CrackFrontPointsProvider 103 : { 104 : public: 105 : /** 106 : * Factory constructor, takes parameters so that all derived classes can be built using the same 107 : * constructor. 108 : */ 109 : static InputParameters validParams(); 110 : 111 : GeometricCutUserObject(const InputParameters & parameters, const bool uses_mesh = false); 112 : 113 : virtual void initialize() override; 114 : virtual void execute() override; 115 : virtual void threadJoin(const UserObject & y) override; 116 : virtual void finalize() override; 117 : 118 0 : virtual unsigned int getNumberOfCrackFrontPoints() const override 119 : { 120 0 : mooseError("getNumberOfCrackFrontPoints() is not implemented for this object."); 121 : } 122 : 123 : /** 124 : * Check to see whether a specified 2D element should be cut based on geometric 125 : * conditions 126 : * @param elem Pointer to the libMesh element to be considered for cutting 127 : * @param cut_edges Data structure filled with information about edges to be cut 128 : * @param cut_nodes Data structure filled with information about nodes to be cut 129 : * @return bool true if element is to be cut 130 : */ 131 : virtual bool cutElementByGeometry(const Elem * elem, 132 : std::vector<Xfem::CutEdge> & cut_edges, 133 : std::vector<Xfem::CutNode> & cut_nodes) const = 0; 134 : 135 : /** 136 : * Check to see whether a specified 3D element should be cut based on geometric 137 : * conditions 138 : * @param elem Pointer to the libMesh element to be considered for cutting 139 : * @param cut_faces Data structure filled with information about edges to be cut 140 : * @return bool true if element is to be cut 141 : */ 142 : virtual bool cutElementByGeometry(const Elem * elem, 143 : std::vector<Xfem::CutFace> & cut_faces) const = 0; 144 : 145 : /** 146 : * Check to see whether a fragment of a 2D element should be cut based on geometric conditions 147 : * @param frag_edges Data structure defining the current fragment to be considered 148 : * @param cut_edges Data structure filled with information about fragment edges to be cut 149 : * @return bool true if fragment is to be cut 150 : */ 151 : virtual bool cutFragmentByGeometry(std::vector<std::vector<Point>> & frag_edges, 152 : std::vector<Xfem::CutEdge> & cut_edges) const = 0; 153 : 154 : /** 155 : * Check to see whether a fragment of a 3D element should be cut based on geometric conditions 156 : * @param frag_faces Data structure defining the current fragment to be considered 157 : * @param cut_faces Data structure filled with information about fragment faces to be cut 158 : * @return bool true if fragment is to be cut 159 : */ 160 : virtual bool cutFragmentByGeometry(std::vector<std::vector<Point>> & frag_faces, 161 : std::vector<Xfem::CutFace> & cut_faces) const = 0; 162 : 163 : /** 164 : * Get the interface ID for this cutting object. 165 : * @return the interface ID 166 : */ 167 54013 : unsigned int getInterfaceID() const { return _interface_id; }; 168 : 169 : /** 170 : * Set the interface ID for this cutting object. 171 : * @param the interface ID 172 : */ 173 501 : void setInterfaceID(unsigned int interface_id) { _interface_id = interface_id; }; 174 : 175 : /** 176 : * Should the elements cut by this cutting object be healed in the current 177 : * time step? 178 : * @return true if the cut element should be healed 179 : */ 180 17495 : bool shouldHealMesh() const { return _heal_always; }; 181 : 182 : /** 183 : * Get CutSubdomainID telling which side the node belongs to relative to the cut. 184 : * The returned ID contains no physical meaning, but should be consistent throughout the 185 : * simulation. 186 : * @param node Pointer to the node 187 : * @return An unsigned int indicating the side 188 : */ 189 0 : virtual CutSubdomainID getCutSubdomainID(const Node * /*node*/) const 190 : { 191 0 : mooseError("Objects that inherit from GeometricCutUserObject should override the " 192 : "getCutSubdomainID method"); 193 : return 0; 194 : } 195 : 196 : /** 197 : * Get the CutSubdomainID for the given element. 198 : * @param node Pointer to the element 199 : * @return The CutSubdomainID 200 : */ 201 : CutSubdomainID getCutSubdomainID(const Elem * elem) const; 202 : 203 : protected: 204 : /// Pointer to the XFEM controller object 205 : std::shared_ptr<XFEM> _xfem; 206 : 207 : /// Associated interface id 208 : unsigned int _interface_id; 209 : 210 : /// Heal the mesh 211 : bool _heal_always; 212 : 213 : /// Time step information needed to advance a 3D crack only at the real beginning of a time step 214 : int _last_step_initialized; 215 : 216 : ///@{Containers with information about all 2D and 3D elements marked for cutting by this object 217 : std::map<unsigned int, std::vector<Xfem::GeomMarkedElemInfo2D>> _marked_elems_2d; 218 : std::map<unsigned int, std::vector<Xfem::GeomMarkedElemInfo3D>> _marked_elems_3d; 219 : ///@} 220 : 221 : ///@{ Methods to pack/unpack the _marked_elems_2d and _marked_elems_3d data into a structure suitable for parallel communication 222 : void serialize(std::string & serialized_buffer); 223 : void deserialize(std::vector<std::string> & serialized_buffers); 224 : ///@} 225 : };