(modules/phase_field/examples/slkks/CrFe.i)
#
# SLKKS two phase example for the BCC and SIGMA phases. The sigma phase contains
# multiple sublattices. Free energy from
# Jacob, Aurelie, Erwin Povoden-Karadeniz, and Ernst Kozeschnik. "Revised thermodynamic
# description of the Fe-Cr system based on an improved sublattice model of the sigma phase."
# Calphad 60 (2018): 16-28.
#
# In this simulation we consider diffusion (Cahn-Hilliard) and phase transformation.
#
# This example requires CrFe_sigma_out_var_0001.csv file, which generated by first
# running the CrFe_sigma.i input file.
[Mesh]
  [gen]
    type = GeneratedMeshGenerator
    dim = 2
    nx = 160
    ny = 1
    nz = 0
    xmin = -25
    xmax = 25
    ymin = -2.5
    ymax = 2.5
    elem_type = QUAD4
  []
[]
[AuxVariables]
  [Fglobal]
    order = CONSTANT
    family = MONOMIAL
  []
[]
[Functions]
  [sigma_cr0]
    type = PiecewiseLinear
    data_file = CrFe_sigma_out_var_0001.csv
    format = columns
    x_index_in_file = 5
    y_index_in_file = 2
    xy_in_file_only = false
  []
  [sigma_cr1]
    type = PiecewiseLinear
    data_file = CrFe_sigma_out_var_0001.csv
    format = columns
    x_index_in_file = 5
    y_index_in_file = 3
    xy_in_file_only = false
  []
  [sigma_cr2]
    type = PiecewiseLinear
    data_file = CrFe_sigma_out_var_0001.csv
    format = columns
    x_index_in_file = 5
    y_index_in_file = 4
    xy_in_file_only = false
  []
[]
[Variables]
  # order parameters
  [eta1]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.5
  []
  [eta2]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.5
  []
  # solute concentration
  [cCr]
    order = FIRST
    family = LAGRANGE
    [InitialCondition]
      type = FunctionIC
      function = '(x+25)/50*0.5+0.1'
    []
  []
  # sublattice concentrations
  [BCC_CR]
    initial_condition = 0.45
  []
  [SIGMA_0CR]
    [InitialCondition]
      type = CoupledValueFunctionIC
      function = sigma_cr0
      v = cCr
      variable = SIGMA_0CR
    []
  []
  [SIGMA_1CR]
    [InitialCondition]
      type = CoupledValueFunctionIC
      function = sigma_cr1
      v = cCr
      variable = SIGMA_1CR
    []
  []
  [SIGMA_2CR]
    [InitialCondition]
      type = CoupledValueFunctionIC
      function = sigma_cr2
      v = cCr
      variable = SIGMA_2CR
    []
  []
  # Lagrange multiplier
  [lambda]
  []
