- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Description:The name of the variable that this residual object operates on
 
MatReaction
Kernel to add -L*v, where L=reaction rate, v=variable
Implements  where  (mob_name) is a reaction rate,  is either a coupled variable (v) or - if not explicitly specified - the non-linear variable the kernel is operating on.
Note the negative sign, which does not appear in Reaction or CoefReaction.
Input Parameters
- argsVector of nonlinear variable arguments this object depends on
C++ Type:std::vector<VariableName>
Options:
Description:Vector of nonlinear variable arguments this object depends on
 - blockThe list of block ids (SubdomainID) that this object will be applied
C++ Type:std::vector<SubdomainName>
Options:
Description:The list of block ids (SubdomainID) that this object will be applied
 - displacementsThe displacements
C++ Type:std::vector<VariableName>
Options:
Description:The displacements
 - mob_nameLThe reaction rate used with the kernel
Default:L
C++ Type:MaterialPropertyName
Options:
Description:The reaction rate used with the kernel
 - vSet this to make v a coupled variable, otherwise it will use the kernel's nonlinear variable for v
C++ Type:std::vector<VariableName>
Options:
Description:Set this to make v a coupled variable, otherwise it will use the kernel's nonlinear variable for v
 
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Options:
Description:Adds user-defined labels for accessing object parameters via control logic.
 - diag_save_inThe name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Options:
Description:The name of auxiliary variables to save this Kernel's diagonal Jacobian contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - enableTrueSet the enabled status of the MooseObject.
Default:True
C++ Type:bool
Options:
Description:Set the enabled status of the MooseObject.
 - implicitTrueDetermines whether this object is calculated using an implicit or explicit form
Default:True
C++ Type:bool
Options:
Description:Determines whether this object is calculated using an implicit or explicit form
 - save_inThe name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
C++ Type:std::vector<AuxVariableName>
Options:
Description:The name of auxiliary variables to save this Kernel's residual contributions to. Everything about that variable must match everything about this variable (the type, what blocks it's on, etc.)
 - seed0The seed for the master random number generator
Default:0
C++ Type:unsigned int
Options:
Description:The seed for the master random number generator
 - use_displaced_meshFalseWhether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
Default:False
C++ Type:bool
Options:
Description:Whether or not this object should use the displaced mesh for computation. Note that in the case this is true but no displacements are provided in the Mesh block the undisplaced mesh will still be used.
 
Advanced Parameters
- extra_matrix_tagsThe extra tags for the matrices this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the matrices this Kernel should fill
 - extra_vector_tagsThe extra tags for the vectors this Kernel should fill
C++ Type:std::vector<TagName>
Options:
Description:The extra tags for the vectors this Kernel should fill
 - matrix_tagssystemThe tag for the matrices this Kernel should fill
Default:system
C++ Type:MultiMooseEnum
Options:nontime, system
Description:The tag for the matrices this Kernel should fill
 - vector_tagsnontimeThe tag for the vectors this Kernel should fill
Default:nontime
C++ Type:MultiMooseEnum
Options:nontime, time
Description:The tag for the vectors this Kernel should fill
 
Tagging Parameters
Input Files
- (modules/phase_field/test/tests/KKS_system/kks_multiphase.i)
 - (modules/fsi/test/tests/fsi_2d/fsi_flat_channel.i)
 - (modules/phase_field/test/tests/phase_field_kernels/MatGradSquareCoupled.i)
 - (modules/combined/examples/publications/rapid_dev/fig6.i)
 - (modules/phase_field/test/tests/phase_field_kernels/CoupledCoefAllenCahn.i)
 - (modules/phase_field/test/tests/phase_field_kernels/CoupledAllenCahn.i)
 - (modules/phase_field/test/tests/mobility_derivative/coupledmatdiffusion.i)
 - (modules/phase_field/test/tests/SimpleACInterface/SimpleCoupledACInterface.i)
 
(modules/phase_field/test/tests/KKS_system/kks_multiphase.i)
#
# This test is for the 3-phase KKS model
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 20
  ny = 20
  nz = 0
  xmin = 0
  xmax = 40
  ymin = 0
  ymax = 40
  zmin = 0
  zmax = 0
  elem_type = QUAD4
[]
[BCs]
  [./Periodic]
    [./all]
      auto_direction = 'x y'
    [../]
  [../]
[]
[AuxVariables]
  [./Energy]
    order = CONSTANT
    family = MONOMIAL
  [../]
[]
[Variables]
  # concentration
  [./c]
    order = FIRST
    family = LAGRANGE
  [../]
  # order parameter 1
  [./eta1]
    order = FIRST
    family = LAGRANGE
  [../]
  # order parameter 2
  [./eta2]
    order = FIRST
    family = LAGRANGE
  [../]
  # order parameter 3
  [./eta3]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.0
  [../]
  # phase concentration 1
  [./c1]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.2
  [../]
  # phase concentration 2
  [./c2]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.5
  [../]
  # phase concentration 3
  [./c3]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.8
  [../]
  # Lagrange multiplier
  [./lambda]
    order = FIRST
    family = LAGRANGE
    initial_condition = 0.0
  [../]
