KKS Model (Kim-Kim-Suzuki)
Overview
The KKS model handles multi-phase, multi-component systems with local thermodynamic equilibrium at the interface. It avoids the artificial interface energy problem by partitioning the composition between phases.
Key feature: Local equilibrium constraint ensures equal diffusion potentials between phases at each point.
Mathematical Formulation
Phase Field Equation
\[\tau \frac{\partial \phi}{\partial t} = W^2 \nabla^2 \phi - g'(\phi) + h'(\phi) \cdot \Delta\omega\]
Concentration Equation
\[\frac{\partial c}{\partial t} = \nabla \cdot \left( M(\phi) \nabla \mu \right)\]
Local Equilibrium Constraints
At each point, the phase compositions (cs, cl) satisfy:
\[c = h(\phi) c_s + (1 - h(\phi)) c_l \quad \text{(mass conservation)}\]
\[\mu_s(c_s) = \mu_l(c_l) \quad \text{(equal chemical potential)}\]
Grand Potential Difference
\[\Delta\omega = f_s(c_s) - f_l(c_l) - \mu (c_s - c_l)\]
Quick Example
Basic Usage with Parabolic Free Energy
using PhaseFields
# Create model
model = KKSModel(τ=1.0, W=1.0, m=1.0, M_s=1.0, M_l=10.0)
# Parabolic free energies for testing
f_s = ParabolicFreeEnergy(A=1.0, c_eq=0.2) # Solid
f_l = ParabolicFreeEnergy(A=1.0, c_eq=0.8) # Liquid
# Partition composition at interface
c_avg = 0.5 # Average composition
φ = 0.5 # At interface
c_s, c_l, μ, converged = kks_partition(c_avg, φ, f_s, f_l)
println("Solid: c_s = $c_s, Liquid: c_l = $c_l")
# Grand potential difference
Δω = kks_grand_potential_diff(f_s, f_l, c_s, c_l, μ)With CALPHAD Thermodynamics
using PhaseFields
using OpenCALPHAD
# Load database
db = read_tdb("AgCu.TDB")
# Create CALPHAD-coupled KKS model
T = 1000.0 # Temperature [K]
model, f_s, f_l = create_calphad_kks_model(db, T, "FCC_A1", "LIQUID")
# Partition with real thermodynamics
c_s, c_l, μ, converged = kks_partition(c_avg, φ, f_s, f_l)API Reference
PhaseFields.KKSModel — Type
KKSModelKim-Kim-Suzuki model for solidification with local equilibrium.
The KKS model introduces phase concentrations (cs, cl) that satisfy:
- Mass conservation: h(φ)·cs + (1-h(φ))·cl = c
- Equal diffusion potential: μs(cs, T) = μl(cl, T)
Evolution equations:
Phase field: τ·∂φ/∂t = W²∇²φ - g'(φ) + h'(φ)·Δω
Concentration: ∂c/∂t = ∇·(M(φ)·∇μ)where Δω = fs(cs) - fl(cl) - μ·(cs - cl) is the grand potential difference.
Fields
τ::Float64: Relaxation time [s]W::Float64: Interface width parameter [m]m::Float64: Driving force scale (typically 1.0 for KKS)M_s::Float64: Mobility in solid phase [m²/(J·s)]M_l::Float64: Mobility in liquid phase [m²/(J·s)]
Example
model = KKSModel(τ=1.0, W=1.0, m=1.0, M_s=1.0, M_l=10.0)Reference
Kim, Kim, Suzuki, Phys. Rev. E 60 (1999) 7186
PhaseFields.ParabolicFreeEnergy — Type
ParabolicFreeEnergySimple parabolic free energy for testing KKS model.
f(c) = A·(c - c_eq)²Fields
A::Float64: Curvature (must be positive for stability)c_eq::Float64: Equilibrium concentration
PhaseFields.free_energy — Function
free_energy(f::ParabolicFreeEnergy, c)Parabolic free energy density: f(c) = A·(c - c_eq)²
PhaseFields.chemical_potential — Function
chemical_potential(f::ParabolicFreeEnergy, c)Chemical potential (first derivative): df/dc = 2A·(c - c_eq)
PhaseFields.d2f_dc2 — Function
d2f_dc2(f::ParabolicFreeEnergy, c)Second derivative of free energy: d²f/dc² = 2A
PhaseFields.kks_partition — Function
kks_partition(c_avg, φ, f_s, f_l; maxiter=50, tol=1e-12, φ_min=1e-8)Solve the KKS local equilibrium problem using Newton iteration.
Given average concentration c and phase fraction φ, find phase concentrations (cs, cl) that satisfy:
- h(φ)·cs + (1-h(φ))·cl = c (mass conservation)
- μs(cs) = μl(cl) (equal chemical potential)
Arguments
c_avg: Average concentration at this pointφ: Phase field value (0=liquid, 1=solid)f_s: Solid phase free energy (ParabolicFreeEnergy or similar)f_l: Liquid phase free energy
Keyword Arguments
maxiter: Maximum Newton iterations (default: 50)tol: Convergence tolerance (default: 1e-12)φ_min: Minimum φ for Newton iteration (default: 1e-8)
Returns
(c_s, c_l, μ, converged): Phase concentrations, chemical potential, convergence flag
Example
f_s = ParabolicFreeEnergy(A=1000.0, c_eq=0.1) # Solid equilibrium at c=0.1
f_l = ParabolicFreeEnergy(A=1000.0, c_eq=0.9) # Liquid equilibrium at c=0.9
c_s, c_l, μ, converged = kks_partition(0.5, 0.5, f_s, f_l)PhaseFields.kks_grand_potential_diff — Function
kks_grand_potential_diff(f_s, f_l, c_s, c_l, μ)Compute the grand potential difference (driving force for phase field).
Δω = f_s(c_s) - f_l(c_l) - μ·(c_s - c_l)Arguments
f_s, f_l: Solid and liquid free energy functionsc_s, c_l: Phase concentrationsμ: Chemical potential (equal in both phases)
Returns
- Grand potential difference Δω [J/mol]
PhaseFields.kks_phase_rhs — Function
kks_phase_rhs(model, φ, ∇²φ, Δω)Compute the right-hand side of the KKS phase field equation.
τ·∂φ/∂t = W²∇²φ - g'(φ) + m·h'(φ)·Δω
⟹ ∂φ/∂t = [W²∇²φ - g'(φ) + m·h'(φ)·Δω] / τArguments
model: KKSModel parametersφ: Phase field value∇²φ: Laplacian of phase fieldΔω: Grand potential difference from kksgrandpotential_diff
Returns
- Time derivative ∂φ/∂t
PhaseFields.kks_concentration_rhs — Function
kks_concentration_rhs(model, φ, ∇²μ)Compute the right-hand side of the KKS concentration equation.
∂c/∂t = ∇·(M(φ)·∇μ) ≈ M(φ)·∇²μ (for constant mobility in each phase)Arguments
model: KKSModel parametersφ: Phase field value (for mobility interpolation)∇²μ: Laplacian of chemical potential
Returns
- Time derivative ∂c/∂t
PhaseFields.kks_mobility — Function
kks_mobility(model, φ)Compute effective mobility using interpolation.
M(φ) = h(φ)·M_s + (1-h(φ))·M_lReturns
- Effective mobility at this phase field value
PhaseFields.kks_interface_width — Function
kks_interface_width(model)Get the interface width parameter W.
See Also
- 301kkssolidification.jl - KKS solidification simulation
- 381calphadcoupling_demo.jl - KKS with CALPHAD thermodynamics