[]
[Materials]
  # CALPHAD free energies
  [F_BCC_A2]
    type = DerivativeParsedMaterial
    property_name = F_BCC_A2
    outputs = exodus
    output_properties = F_BCC_A2
    expression = 'BCC_FE:=1-BCC_CR; G := 8.3145*T*(1.0*if(BCC_CR > 1.0e-15,BCC_CR*log(BCC_CR),0) + '
               '1.0*if(BCC_FE > 1.0e-15,BCC_FE*plog(BCC_FE,eps),0) + 3.0*if(BCC_VA > '
               '1.0e-15,BCC_VA*log(BCC_VA),0))/(BCC_CR + BCC_FE) + 8.3145*T*if(T < '
               '548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - 932.5*BCC_CR*BCC_FE*BCC_VA + '
               '311.5*BCC_CR*BCC_VA - '
               '1043.0*BCC_FE*BCC_VA,-8.13674105561218e-49*T^15/(0.525599232981783*BCC_CR*BCC_FE*BCC_'
               'VA*(BCC_CR - BCC_FE) - 0.894055608820709*BCC_CR*BCC_FE*BCC_VA + '
               '0.298657718120805*BCC_CR*BCC_VA - BCC_FE*BCC_VA + 9.58772770853308e-13)^15 - '
               '4.65558036243985e-30*T^9/(0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA + 0.298657718120805*BCC_CR*BCC_VA - '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^9 - '
               '1.3485349181899e-10*T^3/(0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA + 0.298657718120805*BCC_CR*BCC_VA - '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^3 + 1 - '
               '0.905299382744392*(548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '932.5*BCC_CR*BCC_FE*BCC_VA + 311.5*BCC_CR*BCC_VA - 1043.0*BCC_FE*BCC_VA + '
               '1.0e-9)/T,if(T < -548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '932.5*BCC_CR*BCC_FE*BCC_VA - 311.5*BCC_CR*BCC_VA + '
               '1043.0*BCC_FE*BCC_VA,-8.13674105561218e-49*T^15/(-0.525599232981783*BCC_CR*BCC_FE*BCC'
               '_VA*(BCC_CR - BCC_FE) + 0.894055608820709*BCC_CR*BCC_FE*BCC_VA - '
               '0.298657718120805*BCC_CR*BCC_VA + BCC_FE*BCC_VA + 9.58772770853308e-13)^15 - '
               '4.65558036243985e-30*T^9/(-0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) '
               '+ 0.894055608820709*BCC_CR*BCC_FE*BCC_VA - 0.298657718120805*BCC_CR*BCC_VA + '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^9 - '
               '1.3485349181899e-10*T^3/(-0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA - 0.298657718120805*BCC_CR*BCC_VA + '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^3 + 1 - '
               '0.905299382744392*(-548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '932.5*BCC_CR*BCC_FE*BCC_VA - 311.5*BCC_CR*BCC_VA + 1043.0*BCC_FE*BCC_VA + '
               '1.0e-9)/T,if(T > -548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '932.5*BCC_CR*BCC_FE*BCC_VA - 311.5*BCC_CR*BCC_VA + 1043.0*BCC_FE*BCC_VA & '
               '548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - 932.5*BCC_CR*BCC_FE*BCC_VA + '
               '311.5*BCC_CR*BCC_VA - 1043.0*BCC_FE*BCC_VA < '
               '0,-79209031311018.7*(-0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA - 0.298657718120805*BCC_CR*BCC_VA + '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^5/T^5 - '
               '3.83095660520737e+42*(-0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA - 0.298657718120805*BCC_CR*BCC_VA + '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^15/T^15 - '
               '1.22565886734485e+72*(-0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) + '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA - 0.298657718120805*BCC_CR*BCC_VA + '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^25/T^25,if(T > '
               '548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - 932.5*BCC_CR*BCC_FE*BCC_VA + '
               '311.5*BCC_CR*BCC_VA - 1043.0*BCC_FE*BCC_VA & 548.2*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - '
               'BCC_FE) - 932.5*BCC_CR*BCC_FE*BCC_VA + 311.5*BCC_CR*BCC_VA - 1043.0*BCC_FE*BCC_VA > '
               '0,-79209031311018.7*(0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA + 0.298657718120805*BCC_CR*BCC_VA - '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^5/T^5 - '
               '3.83095660520737e+42*(0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA + 0.298657718120805*BCC_CR*BCC_VA - '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^15/T^15 - '
               '1.22565886734485e+72*(0.525599232981783*BCC_CR*BCC_FE*BCC_VA*(BCC_CR - BCC_FE) - '
               '0.894055608820709*BCC_CR*BCC_FE*BCC_VA + 0.298657718120805*BCC_CR*BCC_VA - '
               'BCC_FE*BCC_VA + 9.58772770853308e-13)^25/T^25,0))))*log((2.15*BCC_CR*BCC_FE*BCC_VA - '
               '0.008*BCC_CR*BCC_VA + 2.22*BCC_FE*BCC_VA)*if(2.15*BCC_CR*BCC_FE*BCC_VA - '
               '0.008*BCC_CR*BCC_VA + 2.22*BCC_FE*BCC_VA <= 0,-1.0,1.0) + 1)/(BCC_CR + BCC_FE) + '
               '1.0*(BCC_CR*BCC_VA*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + '
               'BCC_FE*BCC_VA*if(T >= 298.15 & T < 1811.0,77358.5*1/T - 23.5143*T*log(T) + 124.134*T '
               '- 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= 1811.0 & T < '
               '6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - 25383.581,0)))/(BCC_CR '
               '+ BCC_FE) + 1.0*(BCC_CR*BCC_FE*BCC_VA*(500.0 - 1.5*T)*(BCC_CR - BCC_FE) + '
               'BCC_CR*BCC_FE*BCC_VA*(24600.0 - 14.98*T) + BCC_CR*BCC_FE*BCC_VA*(9.15*T - '
               '14000.0)*(BCC_CR - BCC_FE)^2)/(BCC_CR + BCC_FE); G/100000'
    coupled_variables = 'BCC_CR'
    constant_names = 'BCC_VA T eps'
    constant_expressions = '1 1000 0.01'
  []
  [F_SIGMA]
    type = DerivativeParsedMaterial
    property_name = F_SIGMA
    outputs = exodus
    output_properties = F_SIGMA
    expression = 'SIGMA_0FE := 1-SIGMA_0CR; SIGMA_1FE := 1-SIGMA_1CR; SIGMA_2FE := 1-SIGMA_2CR; G := '
               '8.3145*T*(10.0*if(SIGMA_0CR > 1.0e-15,SIGMA_0CR*plog(SIGMA_0CR,eps),0) + '
               '10.0*if(SIGMA_0FE > 1.0e-15,SIGMA_0FE*plog(SIGMA_0FE,eps),0) + 4.0*if(SIGMA_1CR > '
               '1.0e-15,SIGMA_1CR*plog(SIGMA_1CR,eps),0) + 4.0*if(SIGMA_1FE > '
               '1.0e-15,SIGMA_1FE*plog(SIGMA_1FE,eps),0) + 16.0*if(SIGMA_2CR > '
               '1.0e-15,SIGMA_2CR*plog(SIGMA_2CR,eps),0) + 16.0*if(SIGMA_2FE > '
               '1.0e-15,SIGMA_2FE*plog(SIGMA_2FE,eps),0))/(10.0*SIGMA_0CR + 10.0*SIGMA_0FE + '
               '4.0*SIGMA_1CR + 4.0*SIGMA_1FE + 16.0*SIGMA_2CR + 16.0*SIGMA_2FE) + '
               '(SIGMA_0FE*SIGMA_1CR*SIGMA_2CR*SIGMA_2FE*(-70.0*T - 170400.0) + '
               'SIGMA_0FE*SIGMA_1FE*SIGMA_2CR*SIGMA_2FE*(-10.0*T - 330839.0))/(10.0*SIGMA_0CR + '
               '10.0*SIGMA_0FE + 4.0*SIGMA_1CR + 4.0*SIGMA_1FE + 16.0*SIGMA_2CR + 16.0*SIGMA_2FE) + '
               '(SIGMA_0CR*SIGMA_1CR*SIGMA_2CR*(30.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - '
               '26.908*T*log(T) + 157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= '
               '2180.0 & T < 6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) '
               '+ 132000.0) + SIGMA_0CR*SIGMA_1CR*SIGMA_2FE*(-110.0*T + 16.0*if(T >= 298.15 & T < '
               '1811.0,77358.5*1/T - 23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - '
               '5.89269e-8*T^3.0 + 1225.7,if(T >= 1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - '
               '46.0*T*log(T) + 299.31255*T - 25383.581,0)) + 14.0*if(T >= 298.15 & T < '
               '2180.0,139250.0*1/T - 26.908*T*log(T) + 157.48*T + 0.00189435*T^2.0 - '
               '1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < 6000.0,-2.88526e+32*T^(-9.0) - '
               '50.0*T*log(T) + 344.18*T - 34869.344,0)) + 123500.0) + '
               'SIGMA_0CR*SIGMA_1FE*SIGMA_2CR*(4.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 26.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + 140486.0) '
               '+ SIGMA_0CR*SIGMA_1FE*SIGMA_2FE*(20.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 10.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + 148800.0) '
               '+ SIGMA_0FE*SIGMA_1CR*SIGMA_2CR*(10.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 20.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + 56200.0) + '
               'SIGMA_0FE*SIGMA_1CR*SIGMA_2FE*(26.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 4.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + 152700.0) '
               '+ SIGMA_0FE*SIGMA_1FE*SIGMA_2CR*(14.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 16.0*if(T >= 298.15 & T < 2180.0,139250.0*1/T - 26.908*T*log(T) + '
               '157.48*T + 0.00189435*T^2.0 - 1.47721e-6*T^3.0 - 8856.94,if(T >= 2180.0 & T < '
               '6000.0,-2.88526e+32*T^(-9.0) - 50.0*T*log(T) + 344.18*T - 34869.344,0)) + 46200.0) + '
               'SIGMA_0FE*SIGMA_1FE*SIGMA_2FE*(30.0*if(T >= 298.15 & T < 1811.0,77358.5*1/T - '
               '23.5143*T*log(T) + 124.134*T - 0.00439752*T^2.0 - 5.89269e-8*T^3.0 + 1225.7,if(T >= '
               '1811.0 & T < 6000.0,2.2960305e+31*T^(-9.0) - 46.0*T*log(T) + 299.31255*T - '
               '25383.581,0)) + 173333.0))/(10.0*SIGMA_0CR + 10.0*SIGMA_0FE + 4.0*SIGMA_1CR + '
               '4.0*SIGMA_1FE + 16.0*SIGMA_2CR + 16.0*SIGMA_2FE); G/100000'
    coupled_variables = 'SIGMA_0CR SIGMA_1CR SIGMA_2CR'
    constant_names = 'T eps'
    constant_expressions = '1000 0.01'
  []
  # h(eta)
  [h1]
    type = SwitchingFunctionMaterial
    function_name = h1
    h_order = HIGH
    eta = eta1
  []
  [h2]
    type = SwitchingFunctionMaterial
    function_name = h2
    h_order = HIGH
    eta = eta2
  []
  # g(eta)
  [g1]
    type = BarrierFunctionMaterial
    function_name = g1
    g_order = SIMPLE
    eta = eta1
  []
  [g2]
    type = BarrierFunctionMaterial
    function_name = g2
    g_order = SIMPLE
    eta = eta2
  []
  # constant properties
  [constants]
    type = GenericConstantMaterial
    prop_names = 'D   L   kappa'
    prop_values = '10  1   0.1  '
  []
  # Coefficients for diffusion equation
  [Dh1]
    type = DerivativeParsedMaterial
    material_property_names = 'D h1(eta1)'
    expression = D*h1
    property_name = Dh1
    coupled_variables = eta1
    derivative_order = 1
  []
  [Dh2a]
    type = DerivativeParsedMaterial
    material_property_names = 'D h2(eta2)'
    expression = D*h2*10/30
    property_name = Dh2a
    coupled_variables = eta2
    derivative_order = 1
  []
  [Dh2b]
    type = DerivativeParsedMaterial
    material_property_names = 'D h2(eta2)'
    expression = D*h2*4/30
    property_name = Dh2b
    coupled_variables = eta2
    derivative_order = 1
  []
  [Dh2c]
    type = DerivativeParsedMaterial
    material_property_names = 'D h2(eta2)'
    expression = D*h2*16/30
    property_name = Dh2c
    coupled_variables = eta2
    derivative_order = 1
  []