[]
[ICs]
  [./eta1]
    variable = eta1
    type = SmoothCircleIC
    x1 = 20.0
    y1 = 20.0
    radius = 10
    invalue = 0.9
    outvalue = 0.1
    int_width = 4
  [../]
  [./eta2]
    variable = eta2
    type = SmoothCircleIC
    x1 = 20.0
    y1 = 20.0
    radius = 10
    invalue = 0.1
    outvalue = 0.9
    int_width = 4
  [../]
  [./c]
    variable = c
    type = SmoothCircleIC
    x1 = 20.0
    y1 = 20.0
    radius = 10
    invalue = 0.2
    outvalue = 0.5
    int_width = 2
  [../]
[]
[Materials]
  # simple toy free energies
  [./f1]
    type = DerivativeParsedMaterial
    f_name = F1
    args = 'c1'
    function = '20*(c1-0.2)^2'
  [../]
  [./f2]
    type = DerivativeParsedMaterial
    f_name = F2
    args = 'c2'
    function = '20*(c2-0.5)^2'
  [../]
  [./f3]
    type = DerivativeParsedMaterial
    f_name = F3
    args = 'c3'
    function = '20*(c3-0.8)^2'
  [../]
  # Switching functions for each phase
  # h1(eta1, eta2, eta3)
  [./h1]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta1
    eta_j = eta2
    eta_k = eta3
    f_name = h1
  [../]
  # h2(eta1, eta2, eta3)
  [./h2]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta2
    eta_j = eta3
    eta_k = eta1
    f_name = h2
  [../]
  # h3(eta1, eta2, eta3)
  [./h3]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta3
    eta_j = eta1
    eta_k = eta2
    f_name = h3
  [../]
  # Coefficients for diffusion equation
  [./Dh1]
    type = DerivativeParsedMaterial
    material_property_names = 'D h1'
    function = D*h1
    f_name = Dh1
  [../]
  [./Dh2]
    type = DerivativeParsedMaterial
    material_property_names = 'D h2'
    function = D*h2
    f_name = Dh2
  [../]
  [./Dh3]
    type = DerivativeParsedMaterial
    material_property_names = 'D h3'
    function = D*h3
    f_name = Dh3
  [../]
  # Barrier functions for each phase
  [./g1]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta1
    function_name = g1
  [../]
  [./g2]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta2
    function_name = g2
  [../]
  [./g3]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta3
    function_name = g3
  [../]
  # constant properties
  [./constants]
    type = GenericConstantMaterial
    prop_names  = 'L   kappa  D'
    prop_values = '0.7 1.0    1'
  [../]
[]
[Kernels]
  #Kernels for diffusion equation
  [./diff_time]
    type = TimeDerivative
    variable = c
  [../]
  [./diff_c1]
    type = MatDiffusion
    variable = c
    diffusivity = Dh1
    v = c1
  [../]
  [./diff_c2]
    type = MatDiffusion
    variable = c
    diffusivity = Dh2
    v = c2
  [../]
  [./diff_c3]
    type = MatDiffusion
    variable = c
    diffusivity = Dh3
    v = c3
  [../]
  # Kernels for Allen-Cahn equation for eta1
  [./deta1dt]
    type = TimeDerivative
    variable = eta1
  [../]
  [./ACBulkF1]
    type = KKSMultiACBulkF
    variable  = eta1
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g1
    eta_i     = eta1
    wi        = 1.0
    args      = 'c1 c2 c3 eta2 eta3'
  [../]
  [./ACBulkC1]
    type = KKSMultiACBulkC
    variable  = eta1
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta1
    args      = 'eta2 eta3'
  [../]
  [./ACInterface1]
    type = ACInterface
    variable = eta1
    kappa_name = kappa
  [../]
  [./multipler1]
    type = MatReaction
    variable = eta1
    v = lambda
    mob_name = L
  [../]
  # Kernels for Allen-Cahn equation for eta2
  [./deta2dt]
    type = TimeDerivative
    variable = eta2
  [../]
  [./ACBulkF2]
    type = KKSMultiACBulkF
    variable  = eta2
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g2
    eta_i     = eta2
    wi        = 1.0
    args      = 'c1 c2 c3 eta1 eta3'
  [../]
  [./ACBulkC2]
    type = KKSMultiACBulkC
    variable  = eta2
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta2
    args      = 'eta1 eta3'
  [../]
  [./ACInterface2]
    type = ACInterface
    variable = eta2
    kappa_name = kappa
  [../]
  [./multipler2]
    type = MatReaction
    variable = eta2
    v = lambda
    mob_name = L
  [../]
  # Kernels for the Lagrange multiplier equation
  [./mult_lambda]
    type = MatReaction
    variable = lambda
    mob_name = 3
  [../]
  [./mult_ACBulkF_1]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g1
    eta_i     = eta1
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta2 eta3'
  [../]
  [./mult_ACBulkC_1]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta1
    args      = 'eta2 eta3'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_1]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta1
    kappa_name = kappa
    mob_name = 1
  [../]
  [./mult_ACBulkF_2]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g2
    eta_i     = eta2
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta1 eta3'
  [../]
  [./mult_ACBulkC_2]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta2
    args      = 'eta1 eta3'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_2]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta2
    kappa_name = kappa
    mob_name = 1
  [../]
  [./mult_ACBulkF_3]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g3
    eta_i     = eta3
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta1 eta2'
  [../]
  [./mult_ACBulkC_3]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta3
    args      = 'eta1 eta2'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_3]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta3
    kappa_name = kappa
    mob_name = 1
  [../]
  # Kernels for constraint equation eta1 + eta2 + eta3 = 1
  # eta3 is the nonlinear variable for the constraint equation
  [./eta3reaction]
    type = MatReaction
    variable = eta3
    mob_name = 1
  [../]
  [./eta1reaction]
    type = MatReaction
    variable = eta3
    v = eta1
    mob_name = 1
  [../]
  [./eta2reaction]
    type = MatReaction
    variable = eta3
    v = eta2
    mob_name = 1
  [../]
  [./one]
    type = BodyForce
    variable = eta3
    value = -1.0
  [../]
  # Phase concentration constraints
  [./chempot12]
    type = KKSPhaseChemicalPotential
    variable = c1
    cb       = c2
    fa_name  = F1
    fb_name  = F2
  [../]
  [./chempot23]
    type = KKSPhaseChemicalPotential
    variable = c2
    cb       = c3
    fa_name  = F2
    fb_name  = F3
  [../]
  [./phaseconcentration]
    type = KKSMultiPhaseConcentration
    variable = c3
    cj = 'c1 c2 c3'
    hj_names = 'h1 h2 h3'
    etas = 'eta1 eta2 eta3'
    c = c
  [../]
