Grain Growth Model
Grain boundaries (GBs) migrate to reduce the total free energy of the system. Various sources of free energy drive the GB migration, including stored defect energy, deformation energy, and GB energy. Various modeling approaches have been applied to model GB migration. While all of the various methods predict similar behavior, the phase-field method has emerged as one of the more popular due to its flexibility.
In MOOSE, we have implemented the multiphase grain growth model originally developed in Chen and Yang (1994) and Chen (1995) and further developed in Moelans et al. (2008). The variables evolve to reduce the overall free energy of the system, representing GB migration. The evolution of each grain's order parameter is defined with the Allen-Cahn equation: where is the free energy functional, is the order parameter mobility, and the operator represents a variational derivative. The free energy for this problem is where is the total number of order parameters, represents the local free energy density, and represents any additional sources of energy density. For the grain growth model, where is the free energy weight and for symmetric interfacial profiles.
The model parameters , and are defined in terms of the GB energy , the diffuse GB width , and the GB mobility . In reality, the GB energy and mobility are anisotropic, depending on both the misorientation and inclination of the GB. However, it is common to assume that the GB properties are isotropic. In the case of isotropic GB energy and mobility and symmetric interfacial profiles (), Moelans et al. (2008) defined expressions for the model parameters in terms of , , and :
MOOSE Implementation
The grain growth model implementation in MOOSE can simulate grain growth in 2D and 3D, and can model large polycrystals with thousands of grains. A system of Allen-Cahn equations is solved to define the grain boundary migration. In MOOSE, we solve them in residual weak form: where is the test function. In standard MOOSE fashion, the various terms in the Allen-Cahn equation have been implemented in separate kernels. The MOOSE objects that are used in the grain growth model are summarized, below:
Materials
GBEvolution
- Defines the model parameters , , and using the equations assuming isotropic properties.GBAnistropy
- Defines parameters , , and considering misorientation dependence for the GB energy.
Kernels
TimeDerivative
- Defines the term.ACGrGrPoly
- Defines the term.ACInterface
- Defines the term.
AuxKernels
BndsCalcAux
- To visualize the GB, it is convenient to use .
Simplified MOOSE syntax
The default input file for MOOSE, in which each variable and kernel needed for each variable are added individually, would be very cumbersome for the grain growth model which can have many variables. Therefore, Actions have been created that automate the creation of the various order parameters and their kernels.
GrainGrowthAction
has the simplest syntax. It creates all of the variables and kernels needed for a basic grain growth model.PolycrystalVariablesAction
just creates the variables. It can also be used for other models besides grain growth to create a series of variables with standard naming syntax.PolycrystalKernelAction
just creates the kernels needed for each of the order parameter variables.
Initial Conditions
The initial conditions (ICs) for the polycrystal models are created using the ICs block in the input file. Under the ICs block, there is an additional block called PolycrystalICs
. Under this block, there are various options to create initial conditions for bicrystals, tricrystals and polycrystals. These options are shown in the Table, below.
IC Action name | Description |
---|---|
BicrystalBoundingBoxIC | Two grains defined by a bounding box |
BicrystalCircleGrainIC | Circle grain and matrix grain |
Tricrystal2CircleGrainsIC | Two circle grains and a matrix grain |
PolycrystalRandomIC | Randomly seeds grain order parameters |
PolycrystalColoringIC | Polycrystal initial conditions. Grain shapes are specified by a polycrystal UserObject (see below). |
The polycrystal UserObjects available for use with PolycrystalColoringIC
are:
PolycrystalVoronoi
- Creates polycrystal grain structure using a Voronoi tessellation.PolycrystalHex
- Creates a hexagonal polycrystal grain structure.PolycrystalCircles
- Creates circular/spherical grains from a file defining center coordinates and radii.PolycrystalEBSD
- Reconstructs grain structures from EBSD data in a specific file format. Used in conjunction with theEBSD Reader
.
There is also one tricrystal IC without an associated action: TricrystalTripleJunctionIC
.
Variable Remapping
If each order parameter is used to represent one grain, the computational cost gets very high to model polycrystal grain growth. However, if each order parameter is used to represent more than one grain, unphysical grain coalescence occurs when grains represented by the same variable come in contact. MOOSE uses the remapping algorithm GrainTracker to overcome this issue. Order parameters are used to represent multiple grains and grains are then remapped to other variables when needed to avoid coalescence.
Model Verification
The grain growth model has been verified against analytical solutions for two cases, a circular grain imbedded in a larger grain and a half loop grain imbedded in a larger grain.
Circular Grain
For a circular grain, the grain boundary velocity is calculated according to where is the driving force. In our simple verification cases, the driving force depends only on the curvature of the grain boundary and can be written as where is the GB energy and is the radius of curvature. The change in area with time is equal to Thus, the area at any time can be calculated with where is the initial radius of the circular grain.
This case can be easily implemented using MOOSE, and the change in the area with time is outputted using a Postprocessor, as shown in the below input which leverages data from (Schönfelder et al., 1997).
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
xmin = 0
xmax = 400
ymin = 0
ymax = 400
zmin = 0
zmax = 0
elem_type = QUAD4
uniform_refine = 1
[]
[GlobalParams]
op_num = 2
var_name_base = gr
[]
[Variables]
[./PolycrystalVariables]
[../]
[]
[ICs]
[./PolycrystalICs]
[./BicrystalCircleGrainIC]
radius = 300
x = 400
y = 0
[../]
[../]
[]
[AuxVariables]
[./bnds]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./PolycrystalKernel]
[../]
[]
[AuxKernels]
[./BndsCalc]
type = BndsCalcAux
variable = bnds
[../]
[]
[Materials]
[./Copper]
type = GBEvolution
T = 500 # K
wGB = 60 # nm
GBmob0 = 2.5e-6 #m^4/(Js) from Schoenfelder 1997
Q = 0.23 #Migration energy in eV
GBenergy = 0.708 #GB energy in J/m^2
[../]
[]
[Postprocessors]
[./gr1area]
type = ElementIntegralVariablePostprocessor
variable = gr1
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart'
petsc_options_value = 'hypre boomeramg 31'
l_tol = 1.0e-4
l_max_its = 30
nl_max_its = 20
nl_rel_tol = 1.0e-9
start_time = 0.0
num_steps = 5
dt = 80.0
[./Adaptivity]
initial_adaptivity = 2
refine_fraction = 0.8
coarsen_fraction = 0.05
max_h_level = 2
[../]
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(moose/modules/phase_field/test/tests/grain_growth/test.i)Half Loop Grain
Following the same approach from the circle grain but for the half loop geometry gives the expression
This case can also be easily implemented in MOOSE, with the area outputted using a Postprocessor:
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
nz = 0
xmin = 0
xmax = 1000
ymin = 0
ymax = 1000
zmin = 0
zmax = 0
elem_type = QUAD4
uniform_refine = 2
[]
[GlobalParams]
op_num = 2
var_name_base = gr
v = 'gr0 gr1'
[]
[Variables]
[./gr0]
order = FIRST
family = LAGRANGE
[./InitialCondition]
type = ThumbIC
xcoord = 500.0
height = 600.0
width = 400.0
invalue = 0.0
outvalue = 1.0
[../]
[../]
[./gr1]
order = FIRST
family = LAGRANGE
[./InitialCondition]
type = ThumbIC
xcoord = 500.0
height = 600.0
width = 400.0
invalue = 1.0
outvalue = 0.0
[../]
[../]
[]
[AuxVariables]
[./bnds]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./PolycrystalKernel]
[../]
[]
[AuxKernels]
[./BndsCalc]
type = BndsCalcAux
variable = bnds
[../]
[]
[BCs]
active = ' '
[./Periodic]
[./left_right]
primary = 0
secondary = 2
translation = '0 1000 0'
[../]
[./top_bottom]
primary = 1
secondary = 3
translation = '-1000 0 0'
[../]
[../]
[]
[Materials]
[./Copper]
type = GBEvolution
T = 500 # K
wGB = 60 # nm
GBmob0 = 2.5e-6 #m^4/(Js) from Schoenfelder 1997
Q = 0.23 #Migration energy in eV
GBenergy = 0.708 #GB energy in J/m^2
[../]
[]
[Postprocessors]
[./gr_area]
type = ElementIntegralVariablePostprocessor
variable = gr1
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart'
petsc_options_value = 'hypre boomeramg 31'
l_tol = 1.0e-4
l_max_its = 30
nl_max_its = 20
nl_rel_tol = 1.0e-9
start_time = 0.0
num_steps = 10
dt = 80.0
[./Adaptivity]
initial_adaptivity = 2
refine_fraction = 0.8
coarsen_fraction = 0.05
max_h_level = 2
[../]
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(moose/modules/phase_field/test/tests/grain_growth/thumb.i)Example Input Files
An example 2D input file for a polycrystal simulation using GrainTracker
is available. It should be run on at least 4 processors:
# This simulation predicts GB migration of a 2D copper polycrystal with 100 grains represented with 8 order parameters
# Mesh adaptivity and time step adaptivity are used
# An AuxVariable is used to calculate the grain boundary locations
# Postprocessors are used to record time step and the number of grains
[Mesh]
# Mesh block. Meshes can be read in or automatically generated
type = GeneratedMesh
dim = 2 # Problem dimension
nx = 44 # Number of elements in the x-direction
ny = 44 # Number of elements in the y-direction
xmax = 1000 # maximum x-coordinate of the mesh
ymax = 1000 # maximum y-coordinate of the mesh
elem_type = QUAD4 # Type of elements used in the mesh
uniform_refine = 2 # Initial uniform refinement of the mesh
[]
[GlobalParams]
# Parameters used by several kernels that are defined globally to simplify input file
op_num = 8 # Number of order parameters used
var_name_base = gr # Base name of grains
[]
[Modules]
[PhaseField]
[GrainGrowth]
[]
[]
[]
[UserObjects]
[voronoi]
type = PolycrystalVoronoi
grain_num = 100 # Number of grains
rand_seed = 10
int_width = 7
[]
[grain_tracker]
type = GrainTracker
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = voronoi
[]
[]
[]
[AuxVariables]
# Dependent variables
[unique_grains]
order = CONSTANT
family = MONOMIAL
[]
[var_indices]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
# AuxKernel block, defining the equations used to calculate the auxvars
[bnds_aux]
# AuxKernel that calculates the GB term
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[]
[unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
flood_counter = grain_tracker
field_display = UNIQUE_REGION
execute_on = 'initial timestep_end'
[]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
execute_on = 'initial timestep_end'
[]
[]
[BCs]
# Boundary Condition block
[Periodic]
[All]
auto_direction = 'x y' # Makes problem periodic in the x and y directions
[]
[]
[]
[Materials]
[CuGrGr]
# Material properties
type = GBEvolution
T = 450 # Constant temperature of the simulation (for mobility calculation)
wGB = 14 # Width of the diffuse GB
GBmob0 = 2.5e-6 #m^4(Js) for copper from schonfelder1997molecular bibtex entry
Q = 0.23 #eV for copper from schonfelder1997molecular bibtex entry
GBenergy = 0.708 #J/m^2 from schonfelder1997molecular bibtex entry
[]
[]
[Postprocessors]
# Scalar postprocessors
[dt]
# Outputs the current time step
type = TimestepSize
[]
[]
[Executioner]
type = Transient # Type of executioner, here it is transient with an adaptive time step
scheme = bdf2 # Type of time integration (2nd order backward euler), defaults to 1st order backward euler
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
# Uses newton iteration to solve the problem.
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
l_max_its = 50 # Max number of linear iterations
l_tol = 1e-4 # Relative tolerance for linear solves
nl_max_its = 10 # Max number of nonlinear iterations
end_time = 4000
[TimeStepper]
type = IterationAdaptiveDT
dt = 20 # Initial time step. In this simulation it changes.
optimal_iterations = 6 # Time step will adapt to maintain this number of nonlinear iterations
[]
[Adaptivity]
# Block that turns on mesh adaptivity. Note that mesh will never coarsen beyond initial mesh (before uniform refinement)
initial_adaptivity = 2 # Number of times mesh is adapted to initial condition
refine_fraction = 0.8 # Fraction of high error that will be refined
coarsen_fraction = 0.05 # Fraction of low error that will coarsened
max_h_level = 2 # Max number of refinements used, starting from initial mesh (before uniform refinement)
[]
[]
[Outputs]
exodus = true # Exodus file will be outputted
csv = true
[]
(moose/modules/phase_field/examples/grain_growth/grain_growth_2D_graintracker.i)An example 3D input file is also available. It should be run on at least 2 processors:
# This simulation predicts GB migration of a 3D copper polycrystal with 25 grains represented with 15 order parameters
# Time step adaptivity are used
# An AuxVariable is used to calculate the grain boundary locations
# Postprocessors are used to record time step and the number of grains
[Mesh]
# Mesh block. Meshes can be read in or automatically generated
type = GeneratedMesh
dim = 3 # Problem dimension
nx = 10 # Number of elements in the x-direction
ny = 10 # Number of elements in the y-direction
nz = 10
xmax = 1000 # maximum x-coordinate of the mesh
ymax = 1000 # maximum y-coordinate of the mesh
zmax = 1000
uniform_refine = 1 # Initial uniform refinement of the mesh
parallel_type = distributed
[]
[GlobalParams]
# Parameters used by several kernels that are defined globally to simplify input file
op_num = 15 # Number of order parameters used
var_name_base = gr # Base name of grains
order = CONSTANT
family = MONOMIAL
[]
[Modules]
[PhaseField]
[GrainGrowth]
family = LAGRANGE
order = FIRST
[]
[]
[]
[UserObjects]
[voronoi]
type = PolycrystalVoronoi
grain_num = 25 # Number of grains
rand_seed = 10
coloring_algorithm = jp
[]
[grain_tracker]
type = GrainTracker
threshold = 0.2
connecting_threshold = 0.08
compute_halo_maps = true # Only necessary for displaying HALOS
polycrystal_ic_uo = voronoi
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = voronoi
[]
[]
[]
[AuxVariables]
# Dependent variables
[unique_grains]
[]
[var_indices]
[]
[ghost_regions]
[]
[halos]
[]
[halo0]
[]
[halo1]
[]
[halo2]
[]
[halo3]
[]
[halo4]
[]
[halo5]
[]
[halo6]
[]
[halo7]
[]
[halo8]
[]
[halo9]
[]
[halo10]
[]
[halo11]
[]
[halo12]
[]
[halo13]
[]
[halo14]
[]
[proc]
[]
[]
[AuxKernels]
# AuxKernel block, defining the equations used to calculate the auxvars
[unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
flood_counter = grain_tracker
field_display = UNIQUE_REGION
execute_on = 'initial timestep_end'
[]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
execute_on = 'initial timestep_end'
[]
[ghosted_entities]
type = FeatureFloodCountAux
variable = ghost_regions
flood_counter = grain_tracker
field_display = GHOSTED_ENTITIES
execute_on = 'initial timestep_end'
[]
[halos]
type = FeatureFloodCountAux
variable = halos
flood_counter = voronoi
field_display = HALOS
execute_on = 'initial timestep_end'
[]
[halo0]
type = FeatureFloodCountAux
variable = halo0
map_index = 0
field_display = HALOS
flood_counter = grain_tracker
execute_on = 'initial timestep_end'
[]
[halo1]
type = FeatureFloodCountAux
variable = halo1
map_index = 1
field_display = HALOS
flood_counter = grain_tracker
[]
[halo2]
type = FeatureFloodCountAux
variable = halo2
map_index = 2
field_display = HALOS
flood_counter = grain_tracker
[]
[halo3]
type = FeatureFloodCountAux
variable = halo3
map_index = 3
field_display = HALOS
flood_counter = grain_tracker
[]
[halo4]
type = FeatureFloodCountAux
variable = halo4
map_index = 4
field_display = HALOS
flood_counter = grain_tracker
[]
[halo5]
type = FeatureFloodCountAux
variable = halo5
map_index = 5
field_display = HALOS
flood_counter = grain_tracker
[]
[halo6]
type = FeatureFloodCountAux
variable = halo6
map_index = 6
field_display = HALOS
flood_counter = grain_tracker
[]
[halo7]
type = FeatureFloodCountAux
variable = halo7
map_index = 7
field_display = HALOS
flood_counter = grain_tracker
[]
[halo8]
type = FeatureFloodCountAux
variable = halo8
map_index = 8
field_display = HALOS
flood_counter = grain_tracker
[]
[halo9]
type = FeatureFloodCountAux
variable = halo9
map_index = 9
field_display = HALOS
flood_counter = grain_tracker
[]
[halo10]
type = FeatureFloodCountAux
variable = halo10
map_index = 10
field_display = HALOS
flood_counter = grain_tracker
[]
[halo11]
type = FeatureFloodCountAux
variable = halo11
map_index = 11
field_display = HALOS
flood_counter = grain_tracker
[]
[halo12]
type = FeatureFloodCountAux
variable = halo12
map_index = 12
field_display = HALOS
flood_counter = grain_tracker
[]
[halo13]
type = FeatureFloodCountAux
variable = halo13
map_index = 13
field_display = HALOS
flood_counter = grain_tracker
[]
[halo14]
type = FeatureFloodCountAux
variable = halo14
map_index = 14
field_display = HALOS
flood_counter = grain_tracker
[]
[proc]
type = ProcessorIDAux
variable = proc
execute_on = 'initial timestep_end'
[]
[]
[Materials]
[CuGrGr]
# Material properties
type = GBEvolution
T = 450 # Constant temperature of the simulation (for mobility calculation)
wGB = 125 # Width of the diffuse GB
GBmob0 = 2.5e-6 #m^4(Js) for copper from schonfelder1997molecular bibtex entry
Q = 0.23 #eV for copper from schonfelder1997molecular bibtex entry
GBenergy = 0.708 #J/m^2 from schonfelder1997molecular bibtex entry
[]
[]
[Postprocessors]
# Scalar postprocessors
[dt]
# Outputs the current time step
type = TimestepSize
[]
[]
[Executioner]
type = Transient # Type of executioner, here it is transient with an adaptive time step
scheme = bdf2 # Type of time integration (2nd order backward euler), defaults to 1st order backward euler
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
# Uses newton iteration to solve the problem.
petsc_options_iname = '-pc_type'
petsc_options_value = 'asm'
l_max_its = 30 # Max number of linear iterations
l_tol = 1e-4 # Relative tolerance for linear solves
nl_max_its = 20 # Max number of nonlinear iterations
start_time = 0.0
end_time = 4000
[TimeStepper]
type = IterationAdaptiveDT
dt = 25 # Initial time step. In this simulation it changes.
optimal_iterations = 6 # Time step will adapt to maintain this number of nonlinear iterations
[]
[]
[Outputs]
exodus = true
csv = true
[pg]
type = PerfGraphOutput
execute_on = 'initial final' # Default is "final"
level = 2 # Default is 1
[]
[]
(moose/modules/phase_field/examples/grain_growth/grain_growth_3D.i)There is also an example reconstructing the initial grain structure from EBSD data. It should be run on at least 4 processors:
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = IN100_120x120.txt
pre_refine = 2
[]
[]
[GlobalParams]
op_num = 8
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
enable_var_coloring = true
[]
[grain_tracker]
type = GrainTracker
flood_entity_type = ELEMENTAL
compute_halo_maps = true # For displaying HALO fields
polycrystal_ic_uo = ebsd
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[AuxVariables]
[bnds]
[]
[unique_grains_ic]
order = CONSTANT
family = MONOMIAL
[]
[unique_grains]
order = CONSTANT
family = MONOMIAL
[]
[ghost_elements]
order = CONSTANT
family = MONOMIAL
[]
[halos]
order = CONSTANT
family = MONOMIAL
[]
[var_indices_ic]
order = CONSTANT
family = MONOMIAL
[]
[var_indices]
order = CONSTANT
family = MONOMIAL
[]
[ebsd_grains]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[PolycrystalKernel]
[]
[]
[AuxKernels]
[BndsCalc]
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[]
[ghost_elements]
type = FeatureFloodCountAux
variable = ghost_elements
field_display = GHOSTED_ENTITIES
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
[]
[halos]
type = FeatureFloodCountAux
variable = halos
field_display = HALOS
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
[]
[var_indices_ic]
type = FeatureFloodCountAux
variable = var_indices_ic
execute_on = 'initial'
flood_counter = ebsd
field_display = VARIABLE_COLORING
[]
[unique_grains_ic]
type = FeatureFloodCountAux
variable = unique_grains_ic
execute_on = 'initial'
flood_counter = ebsd
field_display = UNIQUE_REGION
[]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
[]
[unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = UNIQUE_REGION
[]
[grain_aux]
type = EBSDReaderPointDataAux
variable = ebsd_grains
ebsd_reader = ebsd_reader
data_name = 'feature_id'
execute_on = 'initial timestep_end'
[]
[]
[Modules]
[PhaseField]
[EulerAngles2RGB]
crystal_structure = cubic
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
[]
[]
[]
[Materials]
[Copper]
# T = 500 # K
type = GBEvolution
T = 500
wGB = 0.6 # um
GBmob0 = 2.5e-6 # m^4/(Js) from Schoenfelder 1997
Q = 0.23 # Migration energy in eV
GBenergy = 0.708 # GB energy in J/m^2
molar_volume = 7.11e-6 # Molar volume in m^3/mol
length_scale = 1.0e-6
time_scale = 1.0e-6
[]
[]
[Postprocessors]
[dt]
type = TimestepSize
[]
[n_elements]
type = NumElements
execute_on = 'initial timestep_end'
[]
[n_nodes]
type = NumNodes
execute_on = 'initial timestep_end'
[]
[DOFs]
type = NumDOFs
[]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = PJFNK
petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_strong_threshold'
petsc_options_value = 'hypre boomeramg 0.7'
l_tol = 1.0e-4
l_max_its = 20
nl_max_its = 20
nl_rel_tol = 1.0e-8
start_time = 0.0
num_steps = 30
[TimeStepper]
type = IterationAdaptiveDT
cutback_factor = 0.9
dt = 10.0
growth_factor = 1.1
optimal_iterations = 7
[]
[Adaptivity]
initial_adaptivity = 2
refine_fraction = 0.7
coarsen_fraction = 0.1
max_h_level = 2
[]
[]
[Outputs]
exodus = true
checkpoint = true
perf_graph = true
[]
(moose/modules/phase_field/examples/ebsd_reconstruction/IN100-111grn.i)References
- Long-Qing Chen.
A novel computer simulation technique for modeling grain growth.
Scripta Metallurgica et Materialia, 32:115–120, 1995.
URL: https://www.sciencedirect.com/science/article/pii/S0956716X99800223, doi:10.1016/S0956-716X(99)80022-3.[BibTeX]
- Long-Qing Chen and Wei Yang.
Computer simulation of the domain dynamics of a quenched system with a large number of nonconserved order parameters: the grain-growth kinetics.
Phys. Rev. B, 50:15752–15756, Dec 1994.
URL: https://link.aps.org/doi/10.1103/PhysRevB.50.15752, doi:10.1103/PhysRevB.50.15752.[BibTeX]
- N. Moelans, B. Blanpain, and P. Wollants.
Quantitative analysis of grain boundary properties in a generalized phase field model for grain growth in anisotropic systems.
Physical Review B, 78(2):024113, Jul 2008.
URL: http://link.aps.org/doi/10.1103/PhysRevB.78.024113 (visited on 2016-06-02), doi:10.1103/PhysRevB.78.024113.[BibTeX]
- B Schönfelder, D Wolf, SR Phillpot, and M Furtkamp.
Molecular-dynamics method for the simulation of grain-boundary migration.
Interface Science, 5:245–262, 1997.[BibTeX]