[]
[Kernels]
  #Kernels for diffusion equation
  [diff_time]
    type = TimeDerivative
    variable = cCr
  []
  [diff_c1]
    type = MatDiffusion
    variable = cCr
    diffusivity = Dh1
    v = BCC_CR
    args = eta1
  []
  [diff_c2a]
    type = MatDiffusion
    variable = cCr
    diffusivity = Dh2a
    v = SIGMA_0CR
    args = eta2
  []
  [diff_c2b]
    type = MatDiffusion
    variable = cCr
    diffusivity = Dh2b
    v = SIGMA_1CR
    args = eta2
  []
  [diff_c2c]
    type = MatDiffusion
    variable = cCr
    diffusivity = Dh2c
    v = SIGMA_2CR
    args = eta2
  []
  # enforce pointwise equality of chemical potentials
  [chempot1a2a]
    # The BCC phase has only one sublattice
    # we tie it to the first sublattice with site fraction 10/(10+4+16) in the sigma phase
    type = KKSPhaseChemicalPotential
    variable = BCC_CR
    cb = SIGMA_0CR
    kb = '${fparse 10/30}'
    fa_name = F_BCC_A2
    fb_name = F_SIGMA
    args_b = 'SIGMA_1CR SIGMA_2CR'
  []
  [chempot2a2b]
    # This kernel ties the first two sublattices in the sigma phase together
    type = SLKKSChemicalPotential
    variable = SIGMA_0CR
    a = 10
    cs = SIGMA_1CR
    as = 4
    F = F_SIGMA
    coupled_variables = 'SIGMA_2CR'
  []
  [chempot2b2c]
    # This kernel ties the remaining two sublattices in the sigma phase together
    type = SLKKSChemicalPotential
    variable = SIGMA_1CR
    a = 4
    cs = SIGMA_2CR
    as = 16
    F = F_SIGMA
    coupled_variables = 'SIGMA_0CR'
  []
  [phaseconcentration]
    # This kernel ties the sum of the sublattice concentrations to the global concentration cCr
    type = SLKKSMultiPhaseConcentration
    variable = SIGMA_2CR
    c = cCr
    ns = '1      3'
    as = '1      10        4         16'
    cs = 'BCC_CR SIGMA_0CR SIGMA_1CR SIGMA_2CR'
    h_names = 'h1   h2'
    eta = 'eta1 eta2'
  []
  # Kernels for Allen-Cahn equation for eta1
  [deta1dt]
    type = TimeDerivative
    variable = eta1
  []
  [ACBulkF1]
    type = KKSMultiACBulkF
    variable = eta1
    Fj_names = 'F_BCC_A2 F_SIGMA'
    hj_names = 'h1    h2'
    gi_name = g1
    eta_i = eta1
    wi = 0.1
    coupled_variables = 'BCC_CR SIGMA_0CR SIGMA_1CR SIGMA_2CR eta2'
  []
  [ACBulkC1]
    type = SLKKSMultiACBulkC
    variable = eta1
    F = F_BCC_A2
    c = BCC_CR
    ns = '1      3'
    as = '1      10        4         16'
    cs = 'BCC_CR SIGMA_0CR SIGMA_1CR SIGMA_2CR'
    h_names = 'h1   h2'
    eta = 'eta1 eta2'
  []
  [ACInterface1]
    type = ACInterface
    variable = eta1
    kappa_name = kappa
  []
  [lagrange1]
    type = SwitchingFunctionConstraintEta
    variable = eta1
    h_name = h1
    lambda = lambda
    coupled_variables = 'eta2'
  []
  # Kernels for Allen-Cahn equation for eta1
  [deta2dt]
    type = TimeDerivative
    variable = eta2
  []
  [ACBulkF2]
    type = KKSMultiACBulkF
    variable = eta2
    Fj_names = 'F_BCC_A2 F_SIGMA'
    hj_names = 'h1    h2'
    gi_name = g2
    eta_i = eta2
    wi = 0.1
    coupled_variables = 'BCC_CR SIGMA_0CR SIGMA_1CR SIGMA_2CR eta1'
  []
  [ACBulkC2]
    type = SLKKSMultiACBulkC
    variable = eta2
    F = F_BCC_A2
    c = BCC_CR
    ns = '1      3'
    as = '1      10        4         16'
    cs = 'BCC_CR SIGMA_0CR SIGMA_1CR SIGMA_2CR'
    h_names = 'h1   h2'
    eta = 'eta1 eta2'
  []
  [ACInterface2]
    type = ACInterface
    variable = eta2
    kappa_name = kappa
  []
  [lagrange2]
    type = SwitchingFunctionConstraintEta
    variable = eta2
    h_name = h2
    lambda = lambda
    coupled_variables = 'eta1'
  []
  # Lagrange-multiplier constraint kernel for lambda
  [lagrange]
    type = SwitchingFunctionConstraintLagrange
    variable = lambda
    h_names = 'h1   h2'
    etas = 'eta1 eta2'
    epsilon = 1e-6
  []