[]
[AuxKernels]
  [./Energy_total]
    type = KKSMultiFreeEnergy
    Fj_names = 'F1 F2 F3'
    hj_names = 'h1 h2 h3'
    gj_names = 'g1 g2 g3'
    variable = Energy
    w = 1
    interfacial_vars =  'eta1  eta2  eta3'
    kappa_names =       'kappa kappa kappa'
  [../]
[]
[Executioner]
  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 = 30
  nl_max_its = 10
  l_tol = 1.0e-4
  nl_rel_tol = 1.0e-10
  nl_abs_tol = 1.0e-11
  num_steps = 2
  dt = 0.5
[]
[Preconditioning]
  active = 'full'
  [./full]
    type = SMP
    full = true
  [../]
  [./mydebug]
    type = FDP
    full = true
  [../]
[]
[Outputs]
  exodus = true
[]
(modules/fsi/test/tests/fsi_2d/fsi_flat_channel.i)
[GlobalParams]
  gravity = '0 0 0'
  integrate_p_by_parts = true
  laplace = true
  convective_term = true
  transient_term = true
  pspg = true
  displacements = 'disp_x disp_y'
[]
[Mesh]
  [gmg]
  type = GeneratedMeshGenerator
  dim = 2
  xmin = 0
  xmax = 3.0
  ymin = 0
  ymax = 1.0
  nx = 10
  ny = 15
  elem_type = QUAD4
  []
  [subdomain1]
    type = SubdomainBoundingBoxGenerator
    bottom_left = '0.0 0.5 0'
    block_id = 1
    top_right = '3.0 1.0 0'
    input = gmg
  []
  [interface]
    type = SideSetsBetweenSubdomainsGenerator
    primary_block = '0'
    paired_block = '1'
    new_boundary = 'master0_interface'
    input = subdomain1
  []
  [break_boundary]
    type = BreakBoundaryOnSubdomainGenerator
    input = interface
  []
[]
[Variables]
  [./vel_x]
    block = 0
  [../]
  [./vel_y]
    block = 0
  [../]
  [./p]
    block = 0
  [../]
  [./disp_x]
  [../]
  [./disp_y]
  [../]
  [./vel_x_solid]
    block = 1
  [../]
  [./vel_y_solid]
    block = 1
  [../]
[]
[Kernels]
  [./vel_x_time]
    type = INSMomentumTimeDerivative
    variable = vel_x
    block = 0
    use_displaced_mesh = true
  [../]
  [./vel_y_time]
    type = INSMomentumTimeDerivative
    variable = vel_y
    block = 0
    use_displaced_mesh = true
  [../]
  [./mass]
    type = INSMass
    variable = p
    u = vel_x
    v = vel_y
    p = p
    block = 0
    use_displaced_mesh = true
  [../]
  [./x_momentum_space]
    type = INSMomentumLaplaceForm
    variable = vel_x
    u = vel_x
    v = vel_y
    p = p
    component = 0
    block = 0
    use_displaced_mesh = true
  [../]
  [./y_momentum_space]
    type = INSMomentumLaplaceForm
    variable = vel_y
    u = vel_x
    v = vel_y
    p = p
    component = 1
    block = 0
    use_displaced_mesh = true
  [../]
  [./vel_x_mesh]
    type = ConvectedMesh
    disp_x = disp_x
    disp_y = disp_y
    variable = vel_x
    block = 0
    use_displaced_mesh = true
  [../]
  [./vel_y_mesh]
    type = ConvectedMesh
    disp_x = disp_x
    disp_y = disp_y
    variable = vel_y
    block = 0
    use_displaced_mesh = true
  [../]
  [./disp_x_fluid]
    type = Diffusion
    variable = disp_x
    block = 0
  [../]
  [./disp_y_fluid]
    type = Diffusion
    variable = disp_y
    block = 0
  [../]
  [./accel_tensor_x]
    type = CoupledTimeDerivative
    variable = disp_x
    v = vel_x_solid
    block = 1
  [../]
  [./accel_tensor_y]
    type = CoupledTimeDerivative
    variable = disp_y
    v = vel_y_solid
    block = 1
  [../]
  [./vxs_time_derivative_term]
    type = CoupledTimeDerivative
    variable = vel_x_solid
    v = disp_x
    block = 1
  [../]
  [./vys_time_derivative_term]
    type = CoupledTimeDerivative
    variable = vel_y_solid
    v = disp_y
    block = 1
  [../]
  [./source_vxs]
    type = MatReaction
    variable = vel_x_solid
    block = 1
    mob_name = 1
  [../]
  [./source_vys]
    type = MatReaction
    variable = vel_y_solid
    block = 1
    mob_name = 1
  [../]
