- bulk_modulusThe bulk modulus $K$
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The bulk modulus $K$
- phase_fieldName of the phase-field (damage) variable
C++ Type:std::vector<VariableName>
Unit:(no unit assumed)
Controllable:No
Description:Name of the phase-field (damage) variable
- shear_modulusThe shear modulus $G$
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The shear modulus $G$
SmallDeformationIsotropicElasticity
Isotropic elasticity under small strain asumptions.
Overview
The strain energy density is defined as
This class features three most commonly used splits of the strain energy, i.e. no decomposition, spectral decomposition, and volumetric-deviatoric decomposition.
Example Input File Syntax
Input Parameters
- base_nameOptional parameter that allows the user to define multiple mechanics material systems on the same block, i.e. for multiple phases
C++ Type:std::string
Unit:(no unit assumed)
Controllable:No
Description:Optional parameter that allows the user to define multiple mechanics material systems on the same block, i.e. for multiple phases
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Unit:(no unit assumed)
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- boundaryThe list of boundaries (ids or names) from the mesh where this object applies
C++ Type:std::vector<BoundaryName>
Unit:(no unit assumed)
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- constant_onNONEWhen ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
Default:NONE
C++ Type:MooseEnum
Unit:(no unit assumed)
Options:NONE, ELEMENT, SUBDOMAIN
Controllable:No
Description:When ELEMENT, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps.When SUBDOMAIN, MOOSE will only call computeQpProperties() for the 0th quadrature point, and then copy that value to the other qps. Evaluations on element qps will be skipped
- declare_suffixAn optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- decompositionNONEThe decomposition method
Default:NONE
C++ Type:MooseEnum
Unit:(no unit assumed)
Options:NONE, SPECTRAL, VOLDEV
Controllable:No
Description:The decomposition method
- degradation_functiongThe degradation function
Default:g
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:The degradation function
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
- strain_energy_densitypsieName of the strain energy density computed by this material model
Default:psie
C++ Type:MaterialPropertyName
Unit:(no unit assumed)
Controllable:No
Description:Name of the strain energy density computed by this material model
- use_interpolated_stateFalseFor the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:For the old and older state use projected material properties interpolated at the quadrature points. To set up projection use the ProjectedStatefulMaterialStorageAction.
Optional 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.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:Yes
Description:Set the enabled status of the MooseObject.
- implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Determines whether this object is calculated using an implicit or explicit form
- seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:The seed for the master random number generator
- use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Advanced Parameters
- output_propertiesList of material properties, from this material, to output (outputs must also be defined to an output type)
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:List of material properties, from this material, to output (outputs must also be defined to an output type)
- outputsnone Vector of output names where you would like to restrict the output of variables(s) associated with this object
Default:none
C++ Type:std::vector<OutputName>
Unit:(no unit assumed)
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
Input Files
- (tutorials/mode1_cohesive_fracture/elasticity.i)
- (tutorials/small_deformation/elasticity.i)
- (tutorials/surfing_boundary_problem/fracture.i)
- (tutorials/quenching_bibeam/thermoelasticity.i)
- (test/tests/nucleation_models/klbf.i)
- (tutorials/mode1_brittle_fracture/elasticity.i)
- (tutorials/surfing_boundary_problem/elasticity.i)
- (test/tests/nucleation_models/klr.i)
- (tutorials/soil_desiccation/elasticity.i)
- (tutorials/small_deformation/elastoplasticity.i)
- (tutorials/mode1_ductile_fracture/elastoplasticity.i)
- (test/tests/nucleation_models/elasticity.i)
- (tutorials/matrix_fiber/elasticity.i)
- (tutorials/traction_separation/elasticity.i)
- (tutorials/adaptivity/elasticity.i)
- (tutorials/mode2_cohesive_fracture/elasticity.i)
- (tutorials/mode2_brittle_fracture/elasticity.i)
(tutorials/mode1_cohesive_fracture/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
psic = 14.88
l = 0.02
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};psic=${psic};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymax = 0.5
[]
[noncrack]
type = BoundingBoxNodeSetGenerator
input = gen
new_boundary = noncrack
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = BoxMarker
bottom_left = '0.4 0 0'
top_right = '1 0.05 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fy]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
save_in = fy
[]
[]
[BCs]
[ydisp]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = noncrack
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = top
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K} ${G} ${l} ${Gc} ${psic}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'elastic_strain psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fy]
type = NodalSum
variable = fy
boundary = top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 3.5e-3
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/small_deformation/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymax = 0.5
[]
[noncrack]
type = BoundingBoxNodeSetGenerator
input = gen
new_boundary = noncrack
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
construct_side_list_from_node_list = true
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fy]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
save_in = fy
[]
[]
[BCs]
[ydisp]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = noncrack
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = top
value = 0
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'K G'
prop_values = '${K} ${G}'
[]
[no_degradation]
type = NoDegradation
property_name = g
expression = 1
phase_field = d
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
output_properties = 'elastic_strain'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fy]
type = NodalSum
variable = fy
boundary = top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1e-5
end_time = 1e-4
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/surfing_boundary_problem/fracture.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${nx}
ny = ${ny}
xmax = 30
ymin = -5
ymax = 5
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = ${refine}
stop_time = 0
max_h_level = ${refine}
[Markers]
[marker]
type = BoxMarker
bottom_left = '0 -0.7 0'
top_right = '30 0.7 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[d]
[InitialCondition]
type = FunctionIC
function = 'if(y=0&x>=0&x<=5,1,0)'
[]
[]
[]
[AuxVariables]
[bounds_dummy]
[]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[psie_active]
order = CONSTANT
family = MONOMIAL
[]
[]
[Bounds]
[conditional]
type = VariableOldValueBounds
variable = bounds_dummy
bounded_variable = d
bound_type = lower
[]
[upper]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = d
bound_type = upper
bound_value = 1
[]
[]
[Kernels]
[diff]
type = ADPFFDiffusion
variable = d
fracture_toughness = Gc_delta
regularization_length = l
normalization_constant = c0
[]
[source]
type = ADPFFSource
variable = d
free_energy = psi
[]
[nuc_force]
type = ADCoefMatSource
variable = d
prop_names = 'ce'
coefficient = 1.0
[]
[]
[Materials]
[fracture_properties]
type = ADGenericConstantMaterial
prop_names = 'E K G lambda Gc l sigma_ts sigma_hs'
prop_values = '${E} ${K} ${G} ${Lambda} ${Gc} ${l} ${sigma_ts} ${sigma_hs}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 0'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[psi]
type = ADDerivativeParsedMaterial
property_name = psi
expression = 'g*psie_active+(Gc*delta/c0/l)*alpha'
coupled_variables = 'd psie_active'
material_property_names = 'delta alpha(d) g(d) Gc c0 l'
derivative_order = 1
[]
[nucleation_micro_force]
type = LDLNucleationMicroForce
phase_field = d
degradation_function = g
regularization_length = l
normalization_constant = c0
fracture_toughness = Gc
tensile_strength = sigma_ts
hydrostatic_strength = sigma_hs
delta = delta
h_correction = true
external_driving_force_name = ce
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
[]
[Gc_delta]
type = ADParsedMaterial
property_name = Gc_delta
expression = 'Gc*delta'
material_property_names = 'Gc delta'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -snes_type'
petsc_options_value = 'lu superlu_dist vinewtonrsls'
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
[]
[Outputs]
print_linear_residuals = false
[]
(tutorials/quenching_bibeam/thermoelasticity.i)
E_glass = 6.4e4
nu_glass = 0.2
K_glass = '${fparse E_glass/3/(1-2*nu_glass)}'
G_glass = '${fparse E_glass/2/(1+nu_glass)}'
alpha_glass = 3.25e-6
E_steel = 1.93e5
nu_steel = 0.29
K_steel = '${fparse E_steel/3/(1-2*nu_steel)}'
G_steel = '${fparse E_steel/2/(1+nu_steel)}'
alpha_steel = 1.73e-5
Gc = 0.4
psic = 0.05
l = 1
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};psic=${psic};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppGeneralFieldShapeEvaluationTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'gold/plate.msh'
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = BoxMarker
bottom_left = '27 7 0'
top_right = '60 28 0'
inside = REFINE
outside = DO_NOTHING
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[d]
[]
[T]
initial_condition = 1000
[]
[psie_active]
order = CONSTANT
family = MONOMIAL
block = glass
[]
[]
[AuxKernels]
[cooling]
type = FunctionAux
variable = T
function = '1000-t'
[]
[psie_active]
type = ADMaterialRealAux
variable = psie_active
property = psie_active
block = glass
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[]
[BCs]
[xdisp]
type = DirichletBC
variable = 'disp_x'
boundary = 'left'
value = 0
[]
[yfix]
type = DirichletBC
variable = 'disp_y'
boundary = 'left right'
value = 0
[]
[]
[Materials]
# Glass
[bulk_properties_glass]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K_glass} ${G_glass} ${l} ${Gc} ${psic}'
block = glass
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
block = glass
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
block = glass
[]
[eigenstrain_glass]
type = ADComputeThermalExpansionEigenstrain
eigenstrain_name = thermal_eigenstrain
stress_free_temperature = 1000
thermal_expansion_coeff = ${alpha_glass}
temperature = T
block = glass
[]
[strain_glass]
type = ADComputeSmallStrain
eigenstrain_names = thermal_eigenstrain
block = glass
[]
[elasticity_glass]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = SPECTRAL
block = glass
[]
[stress_glass]
type = ComputeSmallDeformationStress
elasticity_model = elasticity_glass
block = glass
[]
# Steel
[eigenstrain_steel]
type = ADComputeThermalExpansionEigenstrain
eigenstrain_name = thermal_eigenstrain
stress_free_temperature = 1000
thermal_expansion_coeff = ${alpha_steel}
temperature = T
block = steel
[]
[strain_steel]
type = ADComputeSmallStrain
eigenstrain_names = thermal_eigenstrain
block = steel
[]
[elasticity_steel]
type = ADComputeIsotropicElasticityTensor
shear_modulus = ${G_steel}
bulk_modulus = ${K_steel}
[]
[stress_steel]
type = ADComputeLinearElasticStress
block = steel
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-08
nl_abs_tol = 1e-10
dt = 0.5
end_time = 200
fixed_point_algorithm = picard
fixed_point_max_its = 20
fixed_point_rel_tol = 1e-08
fixed_point_abs_tol = 1e-10
accept_on_max_fixed_point_iteration = true
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(test/tests/nucleation_models/klbf.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${nx}
ny = ${ny}
xmax = 30
ymin = -5
ymax = 5
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = ${refine}
stop_time = 0
max_h_level = ${refine}
[Markers]
[marker]
type = BoxMarker
bottom_left = '0 -0.7 0'
top_right = '30 0.7 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[d]
[InitialCondition]
type = FunctionIC
function = 'if(y=0&x>=0&x<=5,1,0)'
[]
[]
[]
[AuxVariables]
[bounds_dummy]
[]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[psie_active]
order = CONSTANT
family = MONOMIAL
[]
[ce]
order = CONSTANT
family = MONOMIAL
[]
[]
[Bounds]
[conditional]
type = ConditionalBoundsAux
variable = bounds_dummy
bounded_variable = d
fixed_bound_value = 0
threshold_value = 0.95
[]
[upper]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = d
bound_type = upper
bound_value = 1
[]
[]
[Kernels]
[diff]
type = ADPFFDiffusion
variable = d
fracture_toughness = Gc
regularization_length = l
normalization_constant = c0
[]
[source]
type = ADPFFSource
variable = d
free_energy = psi
[]
[nuc_force]
type = ADCoefMatSource
variable = d
prop_names = 'ce'
[]
[]
[AuxKernels]
[get_ce]
type = ADMaterialRealAux
variable = ce
property = ce
[]
[]
[Materials]
[fracture_properties]
type = ADGenericConstantMaterial
prop_names = 'E K G lambda Gc l sigma_ts sigma_cs delta'
prop_values = '${E} ${K} ${G} ${Lambda} ${Gc} ${l} ${sigma_ts} ${sigma_cs} ${delta}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 0'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[psi]
type = ADDerivativeParsedMaterial
property_name = psi
expression = 'g*psie_active+(Gc/c0/l)*alpha'
coupled_variables = 'd psie_active'
material_property_names = 'alpha(d) g(d) Gc c0 l'
derivative_order = 1
[]
[kumar_material]
type = KLBFNucleationMicroForce
phase_field = d
normalization_constant = c0
tensile_strength = sigma_ts
compressive_strength = sigma_cs
delta = delta
external_driving_force_name = ce
output_properties = 'ce'
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -snes_type'
petsc_options_value = 'lu superlu_dist vinewtonrsls'
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
[]
[Outputs]
print_linear_residuals = false
[]
(tutorials/mode1_brittle_fracture/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
l = 0.02
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymax = 0.5
[]
[noncrack]
type = BoundingBoxNodeSetGenerator
input = gen
new_boundary = noncrack
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = BoxMarker
bottom_left = '0.4 0 0'
top_right = '1 0.05 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fy]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
save_in = fy
[]
[]
[BCs]
[ydisp]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = noncrack
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = top
value = 0
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'K G'
prop_values = '${K} ${G}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'elastic_strain psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fy]
type = NodalSum
variable = fy
boundary = top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 3.5e-3
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/surfing_boundary_problem/elasticity.i)
# This tutorial is based on example of material graphite in section 4.3.4
# See Kumar et. al. https://doi.org/10.1016/j.jmps.2020.104027.
# For regularization length 0.35mm, delta 1.16, mesh size 0.035,
# J intergal over Gc should return a value close to 1
# delta value 0.2 is updated for model 2022
# (see Kumar et. al. https://doi.org/10.1007/s10704-022-00653-z)
E = 9.8e3 #9.8Gpa
nu = 0.13
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Lambda = '${fparse E*nu/(1+nu)/(1-2*nu)}'
Gc = 9.1e-2 # 91N/m
l = 0.35
sigma_ts = 27
sigma_cs = 77
sigma_hs = '${fparse 2/3*sigma_ts*sigma_cs/(sigma_cs - sigma_ts)}'
c1 = '${fparse (1+nu)*sqrt(Gc)/sqrt(2*pi*E)}'
c2 = '${fparse (3-nu)/(1+nu)}'
ahead = 2
V = 20
# mesh size 0.035
nx = 108
ny = 36
refine = 3
[Functions]
[bc_func]
type = ParsedFunction
expression = c1*((x-V*t-ahead)^2+y^2)^(0.25)*(c2-cos(atan2(y,(x-V*t-ahead))))*sin(0.5*atan2(y,(x-V*t-ahead)))
symbol_names = 'c1 c2 V ahead'
symbol_values = '${c1} ${c2} ${V} ${ahead}'
[]
[]
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'E=${E};K=${K};G=${G};Lambda=${Lambda};Gc=${Gc};l=${l};nx=${nx};ny=${ny};refine=${refine};sigma_ts=${sigma_ts};sigma_hs=${sigma_hs}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = fracture
variable = 'd'
source_variable = 'd'
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = fracture
variable = 'disp_x disp_y strain_zz psie_active'
source_variable = 'disp_x disp_y strain_zz psie_active'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${nx}
ny = ${ny}
xmax = 30
ymin = -5
ymax = 5
[]
[gen2]
type = ExtraNodesetGenerator
input = gen
new_boundary = fix_point
coord = '0 -5'
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = ${refine}
stop_time = 0
max_h_level = ${refine}
[Markers]
[marker]
type = BoxMarker
bottom_left = '0 -0.7 0'
top_right = '30 0.7 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[]
[AuxVariables]
[d]
[InitialCondition]
type = FunctionIC
function = 'if(y=0&x>=0&x<=5,1,0)'
[]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[plane_stress]
type = ADWeakPlaneStress
variable = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[fix_x]
type = DirichletBC
variable = disp_x
boundary = fix_point
value = 0
[]
[bottom_y]
type = FunctionDirichletBC
variable = disp_y
boundary = bottom
function = bc_func
[]
[top_y]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = bc_func
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'E K G lambda Gc l'
prop_values = '${E} ${K} ${G} ${Lambda} ${Gc} ${l}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 0'
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'psie_active psie'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
[]
[]
[Postprocessors]
[Jint]
type = PhaseFieldJIntegral
J_direction = '1 0 0'
strain_energy_density = psie
displacements = 'disp_x disp_y'
boundary = 'left top right bottom'
[]
[Jint_over_Gc]
type = ParsedPostprocessor
pp_names = 'Jint'
expression = 'Jint/${Gc}'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1e-2
dtmin = 1e-2
end_time = 0.3
fixed_point_max_its = 500
accept_on_max_fixed_point_iteration = false
fixed_point_rel_tol = 1e-6
fixed_point_abs_tol = 1e-8
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(test/tests/nucleation_models/klr.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${nx}
ny = ${ny}
xmax = 30
ymin = -5
ymax = 5
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = ${refine}
stop_time = 0
max_h_level = ${refine}
[Markers]
[marker]
type = BoxMarker
bottom_left = '0 -0.7 0'
top_right = '30 0.7 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[d]
[InitialCondition]
type = FunctionIC
function = 'if(y=0&x>=0&x<=5,1,0)'
[]
[]
[]
[AuxVariables]
[bounds_dummy]
[]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[psie_active]
order = CONSTANT
family = MONOMIAL
[]
[ce]
order = CONSTANT
family = MONOMIAL
[]
[]
[Bounds]
[conditional]
type = ConditionalBoundsAux
variable = bounds_dummy
bounded_variable = d
fixed_bound_value = 0
threshold_value = 0.95
[]
[upper]
type = ConstantBounds
variable = bounds_dummy
bounded_variable = d
bound_type = upper
bound_value = 1
[]
[]
[AuxKernels]
[get_ce]
type = ADMaterialRealAux
variable = ce
property = ce
[]
[]
[Kernels]
[diff]
type = ADPFFDiffusion
variable = d
fracture_toughness = Gc
regularization_length = l
normalization_constant = c0
[]
[source]
type = ADPFFSource
variable = d
free_energy = psi
[]
[nuc_force]
type = ADCoefMatSource
variable = d
prop_names = 'ce'
[]
[]
[Materials]
[fracture_properties]
type = ADGenericConstantMaterial
prop_names = 'E K G lambda Gc l sigma_ts sigma_cs delta'
prop_values = '${E} ${K} ${G} ${Lambda} ${Gc} ${l} ${sigma_ts} ${sigma_cs} ${delta}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 0'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[psi]
type = ADDerivativeParsedMaterial
property_name = psi
expression = 'g*psie_active+(Gc/c0/l)*alpha'
coupled_variables = 'd psie_active'
material_property_names = 'alpha(d) g(d) Gc c0 l'
derivative_order = 1
[]
[kumar_material]
type = KLRNucleationMicroForce
phase_field = d
normalization_constant = c0
tensile_strength = sigma_ts
compressive_strength = sigma_cs
delta = delta
external_driving_force_name = ce
output_properties = 'ce'
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -snes_type'
petsc_options_value = 'lu superlu_dist vinewtonrsls'
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
[]
[Outputs]
print_linear_residuals = false
[]
(tutorials/soil_desiccation/elasticity.i)
E = 4
nu = 0.2
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
l = 0.5
c = 0.1
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = 'fracture.i'
cli_args = 'l=${l}'
[]
[]
[Transfers]
[send_psie_active]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = fracture
source_variable = psie_active
variable = psie_active
[]
[send_psii_active]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = fracture
source_variable = psii_active
variable = psii_active
[]
[get_d]
type = MultiAppGeneralFieldShapeEvaluationTransfer
from_multi_app = fracture
source_variable = d
variable = d
[]
[]
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'gold/fields_refined.e'
use_for_exodus_restart = true
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[]
[AuxVariables]
[d]
[]
[Gc]
initial_from_file_var = Gc
[]
[psic]
initial_from_file_var = psic
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[plane_stress]
type = ADWeakPlaneStress
variable = strain_zz
[]
[react_x]
type = ADCoefMatReaction
variable = disp_x
coefficient = ${c}
prop_names = 'g'
[]
[react_y]
type = ADCoefMatReaction
variable = disp_y
coefficient = ${c}
prop_names = 'g'
[]
[]
[BCs]
[Periodic]
[x]
variable = disp_x
auto_direction = 'x y'
[]
[y]
variable = disp_y
auto_direction = 'x y'
[]
[]
[]
[Materials]
# fracture
[fracture_toughness]
type = ADParsedMaterial
property_name = Gc
coupled_variables = Gc
expression = 'Gc'
[]
[critial_fracture_energy]
type = ADParsedMaterial
property_name = psic
coupled_variables = psic
expression = 'psic'
[]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l'
prop_values = '${K} ${G} ${l}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
parameter_names = 'p a2 a3 eta '
parameter_values = '2 1 0 1e-6'
material_property_names = 'Gc psic xi c0 l '
[]
# elasticity
[eigen_strain]
type = ComputeEigenstrainFromFunctionInitialStress
initial_stress = 't t 0'
eigenstrain_name = 'es'
bulk_modulus = K
shear_modulus = G
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = strain_zz
eigenstrain_names = 'es'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
[]
[interfacial_energy]
type = ThinFilmInterfaceEnergyDensity
phase_field = d
shear_lag_coef = ${c}
output_properties = 'psii_active'
outputs = exodus
[]
[]
[Postprocessors]
[delta_d]
type = SolutionChangeNorm
variable = d
outputs = none
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options = '-ksp_converged_reason'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package -ksp_type'
petsc_options_value = 'lu superlu_dist gmres'
dt = 0.0005
end_time = 0.22
nl_rel_tol = 1e-06
nl_abs_tol = 1e-08
automatic_scaling = true
fixed_point_max_its = 10
custom_pp = delta_d
custom_abs_tol = 1e-03
custom_rel_tol = 1e-03
disable_fixed_point_residual_norm_check = true
accept_on_max_fixed_point_iteration = true
reuse_preconditioner = true
reuse_preconditioner_max_linear_its = 50
[]
[Outputs]
print_linear_residuals = false
exodus = true
[]
(tutorials/small_deformation/elastoplasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
sigma_y = 50
n = 5
ep0 = 0.001
[GlobalParams]
displacements = 'disp_x disp_y'
volumetric_locking_correction = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymax = 0.5
[]
[noncrack]
type = BoundingBoxNodeSetGenerator
input = gen
new_boundary = noncrack
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
construct_side_list_from_node_list = true
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fy]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
save_in = fy
[]
[]
[BCs]
[ydisp]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = noncrack
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = top
value = 0
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'K G'
prop_values = '${K} ${G}'
[]
[no_degradation]
type = NoDegradation
property_name = g
expression = 1
phase_field = d
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
output_properties = 'elastic_strain'
outputs = exodus
[]
[hardening]
type = PowerLawHardening
yield_stress = ${sigma_y}
exponent = ${n}
reference_plastic_strain = ${ep0}
phase_field = d
degradation_function = g
[]
[plasticity]
type = SmallDeformationJ2Plasticity
hardening_model = hardening
output_properties = 'effective_plastic_strain plastic_strain'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
plasticity_model = plasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fy]
type = NodalSum
variable = fy
boundary = top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1e-5
end_time = 1e-3
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/mode1_ductile_fracture/elastoplasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
psic = 14.88
l = 0.02
sigma_y = 1000
n = 5
ep0 = 0.005
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};psic=${psic};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppGeneralFieldShapeEvaluationTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[to_psip_active]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = 'fracture'
variable = psip_active
source_variable = psip_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
volumetric_locking_correction = true
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymax = 0.5
[]
[noncrack]
type = BoundingBoxNodeSetGenerator
input = gen
new_boundary = noncrack
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = OrientedBoxMarker
center = '0.75 0.25 0'
length = 0.8
width = 0.1
height = 1
length_direction = '1 1 0'
width_direction = '-1 1 0'
inside = REFINE
outside = DO_NOTHING
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fy]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
save_in = fy
[]
[]
[BCs]
[ydisp]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = noncrack
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = top
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K} ${G} ${l} ${Gc} ${psic}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'elastic_strain psie_active'
outputs = exodus
[]
[plasticity]
type = SmallDeformationJ2Plasticity
hardening_model = power_law_hardening
output_properties = 'effective_plastic_strain'
outputs = exodus
[]
[power_law_hardening]
type = PowerLawHardening
degradation_function = g
yield_stress = ${sigma_y}
exponent = ${n}
reference_plastic_strain = ${ep0}
phase_field = d
output_properties = 'psip_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
plasticity_model = plasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fy]
type = NodalSum
variable = fy
boundary = top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 1e-2
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(test/tests/nucleation_models/elasticity.i)
# This tutorial is based on example of material graphite in section 4.3.4
# See Kumar et. al. https://doi.org/10.1016/j.jmps.2020.104027.
# For regularization length 0.35mm, delta 1.16, mesh size 0.035,
# J intergal over Gc should return a value close to 1
# delta value 0.2 is updated for model 2022
# (see Kumar et. al. https://doi.org/10.1007/s10704-022-00653-z)
E = 9.8e3 #9.8Gpa
nu = 0.13
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Lambda = '${fparse E*nu/(1+nu)/(1-2*nu)}'
Gc = 9.1e-2 # 91N/m
l = 0.35
sigma_ts = 27
sigma_cs = 77
# for model 2022 (KLR), select delta=0.2
# for model 2020 (KLBF), select delta=1.16
# delta = 0.2 #1.16
delta = 1.16
c1 = '${fparse (1+nu)*sqrt(Gc)/sqrt(2*pi*E)}'
c2 = '${fparse (3-nu)/(1+nu)}'
ahead = 0
V = 20
# mesh size 0.035
nx = 30
ny = 10
refine = 1
[Functions]
[bc_func]
type = ParsedFunction
expression = c1*((x-V*t-ahead)^2+y^2)^(0.25)*(c2-cos(atan2(y,(x-V*t-ahead))))*sin(0.5*atan2(y,(x-V*t-ahead)))
symbol_names = 'c1 c2 V ahead'
symbol_values = '${c1} ${c2} ${V} ${ahead}'
[]
[]
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = klr.i
cli_args = 'E=${E};K=${K};G=${G};Lambda=${Lambda};Gc=${Gc};l=${l};nx=${nx};ny=${ny};refine=${refine};sigma_ts=${sigma_ts};sigma_cs=${sigma_cs};delta=${delta}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = fracture
variable = 'd ce'
source_variable = 'd ce'
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = fracture
variable = 'disp_x disp_y strain_zz psie_active'
source_variable = 'disp_x disp_y strain_zz psie_active'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = ${nx}
ny = ${ny}
xmax = 30
ymin = -5
ymax = 5
[]
[gen2]
type = ExtraNodesetGenerator
input = gen
new_boundary = fix_point
coord = '0 -5'
[]
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = ${refine}
stop_time = 0
max_h_level = ${refine}
[Markers]
[marker]
type = BoxMarker
bottom_left = '0 -0.7 0'
top_right = '30 0.7 0'
outside = DO_NOTHING
inside = REFINE
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[strain_zz]
[]
[]
[AuxVariables]
[d]
[InitialCondition]
type = FunctionIC
function = 'if(y=0&x>=0&x<=5,1,0)'
[]
[]
[ce]
order = CONSTANT
family = MONOMIAL
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[plane_stress]
type = ADWeakPlaneStress
variable = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[]
[BCs]
[fix_x]
type = DirichletBC
variable = disp_x
boundary = fix_point
value = 0
[]
[bottom_y]
type = FunctionDirichletBC
variable = disp_y
boundary = bottom
function = bc_func
[]
[top_y]
type = FunctionDirichletBC
variable = disp_y
boundary = top
function = bc_func
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'E K G lambda Gc l'
prop_values = '${E} ${K} ${G} ${Lambda} ${Gc} ${l}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 0'
[]
[strain]
type = ADComputePlaneSmallStrain
out_of_plane_strain = 'strain_zz'
displacements = 'disp_x disp_y'
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'psie_active psie'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
[]
[]
[Postprocessors]
[Jint]
type = PhaseFieldJIntegral
J_direction = '1 0 0'
strain_energy_density = psie
displacements = 'disp_x disp_y'
boundary = 'left top right bottom'
[]
[Jint_over_Gc]
type = ScalePostprocessor
value = 'Jint'
scaling_factor = '${fparse 1.0/Gc}'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 1e-2
dtmin = 1e-2
end_time = 0.3
num_steps = 2
fixed_point_max_its = 500
accept_on_max_fixed_point_iteration = false
fixed_point_rel_tol = 1e-6
fixed_point_abs_tol = 1e-8
[]
[Outputs]
file_base = klr_out
exodus = true
print_linear_residuals = false
[]
(tutorials/matrix_fiber/elasticity.i)
E = 4000
nu = 0.2
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
lambda = '${fparse K-2*G/3}'
Gc = 1
psic = 2.64
l = 0.02
k = 1e-6
v = '${fparse sqrt(Gc*3/lambda)}'
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = 'fracture.i'
cli_args = 'Gc=${Gc};l=${l};psic=${psic};k=${k}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
source_variable = d
variable = d
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[Mesh]
[fmg]
type = FileMeshGenerator
file = 'gold/matrix.msh'
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[]
[BCs]
[forcing]
type = FunctionDirichletBC
variable = disp_y
boundary = 'top'
function = '${v}*t'
preset = false
[]
[fixed_hole_x]
type = DirichletBC
variable = disp_x
boundary = 'hole'
value = 0
[]
[fixed_hole_y]
type = DirichletBC
variable = disp_y
boundary = 'hole'
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K} ${G} ${l} ${Gc} ${psic}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = d
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
phase_field = d
parameter_names = 'p a2 a3 eta '
parameter_values = '2 1 0 ${k}'
material_property_names = 'Gc psic xi c0 l '
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = VOLDEV
output_properties = 'psie_active'
outputs = exodus
[]
[strain]
type = ADComputeSmallStrain
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist'
dt = 0.01
end_time = 3
nl_abs_tol = 1e-08
nl_rel_tol = 1e-06
automatic_scaling = true
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/traction_separation/elasticity.i)
E = 3e4
nu = 0.2
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
p = 0.5
Gc = 0.12
sigmac = 3
psic = '${fparse sigmac^2/2/E}'
l = 20
R = 10
h = '${fparse l/R}'
nx = '${fparse ceil(200/h)}'
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'p=${p};Gc=${Gc};psic=${psic};l=${l};nx=${nx}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_fracture]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_fracture]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = 'psie_active disp_x disp_y'
source_variable = 'psie_active disp_x disp_y'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmax = 200
ymax = 1
nx = ${nx}
ny = 1
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fx]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
save_in = fx
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[pressure_x]
type = ADPressurizedCrack
variable = disp_x
pressure = p
phase_field = d
indicator_function = I
component = 0
[]
[pressure_y]
type = ADPressurizedCrack
variable = disp_y
pressure = p
phase_field = d
indicator_function = I
component = 1
[]
[]
[BCs]
[xdisp]
type = FunctionDirichletBC
variable = disp_x
boundary = right
function = 't'
[]
[yfix]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0
[]
[xfix]
type = DirichletBC
variable = disp_x
boundary = left
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic p'
prop_values = '${K} ${G} ${l} ${Gc} ${psic} ${p}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = '2*d-d^2'
phase_field = d
[]
[indicator]
type = ADDerivativeParsedMaterial
property_name = I
coupled_variables = 'd'
expression = 'd^2'
derivative_order = 1
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = NONE
output_properties = 'psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
[]
[pressure_density]
type = PressureDensity
effective_pressure_density = p_density
phase_field = d
pressure = p
indicator_function = I
[]
[]
[Postprocessors]
[Fx_left]
type = NodalSum
variable = fx
boundary = left
[]
[Fx_right]
type = NodalSum
variable = fx
boundary = right
[]
[effective_pressure]
type = ADElementIntegralMaterialProperty
mat_prop = p_density
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 5e-4
end_time = 5e-2
fixed_point_max_its = 200
accept_on_max_fixed_point_iteration = false
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
file_base = 'p_${p}_l_${l}'
exodus = true
csv = true
print_linear_residuals = false
[]
(tutorials/adaptivity/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
psic = 14.88
l = 0.01
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};psic=${psic};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppGeneralFieldShapeEvaluationTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie]
type = MultiAppGeneralFieldShapeEvaluationTransfer
to_multi_app = 'fracture'
variable = 'psie_active psie'
source_variable = 'psie_active psie'
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[top_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = 0
ymax = 0.5
boundary_id_offset = 0
boundary_name_prefix = top_half
[]
[top_stitch]
type = BoundingBoxNodeSetGenerator
input = top_half
new_boundary = top_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[bottom_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = -0.5
ymax = 0
boundary_id_offset = 5
boundary_name_prefix = bottom_half
[]
[bottom_stitch]
type = BoundingBoxNodeSetGenerator
input = bottom_half
new_boundary = bottom_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[stitch]
type = StitchedMeshGenerator
inputs = 'top_stitch bottom_stitch'
stitch_boundaries_pairs = 'top_stitch bottom_stitch'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = combo
max_h_level = 3
[Markers]
[combo]
type = ComboMarker
markers = 'damage_marker strain_energy_marker'
[]
[damage_marker]
type = ValueThresholdMarker
variable = d
refine = 0.5
[]
[strain_energy_marker]
type = ValueThresholdMarker
variable = psie
refine = '${fparse 0.8*psic}'
coarsen = '${fparse 0.8*psic}'
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fx]
[]
[d]
[]
[psie_active]
order = CONSTANT
family = MONOMIAL
[]
[psie]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[psie_active]
type = ADMaterialRealAux
variable = psie_active
property = psie_active
[]
[psie]
type = ADMaterialRealAux
variable = psie
property = psie
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
save_in = fx
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[]
[BCs]
[xdisp]
type = FunctionDirichletBC
variable = 'disp_x'
boundary = 'top_half_top'
function = 't'
preset = false
[]
[yfix]
type = DirichletBC
variable = 'disp_y'
boundary = 'top_half_top bottom_half_bottom top_half_left top_half_right bottom_half_left bottom_half_right'
value = 0
[]
[xfix]
type = DirichletBC
variable = 'disp_x'
boundary = 'bottom_half_bottom'
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K} ${G} ${l} ${Gc} ${psic}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = SPECTRAL
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
[]
[]
[Postprocessors]
[Fx]
type = NodalSum
variable = fx
boundary = top_half_top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 2e-2
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/mode2_cohesive_fracture/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
psic = 14.88
l = 0.02
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};psic=${psic};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[top_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = 0
ymax = 0.5
boundary_id_offset = 0
boundary_name_prefix = top_half
[]
[top_stitch]
type = BoundingBoxNodeSetGenerator
input = top_half
new_boundary = top_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[bottom_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = -0.5
ymax = 0
boundary_id_offset = 5
boundary_name_prefix = bottom_half
[]
[bottom_stitch]
type = BoundingBoxNodeSetGenerator
input = bottom_half
new_boundary = bottom_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[stitch]
type = StitchedMeshGenerator
inputs = 'top_stitch bottom_stitch'
stitch_boundaries_pairs = 'top_stitch bottom_stitch'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = OrientedBoxMarker
center = '0.65 -0.25 0'
length = 0.8
width = 0.2
height = 1
length_direction = '1 -1.5 0'
width_direction = '1.5 1 0'
inside = REFINE
outside = DO_NOTHING
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fx]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
save_in = fx
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[]
[BCs]
[xdisp]
type = FunctionDirichletBC
variable = 'disp_x'
boundary = 'top_half_top'
function = 't'
preset = false
[]
[yfix]
type = DirichletBC
variable = 'disp_y'
boundary = 'top_half_top bottom_half_bottom top_half_left top_half_right bottom_half_left bottom_half_right'
value = 0
[]
[xfix]
type = DirichletBC
variable = 'disp_x'
boundary = 'bottom_half_bottom'
value = 0
[]
[]
[Materials]
[bulk_properties]
type = ADGenericConstantMaterial
prop_names = 'K G l Gc psic'
prop_values = '${K} ${G} ${l} ${Gc} ${psic}'
[]
[crack_geometric]
type = CrackGeometricFunction
property_name = alpha
expression = 'd'
phase_field = d
[]
[degradation]
type = RationalDegradationFunction
property_name = g
expression = (1-d)^p/((1-d)^p+(Gc/psic*xi/c0/l)*d*(1+a2*d+a2*a3*d^2))*(1-eta)+eta
phase_field = d
material_property_names = 'Gc psic xi c0 l '
parameter_names = 'p a2 a3 eta '
parameter_values = '2 -0.5 0 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = SPECTRAL
output_properties = 'elastic_strain psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fx]
type = NodalSum
variable = fx
boundary = top_half_top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 2e-2
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]
(tutorials/mode2_brittle_fracture/elasticity.i)
E = 2.1e5
nu = 0.3
K = '${fparse E/3/(1-2*nu)}'
G = '${fparse E/2/(1+nu)}'
Gc = 2.7
l = 0.02
[MultiApps]
[fracture]
type = TransientMultiApp
input_files = fracture.i
cli_args = 'Gc=${Gc};l=${l}'
execute_on = 'TIMESTEP_END'
[]
[]
[Transfers]
[from_d]
type = MultiAppCopyTransfer
from_multi_app = 'fracture'
variable = d
source_variable = d
[]
[to_psie_active]
type = MultiAppCopyTransfer
to_multi_app = 'fracture'
variable = psie_active
source_variable = psie_active
[]
[]
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
[top_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = 0
ymax = 0.5
boundary_id_offset = 0
boundary_name_prefix = top_half
[]
[top_stitch]
type = BoundingBoxNodeSetGenerator
input = top_half
new_boundary = top_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[bottom_half]
type = GeneratedMeshGenerator
dim = 2
nx = 30
ny = 15
ymin = -0.5
ymax = 0
boundary_id_offset = 5
boundary_name_prefix = bottom_half
[]
[bottom_stitch]
type = BoundingBoxNodeSetGenerator
input = bottom_half
new_boundary = bottom_stitch
bottom_left = '0.5 0 0'
top_right = '1 0 0'
[]
[stitch]
type = StitchedMeshGenerator
inputs = 'top_stitch bottom_stitch'
stitch_boundaries_pairs = 'top_stitch bottom_stitch'
[]
construct_side_list_from_node_list = true
[]
[Adaptivity]
marker = marker
initial_marker = marker
initial_steps = 2
stop_time = 0
max_h_level = 2
[Markers]
[marker]
type = OrientedBoxMarker
center = '0.65 -0.25 0'
length = 0.8
width = 0.2
height = 1
length_direction = '1 -1.5 0'
width_direction = '1.5 1 0'
inside = REFINE
outside = DO_NOTHING
[]
[]
[]
[Variables]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[fx]
[]
[d]
[]
[]
[Kernels]
[solid_x]
type = ADStressDivergenceTensors
variable = disp_x
component = 0
save_in = fx
[]
[solid_y]
type = ADStressDivergenceTensors
variable = disp_y
component = 1
[]
[]
[BCs]
[xdisp]
type = FunctionDirichletBC
variable = 'disp_x'
boundary = 'top_half_top'
function = 't'
preset = false
[]
[yfix]
type = DirichletBC
variable = 'disp_y'
boundary = 'top_half_top bottom_half_bottom'
value = 0
[]
[xfix]
type = DirichletBC
variable = 'disp_x'
boundary = 'bottom_half_bottom'
value = 0
[]
[]
[Materials]
[bulk]
type = ADGenericConstantMaterial
prop_names = 'K G'
prop_values = '${K} ${G}'
[]
[degradation]
type = PowerDegradationFunction
property_name = g
expression = (1-d)^p*(1-eta)+eta
phase_field = d
parameter_names = 'p eta '
parameter_values = '2 1e-6'
[]
[strain]
type = ADComputeSmallStrain
[]
[elasticity]
type = SmallDeformationIsotropicElasticity
bulk_modulus = K
shear_modulus = G
phase_field = d
degradation_function = g
decomposition = SPECTRAL
output_properties = 'elastic_strain psie_active'
outputs = exodus
[]
[stress]
type = ComputeSmallDeformationStress
elasticity_model = elasticity
output_properties = 'stress'
outputs = exodus
[]
[]
[Postprocessors]
[Fx]
type = NodalSum
variable = fx
boundary = top_half_top
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_factor_mat_solver_package'
petsc_options_value = 'lu superlu_dist '
automatic_scaling = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-10
dt = 2e-5
end_time = 2e-2
fixed_point_max_its = 20
accept_on_max_fixed_point_iteration = true
fixed_point_rel_tol = 1e-8
fixed_point_abs_tol = 1e-10
[]
[Outputs]
exodus = true
print_linear_residuals = false
[]