Kim-Kim-Suzuki Example for three or more components
#
# KKS ternary (3 chemical component) system example in the split form
# We track c1 and c2 only, since c1 + c2 + c3 = 1
#
[Mesh<<<{"href": "../../../syntax/Mesh/index.html"}>>>]
type = GeneratedMesh
dim = 2
nx = 150
ny = 15
nz = 0
xmin = -25
xmax = 25
ymin = -2.5
ymax = 2.5
zmin = 0
zmax = 0
elem_type = QUAD4
[]
[AuxVariables<<<{"href": "../../../syntax/AuxVariables/index.html"}>>>]
[./Fglobal]
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
[../]
[]
[Variables<<<{"href": "../../../syntax/Variables/index.html"}>>>]
# order parameter
[./eta]
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
[../]
# solute 1 concentration
[./c1]
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
[../]
# solute 2 concentration
[./c2]
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
[../]
# chemical potential solute 1
[./w1]
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
[../]
# chemical potential solute 2
[./w2]
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
[../]
# Liquid phase solute 1 concentration
[./c1l]
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
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.1
[../]
# Liquid phase solute 2 concentration
[./c2l]
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
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.05
[../]
# Solid phase solute 1 concentration
[./c1s]
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
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.8
[../]
# Solid phase solute 2 concentration
[./c2s]
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
initial_condition<<<{"description": "Specifies a constant initial condition for this variable"}>>> = 0.1
[../]
[]
[Functions<<<{"href": "../../../syntax/Functions/index.html"}>>>]
[./ic_func_eta]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = '0.5*(1.0-tanh((x)/sqrt(2.0)))'
[../]
[./ic_func_c1]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = '0.8*(0.5*(1.0-tanh(x/sqrt(2.0))))^3*(6*(0.5*(1.0-tanh(x/sqrt(2.0))))^2-15*(0.5*(1.0-tanh(x/sqrt(2.0))))+10)+0.1*(1-(0.5*(1.0-tanh(x/sqrt(2.0))))^3*(6*(0.5*(1.0-tanh(x/sqrt(2.0))))^2-15*(0.5*(1.0-tanh(x/sqrt(2.0))))+10))'
[../]
[./ic_func_c2]
type = ParsedFunction<<<{"description": "Function created by parsing a string", "href": "../../../source/functions/MooseParsedFunction.html"}>>>
expression<<<{"description": "The user defined function."}>>> = '0.1*(0.5*(1.0-tanh(x/sqrt(2.0))))^3*(6*(0.5*(1.0-tanh(x/sqrt(2.0))))^2-15*(0.5*(1.0-tanh(x/sqrt(2.0))))+10)+0.05*(1-(0.5*(1.0-tanh(x/sqrt(2.0))))^3*(6*(0.5*(1.0-tanh(x/sqrt(2.0))))^2-15*(0.5*(1.0-tanh(x/sqrt(2.0))))+10))'
[../]
[]
[ICs<<<{"href": "../../../syntax/ICs/index.html"}>>>]
[./eta]
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = eta
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = ic_func_eta
[../]
[./c1]
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = c1
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = ic_func_c1
[../]
[./c2]
variable<<<{"description": "The variable this initial condition is supposed to provide values for."}>>> = c2
type = FunctionIC<<<{"description": "An initial condition that uses a normal function of x, y, z to produce values (and optionally gradients) for a field variable.", "href": "../../../source/ics/FunctionIC.html"}>>>
function<<<{"description": "The initial condition function."}>>> = ic_func_c2
[../]
[]
[Materials<<<{"href": "../../../syntax/Materials/index.html"}>>>]
# Free energy of the liquid
[./fl]
type = DerivativeParsedMaterial<<<{"description": "Parsed Function Material with automatic derivatives.", "href": "../../../source/materials/DerivativeParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = fl
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'c1l c2l'
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '(0.1-c1l)^2+(0.05-c2l)^2'
[../]
# Free energy of the solid
[./fs]
type = DerivativeParsedMaterial<<<{"description": "Parsed Function Material with automatic derivatives.", "href": "../../../source/materials/DerivativeParsedMaterial.html"}>>>
property_name<<<{"description": "Name of the parsed material property"}>>> = fs
coupled_variables<<<{"description": "Vector of variables used in the parsed function"}>>> = 'c1s c2s'
expression<<<{"description": "Parsed function (see FParser) expression for the parsed material"}>>> = '(0.8-c1s)^2+(0.1-c2s)^2'
[../]
# h(eta)
[./h_eta]
type = SwitchingFunctionMaterial<<<{"description": "Helper material to provide $h(\\eta)$ and its derivative in one of two polynomial forms.\nSIMPLE: $3\\eta^2-2\\eta^3$\nHIGH: $\\eta^3(6\\eta^2-15\\eta+10)$", "href": "../../../source/materials/SwitchingFunctionMaterial.html"}>>>
h_order<<<{"description": "Polynomial order of the switching function h(eta)"}>>> = HIGH
eta<<<{"description": "Order parameter variable"}>>> = eta
[../]
# g(eta)
[./g_eta]
type = BarrierFunctionMaterial<<<{"description": "Helper material to provide $g(\\eta)$ and its derivative in a polynomial.\nSIMPLE: $\\eta^2(1-\\eta)^2$\nLOW: $\\eta(1-\\eta)$\nHIGH: $\\eta^2(1-\\eta^2)^2$", "href": "../../../source/materials/BarrierFunctionMaterial.html"}>>>
g_order<<<{"description": "Polynomial order of the barrier function g(eta)"}>>> = SIMPLE
eta<<<{"description": "Order parameter variable"}>>> = eta
[../]
# constant properties
[./constants]
type = GenericConstantMaterial<<<{"description": "Declares material properties based on names and values prescribed by input parameters.", "href": "../../../source/materials/GenericConstantMaterial.html"}>>>
prop_names<<<{"description": "The names of the properties this material will have"}>>> = 'M L eps_sq'
prop_values<<<{"description": "The values associated with the named properties"}>>> = '0.7 0.7 1.0 '
[../]
[]
[Kernels<<<{"href": "../../../syntax/Kernels/index.html"}>>>]
# enforce c1 = (1-h(eta))*c1l + h(eta)*c1s
[./PhaseConc1]
type = KKSPhaseConcentration<<<{"description": "KKS model kernel to enforce the decomposition of concentration into phase concentration $(1-h(\\eta))c_a + h(\\eta)c_b - c = 0$. The non-linear variable of this kernel is $c_b$.", "href": "../../../source/kernels/KKSPhaseConcentration.html"}>>>
ca<<<{"description": "Phase a concentration"}>>> = c1l
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c1s
c<<<{"description": "Real concentration"}>>> = c1
eta<<<{"description": "Phase a/b order parameter"}>>> = eta
[../]
# enforce c2 = (1-h(eta))*c2l + h(eta)*c2s
[./PhaseConc2]
type = KKSPhaseConcentration<<<{"description": "KKS model kernel to enforce the decomposition of concentration into phase concentration $(1-h(\\eta))c_a + h(\\eta)c_b - c = 0$. The non-linear variable of this kernel is $c_b$.", "href": "../../../source/kernels/KKSPhaseConcentration.html"}>>>
ca<<<{"description": "Phase a concentration"}>>> = c2l
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c2s
c<<<{"description": "Real concentration"}>>> = c2
eta<<<{"description": "Phase a/b order parameter"}>>> = eta
[../]
# enforce pointwise equality of chemical potentials
[./ChemPotSolute1]
type = KKSPhaseChemicalPotential<<<{"description": "KKS model kernel to enforce the pointwise equality of phase chemical potentials $dF_a/dc_a = dF_b/dc_b$. The non-linear variable of this kernel is $c_a$.", "href": "../../../source/kernels/KKSPhaseChemicalPotential.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c1l
cb<<<{"description": "Phase b concentration"}>>> = c1s
fa_name<<<{"description": "Base name of the free energy function Fa (f_name in the corresponding derivative function material)"}>>> = fl
fb_name<<<{"description": "Base name of the free energy function Fb (f_name in the corresponding derivative function material)"}>>> = fs
args_a<<<{"description": "Vector of further parameters to Fa (optional, to add in second cross derivatives of Fa)"}>>> = 'c2l'
args_b<<<{"description": "Vector of further parameters to Fb (optional, to add in second cross derivatives of Fb)"}>>> = 'c2s'
[../]
[./ChemPotSolute2]
type = KKSPhaseChemicalPotential<<<{"description": "KKS model kernel to enforce the pointwise equality of phase chemical potentials $dF_a/dc_a = dF_b/dc_b$. The non-linear variable of this kernel is $c_a$.", "href": "../../../source/kernels/KKSPhaseChemicalPotential.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c2l
cb<<<{"description": "Phase b concentration"}>>> = c2s
fa_name<<<{"description": "Base name of the free energy function Fa (f_name in the corresponding derivative function material)"}>>> = fl
fb_name<<<{"description": "Base name of the free energy function Fb (f_name in the corresponding derivative function material)"}>>> = fs
args_a<<<{"description": "Vector of further parameters to Fa (optional, to add in second cross derivatives of Fa)"}>>> = 'c1l'
args_b<<<{"description": "Vector of further parameters to Fb (optional, to add in second cross derivatives of Fb)"}>>> = 'c1s'
[../]
#
# Cahn-Hilliard Equations
#
[./CHBulk1]
type = KKSSplitCHCRes<<<{"description": "KKS model kernel for the split Bulk Cahn-Hilliard term. This kernel operates on the physical concentration 'c' as the non-linear variable", "href": "../../../source/kernels/KKSSplitCHCRes.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c1
ca<<<{"description": "phase concentration corresponding to the non-linear variable of this kernel"}>>> = c1l
fa_name<<<{"description": "Base name of an arbitrary phase free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fl
w<<<{"description": "Chemical potential non-linear helper variable for the split solve"}>>> = w1
args_a<<<{"description": "Vector of additional arguments to Fa"}>>> = 'c2l'
[../]
[./CHBulk2]
type = KKSSplitCHCRes<<<{"description": "KKS model kernel for the split Bulk Cahn-Hilliard term. This kernel operates on the physical concentration 'c' as the non-linear variable", "href": "../../../source/kernels/KKSSplitCHCRes.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = c2
ca<<<{"description": "phase concentration corresponding to the non-linear variable of this kernel"}>>> = c2l
fa_name<<<{"description": "Base name of an arbitrary phase free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fl
w<<<{"description": "Chemical potential non-linear helper variable for the split solve"}>>> = w2
args_a<<<{"description": "Vector of additional arguments to Fa"}>>> = 'c1l'
[../]
[./dc1dt]
type = CoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../../source/kernels/CoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = w1
v<<<{"description": "Coupled variable"}>>> = c1
[../]
[./dc2dt]
type = CoupledTimeDerivative<<<{"description": "Time derivative Kernel that acts on a coupled variable. Weak form: $(\\psi_i, \\frac{\\partial v_h}{\\partial t})$.", "href": "../../../source/kernels/CoupledTimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = w2
v<<<{"description": "Coupled variable"}>>> = c2
[../]
[./w1kernel]
type = SplitCHWRes<<<{"description": "Split formulation Cahn-Hilliard Kernel for the chemical potential variable with a scalar (isotropic) mobility", "href": "../../../source/kernels/SplitCHWRes.html"}>>>
mob_name<<<{"description": "The mobility used with the kernel"}>>> = M
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = w1
[../]
[./w2kernel]
type = SplitCHWRes<<<{"description": "Split formulation Cahn-Hilliard Kernel for the chemical potential variable with a scalar (isotropic) mobility", "href": "../../../source/kernels/SplitCHWRes.html"}>>>
mob_name<<<{"description": "The mobility used with the kernel"}>>> = M
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = w2
[../]
#
# Allen-Cahn Equation
#
[./ACBulkF]
type = KKSACBulkF<<<{"description": "KKS model kernel (part 1 of 2) for the Bulk Allen-Cahn. This includes all terms NOT dependent on chemical potential.", "href": "../../../source/kernels/KKSACBulkF.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = eta
fa_name<<<{"description": "Base name of the free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fl
fb_name<<<{"description": "Base name of the free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fs
w<<<{"description": "Double well height parameter"}>>> = 1.0
coupled_variables<<<{"description": "Vector of nonlinear variable arguments this object depends on"}>>> = 'c1l c1s c2l c2s'
[../]
[./ACBulkC1]
type = KKSACBulkC<<<{"description": "KKS model kernel (part 2 of 2) for the Bulk Allen-Cahn. This includes all terms dependent on chemical potential.", "href": "../../../source/kernels/KKSACBulkC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = eta
ca<<<{"description": "a-phase concentration"}>>> = c1l
cb<<<{"description": "b-phase concentration"}>>> = c1s
fa_name<<<{"description": "Base name of the free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fl
coupled_variables<<<{"description": "Vector of nonlinear variable arguments this object depends on"}>>> = 'c2l'
[../]
[./ACBulkC2]
type = KKSACBulkC<<<{"description": "KKS model kernel (part 2 of 2) for the Bulk Allen-Cahn. This includes all terms dependent on chemical potential.", "href": "../../../source/kernels/KKSACBulkC.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = eta
ca<<<{"description": "a-phase concentration"}>>> = c2l
cb<<<{"description": "b-phase concentration"}>>> = c2s
fa_name<<<{"description": "Base name of the free energy function F (f_base in the corresponding KKSBaseMaterial)"}>>> = fl
coupled_variables<<<{"description": "Vector of nonlinear variable arguments this object depends on"}>>> = 'c1l'
[../]
[./ACInterface]
type = ACInterface<<<{"description": "Gradient energy Allen-Cahn Kernel", "href": "../../../source/kernels/ACInterface.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = eta
kappa_name<<<{"description": "The kappa used with the kernel"}>>> = eps_sq
[../]
[./detadt]
type = TimeDerivative<<<{"description": "The time derivative operator with the weak form of $(\\psi_i, \\frac{\\partial u_h}{\\partial t})$.", "href": "../../../source/kernels/TimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = eta
[../]
[]
[AuxKernels<<<{"href": "../../../syntax/AuxKernels/index.html"}>>>]
[./GlobalFreeEnergy]
variable<<<{"description": "The name of the variable that this object applies to"}>>> = Fglobal
type = KKSGlobalFreeEnergy<<<{"description": "Total free energy in KKS system, including chemical, barrier and gradient terms", "href": "../../../source/auxkernels/KKSGlobalFreeEnergy.html"}>>>
fa_name<<<{"description": "Base name of the free energy function F (f_name in the corresponding derivative function material)"}>>> = fl
fb_name<<<{"description": "Base name of the free energy function F (f_name in the corresponding derivative function material)"}>>> = fs
w<<<{"description": "Double well height parameter"}>>> = 1.0
[../]
[]
[Executioner<<<{"href": "../../../syntax/Executioner/index.html"}>>>]
type = Transient
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm ilu nonzero'
l_max_its = 100
nl_max_its = 100
num_steps = 50
dt = 0.1
[]
#
# Precondition using handcoded off-diagonal terms
#
[Preconditioning<<<{"href": "../../../syntax/Preconditioning/index.html"}>>>]
[./full]
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
[../]
[]
[Outputs<<<{"href": "../../../syntax/Outputs/index.html"}>>>]
exodus<<<{"description": "Output the results using the default settings for Exodus output."}>>> = true
[]
(moose/modules/phase_field/examples/kim-kim-suzuki/kks_example_ternary.i)When additional chemical components are added to the KKS model, a Cahn-Hilliard equation must be added for each additional component. (For components, Cahn-Hilliard equations are required). Each additional Cahn-Hilliard equation requires the kernels:
To enforce the composition and chemical potential constraints, each additional component also requires the kernels
The Allen-Cahn equation is also modified when additional components are added. The residual becomes
where is the number of components. A single KKSACBulkF
kernel is needed as in the 2-component case, and an additional KKSACBulkC
kernel must added for each additional component.