[]
[InterfaceKernels]
  [./penalty_interface_x]
    type = CoupledPenaltyInterfaceDiffusion
    variable = vel_x
    neighbor_var = disp_x
    secondary_coupled_var = vel_x_solid
    boundary = master0_interface
    penalty = 1e6
  [../]
  [./penalty_interface_y]
    type = CoupledPenaltyInterfaceDiffusion
    variable = vel_y
    neighbor_var = disp_y
    secondary_coupled_var = vel_y_solid
    boundary = master0_interface
    penalty = 1e6
  [../]
[]
[Modules/TensorMechanics/Master]
  [./solid_domain]
    strain = SMALL
    incremental = false
    # generate_output = 'strain_xx strain_yy strain_zz' ## Not at all necessary, but nice
    block = '1'
  [../]
[]
[Materials]
  [./elasticity_tensor]
    type = ComputeIsotropicElasticityTensor
    youngs_modulus = 1e2
    poissons_ratio = 0.3
    block = '1'
  [../]
  [./small_stress]
    type = ComputeLinearElasticStress
    block = 1
  [../]
  [./const]
    type = GenericConstantMaterial
    block = 0
    prop_names = 'rho mu'
    prop_values = '1  1'
  [../]
[]
[BCs]
  [./fluid_x_no_slip]
    type = DirichletBC
    variable = vel_x
    boundary = 'bottom'
    value = 0.0
  [../]
  [./fluid_y_no_slip]
    type = DirichletBC
    variable = vel_y
    boundary = 'bottom left_to_0'
    value = 0.0
  [../]
  [./x_inlet]
    type = FunctionDirichletBC
    variable = vel_x
    boundary = 'left_to_0'
    function = 'inlet_func'
  [../]
  [./no_disp_x]
    type = DirichletBC
    variable = disp_x
    boundary = 'bottom top left_to_1 right_to_1 left_to_0 right_to_0'
    value = 0
  [../]
  [./no_disp_y]
    type = DirichletBC
    variable = disp_y
    boundary = 'bottom top left_to_1 right_to_1 left_to_0 right_to_0'
    value = 0
  [../]
  [./solid_x_no_slip]
    type = DirichletBC
    variable = vel_x_solid
    boundary = 'top left_to_1 right_to_1'
    value = 0.0
  [../]
  [./solid_y_no_slip]
    type = DirichletBC
    variable = vel_y_solid
    boundary = 'top left_to_1 right_to_1'
    value = 0.0
  [../]
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  num_steps = 5
  # num_steps = 60
  dt = 0.1
  dtmin = 0.1
  solve_type = 'PJFNK'
  petsc_options_iname = '-pc_type'
  petsc_options_value = 'lu'
  line_search = none
[]
[Outputs]
  [./out]
    type = Exodus
  [../]
[]
[Functions]
  [./inlet_func]
    type = ParsedFunction
    value = '(-16 * (y - 0.25)^2 + 1) * (1 + cos(t))'
  [../]
[]
(modules/phase_field/test/tests/phase_field_kernels/MatGradSquareCoupled.i)
#
# Test the MatGradSquareCoupled kernel
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 15
  ny = 15
  nz = 0
  xmin = 0
  xmax = 50
  ymin = 0
  ymax = 50
  zmin = 0
  zmax = 50
  elem_type = QUAD4
[]
[Variables]
  [./w]
  [../]
  [./eta]
    order = FIRST
    family = LAGRANGE
    [./InitialCondition]
      type = SmoothCircleIC
      x1 = 25.0
      y1 = 25.0
      radius = 6.0
      invalue = 1.0
      outvalue = 0.0
      int_width = 3.0
    [../]
  [../]
[]
[Kernels]
  [./detadt]
    type = TimeDerivative
    variable = eta
  [../]
  [./ACBulk]
    type = CoupledAllenCahn
    variable = w
    v = eta
    f_name = F
    mob_name = 1
  [../]
  [./W]
    type = MatReaction
    variable = w
    mob_name = -1
  [../]
  [./CoupledBulk]
    type = MatReaction
    variable = eta
    v = w
    mob_name = L
  [../]
  [./ACInterface]
    type = ACInterface
    variable = eta
    kappa_name = 1
    mob_name = L
    args = w
  [../]
