- filenameThe name of the file containing the EBSD data
C++ Type:FileName
Unit:(no unit assumed)
Controllable:No
Description:The name of the file containing the EBSD data
EBSDMeshGenerator
Mesh generated from a specified DREAM.3D EBSD data file.
The mesh is generated from the EBSD information, to get an optimal reconstruction of the data. This is accomplished in the mesh block using the EBSDMeshGenerator type. The same data file used with the EBSD reader is used in the EBSDReader UserObject. The mesh is created with one node per data point in the EBSD data file. If you wish to use mesh adaptivity and allow the mesh to get coarser during the simulation, the "uniform_refine" parameter is used to set how many times the mesh can be coarsened. For this to work the number of elements along _each_ dimension has to be divisible by where is the value of the "uniform_refine" parameter.
Contrary to other mesh objects the "uniform_refine" parameter will not affect the resolution of the final mesh. It sets the levels of coarsening that can be applied to the EBSD data.
Input Parameters
- num_cores_for_partitionNumber of cores for partitioning the graph (dafaults to the number of MPI ranks)
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:Number of cores for partitioning the graph (dafaults to the number of MPI ranks)
- pre_refine0Number 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')
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
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')
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Unit:(no unit assumed)
Controllable:No
Description:Adds user-defined labels for accessing object parameters via control logic.
- enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Set the enabled status of the MooseObject.
- save_with_nameKeep the mesh from this mesh generator in memory with the name specified
C++ Type:std::string
Unit:(no unit assumed)
Controllable:No
Description:Keep the mesh from this mesh generator in memory with the name specified
Advanced Parameters
- nemesisFalseWhether or not to output the mesh file in the nemesisformat (only if output = true)
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether or not to output the mesh file in the nemesisformat (only if output = true)
- outputFalseWhether or not to output the mesh file after generating the mesh
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether or not to output the mesh file after generating the mesh
- show_infoFalseWhether or not to show mesh info after generating the mesh (bounding box, element types, sidesets, nodesets, subdomains, etc)
Default:False
C++ Type:bool
Unit:(no unit assumed)
Controllable:No
Description:Whether or not to show mesh info after generating the mesh (bounding box, element types, sidesets, nodesets, subdomains, etc)
Debugging Parameters
Input Files
- (modules/phase_field/test/tests/reconstruction/1phase_evolution.i)
- (modules/phase_field/test/tests/reconstruction/2phase_reconstruction3.i)
- (modules/phase_field/examples/ebsd_reconstruction/IN100-111grn.i)
- (modules/phase_field/test/tests/reconstruction/1phase_reconstruction.i)
- (modules/phase_field/test/tests/GBType/GB_Type_Phase1.i)
- (modules/phase_field/test/tests/grain_tracker_test/grain_tracker_ebsd.i)
- (modules/phase_field/test/tests/reconstruction/EulerAngleVariables2RGBAux.i)
- (modules/phase_field/test/tests/grain_tracker_test/split_grain.i)
- (modules/phase_field/test/tests/EBSDMeshGenerator/test.i)
- (modules/phase_field/test/tests/reconstruction/2phase_reconstruction.i)
- (modules/phase_field/test/tests/reconstruction/euler2rgb_non_uniform_orientation.i)
- (modules/phase_field/test/tests/reconstruction/euler2rgb_no_grain_region.i)
- (modules/phase_field/test/tests/reconstruction/2phase_reconstruction2.i)
- (modules/combined/examples/phase_field-mechanics/EBSD_reconstruction_grain_growth_mech.i)
- (modules/phase_field/test/tests/reconstruction/2phase_reconstruction4.i)
uniform_refine
Default:0
C++ Type:unsigned int
Unit:(no unit assumed)
Controllable:No
Description:Specify the level of uniform refinement applied to the initial mesh
uniform_refine
uniform_refine
(modules/phase_field/test/tests/reconstruction/1phase_evolution.i)
#
# In this test we set the initial condition of a set of order parameters
# by pulling out the grain data from given EBSD data file ignoring the phase completely.
#
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = IN100_001_28x28_Marmot.txt
[]
[]
[GlobalParams]
op_num = 5
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
polycrystal_ic_uo = ebsd
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[Kernels]
[PolycrystalKernel]
[]
[]
[AuxVariables]
[feature]
family = MONOMIAL
order = CONSTANT
[]
[bnds]
[]
[]
[AuxKernels]
[feature]
type = EBSDReaderAvgDataAux
variable = feature
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
data_name = feature_id
execute_on = 'initial timestep_end'
[]
[bnds]
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[]
[]
[Materials]
[CuGrGr]
# Material properties
type = GBEvolution # Quantitative material properties for copper grain growth. Dimensions are nm and ns
block = 0 # Block ID (only one block in this problem)
GBmob0 = 2.5e-6 #Mobility prefactor for Cu from schonfelder1997molecular bibtex entry
GBenergy = 0.708 # GB energy in J/m^2
Q = 0.23 #Activation energy for grain growth from Schonfelder 1997
T = 500 # K #Constant temperature of the simulation (for mobility calculation)
wGB = 1 # nm #Width of the diffuse GB
#outputs = exodus
length_scale = 1e-06
time_scale = 1e-6
[]
[]
[Executioner]
type = Transient
num_steps = 2
dt = 10
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/reconstruction/2phase_reconstruction3.i)
#
# In this test , which is set up similarly to 2phase_reconstruction_test2.i
# we demonstrate that the feature numbers in the EBSD file can be chosen arbitrarily.
# There is no need for then to start at a certain index or even to be contiguous!
# The EBSDReaderPointDataAux AuxKernel outputs the original feature IDs (grain numbers)
#
[Problem]
type = FEProblem
solve = false
kernel_coverage_check = false
[]
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = Renumbered.txt
[]
[]
[GlobalParams]
op_num = 2
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
phase = 1
output_adjacency_matrix = true
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[AuxVariables]
[GRAIN]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[grain_aux]
type = EBSDReaderPointDataAux
variable = GRAIN
ebsd_reader = ebsd_reader
data_name = 'feature_id'
execute_on = 'initial'
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[Executioner]
type = Transient
num_steps = 0
[]
[Outputs]
exodus = true
[]
(modules/phase_field/examples/ebsd_reconstruction/IN100-111grn.i)
[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
[]
(modules/phase_field/test/tests/reconstruction/1phase_reconstruction.i)
#
# In this test we set the initial condition of a set of order parameters
# by pulling out the grain data from given EBSD data file ignoring the phase completely.
#
[Problem]
type = FEProblem
solve = false
kernel_coverage_check = false
[]
# The following sections are extracted in the documentation in
# moose/docs/content/modules/phase_field/ICs/EBSD.md
[Mesh]
# Create a mesh representing the EBSD data
[ebsd_mesh]
type = EBSDMeshGenerator
filename = IN100_001_28x28_Marmot.txt
[]
[]
[GlobalParams]
# Define the number and names of the order parameters used to represent the grains
op_num = 4
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
# Read in the EBSD data. Uses the filename given in the mesh block.
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
# For displaying HALO fields
compute_halo_maps = true
# Link in the ebsd userobject here so that grain tracker can extract info from it
polycrystal_ic_uo = ebsd
[]
[]
[Variables]
[PolycrystalVariables]
# Create all the order parameters
order = FIRST
family = LAGRANGE
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
# Uses the data from the user object 'ebsd' to initialize the variables for all the order parameters.
polycrystal_ic_uo = ebsd
[]
[]
[]
#ENDDOC - End of the file section that is included in the documentation. Do not change this line!
[GlobalParams]
execute_on = 'initial'
family = MONOMIAL
order = CONSTANT
[]
[AuxVariables]
[PHI1]
[]
[PHI]
[]
[PHI2]
[]
[GRAIN]
[]
[unique_grains]
[]
[var_indices]
[]
[halo0]
[]
[halo1]
[]
[halo2]
[]
[halo3]
[]
[]
[AuxKernels]
[phi1_aux]
type = EBSDReaderPointDataAux
variable = PHI1
ebsd_reader = ebsd_reader
data_name = 'phi1'
[]
[phi_aux]
type = EBSDReaderPointDataAux
variable = PHI
ebsd_reader = ebsd_reader
data_name = 'phi'
[]
[phi2_aux]
type = EBSDReaderPointDataAux
variable = PHI2
ebsd_reader = ebsd_reader
data_name = 'phi2'
[]
[grain_aux]
type = EBSDReaderPointDataAux
variable = GRAIN
ebsd_reader = ebsd_reader
data_name = 'feature_id'
[]
[unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
flood_counter = grain_tracker
field_display = UNIQUE_REGION
[]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
[]
[halo0]
type = FeatureFloodCountAux
variable = halo0
map_index = 0
field_display = HALOS
flood_counter = grain_tracker
[]
[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
[]
[]
[Executioner]
type = Steady
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/GBType/GB_Type_Phase1.i)
# MOOSE input file
# Written by Pierre-Clement Simon - Idaho National Laboratory
#
# Project:
# TRISO fuel fission gas transport: Silver diffusion in silicon carbide
#
# Published with:
# ---
#
# Phase Field Model: Isotropic diffusion equation
# type: Transient
# Grain structure: Single grain
# BCs: Fixed value on the right, flux on the left
#
#
# Info:
# - Input file used to generate polycrystals for SiC
#
# Updates from previous file:
# -
#
# Units
# length: --
# time: --
# energy: --
# quantity: --
# This simulation predicts GB migration of a 2D copper polycrystal with 15 grains
# Mesh adaptivity (new system) 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
# We are not using the GrainTracker in this example so the number
# of order paramaters must match the number of grains.
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
# Two Parallel Grains
filename = 'EBSD_ThreeGrains.txt'
[]
[]
[GlobalParams]
# Parameters used by several kernels that are defined globally to simplify input file
op_num = 6 # Number of grains
var_name_base = gr # Base name of grains
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
enable_var_coloring = true
# output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
threshold = 0.001
connecting_threshold = 0.008
compute_var_to_feature_map = true
compute_halo_maps = true # For displaying HALO fields
remap_grains = true
polycrystal_ic_uo = ebsd
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
# Variable block, where all variables in the simulation are declared
[./PolycrystalVariables]
# Custom action that created all of the grain variables and sets their initial condition
[../]
[]
[AuxVariables]
# Dependent variables
[./bnds]
# Variable used to visualize the grain boundaries in the simulation
[../]
[./unique_grains]
order = CONSTANT
family = MONOMIAL
[../]
[./aphi1]
order = CONSTANT
family = MONOMIAL
[../]
[./bPhi]
order = CONSTANT
family = MONOMIAL
[../]
[./cphi2]
order = CONSTANT
family = MONOMIAL
[../]
[./ebsd_numbers]
order = CONSTANT
family = MONOMIAL
[../]
[]
[Kernels]
# Kernel block, where the kernels defining the residual equations are set up.
[./PolycrystalKernel]
# Custom action creating all necessary kernels for grain growth. All input parameters are up in GlobalParams
[../]
[]
[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'
[../]
# generate the unique ID from grain_tracker
[./unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = UNIQUE_REGION
[../]
# The phi will output the Euler angle from EBSD data, and the data structure
# will change with the guide from grain_tracker
[./aphi1]
type = OutputEulerAngles
variable = aphi1
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'phi1'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[./bPhi]
type = OutputEulerAngles
variable = bPhi
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'Phi'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[./cphi2]
type = OutputEulerAngles
variable = cphi2
euler_angle_provider = ebsd_reader
grain_tracker = grain_tracker
output_euler_angle = 'phi2'
execute_on = 'INITIAL TIMESTEP_END'
[../]
# Import the unique grain ID from ebsd data, and the data structure
# will change with the guide from grain_tracker
[ebsd_numbers]
type = EBSDReaderAvgDataAux
data_name = feature_id
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
variable = ebsd_numbers
execute_on = 'initial timestep_end'
[../]
[]
[BCs]
# Boundary Condition block
[./Periodic]
[./top_bottom]
auto_direction = 'x y' # Makes problem periodic in the x and y directions
[../]
[../]
[]
[Materials]
[./CuGrGr]
# Material properties
type = GBEvolution # Quantitative material properties for copper grain growth. Dimensions are nm and ns
GBmob0 = 2.5e-6 # Mobility prefactor for Cu from schonfelder1997molecular bibtex entry
GBenergy = 0.708 # GB energy for Cu from schonfelder1997molecular bibtex entry
Q = 0.23 # Activation energy for grain growth from Schonfelder 1997
T = 450 # Constant temperature of the simulation (for mobility calculation)
wGB = 6 # Width of the diffuse GB
[../]
[./GB_type]
# The new developed Miso Bnds Aux Kernel
type = ComputeGBMisorientationType
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
output_properties = 'gb_type'
outputs = exodus
[../]
[]
[Postprocessors]
# Scalar postprocessors
[./dt]
# Outputs the current time step
type = TimestepSize
[../]
[n_elements]
type = NumElements
execute_on = 'initial timestep_end'
[]
[n_nodes]
type = NumNodes
execute_on = 'initial timestep_end'
[]
[DOFs]
type = NumDOFs
[]
[]
[Adaptivity]
initial_steps = 1
max_h_level = 1
marker = combined
[./Indicators]
[./error]
type = GradientJumpIndicator
variable = bnds
[../]
[../]
[./Markers]
[./bound_adapt]
type = ValueThresholdMarker
third_state = DO_NOTHING
coarsen = 0.999 #1.0
refine = 0.95 #0.95
variable = bnds
invert = true
[../]
[./errorfrac]
type = ErrorFractionMarker
coarsen = 0.1
indicator = error
refine = 0.7
[../]
[./combined]
type = ComboMarker
markers = 'bound_adapt errorfrac'
[../]
[../]
[]
[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'
petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_strong_threshold'
petsc_options_value = ' hypre boomeramg 0.7'
l_max_its = 30 # Max number of linear iterations
l_tol = 1e-4 # Relative tolerance for linear solves
nl_max_its = 40 # Max number of nonlinear iterations
nl_abs_tol = 1e-11 # Relative tolerance for nonlienar solves
nl_rel_tol = 1e-10 # Absolute tolerance for nonlienar solves
[TimeStepper]
type = IterationAdaptiveDT
cutback_factor = 0.9
dt = 1
growth_factor = 1.1
optimal_iterations = 7
[]
start_time = 0.0
num_steps = 2
[]
[Outputs]
perf_graph = true
exodus = true
[./console]
type = Console
max_rows = 10
[../]
[]
(modules/phase_field/test/tests/grain_tracker_test/grain_tracker_ebsd.i)
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = 'ebsd_9.txt'
[]
[]
[GlobalParams]
op_num = 4
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
threshold = 0.2
connecting_threshold = 0.08
flood_entity_type = ELEMENTAL
compute_halo_maps = true # For displaying HALO fields
polycrystal_ic_uo = ebsd
execute_on = 'initial timestep_end'
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[AuxVariables]
[bnds]
[]
[unique_grains]
family = MONOMIAL
order = CONSTANT
[]
[var_indices]
family = MONOMIAL
order = CONSTANT
[]
[ebsd_grains]
family = MONOMIAL
order = CONSTANT
[]
[phi1]
family = MONOMIAL
order = CONSTANT
[]
[halo0]
order = CONSTANT
family = MONOMIAL
[]
[halo1]
order = CONSTANT
family = MONOMIAL
[]
[halo2]
order = CONSTANT
family = MONOMIAL
[]
[halo3]
order = CONSTANT
family = MONOMIAL
[]
[]
[Kernels]
[PolycrystalKernel]
[]
[]
[AuxKernels]
[BndsCalc]
type = BndsCalcAux
variable = bnds
execute_on = timestep_end
[]
[unique_grains]
type = FeatureFloodCountAux
variable = unique_grains
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = UNIQUE_REGION
[]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
execute_on = 'initial timestep_end'
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
[]
[grain_aux]
type = EBSDReaderPointDataAux
variable = ebsd_grains
ebsd_reader = ebsd_reader
data_name = 'feature_id'
execute_on = 'initial timestep_end'
[]
[phi1]
type = OutputEulerAngles
euler_angle_provider = ebsd_reader
output_euler_angle = phi1
grain_tracker = grain_tracker
variable = phi1
[]
[halo0]
type = FeatureFloodCountAux
variable = halo0
map_index = 0
field_display = HALOS
flood_counter = grain_tracker
[]
[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
[]
[]
[Materials]
[CuGrGr]
type = GBEvolution
T = 500 #K
wGB = 0.75 #micron
length_scale = 1.0e-6
time_scale = 1.0e-4
GBmob0 = 2.5e-6
Q = 0.23
GBenergy = 0.708
molar_volume = 7.11e-6
[]
[]
[Postprocessors]
[n_nodes]
type = NumNodes
execute_on = timestep_end
[]
[DOFs]
type = NumDOFs
[]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type -ksp_gmres_restart '
'-pc_hypre_boomeramg_strong_threshold'
petsc_options_value = 'hypre boomeramg 31 0.7'
l_tol = 1.0e-4
l_max_its = 20
nl_rel_tol = 1.0e-9
nl_max_its = 20
start_time = 0.0
num_steps = 1
dt = 0.05
[]
[Outputs]
execute_on = 'initial'
exodus = true
perf_graph = true
[]
(modules/phase_field/test/tests/reconstruction/EulerAngleVariables2RGBAux.i)
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = IN100_001_28x28_Clean_Marmot.txt
[]
[]
[UserObjects]
[ebsd]
type = EBSDReader
[]
[]
[AuxVariables]
[phi1]
family = MONOMIAL
order = CONSTANT
[]
[phi]
family = MONOMIAL
order = CONSTANT
[]
[phi2]
family = MONOMIAL
order = CONSTANT
[]
[phase]
order = CONSTANT
family = MONOMIAL
[]
[symmetry]
order = CONSTANT
family = MONOMIAL
[]
[rgb]
order = CONSTANT
family = MONOMIAL
[]
[]
[AuxKernels]
[phi1_aux]
type = EBSDReaderPointDataAux
variable = phi1
ebsd_reader = ebsd
data_name = phi1
execute_on = initial
[]
[phi_aux]
type = EBSDReaderPointDataAux
variable = phi
ebsd_reader = ebsd
data_name = phi
execute_on = initial
[]
[phi2_aux]
type = EBSDReaderPointDataAux
variable = phi2
ebsd_reader = ebsd
data_name = phi2
execute_on = initial
[]
[phase_aux]
type = EBSDReaderPointDataAux
variable = phase
ebsd_reader = ebsd
data_name = phase
execute_on = initial
[]
[symmetry_aux]
type = EBSDReaderPointDataAux
variable = symmetry
ebsd_reader = ebsd
data_name = symmetry
execute_on = initial
[]
[rgb]
type = EulerAngleVariables2RGBAux
variable = rgb
phi1 = phi1
phi = phi
phi2 = phi2
phase = phase
symmetry = symmetry
execute_on = initial
[]
[]
[Executioner]
type = Steady
[]
[Problem]
solve = false
kernel_coverage_check = false
[]
[Outputs]
exodus = true
show = 'rgb'
[]
(modules/phase_field/test/tests/grain_tracker_test/split_grain.i)
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = EBSD_split_grain.txt
[]
[]
[GlobalParams]
op_num = 4
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
enable_var_coloring = true
output_adjacency_matrix = 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]
order = CONSTANT
family = MONOMIAL
[]
[ghost_elements]
order = CONSTANT
family = MONOMIAL
[]
[halos]
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]
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]
[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 = 2
[TimeStepper]
type = IterationAdaptiveDT
cutback_factor = 0.9
dt = 10.0
growth_factor = 1.1
optimal_iterations = 7
[]
[]
[Outputs]
exodus = true
csv = true
perf_graph = true
[]
(modules/phase_field/test/tests/EBSDMeshGenerator/test.i)
[Mesh]
[ebsd]
type = EBSDMeshGenerator
[]
[]
(modules/phase_field/test/tests/reconstruction/2phase_reconstruction.i)
#
# In this test we set the initial condition of two variables
# based on solely the phase information in a given EBSD data file,
# ignoring the feature IDs entirely
#
[Problem]
type = FEProblem
solve = false
kernel_coverage_check = false
[]
# The following sections are extracted in the documentation in
# moose/docs/content/modules/phase_field/ICs/EBSD.md
[Mesh]
# Create a mesh representing the EBSD data
[ebsd_mesh]
type = EBSDMeshGenerator
filename = 'Ti_2Phase_28x28_ebsd.txt'
[]
[]
[UserObjects]
[ebsd]
# Read in the EBSD data. Uses the filename given in the mesh block.
type = EBSDReader
[]
[]
[Variables]
# Creates the two variables being initialized
[c1]
[]
[c2]
[]
[]
[ICs]
[phase1_recon]
# Initializes the variable info from the ebsd data
type = ReconPhaseVarIC
ebsd_reader = ebsd
phase = 1
variable = c1
[]
[phase2_recon]
type = ReconPhaseVarIC
ebsd_reader = ebsd
phase = 2
variable = c2
[]
[]
#ENDDOC - End of the file section that is included in the documentation. Do not change this line!
[AuxVariables]
[PHI1]
family = MONOMIAL
order = CONSTANT
[]
[PHI]
family = MONOMIAL
order = CONSTANT
[]
[APHI2]
family = MONOMIAL
order = CONSTANT
[]
[PHI2]
family = MONOMIAL
order = CONSTANT
[]
[PHASE]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[phi1_aux]
type = EBSDReaderPointDataAux
variable = PHI1
ebsd_reader = ebsd
data_name = 'phi1'
execute_on = 'initial'
[]
[phi_aux]
type = EBSDReaderPointDataAux
variable = PHI
ebsd_reader = ebsd
data_name = 'phi'
execute_on = 'initial'
[]
[phi2_aux]
type = EBSDReaderPointDataAux
variable = PHI2
ebsd_reader = ebsd
data_name = 'phi2'
execute_on = 'initial'
[]
[phase_aux]
type = EBSDReaderPointDataAux
variable = PHASE
ebsd_reader = ebsd
data_name = 'phase'
execute_on = 'initial'
[]
[]
[Executioner]
type = Transient
num_steps = 0
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/reconstruction/euler2rgb_non_uniform_orientation.i)
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = ebsd_scan.txt
[]
[]
[GlobalParams]
op_num = 10
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
bins = 40
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = jp
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]
[]
[]
[Kernels]
[PolycrystalKernel]
[]
[]
[AuxKernels]
[BndsCalc]
type = BndsCalcAux
variable = bnds
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 = Newton
petsc_options_iname = '-pc_type -pc_hypre_type -pc_hypre_boomeramg_strong_threshold'
petsc_options_value = 'hypre boomeramg 0.7'
l_tol = 1.0e-6
l_max_its = 100
nl_max_its = 20
nl_rel_tol = 1.0e-8
start_time = 0.0
num_steps = 0
[]
[Outputs]
exodus = true
perf_graph = true
[]
(modules/phase_field/test/tests/reconstruction/euler2rgb_no_grain_region.i)
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = ebsd_small.txt
[]
[]
[GlobalParams]
op_num = 8
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
execute_on = initial
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
phase = 2
output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
polycrystal_ic_uo = ebsd
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[void_phase]
type = ReconPhaseVarIC
variable = c
ebsd_reader = ebsd_reader
phase = 1
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[AuxVariables]
# active = 'c bnds'
[c]
[]
[bnds]
[]
[ebsd_numbers]
family = MONOMIAL
order = CONSTANT
[]
# Note: Not active
[unique_grains]
family = MONOMIAL
order = CONSTANT
[]
[var_indices]
family = MONOMIAL
order = CONSTANT
[]
[]
[Kernels]
[PolycrystalKernel]
c = c
[]
[]
[AuxKernels]
# active = 'BndsCalc'
[BndsCalc]
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[]
[ebsd_numbers]
type = EBSDReaderAvgDataAux
data_name = feature_id
ebsd_reader = ebsd_reader
grain_tracker = grain_tracker
variable = ebsd_numbers
phase = 2
execute_on = 'initial timestep_end'
[]
# Note: Not active
[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'
[]
[]
[Modules]
[PhaseField]
[EulerAngles2RGB]
crystal_structure = cubic
grain_tracker = grain_tracker
euler_angle_provider = ebsd_reader
no_grain_color = '.1 .1 .1'
phase = 2
[]
[]
[]
[Materials]
[bulk]
type = GBEvolution
block = 0
T = 2273
wGB = 10.0
GBenergy = 1.58
GBmob0 = 9.2124e-9
Q = 2.77
length_scale = 1.0e-6
time_scale = 60.0
[]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = PJFNK
petsc_options_iname = '-pc_type -sub_pc_type -pc_asm_overlap -ksp_grmres_restart '
petsc_options_value = ' asm lu 1 21'
start_time = 0.0
dt = 0.2
num_steps = 1
[]
[Outputs]
csv = true
exodus = true
execute_on = 'INITIAL TIMESTEP_END'
perf_graph = true
[]
(modules/phase_field/test/tests/reconstruction/2phase_reconstruction2.i)
#
# In this test we set the initial condition of a set of order parameters
# by pulling out the only grains from given EBSD data file that belong to a specified phase
#
[Problem]
type = FEProblem
solve = false
kernel_coverage_check = false
[]
# The following sections are extracted in the documentation in
# moose/docs/content/modules/phase_field/ICs/EBSD.md
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = Ti_2Phase_28x28_ebsd.txt
[]
[]
[GlobalParams]
op_num = 2
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
phase = 1
output_adjacency_matrix = true
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
# select only data for phase 1 from the EBSD file
polycrystal_ic_uo = ebsd
[]
[]
[]
#ENDDOC - End of the file section that is included in the documentation. Do not change this line!
[Executioner]
type = Transient
num_steps = 0
[]
[Outputs]
exodus = true
[]
(modules/combined/examples/phase_field-mechanics/EBSD_reconstruction_grain_growth_mech.i)
# This example reconstructs the grain structure from an EBSD data file
# Then, an isotropic grain model is run with linear elasticity and an anisotropic
# elasticity tensor that uses the measured EBSD angles.
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
uniform_refine = 2 #Mesh can go two levels coarser than the EBSD grid
filename = IN100_128x128.txt
[]
[]
[GlobalParams]
op_num = 8
var_name_base = gr
displacements = 'disp_x disp_y'
[]
[Variables]
[PolycrystalVariables] #Polycrystal variable generation (30 order parameters)
[]
[disp_x]
[]
[disp_y]
[]
[]
[AuxVariables]
[bnds]
[]
[gt_indices]
order = CONSTANT
family = MONOMIAL
[]
[unique_grains]
order = CONSTANT
family = MONOMIAL
[]
[vonmises_stress]
order = CONSTANT
family = MONOMIAL
[]
[C1111]
order = CONSTANT
family = MONOMIAL
[]
[phi1]
order = CONSTANT
family = MONOMIAL
[]
[Phi]
order = CONSTANT
family = MONOMIAL
[]
[phi2]
order = CONSTANT
family = MONOMIAL
[]
[EBSD_grain]
family = MONOMIAL
order = CONSTANT
[]
[]
[ICs]
[PolycrystalICs]
[ReconVarIC]
ebsd_reader = ebsd
coloring_algorithm = bt
[]
[]
[]
[Kernels]
[PolycrystalKernel]
[]
[PolycrystalElasticDrivingForce]
[]
[TensorMechanics]
[]
[]
[AuxKernels]
[BndsCalc]
type = BndsCalcAux
variable = bnds
execute_on = 'initial timestep_end'
[]
[gt_indices]
type = FeatureFloodCountAux
variable = gt_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
[]
[C1111]
type = RankFourAux
variable = C1111
rank_four_tensor = elasticity_tensor
index_l = 0
index_j = 0
index_k = 0
index_i = 0
execute_on = timestep_end
[]
[vonmises_stress]
type = RankTwoScalarAux
variable = vonmises_stress
rank_two_tensor = stress
scalar_type = VonMisesStress
execute_on = timestep_end
[]
[phi1]
type = OutputEulerAngles
variable = phi1
euler_angle_provider = ebsd
grain_tracker = grain_tracker
output_euler_angle = 'phi1'
execute_on = 'initial'
[]
[Phi]
type = OutputEulerAngles
variable = Phi
euler_angle_provider = ebsd
grain_tracker = grain_tracker
output_euler_angle = 'Phi'
execute_on = 'initial'
[]
[phi2]
type = OutputEulerAngles
variable = phi2
euler_angle_provider = ebsd
grain_tracker = grain_tracker
output_euler_angle = 'phi2'
execute_on = 'initial'
[]
[grain_aux]
type = EBSDReaderPointDataAux
variable = EBSD_grain
ebsd_reader = ebsd
data_name = 'feature_id'
execute_on = 'initial'
[]
[]
[BCs]
[top_displacement]
type = DirichletBC
variable = disp_y
boundary = top
value = -2.0
[]
[x_anchor]
type = DirichletBC
variable = disp_x
boundary = 'left right'
value = 0.0
[]
[y_anchor]
type = DirichletBC
variable = disp_y
boundary = bottom
value = 0.0
[]
[]
[Modules]
[PhaseField]
[EulerAngles2RGB]
crystal_structure = cubic
euler_angle_provider = ebsd
grain_tracker = grain_tracker
[]
[]
[]
[Materials]
[Copper]
# T = 500 # K
type = GBEvolution
block = 0
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
[]
[ElasticityTensor]
type = ComputePolycrystalElasticityTensor
grain_tracker = grain_tracker
[]
[strain]
type = ComputeSmallStrain
block = 0
displacements = 'disp_x disp_y'
[]
[stress]
type = ComputeLinearElasticStress
block = 0
[]
[]
[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
[]
[]
[UserObjects]
[ebsd]
type = EBSDReader
[]
[grain_tracker]
type = GrainTrackerElasticity
compute_var_to_feature_map = true
ebsd_reader = ebsd
fill_method = symmetric9
C_ijkl = '1.27e5 0.708e5 0.708e5 1.27e5 0.708e5 1.27e5 0.7355e5 0.7355e5 0.7355e5'
euler_angle_provider = ebsd
[]
[]
[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
dt = 10
[Adaptivity]
initial_adaptivity = 0
refine_fraction = 0.7
coarsen_fraction = 0.1
max_h_level = 2
[]
[TimeStepper]
type = IterationAdaptiveDT
cutback_factor = 0.9
dt = 10.0
growth_factor = 1.1
optimal_iterations = 7
[]
[]
[Outputs]
csv = true
exodus = true
[]
(modules/phase_field/test/tests/reconstruction/2phase_reconstruction4.i)
#
# In this test we set the initial condition of a set of order parameters
# by pulling out the only grains from given EBSD data file that belong to a specified phase
#
[Problem]
type = FEProblem
solve = false
kernel_coverage_check = false
[]
[Mesh]
[ebsd_mesh]
type = EBSDMeshGenerator
filename = ebsd_40x40_2_phase.txt
[]
[]
[GlobalParams]
op_num = 6
var_name_base = gr
[]
[UserObjects]
[ebsd_reader]
type = EBSDReader
[]
[ebsd]
type = PolycrystalEBSD
coloring_algorithm = bt
ebsd_reader = ebsd_reader
phase = 2
output_adjacency_matrix = true
[]
[grain_tracker]
type = GrainTracker
polycrystal_ic_uo = ebsd
remap_grains = false
[]
[]
[AuxVariables]
[var_indices]
family = MONOMIAL
order = CONSTANT
[]
[]
[AuxKernels]
[var_indices]
type = FeatureFloodCountAux
variable = var_indices
flood_counter = grain_tracker
field_display = VARIABLE_COLORING
[]
[]
[ICs]
[PolycrystalICs]
[PolycrystalColoringIC]
polycrystal_ic_uo = ebsd
[]
[]
[]
[Variables]
[PolycrystalVariables]
[]
[]
[Executioner]
type = Transient
num_steps = 0
[]
[Outputs]
exodus = true
[]