[]
[AuxKernels]
  [GlobalFreeEnergy]
    type = KKSMultiFreeEnergy
    variable = Fglobal
    Fj_names = 'F_BCC_A2 F_SIGMA'
    hj_names = 'h1 h2'
    gj_names = 'g1 g2'
    interfacial_vars = 'eta1 eta2'
    kappa_names = 'kappa kappa'
    w = 0.1
  []
[]
[Executioner]
  type = Transient
  solve_type = 'NEWTON'
  line_search = none
  petsc_options_iname = '-pc_type -sub_pc_type -sub_pc_factor_shift_type -ksp_gmres_restart'
  petsc_options_value = 'asm      lu          nonzero                    30'
  l_max_its = 100
  nl_max_its = 20
  nl_abs_tol = 1e-10
  end_time = 10000
  [TimeStepper]
    type = IterationAdaptiveDT
    optimal_iterations = 12
    iteration_window = 2
    growth_factor = 1.5
    cutback_factor = 0.7
    dt = 0.1
  []
[]
[VectorPostprocessors]
  [var]
    type = LineValueSampler
    start_point = '-25 0 0'
    end_point = '25 0 0'
    variable = 'cCr eta1 eta2 SIGMA_0CR SIGMA_1CR SIGMA_2CR'
    num_points = 151
    sort_by = id
    execute_on = 'initial timestep_end'
  []
  [mat]
    type = LineMaterialRealSampler
    start = '-25 0 0'
    end = '25 0 0'
    property = 'F_BCC_A2 F_SIGMA'
    sort_by = id
    execute_on = 'initial timestep_end'
  []
[]
[Postprocessors]
  [F]
    type = ElementIntegralVariablePostprocessor
    variable = Fglobal
    execute_on = 'initial timestep_end'
  []
  [cmin]
    type = NodalExtremeValue
    value_type = min
    variable = cCr
    execute_on = 'initial timestep_end'
  []
  [cmax]
    type = NodalExtremeValue
    value_type = max
    variable = cCr
    execute_on = 'initial timestep_end'
  []
  [ctotal]
    type = ElementIntegralVariablePostprocessor
    variable = cCr
    execute_on = 'initial timestep_end'
  []
[]
[Outputs]
  exodus = true
  print_linear_residuals = false
  csv = true
  perf_graph = true
[]