# MatGradSquareCoupled kernel
  [./nabla_eta]
    type = MatGradSquareCoupled
    variable = w
    elec_potential = eta
    prefactor = 0.5
  [../]
[]
[Materials]
  [./mobility]
    type = DerivativeParsedMaterial
    f_name  = L
    args = 'eta w'
    function = '(1.5-eta)^2+(1.5-w)^2'
    derivative_order = 2
  [../]
  [./free_energy]
    type = DerivativeParsedMaterial
    f_name = F
    args = 'eta'
    function = 'eta^2 * (1-eta)^2'
    derivative_order = 2
  [../]
[]
[Preconditioning]
  [./smp]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  scheme = 'bdf2'
  solve_type = 'PJFNK'
  l_max_its = 15
  l_tol = 1.0e-4
  nl_max_its = 10
  nl_rel_tol = 1.0e-11
  start_time = 0.0
  num_steps = 2
  dt = 0.5
[]
[Outputs]
  hide = w
  exodus = true
  console = true
[]
(modules/combined/examples/publications/rapid_dev/fig6.i)
#
# Fig. 6 input for 10.1016/j.commatsci.2017.02.017
# D. Schwen et al./Computational Materials Science 132 (2017) 36-45
# Three phase interface simulation demonstrating the interfacial stability
# w.r.t. formation of a tspurious third phase
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 120
  ny = 120
  nz = 0
  xmin = 0
  xmax = 40
  ymin = 0
  ymax = 40
  zmin = 0
  zmax = 0
  elem_type = QUAD4
[]
[Variables]
  # concentration
  [./c]
  [../]
  # order parameter 1
  [./eta1]
  [../]
  # order parameter 2
  [./eta2]
  [../]
  # order parameter 3
  [./eta3]
  [../]
  # phase concentration 1
  [./c1]
    initial_condition = 0.4
  [../]
  # phase concentration 2
  [./c2]
    initial_condition = 0.5
  [../]
  # phase concentration 3
  [./c3]
    initial_condition = 0.8
  [../]
  # Lagrange multiplier
  [./lambda]
    initial_condition = 0.0
  [../]
[]
[AuxVariables]
  [./T]
    [./InitialCondition]
      type = FunctionIC
      function = 'x-10'
    [../]
  [../]
[]
[Functions]
  [./ic_func_eta1]
    type = ParsedFunction
    value = '0.5*(1.0+tanh((x-10)/sqrt(2.0))) * 0.5*(1.0+tanh((y-10)/sqrt(2.0)))'
  [../]
  [./ic_func_eta2]
    type = ParsedFunction
    value = '0.5*(1.0-tanh((x-10)/sqrt(2.0)))'
  [../]
  [./ic_func_eta3]
    type = ParsedFunction
    value = '1 - 0.5*(1.0-tanh((x-10)/sqrt(2.0)))
              - 0.5*(1.0+tanh((x-10)/sqrt(2.0))) * 0.5*(1.0+tanh((y-10)/sqrt(2.0)))'
  [../]
  [./ic_func_c]
    type = ParsedFunction
    value = '0.5 * 0.5*(1.0-tanh((x-10)/sqrt(2.0)))
              + 0.4 * 0.5*(1.0+tanh((x-10)/sqrt(2.0))) * 0.5*(1.0+tanh((y-10)/sqrt(2.0)))
              + 0.8 * (1 - 0.5*(1.0-tanh((x-10)/sqrt(2.0)))
                        - 0.5*(1.0+tanh((x-10)/sqrt(2.0))) * 0.5*(1.0+tanh((y-10)/sqrt(2.0))))'
  [../]
[]
[ICs]
  [./eta1]
    variable = eta1
    type = FunctionIC
    function = ic_func_eta1
  [../]
  [./eta2]
    variable = eta2
    type = FunctionIC
    function = ic_func_eta2
  [../]
  [./eta3]
    variable = eta3
    type = FunctionIC
    function = ic_func_eta3
  [../]
  [./c]
    variable = c
    type = FunctionIC
    function = ic_func_c
  [../]
[]
[Materials]
  # simple toy free energies
  [./f1]
    type = DerivativeParsedMaterial
    f_name = F1
    args = 'c1'
    function = '20*(c1-0.4)^2'
  [../]
  [./f2]
    type = DerivativeParsedMaterial
    f_name = F2
    args = 'c2 T'
    function = '20*(c2-0.5)^2 + 0.01*T'
  [../]
  [./f3]
    type = DerivativeParsedMaterial
    f_name = F3
    args = 'c3'
    function = '20*(c3-0.8)^2'
  [../]
  # Switching functions for each phase
  # h1(eta1, eta2, eta3)
  [./h1]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta1
    eta_j = eta2
    eta_k = eta3
    f_name = h1
  [../]
  # h2(eta1, eta2, eta3)
  [./h2]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta2
    eta_j = eta3
    eta_k = eta1
    f_name = h2
  [../]
  # h3(eta1, eta2, eta3)
  [./h3]
    type = SwitchingFunction3PhaseMaterial
    eta_i = eta3
    eta_j = eta1
    eta_k = eta2
    f_name = h3
  [../]
  # Coefficients for diffusion equation
  [./Dh1]
    type = DerivativeParsedMaterial
    material_property_names = 'D h1'
    function = D*h1
    f_name = Dh1
  [../]
  [./Dh2]
    type = DerivativeParsedMaterial
    material_property_names = 'D h2'
    function = D*h2
    f_name = Dh2
  [../]
  [./Dh3]
    type = DerivativeParsedMaterial
    material_property_names = 'D h3'
    function = D*h3
    f_name = Dh3
  [../]
  # Barrier functions for each phase
  [./g1]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta1
    function_name = g1
  [../]
  [./g2]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta2
    function_name = g2
  [../]
  [./g3]
    type = BarrierFunctionMaterial
    g_order = SIMPLE
    eta = eta3
    function_name = g3
  [../]
  # constant properties
  [./constants]
    type = GenericConstantMaterial
    prop_names  = 'L   kappa  D'
    prop_values = '1.0 1.0    1'
  [../]
