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 "CohesiveZoneModelBase.h" 13 : #include "TwoVector.h" 14 : 15 : /** 16 : * User object that computes bilinear mixed mode traction separation law. 17 : */ 18 : class BilinearMixedModeCohesiveZoneModel : virtual public CohesiveZoneModelBase 19 : { 20 : public: 21 : static InputParameters validParams(); 22 : 23 : BilinearMixedModeCohesiveZoneModel(const InputParameters & parameters); 24 : 25 : virtual void initialize() override; 26 : virtual void finalize() override; 27 : 28 : // Getters for analysis output 29 : Real getModeMixityRatio(const Node * const node) const; 30 : Real getCohesiveDamage(const Node * const node) const; 31 : Real getLocalDisplacementNormal(const Node * const node) const; 32 : Real getLocalDisplacementTangential(const Node * const node) const; 33 : 34 : protected: 35 : virtual void computeQpProperties() override; 36 : virtual void computeQpIProperties() override; 37 : 38 28348 : virtual bool constrainedByOwner() const override { return false; } 39 : 40 : // @{ 41 : // Compute CZM bilinear traction law. 42 : virtual void computeModeMixity(const Node * const node); 43 : virtual void computeCriticalDisplacementJump(const Node * const node); 44 : virtual void computeFinalDisplacementJump(const Node * const node); 45 : virtual void computeEffectiveDisplacementJump(const Node * const node); 46 : virtual void computeDamage(const Node * const node) override; 47 : // @} 48 : 49 : /// Encapsulate the CZM constitutive behavior. 50 : virtual void computeCZMTraction(const Node * const node) override; 51 : 52 : /// Zero compressive traction 53 : const bool _set_compressive_traction_to_zero; 54 : 55 : /// Map from degree of freedom to mode mixity ratio (AD needed?) 56 : std::unordered_map<const DofObject *, ADReal> _dof_to_mode_mixity_ratio; 57 : 58 : /// The normal strength material property 59 : const MaterialProperty<Real> & _normal_strength; 60 : 61 : /// The shear strength material property 62 : const MaterialProperty<Real> & _shear_strength; 63 : 64 : /// Fracture parameter mode I 65 : const MaterialProperty<Real> & _GI_c; 66 : 67 : /// Fracture parameter mode II 68 : const MaterialProperty<Real> & _GII_c; 69 : 70 : /// Interpolated value of normal_strength 71 : ADReal _normal_strength_interpolation; 72 : 73 : /// Interpolated value of shear_strength 74 : ADReal _shear_strength_interpolation; 75 : 76 : /// Interpolated value of fracture paramter mode I 77 : ADReal _GI_c_interpolation; 78 : 79 : /// Interpolated value of fracture paramter mode II 80 : ADReal _GII_c_interpolation; 81 : 82 : // Penalty stiffness for bilinear traction model 83 : const Real _penalty_stiffness_czm; 84 : 85 : /// Mixed-mode propagation criterion 86 : enum class MixedModeCriterion 87 : { 88 : POWER_LAW, 89 : BK 90 : } _mix_mode_criterion; 91 : 92 : /// Power law parameter for bilinear traction model 93 : const Real _power_law_parameter; 94 : 95 : /// Viscosity for damage model 96 : const Real _viscosity; 97 : 98 : /// Parameter for the regularization of the Macaulay bracket 99 : const Real _regularization_alpha; 100 : 101 : // @{ 102 : // Strength material properties at the nodes 103 : std::unordered_map<const DofObject *, ADReal> _dof_to_normal_strength; 104 : std::unordered_map<const DofObject *, ADReal> _dof_to_shear_strength; 105 : // @} 106 : 107 : // @{ 108 : // Fracture properties at the nodes 109 : std::unordered_map<const DofObject *, ADReal> _dof_to_GI_c; 110 : std::unordered_map<const DofObject *, ADReal> _dof_to_GII_c; 111 : // @} 112 : 113 : // @{ 114 : // The parameters in the damage evolution law: Maps node to damage 115 : std::unordered_map<const DofObject *, ADReal> _dof_to_delta_initial; 116 : std::unordered_map<const DofObject *, ADReal> _dof_to_delta_final; 117 : std::unordered_map<const DofObject *, ADReal> _dof_to_delta_max; 118 : // @} 119 : };