- grain_sizeValue of the crystal grain size, in mm
C++ Type:double
Controllable:No
Description:Value of the crystal grain size, in mm
- initial_forest_dislocation_densityThe initial density of the forest dislocations, in 1/mm^2, assumed to be split evenly among all slip systems
C++ Type:double
Controllable:No
Description:The initial density of the forest dislocations, in 1/mm^2, assumed to be split evenly among all slip systems
- initial_substructure_densityThe initial total density of the sessile dislocations, in 1/mm^2
C++ Type:double
Controllable:No
Description:The initial total density of the sessile dislocations, in 1/mm^2
- number_slip_systemsThe total number of possible active slip systems for the crystalline material
C++ Type:unsigned int
Controllable:No
Description:The total number of possible active slip systems for the crystalline material
- slip_sys_file_nameName of the file containing the slip systems, one slip system per row, with the slip plane normal given before the slip plane direction.
C++ Type:FileName
Controllable:No
Description:Name of the file containing the slip systems, one slip system per row, with the slip plane normal given before the slip plane direction.
Crystal Plasticity HCP Dislocation Slip Beyerlein Update
Two-term dislocation slip model for hexagonal close packed crystals from Beyerline and Tome
Description
CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
is designed to be used in conjunction with the ComputeMultipleCrystalPlasticityStress class to calculate the response of a FCC crystalline solid. Details about the algorithm and specific stress and strain measures used in the CrystalPlasticityUpdate
base class are given on the documentation page for ComputeMultipleCrystalPlasticityStress.
As in other crystal plasticity consitutive models compatible with ComputeMultipleCrystalPlasticityStress, the user must supply the slip plane normal and slip direction information, and the lengths of the unit cell lattice parameters. Unlike other crystal plasticity consistutive models, CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
is set for use with the HCP lattice type only. Because of this lattice type requirement, four entries each are expected for both the slip plane normal and the slip direction vector entries.
Constitutive Model Definition
A constitutive model for only the glide and evolution of forest dislocations within a Hexagonal Close-Packed (HCP) crystal lattice is implemented in this class. The constitutive model is taken from the work of Beyerlein and Tomé (2008), Capolungo et al. (2009), and Beyerlein and Tomé (2010), which was developed for a viscoplastic self-consistent application. Here the constitutive model is employed within a crystal plasticity implementation.
CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
allows for the evaluation of the slip resistance, forest dislocation density, and substructure density on different types of slip systems, such as basal, prismatic, pyramidal, and so on. Within this class these different types of slip systems are termed modes
, following Beyerlein and Tomé (2008). Vectors of values for each slip mode dependent parameter may be supplied within the input file. The number of slip systems, value of the burgers vector, and the value of the initial lattice friction are examples of slip mode dependent input file parameters. The specification of modes are not explicityly identified within the code, and no particular order is enforced.
Users are responsible for ensuring that the order of mode-specific parameters is consistent throughout the input file, and that the order of slip system modes corresponds correctly to the order of slip plane normals and slip directions given in the slip system text file, slip_sys_file_name
.
Slip System Resistance Calculation
The resistance to dislocation glide is represented in this constitutive model as the strength value that the applied shear stress must overcome to produce plastic slip. The slip resistance is calculated individually for each of the user-supplied slip systems. The resistance is considered here as the sum of four terms: (1) where g is a user-defined constant initial lattice friction, g represents a Hall-Petch type treatment of the slip resistance dependence on grain size, and g and g represent the hardening contributions from the two dislocation density populations, forest and substructure. (Beyerlein and Tomé, 2008).The grain-size dependence hardening term employs a Hall-Petch type dependence on grain size, following Beyerlein and Tomé (2008): (2) where HP is a user-defined constant, is the shear modulus for the slip system, b is the slip system Burgers vector, and d is the average grain size. This value is computed once during the initalization stages of the simulation.
The hardening contribution from forest dislocations, see forest dislocation evolution discussion, is calculated on a per slip basis(Capolungo et al., 2009). (3) where is a user-defined interaction coefficient and is the density of forest dislocations on the slip system . This calculation is completed at every timestep to capture the dislocation density evolution.
Similarly, the hardening contribution due to the substructure density is computed once per timestep; however, only a single total value of substructure density from all slip systems is tracked, as described in the substructure evolution section. Characteristics of the slip systems, such as the Burgers vector, do affect the hardening contribution from substructure density on each slip system: (4) where k is a user-defined interaction coefficient and is the substructure density.
Slip Rate Calculation
The slip rate due to dislocation glide is calculated with a power law relation: (5) where is the reference strain rate, is the resolved applied shear stress on the slip system , m is the slip rate exponent, and t is the current timestep size. As in the work of Beyerlein and Tomé (2008) a low slip rate exponent is recommended here to produce slip only when the applied shear stress is near the value fo the slip system resistance, Eq. (1). Note that lower values of the strain rate exponent can cause strain rate sensitivity; following Capolungo et al. (2009) and Beyerlein and Tomé (2010), CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
sets the reference strain rate, equal to the macroscopic strain rate through a user-supplied parameter. To reduce potential issues arising from very low slip increments, the requirement is implemented in this class, where the tolerance is a user-defined parameter. For slip systems which do not meet this threshold requirement, the slip increment value is set to zero for that timestep.
Forest Dislocation Evolution
The forest dislocations evolve in response to the resolved applied shear stress, which provides the driving force for dicloscation glide. The forest dislocation density evolution model is a direct function of the slip increment (Beyerlein and Tomé, 2008): where k is the slip mode dependent forest dislocation generation coefficient and is the removed forest dislocation density increment for the current timestep. These forest dislocations are either annihilated or generated substructure debris. The evolution of removed forest dislocations has the form (Beyerlein and Tomé, 2008): (6) where G is a slip-mode specific normalized activation energy, k is the Boltzman constant, T is the temperature, D is the slip mode dependent propotionality factor, is the applied strain rate, bounded within the range 10 to 10 1/s, and is a material value estimated to be 10 (Beyerlein and Tomé, 2008).
The final forest dislocation density per slip system is updated at the end of each time step as where n is the current timestep and n-1 is the previous timestep, provided the increment of the forest dislocation density is within the user-set tolerance. If the forest density increment is out of tolerance, CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
triggers the use of a smaller timestep or substep.
Substructure Dislocation Evolution
The substructure density is calculated as a function of the removed forest dislocation density, Eq. (6), following Capolungo et al. (2009): (7) where f is the slip-mode dependent substructure generation rate coefficient. The substructure density is considered to affect the forest dislocation evolution on all slip systesms eqvenly; thus, only the total sum of substructure density across all slip systems is modeled here. As noted by (Capolungo et al., 2009), the constitutive model does nto allow for recovery of dislocations from the substructure density.
Example Input File Syntax
[trial_xtalpl]
type = CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
number_slip_systems = 15
slip_sys_file_name = hcp_aprismatic_capyramidal_slip_sys.txt
unit_cell_dimension = '2.934e-7 2.934e-7 4.657e-7' #Ti, in mm, https://materialsproject.org/materials/mp-46/
temperature = temperature
initial_forest_dislocation_density = 15.0e5
initial_substructure_density = 1.0e3
slip_system_modes = 2
number_slip_systems_per_mode = '3 12'
lattice_friction_per_mode = '98 224' #Knezevic et al MSEA 654 (2013)
effective_shear_modulus_per_mode = '4.7e4 4.7e4' #Ti, in MPa, https://materialsproject.org/materials/mp-46/
burgers_vector_per_mode = '2.934e-7 6.586e-7' #Ti, in mm, https://materialsproject.org/materials/mp-46/
slip_generation_coefficient_per_mode = '1.25e5 2.25e7' #from Beyerlein and Tome 2008 IJP
normalized_slip_activiation_energy_per_mode = '3.73e-3 3.2e-2' #from Beyerlein and Tome 2008 IJP
slip_energy_proportionality_factor_per_mode = '330 100' #from Beyerlein and Tome 2008 IJP
substructure_rate_coefficient_per_mode = '355 0.4' #from Capolungo et al MSEA (2009)
applied_strain_rate = 0.001
gamma_o = 1.0e-3
Hall_Petch_like_constant_per_mode = '0.2 0.2' #Estimated to match graph in Capolungo et al MSEA (2009), Figure 2
grain_size = 20.0e-3 #20 microns, Beyerlein and Tome IJP (2008)
[]
(moose/modules/tensor_mechanics/test/tests/crystal_plasticity/hcp_single_crystal/update_method_hcp_aprismatic_capyramidal.i)CrystalPlasticityHCPDislocationSlipBeyerleinUpdate
must be run in conjunction with the crystal plasticity specific stress calculator as shown below:
[stress]
type = ComputeMultipleCrystalPlasticityStress
crystal_plasticity_models = 'trial_xtalpl'
tan_mod_type = exact
[]
(moose/modules/tensor_mechanics/test/tests/crystal_plasticity/hcp_single_crystal/update_method_hcp_aprismatic_capyramidal.i)Input Parameters
- Boltzman_constant1.38065e-20Boltzman constant, in MPa-mm^3/K
Default:1.38065e-20
C++ Type:double
Controllable:No
Description:Boltzman constant, in MPa-mm^3/K
- Hall_Petch_like_constant_per_modeThe microstructure Hall-Petch like coefficient value used to capture the influence of grain size on slip system resistance in the absence of twin dislocations, dimensionless
C++ Type:std::vector<double>
Controllable:No
Description:The microstructure Hall-Petch like coefficient value used to capture the influence of grain size on slip system resistance in the absence of twin dislocations, dimensionless
- applied_strain_rate0.0001Value of the applied macroscopic strain rate and should correspond to the simulation loading conditions, in 1/s.
Default:0.0001
C++ Type:double
Controllable:No
Description:Value of the applied macroscopic strain rate and should correspond to the simulation loading conditions, in 1/s.
- base_nameOptional parameter that allows the user to define multiple crystal plasticity mechanisms
C++ Type:std::string
Controllable:No
Description:Optional parameter that allows the user to define multiple crystal plasticity mechanisms
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
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>
Controllable:No
Description:The list of boundaries (ids or names) from the mesh where this object applies
- burgers_vector_per_modeValue of the Burgers vector, b, for each type of the slip system, units of mm. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Value of the Burgers vector, b, for each type of the slip system, units of mm. The order must be consistent with the number of slip systems per type vector.
- 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
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
Controllable:No
Description:An optional suffix parameter that can be appended to any declared properties. The suffix will be prepended with a '_' character.
- effective_shear_modulus_per_modeEffective isotropic shear modulus value, mu, in MPa. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Effective isotropic shear modulus value, mu, in MPa. The order must be consistent with the number of slip systems per type vector.
- forest_interaction_parameter0.9Forest dislocation interaction parameter, Chi, dimensionless.
Default:0.9
C++ Type:double
Controllable:No
Description:Forest dislocation interaction parameter, Chi, dimensionless.
- gamma_o0.001Reference strain rate on each slip system, in 1/s
Default:0.001
C++ Type:double
Controllable:No
Description:Reference strain rate on each slip system, in 1/s
- lattice_friction_per_modeValue of the lattice friction for each type of the slip system, units of MPa. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Value of the lattice friction for each type of the slip system, units of MPa. The order must be consistent with the number of slip systems per type vector.
- normalized_slip_activiation_energy_per_modeValue of the slip dislocation attraction activation energy for each type of the slip system, g, dimensionless. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Value of the slip dislocation attraction activation energy for each type of the slip system, g, dimensionless. The order must be consistent with the number of slip systems per type vector.
- number_cross_slip_directions0Quanity of unique slip directions, used to determine cross slip familes
Default:0
C++ Type:double
Controllable:No
Description:Quanity of unique slip directions, used to determine cross slip familes
- number_cross_slip_planes0Quanity of slip planes belonging to a single cross slip direction; used to determine cross slip families
Default:0
C++ Type:double
Controllable:No
Description:Quanity of slip planes belonging to a single cross slip direction; used to determine cross slip families
- number_slip_systems_per_modeThe number of slip systems per each slip system type. The sum of the entries of the vector given here must equal the value given for the total number of slip systems.
C++ Type:std::vector<unsigned int>
Controllable:No
Description:The number of slip systems per each slip system type. The sum of the entries of the vector given here must equal the value given for the total number of slip systems.
- print_state_variable_convergence_error_messagesFalseWhether or not to print warning messages from the crystal plasticity specific convergence checks on both the constiutive model internal state variables.
Default:False
C++ Type:bool
Controllable:No
Description:Whether or not to print warning messages from the crystal plasticity specific convergence checks on both the constiutive model internal state variables.
- 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
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.
- reference_macroscopic_strain_rate1e+07Value of the reference macroscopic strain rate for the thermal dislocation attraction, in 1/s.
Default:1e+07
C++ Type:double
Controllable:No
Description:Value of the reference macroscopic strain rate for the thermal dislocation attraction, in 1/s.
- resistance_tol0.01Constitutive slip system resistance relative residual tolerance for each individual constitutive model
Default:0.01
C++ Type:double
Controllable:No
Description:Constitutive slip system resistance relative residual tolerance for each individual constitutive model
- slip_energy_proportionality_factor_per_modeValue of the the dislocation slip attraction energy proportionality factor for each type of the slip system, D, units of MPa. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Value of the the dislocation slip attraction energy proportionality factor for each type of the slip system, D, units of MPa. The order must be consistent with the number of slip systems per type vector.
- slip_generation_coefficient_per_modeSlip dislocation generation coefficient value for each type of the slip system, k_1, units of 1/mm. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Slip dislocation generation coefficient value for each type of the slip system, k_1, units of 1/mm. The order must be consistent with the number of slip systems per type vector.
- slip_increment_tolerance0.02Maximum allowable slip in an increment for each individual constitutive model
Default:0.02
C++ Type:double
Controllable:No
Description:Maximum allowable slip in an increment for each individual constitutive model
- slip_system_modes1Number of different types of slip systems in this HCP crystal, e.g. for a material with basal, prismatic, and pyramidal active slip systems, this number would be 3
Default:1
C++ Type:unsigned int
Controllable:No
Description:Number of different types of slip systems in this HCP crystal, e.g. for a material with basal, prismatic, and pyramidal active slip systems, this number would be 3
- stol0.01Constitutive internal state variable relative change tolerance
Default:0.01
C++ Type:double
Controllable:No
Description:Constitutive internal state variable relative change tolerance
- strain_rate_sensitivity_exponent0.05The strain rate sensitivity exponent for the power law relationship of resolved shear stress
Default:0.05
C++ Type:double
Controllable:No
Description:The strain rate sensitivity exponent for the power law relationship of resolved shear stress
- substructure_hardening_coefficient0.086Value of the coefficient for the expanded Taylor hardening substructure hardening relation, set to recover the Taylor hardening law for low substructure densities, k_{sub}, dimensionless.
Default:0.086
C++ Type:double
Controllable:No
Description:Value of the coefficient for the expanded Taylor hardening substructure hardening relation, set to recover the Taylor hardening law for low substructure densities, k_{sub}, dimensionless.
- substructure_rate_coefficient_per_modeMaterial-independent rate constant that accounts for locking of slip dislocations in sessile substructure dislocation segments, q, dimensionless. This value is often determined through dislocation dynamics calculations. The order must be consistent with the number of slip systems per type vector.
C++ Type:std::vector<double>
Controllable:No
Description:Material-independent rate constant that accounts for locking of slip dislocations in sessile substructure dislocation segments, q, dimensionless. This value is often determined through dislocation dynamics calculations. The order must be consistent with the number of slip systems per type vector.
- temperatureThe name of the temperature variable
C++ Type:std::vector<VariableName>
Controllable:No
Description:The name of the temperature variable
- total_twin_volume_fractionTotal twin volume fraction, if twinning is considered in the simulation
C++ Type:MaterialPropertyName
Controllable:No
Description:Total twin volume fraction, if twinning is considered in the simulation
- unit_cell_dimension1 1 1 The dimension of the unit cell along three directions, where a cubic unit cell is assumed for cubic crystals and a hexagonal unit cell (a, a, c) is assumed for HCP crystals. These dimensions will be taken into account while computing the slip systems. Default size is 1.0 along all three directions.
Default:1 1 1
C++ Type:std::vector<double>
Controllable:No
Description:The dimension of the unit cell along three directions, where a cubic unit cell is assumed for cubic crystals and a hexagonal unit cell (a, a, c) is assumed for HCP crystals. These dimensions will be taken into account while computing the slip systems. Default size is 1.0 along all three directions.
- zero_tol1e-12Tolerance for residual check when variable value is zero for each individual constitutive model
Default:1e-12
C++ Type:double
Controllable:No
Description:Tolerance for residual check when variable value is zero for each individual constitutive model
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
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
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
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
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
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>
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>
Controllable:No
Description:Vector of output names where you would like to restrict the output of variables(s) associated with this object
Outputs Parameters
References
- I.J. Beyerlein and C.N. Tomé.
A dislocation-based constitutive law for pure zr including temperature effects.
International Journal of Plasticity, 24(5):867–895, 2008.
URL: https://www.sciencedirect.com/science/article/pii/S074964190700109X, doi:https://doi.org/10.1016/j.ijplas.2007.07.017.[BibTeX]
- I.J. Beyerlein and C.N. Tomé.
A probabilistic twin nucleatin model for HCP polycrystalline metals.
Proceedings of the Royal Society A, 466:2517–2544, 2010.
doi:https://doi.org/10.1098/rspa.2009.0661.[BibTeX]
- L. Capolungo, I.J. Beyerlein, G.C. Kaschner, and C.N. Tomé.
On the interaction between slip dislocations and twins in HCP Zr.
Materials Science and Engineering A, 513-514:42–51, 2009.
URL: https://www.sciencedirect.com/science/article/pii/S0921509309000616, doi:https://doi.org/10.1016/j.msea.2009.01.035.[BibTeX]