[]
[Kernels]
  #Kernels for diffusion equation
  [./diff_time]
    type = TimeDerivative
    variable = c
  [../]
  [./diff_c1]
    type = MatDiffusion
    variable = c
    diffusivity = Dh1
    v = c1
  [../]
  [./diff_c2]
    type = MatDiffusion
    variable = c
    diffusivity = Dh2
    v = c2
  [../]
  [./diff_c3]
    type = MatDiffusion
    variable = c
    diffusivity = Dh3
    v = c3
  [../]
  # Kernels for Allen-Cahn equation for eta1
  [./deta1dt]
    type = TimeDerivative
    variable = eta1
  [../]
  [./ACBulkF1]
    type = KKSMultiACBulkF
    variable  = eta1
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g1
    eta_i     = eta1
    wi        = 1.0
    args      = 'c1 c2 c3 eta2 eta3'
  [../]
  [./ACBulkC1]
    type = KKSMultiACBulkC
    variable  = eta1
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta1
    args      = 'eta2 eta3'
  [../]
  [./ACInterface1]
    type = ACInterface
    variable = eta1
    kappa_name = kappa
  [../]
  [./multipler1]
    type = MatReaction
    variable = eta1
    v = lambda
    mob_name = L
  [../]
  # Kernels for Allen-Cahn equation for eta2
  [./deta2dt]
    type = TimeDerivative
    variable = eta2
  [../]
  [./ACBulkF2]
    type = KKSMultiACBulkF
    variable  = eta2
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g2
    eta_i     = eta2
    wi        = 1.0
    args      = 'c1 c2 c3 eta1 eta3'
  [../]
  [./ACBulkC2]
    type = KKSMultiACBulkC
    variable  = eta2
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta2
    args      = 'eta1 eta3'
  [../]
  [./ACInterface2]
    type = ACInterface
    variable = eta2
    kappa_name = kappa
  [../]
  [./multipler2]
    type = MatReaction
    variable = eta2
    v = lambda
    mob_name = L
  [../]
  # Kernels for the Lagrange multiplier equation
  [./mult_lambda]
    type = MatReaction
    variable = lambda
    mob_name = 3
  [../]
  [./mult_ACBulkF_1]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g1
    eta_i     = eta1
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta2 eta3'
  [../]
  [./mult_ACBulkC_1]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta1
    args      = 'eta2 eta3'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_1]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta1
    kappa_name = kappa
    mob_name = 1
  [../]
  [./mult_ACBulkF_2]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g2
    eta_i     = eta2
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta1 eta3'
  [../]
  [./mult_ACBulkC_2]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta2
    args      = 'eta1 eta3'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_2]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta2
    kappa_name = kappa
    mob_name = 1
  [../]
  [./mult_ACBulkF_3]
    type = KKSMultiACBulkF
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    gi_name   = g3
    eta_i     = eta3
    wi        = 1.0
    mob_name  = 1
    args      = 'c1 c2 c3 eta1 eta2'
  [../]
  [./mult_ACBulkC_3]
    type = KKSMultiACBulkC
    variable  = lambda
    Fj_names  = 'F1 F2 F3'
    hj_names  = 'h1 h2 h3'
    cj_names  = 'c1 c2 c3'
    eta_i     = eta3
    args      = 'eta1 eta2'
    mob_name  = 1
  [../]
  [./mult_CoupledACint_3]
    type = SimpleCoupledACInterface
    variable = lambda
    v = eta3
    kappa_name = kappa
    mob_name = 1
  [../]
  # Kernels for constraint equation eta1 + eta2 + eta3 = 1
  # eta3 is the nonlinear variable for the constraint equation
  [./eta3reaction]
    type = MatReaction
    variable = eta3
    mob_name = 1
  [../]
  [./eta1reaction]
    type = MatReaction
    variable = eta3
    v = eta1
    mob_name = 1
  [../]
  [./eta2reaction]
    type = MatReaction
    variable = eta3
    v = eta2
    mob_name = 1
  [../]
  [./one]
    type = BodyForce
    variable = eta3
    value = -1.0
  [../]
  # Phase concentration constraints
  [./chempot12]
    type = KKSPhaseChemicalPotential
    variable = c1
    cb       = c2
    fa_name  = F1
    fb_name  = F2
  [../]
  [./chempot23]
    type = KKSPhaseChemicalPotential
    variable = c2
    cb       = c3
    fa_name  = F2
    fb_name  = F3
  [../]
  [./phaseconcentration]
    type = KKSMultiPhaseConcentration
    variable = c3
    cj = 'c1 c2 c3'
    hj_names = 'h1 h2 h3'
    etas = 'eta1 eta2 eta3'
    c = c
  [../]
