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 nameDescription
BicrystalBoundingBoxICTwo grains defined by a bounding box
BicrystalCircleGrainICCircle grain and matrix grain
Tricrystal2CircleGrainsICTwo circle grains and a matrix grain
PolycrystalRandomICRandomly seeds grain order parameters
PolycrystalColoringICPolycrystal 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 the EBSD 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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  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<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 1
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  op_num = 2
  var_name_base = gr
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [./PolycrystalVariables<<<{"href": "../../syntax/Variables/PolycrystalVariables/index.html"}>>>]
  [../]
[]

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  [./PolycrystalICs<<<{"href": "../../syntax/ICs/PolycrystalICs/index.html"}>>>]
    [./BicrystalCircleGrainIC<<<{"href": "../../syntax/ICs/PolycrystalICs/BicrystalCircleGrainIC/index.html"}>>>]
      radius<<<{"description": "Void radius"}>>> = 300
      x<<<{"description": "The x coordinate of the circle grain center"}>>> = 400
      y<<<{"description": "The y coordinate of the circle grain center"}>>> = 0
    [../]
  [../]
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [./bnds]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
  [../]
[]

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [./PolycrystalKernel<<<{"href": "../../syntax/Kernels/PolycrystalKernel/index.html"}>>>]
  [../]
[]

