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 : #include "GeneralUserObject.h" 13 : #include "CrackFrontPointsProvider.h" 14 : #include "BoundaryRestrictable.h" 15 : #include <set> 16 : #include "ADRankTwoTensorForward.h" 17 : 18 : class AuxiliarySystem; 19 : 20 : // libMesh forward declarations 21 : namespace libMesh 22 : { 23 : class QBase; 24 : } 25 : 26 : // fixme lynn this needs to be removed from the globalname space in a seperate commit after Grizzly 27 : // is fixed to use includeCrackFrontDefinitionParams 28 : void addCrackFrontDefinitionParams(InputParameters & params); 29 : /** 30 : * Class used in fracture integrals to define geometric characteristics of the crack front 31 : */ 32 : class CrackFrontDefinition : public GeneralUserObject, public BoundaryRestrictable 33 : { 34 : public: 35 : static InputParameters validParams(); 36 : 37 : CrackFrontDefinition(const InputParameters & parameters); 38 : 39 : virtual void initialSetup() override; 40 : virtual void initialize() override; 41 : virtual void finalize() override; 42 : virtual void execute() override; 43 : 44 : /// used by Actions to add CrackFrontDefinitionParams 45 : static void includeCrackFrontDefinitionParams(InputParameters & params); 46 : 47 : /** 48 : * Change the number of crack front nodes. As the crack grows, the number of crack fronts nodes 49 : * may keep increasing in many cases. 50 : */ 51 : void updateNumberOfCrackFrontPoints(const std::size_t num_points); 52 : 53 : /** 54 : * Get the node pointer for a specified node on the crack front 55 : * @param node_index Index of the node 56 : * @return Pointer to node 57 : */ 58 : const Node * getCrackFrontNodePtr(const std::size_t node_index) const; 59 : 60 : /** 61 : * Get a Point object for a specified point on the crack front 62 : * @param point_index Index of the point 63 : * @return Point object 64 : */ 65 : const Point * getCrackFrontPoint(const std::size_t point_index) const; 66 : 67 : /** 68 : * Get the vector tangent to the crack front at a specified position 69 : * @param point_index Index of the point 70 : * @return tangent vector 71 : */ 72 : const RealVectorValue & getCrackFrontTangent(const std::size_t point_index) const; 73 : 74 : /** 75 : * Get the vector normal to the crack front at a specified position 76 : * @param point_index Index of the point 77 : * @return normal vector 78 : */ 79 : const RealVectorValue & getCrackFrontNormal(const std::size_t point_index) const; 80 : 81 : /** 82 : * Get the length of the line segment on the crack front ahead of the specified position 83 : * @param point_index Index of the point 84 : * @return Line segment length 85 : */ 86 : Real getCrackFrontForwardSegmentLength(const std::size_t point_index) const; 87 : 88 : /** 89 : * Get the length of the line segment on the crack front behind the specified position 90 : * @param point_index Index of the point 91 : * @return Line segment length 92 : */ 93 : Real getCrackFrontBackwardSegmentLength(const std::size_t point_index) const; 94 : 95 : /** 96 : * Get the unit vector of the crack extension direction at the specified position 97 : * @param point_index Index of the point 98 : * @return Crack extension direction vector 99 : */ 100 : const RealVectorValue & getCrackDirection(const std::size_t point_index) const; 101 : 102 : /** 103 : * Get the distance along the crack front from the beginning of the crack to the 104 : * specified position 105 : * @param point_index Index of the point 106 : * @return Distance along crack 107 : */ 108 : Real getDistanceAlongFront(const std::size_t point_index) const; 109 : 110 : /** 111 : * Whether the distance along the crack front is available as an angle 112 : * @return true if it is available as an angle 113 : */ 114 : bool hasAngleAlongFront() const; 115 : 116 : /** 117 : * Get the angle along the crack front from the beginning of the crack to the 118 : * specified position 119 : * @param point_index Index of the point 120 : * @return Angle along crack 121 : */ 122 : Real getAngleAlongFront(const std::size_t point_index) const; 123 : 124 : /** 125 : * Get the number of points defining the crack front as a set of line segments 126 : * @return Number of points 127 : */ 128 : std::size_t getNumCrackFrontPoints() const; 129 : 130 : /** 131 : * Whether the fracture computations are treated as 2D for the model 132 : * @return true if treated as 2D 133 : */ 134 3655291 : bool treatAs2D() const { return _treat_as_2d; } 135 : 136 : /** 137 : * Is the crack defined by a mesh cutter object 138 : * @return true if using a mesh cutter 139 : */ 140 2459 : bool usingMeshCutter() const { return _use_mesh_cutter; } 141 : 142 : /** 143 : * Rotate a vector in the global coordinate coordinate system to the crack 144 : * front local coordinate system at a specified point on the crack. 145 : * @param vector Vector in global coordinates 146 : * @param point_index Index of the point 147 : * @return Vector in crack local coordinates 148 : */ 149 : RealVectorValue rotateToCrackFrontCoords(const RealVectorValue vector, 150 : const std::size_t point_index) const; 151 : /** 152 : * Rotate a RankTwoTensor in the global coordinate coordinate system to the crack 153 : * front local coordinate system at a specified point on the crack. 154 : * @param tensor Tensor in global coordinates 155 : * @param point_index Index of the point 156 : * @return Tensor in crack local coordinates 157 : */ 158 : RankTwoTensor rotateToCrackFrontCoords(const RankTwoTensor tensor, 159 : const std::size_t point_index) const; 160 : 161 : /** 162 : * Rotate a vector from crack front cartesian coordinate to global cartesian coordinate 163 : * @param vector Vector in crack local coordinates 164 : * @param point_index Index of the point 165 : * @return Vector in global coordinates 166 : */ 167 : RealVectorValue rotateFromCrackFrontCoordsToGlobal(const RealVectorValue vector, 168 : const std::size_t point_index) const; 169 : 170 : /** 171 : * Calculate r and theta of a point in the crack front polar coordinates for a given 172 : * crack point index. 173 : * @param qp The Point for which coordinates are evaluated 174 : * @param point_index the crack front point index 175 : * @param r Value of the radial coordinate computed in this function 176 : * @param theta Value of the theta coordinate computed in this function 177 : */ 178 : void calculateRThetaToCrackFront(const Point qp, 179 : const std::size_t point_index, 180 : Real & r, 181 : Real & theta) const; 182 : 183 : /** 184 : * Calculate r and theta of a point in the crack front polar coordinate relative to the 185 : * closest crack front point. This function loops over all crack front points 186 : * to find the one closest to the specified point 187 : * @param qp The Point for which coordinates are evaluated 188 : * @param r Value of the radial coordinate computed in this function 189 : * @param theta Value of the theta coordinate computed in this function 190 : * @return Index of the closest crack front point 191 : */ 192 : std::size_t calculateRThetaToCrackFront(const Point qp, Real & r, Real & theta) const; 193 : 194 : /** 195 : * Determine whether a given node is on one of the boundaries that intersects an end 196 : * of the crack front 197 : * @param node Pointer to node 198 : * @return true if the node is on an intersecting boundary 199 : */ 200 : bool isNodeOnIntersectingBoundary(const Node * const node) const; 201 : 202 : /** 203 : * Determine whether a given crack front point is on one of the boundaries that 204 : * intersects an end of the crack front 205 : * @param point_index the crack front point index 206 : * @return true if the point is on an intersecting boundary 207 : */ 208 : bool isPointWithIndexOnIntersectingBoundary(const std::size_t point_index) const; 209 : 210 : /** 211 : * Get the strain in the direction tangent to the crack front at a given point 212 : * @param node_index the crack front node index 213 : * @return Tangential strain 214 : */ 215 : Real getCrackFrontTangentialStrain(const std::size_t node_index) const; 216 : 217 : /** 218 : * Determine whether the crack front was defined using nodes 219 : * @return true if it was defined using nodes 220 : */ 221 : bool hasCrackFrontNodes() const 222 : { 223 160 : return _geom_definition_method == CRACK_GEOM_DEFINITION::CRACK_FRONT_NODES; 224 : } 225 : 226 : /** 227 : * Determine whether a node is contained within a specified volume integral element ring 228 : * for a given node on the crack front 229 : * @param ring_index Index of the ring 230 : * @param connected_node_id ID of the node 231 : * @param node_index Index of the crack front node 232 : * @return true if the node is in the ring 233 : */ 234 : bool isNodeInRing(const std::size_t ring_index, 235 : const dof_id_type connected_node_id, 236 : const std::size_t node_index) const; 237 : 238 : /** 239 : * Compute the q function for the case where it is defined geometrically 240 : * @param crack_front_point_index Index of the point on the crack front 241 : * @param ring_index Index of the volume integral ring 242 : * @param current_node Node at which q is evaluated 243 : * @return q 244 : */ 245 : Real DomainIntegralQFunction(std::size_t crack_front_point_index, 246 : std::size_t ring_index, 247 : const Node * const current_node) const; 248 : 249 : /** 250 : * Compute the q function for the case where it is defined through element connectivity 251 : * @param crack_front_point_index Index of the point on the crack front 252 : * @param ring_index Index of the volume integral ring 253 : * @param current_node Node at which q is evaluated 254 : * @return q 255 : */ 256 : Real DomainIntegralTopologicalQFunction(std::size_t crack_front_point_index, 257 : std::size_t ring_index, 258 : const Node * const current_node) const; 259 : 260 : /** 261 : Set the value of _is_cutter_modified 262 : */ 263 : void isCutterModified(const bool is_cutter_modified); 264 : 265 : protected: 266 : /// Enum used to define the method for computing the crack extension direction 267 : const enum class DIRECTION_METHOD { 268 : CRACK_DIRECTION_VECTOR, 269 : CRACK_MOUTH, 270 : CURVED_CRACK_FRONT 271 : } _direction_method; 272 : 273 : /// Enum used to define the method for computing the crack extension direction 274 : /// at the ends of the crack 275 : const enum class END_DIRECTION_METHOD { 276 : NO_SPECIAL_TREATMENT, 277 : END_CRACK_DIRECTION_VECTOR, 278 : END_CRACK_TANGENT_VECTOR 279 : } _end_direction_method; 280 : 281 : /// Enum used to define the type of the nodes on the crack front (end or middle) 282 : enum CRACK_NODE_TYPE 283 : { 284 : MIDDLE_NODE, 285 : END_1_NODE, 286 : END_2_NODE 287 : }; 288 : 289 : /// Enum used to define whether the crack front is defined using nodes or points 290 : enum class CRACK_GEOM_DEFINITION 291 : { 292 : CRACK_FRONT_NODES, 293 : CRACK_FRONT_POINTS 294 : } _geom_definition_method; 295 : 296 : /// Reference to the auxiliary system 297 : AuxiliarySystem & _aux; 298 : /// Reference to the mesh 299 : MooseMesh & _mesh; 300 : 301 : /// Crack front nodes ordered from the start to end of the crack front 302 : std::vector<dof_id_type> _ordered_crack_front_nodes; 303 : /// Vector of points along the crack front 304 : std::vector<Point> _crack_front_points; 305 : /// Vector of tangent directions along the crack front 306 : std::vector<RealVectorValue> _tangent_directions; 307 : /// Vector of crack extension directions along the crack front 308 : std::vector<RealVectorValue> _crack_directions; 309 : /// Vector of segment lengths along the crack front 310 : std::vector<std::pair<Real, Real>> _segment_lengths; 311 : /// Vector of distances along the crack front 312 : std::vector<Real> _distances_along_front; 313 : /// Vector of angles along the crack front 314 : std::vector<Real> _angles_along_front; 315 : /// Vector of tangential strain along the crack front 316 : std::vector<Real> _strain_along_front; 317 : /// Vector of rotation matrices along the crack front 318 : std::vector<RankTwoTensor> _rot_matrix; 319 : /// Overall length of the crack 320 : Real _overall_length; 321 : /// Fixed vector optionally used to define crack extension direction 322 : RealVectorValue _crack_direction_vector; 323 : /// Fixed vector optionally used to define crack extension direction at end 1 of crack front 324 : RealVectorValue _crack_direction_vector_end_1; 325 : /// Fixed vector optionally used to define crack extension direction at end 2 of crack front 326 : RealVectorValue _crack_direction_vector_end_2; 327 : /// Fixed vector optionally used to define crack tangent direction at end 1 of crack front 328 : RealVectorValue _crack_tangent_vector_end_1; 329 : /// Fixed vector optionally used to define crack tangent direction at end 2 of crack front 330 : RealVectorValue _crack_tangent_vector_end_2; 331 : /// Names of boundaries used to define location of crack mouth 332 : std::vector<BoundaryName> _crack_mouth_boundary_names; 333 : /// IDs of boundaries used to define location of crack mouth 334 : std::vector<BoundaryID> _crack_mouth_boundary_ids; 335 : /// Names of boundaries that intersect crack at its ends 336 : std::vector<BoundaryName> _intersecting_boundary_names; 337 : /// IDs of boundaries that intersect crack at its ends 338 : std::vector<BoundaryID> _intersecting_boundary_ids; 339 : /// Coordinates of crack mouth 340 : RealVectorValue _crack_mouth_coordinates; 341 : /// Vector normals to a nonplanar crack 342 : std::vector<RealVectorValue> _crack_plane_normals; 343 : /// Whether to treat the model as 2D for computation of fracture integrals 344 : bool _treat_as_2d; 345 : /// Whether to describe the crack as a mesh cutter 346 : bool _use_mesh_cutter; 347 : /// Indicator that shows if the cutter mesh is modified or not in the calculation step 348 : bool _is_cutter_modified; 349 : /// Whether the crack forms a closed loop 350 : bool _closed_loop; 351 : /// Out of plane axis when crack is treated as 2D 352 : unsigned int _axis_2d; 353 : /// Whether the crack plane is also a symmetry plane in the model 354 : bool _has_symmetry_plane; 355 : /// Which plane is the symmetry plane 356 : unsigned int _symmetry_plane; 357 : /// Names of the x, y, and z displacement variables 358 : ///@{ 359 : std::string _disp_x_var_name; 360 : std::string _disp_y_var_name; 361 : std::string _disp_z_var_name; 362 : ///@} 363 : /// Whether the T-stress is being computed 364 : bool _t_stress; 365 : /// Whether topological rings are used to define the q functions 366 : bool _q_function_rings; 367 : /// Numer of elements from crack tip to last topological ring 368 : std::size_t _last_ring; 369 : /// Numer of elements from crack tip to first topological ring 370 : std::size_t _first_ring; 371 : /// Data structure used to store information about topological rings 372 : /// Key is a pair of the crack front node index and ring id 373 : /// Data is a set of the IDs of the nodes in the ring for that crack front node 374 : std::map<std::pair<dof_id_type, std::size_t>, std::set<dof_id_type>> 375 : _crack_front_node_to_node_map; 376 : /// Method used to define the q function 377 : MooseEnum _q_function_type; 378 : /// Vector of bools indicating whether individual crack front points are on 379 : /// an intersecting boundary 380 : std::vector<bool> _is_point_on_intersecting_boundary; 381 : /// Vector of inner radii of the rings used for geometric q functions 382 : std::vector<Real> _j_integral_radius_inner; 383 : /// Vector of outer radii of the rings used for geometric q functions 384 : std::vector<Real> _j_integral_radius_outer; 385 : /// Pointer to a CrackFrontPointsProvider object optionally used to define 386 : /// the crack front points 387 : const CrackFrontPointsProvider * _crack_front_points_provider; 388 : /// Number of points coming from the CrackFrontPointsProvider 389 : std::size_t _num_points_from_provider; 390 : /// tolerance for matching nodes at crack front 391 : const Real _tol; 392 : 393 : /** 394 : * Get the set of all crack front nodes 395 : * @param nodes Set of nodes -- populated by this method 396 : */ 397 : void getCrackFrontNodes(std::set<dof_id_type> & nodes); 398 : 399 : /** 400 : * Arrange the crack front nodes by their position along the crack front, 401 : * and put them in the _ordered_crack_front_nodes member variable. 402 : * @param nodes Set of nodes to be ordered 403 : */ 404 : void orderCrackFrontNodes(std::set<dof_id_type> & nodes); 405 : 406 : /** 407 : * Determine which of the end nodes should be the starting point of the 408 : * crack front. 409 : * @param end_nodes Vector containing two end nodes. The order of this is 410 : * rearranged so that the first end node is the start of 411 : * the crack front, and the second is at the end. 412 : */ 413 : void orderEndNodes(std::vector<dof_id_type> & end_nodes); 414 : 415 : /** 416 : * For the case of a crack that is a complete loop, determine which of the 417 : * nodes should be the start and end nodes in a repeatable way. 418 : * @param end_nodes Vector containing two end nodes -- populated by this method 419 : * @param nodes Set of all nodes on the crack front 420 : * @param node_to_line_elem_map Map from crack front nodes to line elements on crack 421 : * front (see line_elems param) 422 : * @param line_elems Line elements on crack front defined as vectors of the 423 : * nodes on each end of the line elements. 424 : */ 425 : void 426 : pickLoopCrackEndNodes(std::vector<dof_id_type> & end_nodes, 427 : std::set<dof_id_type> & nodes, 428 : std::map<dof_id_type, std::vector<dof_id_type>> & node_to_line_elem_map, 429 : std::vector<std::vector<dof_id_type>> & line_elems); 430 : /** 431 : * Find the node with the maximum value of its coordinate. This is used to 432 : * deterministically find the first node on the crack front. First, the nodes 433 : * with the maximum coordinate in one direction are found, and then if there 434 : * are duplicates with that same coordinate, search through the other two 435 : * coordinates to find the node with the maximum coordinate in all 3 directions. 436 : * @param nodes Set of all nodes on the crack front 437 : * @param dir0 First coordinate direction in which to order the coordinates 438 : */ 439 : dof_id_type maxNodeCoor(std::vector<Node *> & nodes, unsigned int dir0 = 0); 440 : 441 : /** 442 : * Update the data structures defining the crack front geometry such as 443 : * the ordered crack front nodes/points and other auxiliary data. 444 : */ 445 : void updateCrackFrontGeometry(); 446 : 447 : /** 448 : * compute node and coordinate data for crack fronts defined by crack_mouth_boundary_ids sidesets 449 : */ 450 : void computeCrackMouthNodes(); 451 : 452 : /** 453 : * Compute crack plane face normals for cracks that have a curved crack front but do not use a 454 : * mesh cutter. 455 : */ 456 : void computeCurvedCrackFrontCrackPlaneNormals(); 457 : 458 : /** 459 : * Compute the direction of crack extension for a given point on the crack front. 460 : * @param crack_front_point Point on the crack front 461 : * @param tangent_direction Tangent direction vector for the crack front point 462 : * @param ntype Node type such as MIDDLE_NODE, END_1_NODE, END_2_NODE 463 : * @param crack_front_point_index Index of the point on the crack front 464 : */ 465 : RealVectorValue calculateCrackFrontDirection(const Point & crack_front_point, 466 : const RealVectorValue & tangent_direction, 467 : const CRACK_NODE_TYPE ntype, 468 : const std::size_t crack_front_point_index = 0) const; 469 : 470 : /** 471 : * Compute the strain in the direction tangent to the crack at all points on the 472 : * crack front. 473 : */ 474 : void calculateTangentialStrainAlongFront(); 475 : 476 : /** 477 : * Create the data defining the rings used to define the q function when the topological 478 : * option is used to define the q function. 479 : */ 480 : void createQFunctionRings(); 481 : 482 : /** 483 : * Find nodes that are connected through elements to the nodes in the previous 484 : * node ring. 485 : * @param nodes_new_ring Nodes in the new ring -- populated by this method 486 : * @param nodes_old_ring Nodes in the previous ring 487 : * @param nodes_all_rings Nodes in all other rings to be excluded from the new ring 488 : * @param nodes_neighbor1 Nodes in the neighboring ring to one side to be excluded from the new 489 : * ring 490 : * @param nodes_neighbor2 Nodes in the neighboring ring to the other side to be excluded from the 491 : * new ring 492 : * @param nodes_to_elem_map Map of nodes to connected elements 493 : */ 494 : void addNodesToQFunctionRing(std::set<dof_id_type> & nodes_new_ring, 495 : const std::set<dof_id_type> & nodes_old_ring, 496 : const std::set<dof_id_type> & nodes_all_rings, 497 : const std::set<dof_id_type> & nodes_neighbor1, 498 : const std::set<dof_id_type> & nodes_neighbor2, 499 : std::vector<std::vector<const Elem *>> & nodes_to_elem_map); 500 : 501 : /** 502 : * Project a point to a specified point along the crack front and compute the 503 : * projected normal and tangential distance to the front 504 : * @param dist_to_front Projected normal distance to the front -- computed by this method 505 : * @param dist_along_tangent Project tangent distance to the front -- computed by this method 506 : * @param crack_front_point_index Index of the point on the crack front 507 : * @param current_node Node to be projected to the front 508 : */ 509 : void projectToFrontAtPoint(Real & dist_to_front, 510 : Real & dist_along_tangent, 511 : std::size_t crack_front_point_index, 512 : const Node * const current_node) const; 513 : };