[]
[Executioner]
  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 = 30
  nl_max_its = 10
  l_tol = 1.0e-4
  nl_rel_tol = 1.0e-10
  nl_abs_tol = 1.0e-11
  num_steps = 1000
  [./TimeStepper]
    type = IterationAdaptiveDT
    dt = 0.2
    optimal_iterations = 10
    iteration_window = 2
  [../]
[]
[Preconditioning]
  active = 'full'
  [./full]
    type = SMP
    full = true
  [../]
  [./mydebug]
    type = FDP
    full = true
  [../]
[]
[Outputs]
  exodus = true
  checkpoint = true
  print_linear_residuals = false
  [./csv]
    type = CSV
    execute_on = 'final'
  [../]
[]
#[VectorPostprocessors]
#  [./c]
#    type =  LineValueSampler
#    start_point = '-25 0 0'
#    end_point = '25 0 0'
#    variable = c
#    num_points = 151
#    sort_by =  id
#    execute_on = timestep_end
#  [../]
#  [./eta1]
#    type =  LineValueSampler
#    start_point = '-25 0 0'
#    end_point = '25 0 0'
#    variable = eta1
#    num_points = 151
#    sort_by =  id
#    execute_on = timestep_end
#  [../]
#  [./eta2]
#    type =  LineValueSampler
#    start_point = '-25 0 0'
#    end_point = '25 0 0'
#    variable = eta2
#    num_points = 151
#    sort_by =  id
#    execute_on = timestep_end
#  [../]
#  [./eta3]
#    type =  LineValueSampler
#    start_point = '-25 0 0'
#    end_point = '25 0 0'
#    variable = eta3
#    num_points = 151
#    sort_by =  id
#    execute_on = timestep_end
#  [../]
#[]
(modules/phase_field/test/tests/phase_field_kernels/CoupledCoefAllenCahn.i)
#
# Test the CoefReaction kernel (which adds -L*v to the residual) for the case
# where v is a coupled variable
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 15
  ny = 15
  nz = 0
  xmin = 0
  xmax = 50
  ymin = 0
  ymax = 50
  zmin = 0
  zmax = 50
  elem_type = QUAD4
[]
[Variables]
  [./w]
  [../]
  [./eta]
    order = FIRST
    family = LAGRANGE
    [./InitialCondition]
      type = SmoothCircleIC
      x1 = 25.0
      y1 = 25.0
      radius = 6.0
      invalue = 1.0
      outvalue = 0.0
      int_width = 3.0
    [../]
  [../]
[]
[Kernels]
  [./detadt]
    type = TimeDerivative
    variable = eta
  [../]
  [./ACBulk]
    type = CoupledAllenCahn
    variable = w
    v = eta
    f_name = F
    mob_name = 1
  [../]
  [./W]
    type = MatReaction
    variable = w
    mob_name = -1
  [../]
  [./CoupledBulk]
    type = MatReaction
    variable = eta
    v = w
    mob_name = L
  [../]
  [./ACInterface]
    type = ACInterface
    variable = eta
    kappa_name = 1
    mob_name = L
    args = w
  [../]
[]
[Materials]
  [./mobility]
    type = DerivativeParsedMaterial
    f_name  = L
    args = 'eta w'
    function = '(1.5-eta)^2+(1.5-w)^2'
    derivative_order = 2
  [../]
  [./free_energy]
    type = DerivativeParsedMaterial
    f_name = F
    args = 'eta'
    function = 'eta^2 * (1-eta)^2'
    derivative_order = 2
  [../]
[]
[Preconditioning]
  [./smp]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  scheme = 'bdf2'
  solve_type = 'PJFNK'
  l_max_its = 15
  l_tol = 1.0e-4
  nl_max_its = 10
  nl_rel_tol = 1.0e-11
  start_time = 0.0
  num_steps = 2
  dt = 0.5
[]
[Outputs]
  hide = w
  exodus = true
[]
(modules/phase_field/test/tests/phase_field_kernels/CoupledAllenCahn.i)
#
# Test the coupled Allen-Cahn Bulk kernel
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 10
  ny = 10
  xmax = 12
  ymax = 12
  elem_type = QUAD4
[]
[Variables]
  [./w]
  [../]
  [./eta]
    order = FIRST
    family = LAGRANGE
    [./InitialCondition]
      type = SmoothCircleIC
      x1 = 0.0
      y1 = 0.0
      radius = 6.0
      invalue = 0.9
      outvalue = 0.1
      int_width = 3.0
    [../]
  [../]
[]
[Kernels]
  [./detadt]
    type = TimeDerivative
    variable = eta
  [../]
  [./ACBulk]
    type = CoupledAllenCahn
    variable = w
    v = eta
    f_name = F
  [../]
  [./W]
    type = Reaction
    variable = w
  [../]
  [./CoupledBulk]
    type = MatReaction
    variable = eta
    v = w
  [../]
  [./ACInterface]
    type = ACInterface
    variable = eta
    kappa_name = 1
  [../]