[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  [./BndsCalc]
    type = BndsCalcAux<<<{"description": "Calculate location of grain boundaries in a polycrystalline sample", "href": "../../source/auxkernels/BndsCalcAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = bnds
  [../]
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [./Copper]
    type = GBEvolution<<<{"description": "Computes necessary material properties for the isotropic grain growth model", "href": "../../source/materials/GBEvolution.html"}>>>
    T<<<{"description": "Temperature in Kelvin"}>>> = 500 # K
    wGB<<<{"description": "Diffuse GB width in the length scale of the model"}>>> = 60 # nm
    GBmob0<<<{"description": "Grain boundary mobility prefactor in m^4/(J*s)"}>>> = 2.5e-6 #m^4/(Js) from Schoenfelder 1997
    Q<<<{"description": "Grain boundary migration activation energy in eV"}>>> = 0.23 #Migration energy in eV
    GBenergy<<<{"description": "Grain boundary energy in J/m^2"}>>> = 0.708 #GB energy in J/m^2
  [../]
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [./gr1area]
    type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = gr1
  [../]
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  [./SMP]
    type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
    full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
  [../]
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Executioner/Adaptivity/index.html"}>>>]
    initial_adaptivity<<<{"description": "The number of adaptivity steps to perform using the initial conditions"}>>> = 2
    refine_fraction<<<{"description": "The fraction of elements or error to refine. Should be between 0 and 1."}>>> = 0.8
    coarsen_fraction<<<{"description": "The fraction of elements or error to coarsen. Should be between 0 and 1."}>>> = 0.05
    max_h_level<<<{"description": "Maximum number of times a single element can be refined. If 0 then infinite."}>>> = 2
  [../]
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = 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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  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<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 2
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  op_num = 2
  var_name_base = gr
  v = 'gr0 gr1'
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [./gr0]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
    [./InitialCondition]
      type = ThumbIC
      xcoord = 500.0
      height = 600.0
      width = 400.0
      invalue = 0.0
      outvalue = 1.0
    [../]
  [../]

  [./gr1]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
    [./InitialCondition]
      type = ThumbIC
      xcoord = 500.0
      height = 600.0
      width = 400.0
      invalue = 1.0
      outvalue = 0.0
    [../]
  [../]
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [./bnds]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = FIRST
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = LAGRANGE
  [../]
[]

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [./PolycrystalKernel<<<{"href": "../../syntax/Kernels/PolycrystalKernel/index.html"}>>>]
  [../]
[]

[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  [./BndsCalc]
    type = BndsCalcAux<<<{"description": "Calculate location of grain boundaries in a polycrystalline sample", "href": "../../source/auxkernels/BndsCalcAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = bnds
  [../]
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  active<<<{"description": "If specified only the blocks named will be visited and made active"}>>> = ' '

  [./Periodic<<<{"href": "../../syntax/BCs/Periodic/index.html"}>>>]
    [./left_right]
      primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = 0
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = 2
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '0 1000 0'
    [../]

    [./top_bottom]
      primary<<<{"description": "Boundary ID associated with the primary boundary."}>>> = 1
      secondary<<<{"description": "Boundary ID associated with the secondary boundary."}>>> = 3
      translation<<<{"description": "Vector that translates coordinates on the primary boundary to coordinates on the secondary boundary."}>>> = '-1000 0 0'
    [../]
  [../]
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [./Copper]
    type = GBEvolution<<<{"description": "Computes necessary material properties for the isotropic grain growth model", "href": "../../source/materials/GBEvolution.html"}>>>
    T<<<{"description": "Temperature in Kelvin"}>>> = 500 # K
    wGB<<<{"description": "Diffuse GB width in the length scale of the model"}>>> = 60 # nm
    GBmob0<<<{"description": "Grain boundary mobility prefactor in m^4/(J*s)"}>>> = 2.5e-6 #m^4/(Js) from Schoenfelder 1997
    Q<<<{"description": "Grain boundary migration activation energy in eV"}>>> = 0.23 #Migration energy in eV
    GBenergy<<<{"description": "Grain boundary energy in J/m^2"}>>> = 0.708 #GB energy in J/m^2
  [../]
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [./gr_area]
    type = ElementIntegralVariablePostprocessor<<<{"description": "Computes a volume integral of the specified variable", "href": "../../source/postprocessors/ElementIntegralVariablePostprocessor.html"}>>>
    variable<<<{"description": "The name of the variable that this object operates on"}>>> = gr1
  [../]
[]

[Preconditioning<<<{"href": "../../syntax/Preconditioning/index.html"}>>>]
  [./SMP]
   type = SMP<<<{"description": "Single matrix preconditioner (SMP) builds a preconditioner using user defined off-diagonal parts of the Jacobian.", "href": "../../source/preconditioners/SingleMatrixPreconditioner.html"}>>>
   full<<<{"description": "Set to true if you want the full set of couplings between variables simply for convenience so you don't have to set every off_diag_row and off_diag_column combination."}>>> = true
  [../]
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Executioner/Adaptivity/index.html"}>>>]
    initial_adaptivity<<<{"description": "The number of adaptivity steps to perform using the initial conditions"}>>> = 2
    refine_fraction<<<{"description": "The fraction of elements or error to refine. Should be between 0 and 1."}>>> = 0.8
    coarsen_fraction<<<{"description": "The fraction of elements or error to coarsen. Should be between 0 and 1."}>>> = 0.05
    max_h_level<<<{"description": "Maximum number of times a single element can be refined. If 0 then infinite."}>>> = 2
  [../]
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'timestep_end'
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = 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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  # 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<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 2 # Initial uniform refinement of the mesh
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  # 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<<<{"href": "../../syntax/Modules/index.html"}>>>]
  [PhaseField<<<{"href": "../../syntax/Modules/PhaseField/index.html"}>>>]
    [GrainGrowth<<<{"href": "../../syntax/Modules/PhaseField/GrainGrowth/index.html"}>>>]
    []
  []
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [voronoi]
    type = PolycrystalVoronoi<<<{"description": "Random Voronoi tessellation polycrystal (used by PolycrystalVoronoiAction)", "href": "../../source/userobjects/PolycrystalVoronoi.html"}>>>
    grain_num<<<{"description": "Number of grains being represented by the order parameters"}>>> = 100 # Number of grains
    rand_seed<<<{"description": "The random seed"}>>> = 10
    int_width<<<{"description": "Width of diffuse interfaces"}>>> = 7
  []
  [grain_tracker]
    type = GrainTracker<<<{"description": "Grain Tracker object for running reduced order parameter simulations without grain coalescence.", "href": "../../source/postprocessors/GrainTracker.html"}>>>
  []
[]

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  [PolycrystalICs<<<{"href": "../../syntax/ICs/PolycrystalICs/index.html"}>>>]
    [PolycrystalColoringIC<<<{"href": "../../syntax/ICs/PolycrystalICs/PolycrystalColoringIC/index.html"}>>>]
      polycrystal_ic_uo<<<{"description": "Optional: TODO"}>>> = voronoi
    []
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  # Dependent variables
  [unique_grains]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [var_indices]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
[]

[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  # AuxKernel block, defining the equations used to calculate the auxvars
  [bnds_aux]
    # AuxKernel that calculates the GB term
    type = BndsCalcAux<<<{"description": "Calculate location of grain boundaries in a polycrystalline sample", "href": "../../source/auxkernels/BndsCalcAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = bnds
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [unique_grains]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = unique_grains
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = UNIQUE_REGION
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [var_indices]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = var_indices
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = VARIABLE_COLORING
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
[]

[BCs<<<{"href": "../../syntax/BCs/index.html"}>>>]
  # Boundary Condition block
  [Periodic<<<{"href": "../../syntax/BCs/Periodic/index.html"}>>>]
    [All]
      auto_direction<<<{"description": "If using a generated mesh, you can specify just the dimension(s) you want to mark as periodic"}>>> = 'x y' # Makes problem periodic in the x and y directions
    []
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [CuGrGr]
    # Material properties
    type = GBEvolution<<<{"description": "Computes necessary material properties for the isotropic grain growth model", "href": "../../source/materials/GBEvolution.html"}>>>
    T<<<{"description": "Temperature in Kelvin"}>>> = 450 # Constant temperature of the simulation (for mobility calculation)
    wGB<<<{"description": "Diffuse GB width in the length scale of the model"}>>> = 14 # Width of the diffuse GB
    GBmob0<<<{"description": "Grain boundary mobility prefactor in m^4/(J*s)"}>>> = 2.5e-6 #m^4(Js) for copper from schonfelder1997molecular bibtex entry
    Q<<<{"description": "Grain boundary migration activation energy in eV"}>>> = 0.23 #eV for copper from schonfelder1997molecular bibtex entry
    GBenergy<<<{"description": "Grain boundary energy in J/m^2"}>>> = 0.708 #J/m^2 from schonfelder1997molecular bibtex entry
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  # Scalar postprocessors
  [dt]
    # Outputs the current time step
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../../source/postprocessors/TimestepSize.html"}>>>
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
    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<<<{"href": "../../syntax/Executioner/Adaptivity/index.html"}>>>]
    # Block that turns on mesh adaptivity. Note that mesh will never coarsen beyond initial mesh (before uniform refinement)
    initial_adaptivity<<<{"description": "The number of adaptivity steps to perform using the initial conditions"}>>> = 2 # Number of times mesh is adapted to initial condition
    refine_fraction<<<{"description": "The fraction of elements or error to refine. Should be between 0 and 1."}>>> = 0.8 # Fraction of high error that will be refined
    coarsen_fraction<<<{"description": "The fraction of elements or error to coarsen. Should be between 0 and 1."}>>> = 0.05 # Fraction of low error that will coarsened
    max_h_level<<<{"description": "Maximum number of times a single element can be refined. If 0 then infinite."}>>> = 2 # Max number of refinements used, starting from initial mesh (before uniform refinement)
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true # Exodus file will be outputted
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = 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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  # 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<<<{"description": "Specify the level of uniform refinement applied to the initial mesh"}>>> = 1 # Initial uniform refinement of the mesh

  parallel_type = distributed
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  # 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<<<{"href": "../../syntax/Modules/index.html"}>>>]
  [PhaseField<<<{"href": "../../syntax/Modules/PhaseField/index.html"}>>>]
    [GrainGrowth<<<{"href": "../../syntax/Modules/PhaseField/GrainGrowth/index.html"}>>>]
      family<<<{"description": "Specifies the family of FE shape function to use for the order parameters"}>>> = LAGRANGE
      order<<<{"description": "Specifies the order of the FE shape function to use for the order parameters"}>>> = FIRST
    []
  []
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [voronoi]
    type = PolycrystalVoronoi<<<{"description": "Random Voronoi tessellation polycrystal (used by PolycrystalVoronoiAction)", "href": "../../source/userobjects/PolycrystalVoronoi.html"}>>>
    grain_num<<<{"description": "Number of grains being represented by the order parameters"}>>> = 25 # Number of grains
    rand_seed<<<{"description": "The random seed"}>>> = 10
    coloring_algorithm<<<{"description": "The grain neighbor graph coloring algorithm to use: \"jp\" (DEFAULT) Jones and Plassmann, an efficient coloring algorithm, \"power\" an alternative stochastic algorithm, \"greedy\", a greedy assignment algorithm with stochastic updates to guarantee a valid coloring, \"bt\", a back tracking algorithm that produces good distributions but may experience exponential run time in the worst case scenario (works well on medium to large 2D problems)"}>>> = jp
  []
  [grain_tracker]
    type = GrainTracker<<<{"description": "Grain Tracker object for running reduced order parameter simulations without grain coalescence.", "href": "../../source/postprocessors/GrainTracker.html"}>>>
    threshold<<<{"description": "The threshold value for which a new feature may be started"}>>> = 0.2
    connecting_threshold<<<{"description": "The threshold for which an existing feature may be extended (defaults to \"threshold\")"}>>> = 0.08
    compute_halo_maps<<<{"description": "Instruct the Postprocessor to communicate proper halo information to all ranks"}>>> = true # Only necessary for displaying HALOS
    polycrystal_ic_uo<<<{"description": "Optional: Polycrystal IC object"}>>> = voronoi
  []
[]

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  [PolycrystalICs<<<{"href": "../../syntax/ICs/PolycrystalICs/index.html"}>>>]
    [PolycrystalColoringIC<<<{"href": "../../syntax/ICs/PolycrystalICs/PolycrystalColoringIC/index.html"}>>>]
      polycrystal_ic_uo<<<{"description": "Optional: TODO"}>>> = voronoi
    []
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  # 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<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  # AuxKernel block, defining the equations used to calculate the auxvars
  [unique_grains]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = unique_grains
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = UNIQUE_REGION
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [var_indices]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = var_indices
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = VARIABLE_COLORING
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [ghosted_entities]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = ghost_regions
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = GHOSTED_ENTITIES
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [halos]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halos
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = voronoi
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [halo0]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo0
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 0
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [halo1]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo1
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 1
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo2]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo2
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 2
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo3]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo3
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 3
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo4]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo4
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 4
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo5]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo5
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 5
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo6]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo6
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 6
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo7]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo7
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 7
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo8]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo8
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 8
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo9]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo9
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 9
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo10]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo10
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 10
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo11]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo11
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 11
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo12]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo12
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 12
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo13]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo13
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 13
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halo14]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halo14
    map_index<<<{"description": "The index of which map to retrieve values from when using FeatureFloodCount with multiple maps."}>>> = 14
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [proc]
    type = ProcessorIDAux<<<{"description": "Creates a field showing the processors and partitioning.", "href": "../../source/auxkernels/ProcessorIDAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = proc
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []

[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [CuGrGr]
    # Material properties
    type = GBEvolution<<<{"description": "Computes necessary material properties for the isotropic grain growth model", "href": "../../source/materials/GBEvolution.html"}>>>
    T<<<{"description": "Temperature in Kelvin"}>>> = 450 # Constant temperature of the simulation (for mobility calculation)
    wGB<<<{"description": "Diffuse GB width in the length scale of the model"}>>> = 125 # Width of the diffuse GB
    GBmob0<<<{"description": "Grain boundary mobility prefactor in m^4/(J*s)"}>>> = 2.5e-6 #m^4(Js) for copper from schonfelder1997molecular bibtex entry
    Q<<<{"description": "Grain boundary migration activation energy in eV"}>>> = 0.23 #eV for copper from schonfelder1997molecular bibtex entry
    GBenergy<<<{"description": "Grain boundary energy in J/m^2"}>>> = 0.708 #J/m^2 from schonfelder1997molecular bibtex entry
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  # Scalar postprocessors
  [dt]
    # Outputs the current time step
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../../source/postprocessors/TimestepSize.html"}>>>
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
    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<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  csv<<<{"description": "Output the scalar variable and postprocessors to a *.csv file using the default CSV output."}>>> = true
  [pg]
    type = PerfGraphOutput<<<{"description": "Controls output of the PerfGraph: the performance log for MOOSE", "href": "../../source/outputs/PerfGraphOutput.html"}>>>
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial final' # Default is "final"
    level<<<{"description": "The level of detail to output.  Higher levels will yield more detail."}>>> = 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<<<{"href": "../../syntax/Mesh/index.html"}>>>]
  [ebsd_mesh]
    type = EBSDMeshGenerator<<<{"description": "Mesh generated from a specified DREAM.3D EBSD data file.", "href": "../../source/meshgenerators/EBSDMeshGenerator.html"}>>>
    filename<<<{"description": "The name of the file containing the EBSD data"}>>> = IN100_120x120.txt
    pre_refine<<<{"description": "Number of coarsening levels available in adaptive mesh refinement. The resulting mesh will have one mesh element per EBSD data cell, but will be based on a refined coarser mesh with 'pre_refine' levels of refinement. This requires all dimension of the EBSD data to be divisible by 2^pre_refine.(this parameter was formerly called 'uniform_refine')"}>>> = 2
  []
[]

[GlobalParams<<<{"href": "../../syntax/GlobalParams/index.html"}>>>]
  op_num = 8
  var_name_base = gr
[]

[UserObjects<<<{"href": "../../syntax/UserObjects/index.html"}>>>]
  [ebsd_reader]
    type = EBSDReader<<<{"description": "Load and manage DREAM.3D EBSD data files for running simulations on reconstructed microstructures.", "href": "../../source/userobjects/EBSDReader.html"}>>>
  []
  [ebsd]
    type = PolycrystalEBSD<<<{"description": "Object for setting up a polycrystal structure from an EBSD Datafile", "href": "../../source/userobjects/PolycrystalEBSD.html"}>>>
    coloring_algorithm<<<{"description": "The grain neighbor graph coloring algorithm to use: \"jp\" (DEFAULT) Jones and Plassmann, an efficient coloring algorithm, \"power\" an alternative stochastic algorithm, \"greedy\", a greedy assignment algorithm with stochastic updates to guarantee a valid coloring, \"bt\", a back tracking algorithm that produces good distributions but may experience exponential run time in the worst case scenario (works well on medium to large 2D problems)"}>>> = bt
    ebsd_reader<<<{"description": "EBSD Reader for initial condition"}>>> = ebsd_reader
    enable_var_coloring<<<{"description": "Instruct the Postprocessor to populate the variable index map."}>>> = true
  []
  [grain_tracker]
    type = GrainTracker<<<{"description": "Grain Tracker object for running reduced order parameter simulations without grain coalescence.", "href": "../../source/postprocessors/GrainTracker.html"}>>>
    flood_entity_type<<<{"description": "Determines whether the flood algorithm runs on nodes or elements"}>>> = ELEMENTAL
    compute_halo_maps<<<{"description": "Instruct the Postprocessor to communicate proper halo information to all ranks"}>>> = true # For displaying HALO fields
    polycrystal_ic_uo<<<{"description": "Optional: Polycrystal IC object"}>>> = ebsd
  []
[]

[ICs<<<{"href": "../../syntax/ICs/index.html"}>>>]
  [PolycrystalICs<<<{"href": "../../syntax/ICs/PolycrystalICs/index.html"}>>>]
    [PolycrystalColoringIC<<<{"href": "../../syntax/ICs/PolycrystalICs/PolycrystalColoringIC/index.html"}>>>]
      polycrystal_ic_uo<<<{"description": "Optional: TODO"}>>> = ebsd
    []
  []
[]

[Variables<<<{"href": "../../syntax/Variables/index.html"}>>>]
  [PolycrystalVariables<<<{"href": "../../syntax/Variables/PolycrystalVariables/index.html"}>>>]
  []
[]

[AuxVariables<<<{"href": "../../syntax/AuxVariables/index.html"}>>>]
  [bnds]
  []
  [unique_grains_ic]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [unique_grains]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [ghost_elements]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [halos]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [var_indices_ic]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [var_indices]
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
  []
  [ebsd_grains]
    family<<<{"description": "Specifies the family of FE shape functions to use for this variable"}>>> = MONOMIAL
    order<<<{"description": "Specifies the order of the FE shape function to use for this variable (additional orders not listed are allowed)"}>>> = CONSTANT
  []
[]

[Kernels<<<{"href": "../../syntax/Kernels/index.html"}>>>]
  [PolycrystalKernel<<<{"href": "../../syntax/Kernels/PolycrystalKernel/index.html"}>>>]
  []
[]

[AuxKernels<<<{"href": "../../syntax/AuxKernels/index.html"}>>>]
  [BndsCalc]
    type = BndsCalcAux<<<{"description": "Calculate location of grain boundaries in a polycrystalline sample", "href": "../../source/auxkernels/BndsCalcAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = bnds
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [ghost_elements]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = ghost_elements
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = GHOSTED_ENTITIES
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [halos]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = halos
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = HALOS
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
  []
  [var_indices_ic]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = var_indices_ic
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = ebsd
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = VARIABLE_COLORING
  []
  [unique_grains_ic]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = unique_grains_ic
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = ebsd
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = UNIQUE_REGION
  []
  [var_indices]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = var_indices
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = VARIABLE_COLORING
  []
  [unique_grains]
    type = FeatureFloodCountAux<<<{"description": "Feature detection by connectivity analysis", "href": "../../source/auxkernels/FeatureFloodCountAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = unique_grains
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
    flood_counter<<<{"description": "The FeatureFloodCount UserObject to get values from."}>>> = grain_tracker
    field_display<<<{"description": "Determines how the auxilary field should be colored. (UNIQUE_REGION and VARIABLE_COLORING are nodal, CENTROID is elemental, default: UNIQUE_REGION)"}>>> = UNIQUE_REGION
  []
  [grain_aux]
    type = EBSDReaderPointDataAux<<<{"description": "Outputs the requested EBSD reader point data.", "href": "../../source/auxkernels/EBSDReaderPointDataAux.html"}>>>
    variable<<<{"description": "The name of the variable that this object applies to"}>>> = ebsd_grains
    ebsd_reader<<<{"description": "The EBSDReader GeneralUserObject"}>>> = ebsd_reader
    data_name<<<{"description": "The data to be extracted from the EBSD data by this AuxKernel"}>>> = 'feature_id'
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
[]

[Modules<<<{"href": "../../syntax/Modules/index.html"}>>>]
  [PhaseField<<<{"href": "../../syntax/Modules/PhaseField/index.html"}>>>]
    [EulerAngles2RGB<<<{"href": "../../syntax/Modules/PhaseField/EulerAngles2RGB/index.html"}>>>]
      crystal_structure<<<{"description": "Crystal structure of the material"}>>> = cubic
      euler_angle_provider<<<{"description": "Name of Euler angle provider user object"}>>> = ebsd_reader
      grain_tracker<<<{"description": "The GrainTracker UserObject to get values from."}>>> = grain_tracker
    []
  []
[]

[Materials<<<{"href": "../../syntax/Materials/index.html"}>>>]
  [Copper]
    # T = 500 # K
    type = GBEvolution<<<{"description": "Computes necessary material properties for the isotropic grain growth model", "href": "../../source/materials/GBEvolution.html"}>>>
    T<<<{"description": "Temperature in Kelvin"}>>> = 500
    wGB<<<{"description": "Diffuse GB width in the length scale of the model"}>>> = 0.6 # um
    GBmob0<<<{"description": "Grain boundary mobility prefactor in m^4/(J*s)"}>>> = 2.5e-6 # m^4/(Js) from Schoenfelder 1997
    Q<<<{"description": "Grain boundary migration activation energy in eV"}>>> = 0.23 # Migration energy in eV
    GBenergy<<<{"description": "Grain boundary energy in J/m^2"}>>> = 0.708 # GB energy in J/m^2
    molar_volume<<<{"description": "Molar volume in m^3/mol, needed for temperature gradient driving force"}>>> = 7.11e-6 # Molar volume in m^3/mol
    length_scale<<<{"description": "Length scale in m, where default is nm"}>>> = 1.0e-6
    time_scale<<<{"description": "Time scale in s, where default is ns"}>>> = 1.0e-6
  []
[]

[Postprocessors<<<{"href": "../../syntax/Postprocessors/index.html"}>>>]
  [dt]
    type = TimestepSize<<<{"description": "Reports the timestep size", "href": "../../source/postprocessors/TimestepSize.html"}>>>
  []
  [n_elements]
    type = NumElements<<<{"description": "Return the number of active or total elements in the simulation.", "href": "../../source/postprocessors/NumElements.html"}>>>
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [n_nodes]
    type = NumNodes<<<{"description": "Returns the total number of nodes in a simulation (works with DistributedMesh)", "href": "../../source/postprocessors/NumNodes.html"}>>>
    execute_on<<<{"description": "The list of flag(s) indicating when this object should be executed. For a description of each flag, see https://mooseframework.inl.gov/source/interfaces/SetupInterface.html."}>>> = 'initial timestep_end'
  []
  [DOFs]
    type = NumDOFs<<<{"description": "Return the number of Degrees of freedom from either the NL, Aux or both systems.", "href": "../../source/postprocessors/NumDOFs.html"}>>>
  []
[]

[Executioner<<<{"href": "../../syntax/Executioner/index.html"}>>>]
  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<<<{"href": "../../syntax/Executioner/TimeStepper/index.html"}>>>]
    type = IterationAdaptiveDT
    cutback_factor = 0.9
    dt = 10.0
    growth_factor = 1.1
    optimal_iterations = 7
  []

  [Adaptivity<<<{"href": "../../syntax/Executioner/Adaptivity/index.html"}>>>]
    initial_adaptivity<<<{"description": "The number of adaptivity steps to perform using the initial conditions"}>>> = 2
    refine_fraction<<<{"description": "The fraction of elements or error to refine. Should be between 0 and 1."}>>> = 0.7
    coarsen_fraction<<<{"description": "The fraction of elements or error to coarsen. Should be between 0 and 1."}>>> = 0.1
    max_h_level<<<{"description": "Maximum number of times a single element can be refined. If 0 then infinite."}>>> = 2
  []
[]

[Outputs<<<{"href": "../../syntax/Outputs/index.html"}>>>]
  exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
  checkpoint<<<{"description": "Create checkpoint files using the default options."}>>> = true
  perf_graph<<<{"description": "Enable printing of the performance graph to the screen (Console)"}>>> = true
[]
(moose/modules/phase_field/examples/ebsd_reconstruction/IN100-111grn.i)

References

  1. 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]
  2. 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]
  3. 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]
  4. 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]