- boundaryThe list of boundary IDs from the mesh where the cohesive zone will be applied
C++ Type:std::vector<BoundaryName>
Unit:(no unit assumed)
Controllable:No
Description:The list of boundary IDs from the mesh where the cohesive zone will be applied
- displacementsThe nonlinear displacement variables for the problem
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:The nonlinear displacement variables for the problem
Cohesive Zone Action System
Action to create an instance of the cohesive zone model kernel for each displacement component
Description
The TensorMechanics
system provides a cohesive zone modeling capability that can be used to introduce traction-separation models running at the interfaces between regions modeled with continuum finite elements. The implemented cohesive zone model (CZM) is based on a Discrete Galerkin approach and therefore does not require cohesive elements. The CZM is formulated in terms of traction separation laws and requires five ingredients:
BreakMeshByBlockGenerator
ComputeDisplacementJump
ComputeLocalTraction
ComputeGlobalTraction
CZMInterfaceKernel
The BreakMeshByBlockGenerator is utilized to create the cohesive zone interface by splitting a monolithic mesh into blocks by adding the required nodes and boundaries between each block pair. The split mesh allows to compute a displacement jump at each quadrature point on the interface. The schematic below is an example of using BreakMeshByBlockGenerator
on a 3-blocks, 2-dimensional mesh. The generated interfaces are highlighted in yellow.
The ComputeDisplacementJump
object computes the displacement jump across the cohesive zone according to the selected formulation. The ComputeLocalTraction
provides the cohesive zone response in the natural interface coordinate system. The ComputeGlobalTraction
object computes the traction in global coordinates and its derivative w.r.t. the displacement jump in global coordinates, . The CZMInterfaceKernel
imposes equilibrium by adding the proper residual to the system, and provides the analytic Jacobian.
The provided CZMInterfaceKernel
assume the ComputeLocalTraction
is only function of the displacement . If one wants to implement a traction separation law depending upon other variables, such as bulk stress or temperature, it is responsibility of the user to implement the proper Jacobian terms.
The CohesiveZone
action automatically adds the proper ComputeDisplacementJump
, ComputeGlobalTraction
, CZMInterfaceKernel
based on the kinematic_type
parameter value (see inputs). The flowchart below summarizes the flow of information of the cohesive zone modeling frameworks, and highlights the objects automatically added by the CohesizeZoneMaster
action.
Even when using the CohesiveZone
action it is the responsibility of the user to add the appropriate ComputeLocalTraction
constitutive model and BreakMeshByBlockGenerator
in the input file.
Supported Kinematic Formulations
The system supports two different kinematic formulations:
Small Strain
Total Lagrangian
Each type of formulations requires adding a different ComputeDisplacementJump
, ComputeGlobalTraction
, CZMInterfaceKernel
.
The Small Strain
formulations requires adding:
The Total Lagrangian
formulations requires adding:
As mentioned in previous section, the CohesiveZone
action adds the appropriate objects depending on the selected kinematic formulation.
To obtain a traction consistent with the bulk stress, the Small Strain
kinematic should be used together with a bulk Small Strain
formulation, and the Total Lagrangian
kinematics should be used together with a Finite Strain
formulation.
Cohesive Zone Constitutive Models
The implemented framework allows to implement two different types of cohesive zone constitutive models (e.g. traction separation laws): i) path independent, and ii) path dependent. The first type of models assumes the traction is only function of the total interface displacement jump, i.e. . This kind of models do not have history variables and can be implemented overriding the ComputeLocalTractionTotalBase
class. Instead, path dependent models have history variables and assume the traction increment is a function of the interface displacement jump increment, and its old value, . Path dependent models can be implemented overriding the ComputeLocalTractionIncrementalBase
class. In both cases, the traction separation law should always be formulated in terms of the local interface response and assuming Small strain
. Furthermore, all constitutive laws should be implemented for the 3D general case, even if used for 1D or 2D simulations. The kinematic is already embedded in the ComputeDisplacementJump
and ComputeGlobalTraction
objects. Both types of ComputeLocalTraction
objects allow using either the Small Strain
and the Total Lagrangian
kinematics.
Input Examples
The following example show how to use the CohesiveZone
action.
[Physics]
[SolidMechanics]
[CohesiveZone]
[./czm_ik]
boundary = 'interface'
strain = FINITE
generate_output = 'traction_x traction_y traction_z jump_x jump_y jump_z normal_traction tangent_traction normal_jump tangent_jump pk1_traction_x pk1_traction_y pk1_traction_z'
[../]
[]
[]
[]
(moose/modules/solid_mechanics/test/tests/cohesive_zone_model/stretch_rotate_large_deformation.i)If necessary, multiple instances of the CohesiveZone
action can be added, for instance when different material properties base_name
are needed for different boundaries. The base_name
parameter used in the action should also be provided to the associated materials. The generate_output
parameter adds scalar quantities of the traction and displacement jump to the outputs. Available options are: traction_x traction_y traction_z normal_traction tangent_traction jump_x jump_y jump_z normal_jump tangent_jump pk1_traction_x pk1_traction_y pk1_traction_z `.
The name `traction
refers to the Cauchy traction, pk1_traction
refers to the first Piola-Kirchhoff traction, and jump
refers to the displacement jump across the cohesive interface. All the above vectors are defined in the global coordinate system. The normal_traction
and tangent_traction
are scalar values compute using CZMRealVectorScalar (the same is true for the normal_jump
and tangent_jump
).
[Physics/SolidMechanics/CohesiveZone]
[./czm_ik_012]
boundary = 'Block0_Block1 Block1_Block2'
base_name = 'czm_b012'
[../]
[./czm_ik_23]
boundary = 'Block2_Block3'
base_name = 'czm_b23'
[../]
[]
[Materials]
# cohesive materials
[./czm_3dc]
type = SalehaniIrani3DCTraction
boundary = 'Block0_Block1 Block1_Block2'
normal_gap_at_maximum_normal_traction = 1
tangential_gap_at_maximum_shear_traction = 0.5
maximum_normal_traction = 500
maximum_shear_traction = 300
base_name = 'czm_b012'
[../]
[./czm_elastic_incremental]
type = PureElasticTractionSeparationIncremental
boundary = 'Block2_Block3'
normal_stiffness = 500
tangent_stiffness = 300
base_name = 'czm_b23'
[../]
# bulk materials
[./stress]
type = ADComputeFiniteStrainElasticStress
[../]
[./elasticity_tensor]
type = ADComputeIsotropicElasticityTensor
youngs_modulus = 200e4
poissons_ratio = 0.3
[../]
[]
(moose/modules/solid_mechanics/test/tests/cohesive_zone_model/czm_multiple_action_and_materials.i)Input Parameters
- active__all__ If specified only the blocks named will be visited and made active
Default:__all__
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:If specified only the blocks named will be visited and made active
- base_nameMaterial property base name
C++ Type:std::string
Unit:(no unit assumed)
Controllable:No
Description:Material property base name
- inactiveIf specified blocks matching these identifiers will be skipped.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:If specified blocks matching these identifiers will be skipped.
- strainSMALLStrain formulation
Default:SMALL
C++ Type:MooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Strain formulation
- use_automatic_differentiationFalseWhether to use automatic differentiation to compute the Jacobian
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether to use automatic differentiation to compute the Jacobian
- verboseFalseDisplay extra information.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Display extra information.
Optional Parameters
- additional_generate_outputAdd scalar quantity output for stress and/or strain (will be appended to the list in `generate_output`)
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Add scalar quantity output for stress and/or strain (will be appended to the list in `generate_output`)
- additional_material_output_familySpecifies the family of FE shape functions to use for this variable.
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Specifies the family of FE shape functions to use for this variable.
- additional_material_output_orderSpecifies the order of the FE shape function to use for this variable.
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Specifies the order of the FE shape function to use for this variable.
- generate_outputAdd scalar quantity output for stress and/or strain
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Add scalar quantity output for stress and/or strain
- material_output_familySpecifies the family of FE shape functions to use for this variable.
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Specifies the family of FE shape functions to use for this variable.
- material_output_orderSpecifies the order of the FE shape function to use for this variable.
C++ Type:MultiMooseEnum
Unit:(no unit assumed)
Controllable:No
Description:Specifies the order of the FE shape function to use for this variable.
Output Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- diag_save_in_masterThe displacement diagonal preconditioner terms on the master side
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacement diagonal preconditioner terms on the master side
- diag_save_in_slaveThe displacement diagonal preconditioner terms on the slave side
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacement diagonal preconditioner terms on the slave side
- save_in_masterThe displacement residuals on the master side
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacement residuals on the master side
- save_in_slaveThe displacement residuals on the slave side
C++ Type:std::vector<AuxVariableName>
Unit:(no unit assumed)
Controllable:No
Description:The displacement residuals on the slave side
Advanced Parameters
Associated Actions
Available Actions
- Solid Mechanics App
- CommonCohesiveZoneActionStore common cohesive zone paramters
- CohesiveZoneActionAction to create an instance of the cohesive zone model kernel for each displacement component