[]
[Materials]
  [./consts]
    type = GenericConstantMaterial
    prop_names  = 'L'
    prop_values = '1'
  [../]
  [./free_energy]
    type = DerivativeParsedMaterial
    f_name = F
    args = 'eta'
    function = '2 * eta^2 * (1-eta)^2 - 0.2*eta'
    derivative_order = 2
  [../]
[]
[Preconditioning]
  [./smp]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  scheme = 'bdf2'
  solve_type = 'PJFNK'
  l_max_its = 15
  l_tol = 1.0e-4
  nl_max_its = 10
  nl_rel_tol = 1.0e-11
  start_time = 0.0
  num_steps = 2
  dt = 0.5
[]
[Debug]
  show_var_residual_norms = true
[]
[Outputs]
  hide = w
  file_base = AllenCahn_out
  exodus = true
[]
(modules/phase_field/test/tests/mobility_derivative/coupledmatdiffusion.i)
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 15
  ny = 15
  xmax = 15.0
  ymax = 15.0
  elem_type = QUAD4
[]
[Variables]
  [./c]
    [./InitialCondition]
      type = CrossIC
      x1 = 0.0
      x2 = 30.0
      y1 = 0.0
      y2 = 30.0
    [../]
  [../]
  [./d]
    [./InitialCondition]
      type = SmoothCircleIC
      x1 = 15
      y1 = 15
      radius = 8
      int_width = 3
      invalue = 2
      outvalue = 0
    [../]
  [../]
  [./u]
  [../]
  [./w]
  [../]
[]
[Kernels]
  [./ctime]
    type = TimeDerivative
    variable = c
  [../]
  [./umat]
    type = MatReaction
    variable = c
    v = u
    mob_name = 1
  [../]
  [./urxn]
    type = Reaction
    variable = u
  [../]
  [./cres]
    type = MatDiffusion
    variable = u
    diffusivity = Dc
    args = d
    v = c
  [../]
  [./dtime]
    type = TimeDerivative
    variable = d
  [../]
  [./wmat]
    type = MatReaction
    variable = d
    v = w
    mob_name = 1
  [../]
  [./wrxn]
    type = Reaction
    variable = w
  [../]
  [./dres]
    type = MatDiffusion
    variable = w
    diffusivity = Dd
    args = c
    v = d
  [../]
[]
[Materials]
  [./Dc]
    type = DerivativeParsedMaterial
    f_name = Dc
    function = '0.01+c^2+d'
    args = 'c d'
    derivative_order = 1
  [../]
  [./Dd]
    type = DerivativeParsedMaterial
    f_name = Dd
    function = 'd^2+c+1.5'
    args = 'c d'
    derivative_order = 1
  [../]
[]
[Preconditioning]
  [./SMP]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  scheme = 'BDF2'
  solve_type = 'NEWTON'
  petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -pc_asm_overlap'
  petsc_options_value = 'asm      31                  lu           1'
  dt = 1
  num_steps = 2
[]
[Outputs]
  exodus = true
[]
(modules/phase_field/test/tests/SimpleACInterface/SimpleCoupledACInterface.i)
#
# Test the coupled Allen-Cahn Bulk kernel
#
[Mesh]
  type = GeneratedMesh
  dim = 2
  nx = 20
  ny = 20
  nz = 0
  xmin = 0
  xmax = 50
  ymin = 0
  ymax = 50
  zmin = 0
  zmax = 50
  elem_type = QUAD4
  uniform_refine = 1
[]
[Variables]
  [./w]
  [../]
  [./eta]
    order = FIRST
    family = LAGRANGE
    [./InitialCondition]
      type = SmoothCircleIC
      x1 = 25.0
      y1 = 25.0
      radius = 6.0
      invalue = 1.0
      outvalue = 0.0
      int_width = 5.0
    [../]
  [../]
[]
[Kernels]
  [./detadt]
    type = TimeDerivative
    variable = eta
  [../]
  [./ACBulk]
    type = AllenCahn
    variable = eta
    f_name = F
  [../]
  [./CoupledBulk]
    type = MatReaction
    variable = eta
    v = w
  [../]
  [./W]
    type = Reaction
    variable = w
  [../]
  [./CoupledACInterface]
    type = SimpleCoupledACInterface
    variable = w
    v = eta
    kappa_name = 1
  [../]
[]
[Materials]
  [./consts]
    type = GenericConstantMaterial
    prop_names  = 'L'
    prop_values = '1'
  [../]
  [./free_energy]
    type = DerivativeParsedMaterial
    f_name = F
    args = 'eta'
    function = 'eta^2 * (1-eta)^2'
    derivative_order = 2
  [../]
[]
[Preconditioning]
  [./smp]
    type = SMP
    full = true
  [../]
[]
[Executioner]
  type = Transient
  scheme = 'bdf2'
  solve_type = 'PJFNK'
  l_max_its = 15
  l_tol = 1.0e-4
  nl_max_its = 10
  nl_rel_tol = 1.0e-11
  start_time = 0.0
  num_steps = 2
  dt = 2
[]
[Debug]
  show_var_residual_norms = true
[]
[Outputs]
  hide = w
  exodus = true
[]