- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Controllable:No
Description:The name of the variable that this residual object operates on
ADTimeDerivative
Description
The ADTimeDerivative
kernel implements a simple time derivative for the domain given by
where the second term on the left hand side corresponds to the strong forms of other kernels. The corresponding ADTimeDerivative
weak form using inner-product notation is
where is the approximate solution and is a finite element test function.
The Jacobian is computed through automatic differentiation. More information about time kernels can be found on the Kernels description page.
Example Syntax
Time derivative terms are ubiquitous in any transient simulation. The kernel block for a transient diffusion problem that demonstrates the ADTimeDerivative
syntax is shown below:
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[]
(test/tests/kernels/ad_transient_diffusion/ad_transient_diffusion.i)Input Parameters
- blockThe list of blocks (ids or names) that this object will be applied
C++ Type:std::vector<SubdomainName>
Controllable:No
Description:The list of blocks (ids or names) that this object will be applied
- displacementsThe displacements
C++ Type:std::vector<VariableName>
Controllable:No
Description:The displacements
- prop_getter_suffixAn optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
C++ Type:MaterialPropertyName
Controllable:No
Description:An optional suffix parameter that can be appended to any attempt to retrieve/get material properties. The suffix will be prepended with a '_' character.
Optional Parameters
- control_tagsAdds user-defined labels for accessing object parameters via control logic.
C++ Type:std::vector<std::string>
Controllable:No
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>
Controllable:No
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
Controllable:Yes
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
Controllable:No
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>
Controllable:No
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
Controllable:No
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
Controllable:No
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>
Controllable:No
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>
Controllable:No
Description:The extra tags for the vectors this Kernel should fill
- matrix_tagssystem timeThe tag for the matrices this Kernel should fill
Default:system time
C++ Type:MultiMooseEnum
Options:nontime, system, time
Controllable:No
Description:The tag for the matrices this Kernel should fill
- vector_tagstimeThe tag for the vectors this Kernel should fill
Default:time
C++ Type:MultiMooseEnum
Options:nontime, time
Controllable:No
Description:The tag for the vectors this Kernel should fill
Tagging Parameters
Input Files
- (test/tests/misc/ad_robustness/ad_two_nl_var_transient_diffusion_jac.i)
- (modules/phase_field/test/tests/ad_coupled_gradient_dot/diffusion_rate.i)
- (modules/stochastic_tools/examples/parameter_study/diffusion_time.i)
- (test/tests/kernels/ad_jacobians/adfunction.i)
- (modules/phase_field/test/tests/phase_field_kernels/ADAllenCahn.i)
- (modules/phase_field/test/tests/ADCHSoretDiffusion/simple_transient_diffusion_with_soret.i)
- (test/tests/nodalkernels/constraint_enforcement/ad-upper-and-lower-bound.i)
- (modules/tensor_mechanics/test/tests/ad_viscoplasticity_stress_update/lps_dual.i)
- (test/tests/materials/derivative_material_interface/ad_construction_order.i)
- (test/tests/materials/derivative_sum_material/ad_random_ic.i)
- (test/tests/misc/ad_robustness/ad_two_nl_var_transient_diffusion.i)
- (modules/heat_conduction/test/tests/radiative_bcs/ad_function_radiative_bc.i)
- (modules/stochastic_tools/examples/parameter_study/diffusion.i)
- (modules/stochastic_tools/test/tests/transfers/batch_sampler_transfer/sub.i)
- (test/tests/materials/coupled_value_function/adjac.i)
- (modules/stochastic_tools/examples/parameter_study/diffusion_vector.i)
- (modules/stochastic_tools/test/tests/multiapps/batch_full_solve_multiapp/sub.i)
- (python/mms/test/mms_temporal.i)
- (test/tests/time_integrators/newmark-beta/ad_newmark_beta_dotdot.i)
- (test/tests/bcs/mat_neumann_bc/ad_mat_neumann.i)
- (test/tests/reporters/iteration_info/iteration_info.i)
- (modules/phase_field/test/tests/phase_field_kernels/ADAllenCahnVariableL.i)
- (modules/stochastic_tools/examples/sobol/diffusion.i)
- (test/tests/kernels/ad_transient_diffusion/ad_transient_diffusion.i)
- (test/tests/time_integrators/central-difference/ad_central_difference_dotdot.i)
- (modules/phase_field/test/tests/ADCHSplitChemicalPotential/simple_transient_diffusion.i)
- (test/tests/kernels/ad_jacobians/test.i)
- (modules/phase_field/test/tests/anisotropic_mobility/ad_diffusion.i)
- (modules/phase_field/examples/anisotropic_interfaces/ad_snow.i)
- (test/tests/outputs/xml/xml_iterations.i)
- (test/tests/multiapps/full_solve_multiapp_reset/sub.i)
- (modules/xfem/test/tests/moving_interface/ad_phase_transition_2d.i)
- (modules/heat_conduction/test/tests/ad_convective_heat_flux/coupled.i)
- (modules/phase_field/test/tests/anisotropic_interfaces/adkobayashi.i)
- (modules/level_set/test/tests/functions/olsson_bubble/olsson_bubble_adjac.i)
- (test/tests/misc/ad_robustness/ad_two_var_transient_diffusion.i)
- (modules/thermal_hydraulics/test/tests/components/hs_boundary_external_app_temperature/phy.master.i)
- (modules/phase_field/test/tests/automatic_differentiation/admatreaction.i)
- (modules/thermal_hydraulics/test/tests/components/hs_boundary_external_app_convection/plate.master.i)
- (modules/xfem/test/tests/moving_interface/moving_ad_diffusion.i)
- (modules/stochastic_tools/examples/batch/sub.i)
- (modules/heat_conduction/test/tests/ad_convective_heat_flux/equilibrium.i)
Child Objects
(test/tests/kernels/ad_transient_diffusion/ad_transient_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/misc/ad_robustness/ad_two_nl_var_transient_diffusion_jac.i)
penalty=1
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2
[]
[Variables]
[./u]
family = MONOMIAL
order = FIRST
[../]
[v]
family = MONOMIAL
order = FIRST
[]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[coupled]
type = ADCoupledValueTest
variable = u
v = v
[]
[v_diff]
type = Diffusion
variable = v
[]
[]
[DGKernels]
[dummy]
type = ADDGCoupledTest
variable = u
v = v
[]
[]
[BCs]
[./left]
type = PenaltyDirichletBC
variable = u
boundary = left
value = 0
penalty = ${penalty}
[../]
[./right]
type = PenaltyDirichletBC
variable = u
boundary = right
value = 1
penalty = ${penalty}
[../]
[./left_v]
type = PenaltyDirichletBC
variable = v
boundary = left
value = 0
penalty = ${penalty}
[../]
[./right_v]
type = PenaltyDirichletBC
variable = v
boundary = right
value = 1
penalty = ${penalty}
[../]
[]
[Executioner]
type = Transient
num_steps = 2
dt = 0.1
dtmin = 0.1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[dof_map]
type = DOFMap
execute_on = 'initial'
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
(modules/phase_field/test/tests/ad_coupled_gradient_dot/diffusion_rate.i)
# Solves the problem
# -mu * Lap(u_dot) + u_dot = alpha * Lap(u) - 2*u*(1-3*u+2*u^2)
# for mu = 1 and alpha = 0.01
# (see appendix B of A. Guevel et al. "Viscous phase-field modeling for chemo-mechanical microstructural evolution: application to geomaterials and pressure solution." In print.)
n_elem = 100
alpha = 0.01
mu = 1
[Mesh]
type = GeneratedMesh
dim = 1
xmin = 0
xmax = 1
nx = '${n_elem}'
elem_type = EDGE2
[]
[Variables]
[u]
[InitialCondition]
type = ConstantIC
value = 0.51
[]
[]
[]
[Kernels]
[Lap]
type = ADMatDiffusion
variable = u
diffusivity = '${alpha}'
[]
[LapDot]
type = ADDiffusionRate
variable = u
mu = '${mu}'
[]
[Reac]
type = ADMatReaction
variable = u
mob_name = L
[]
[Visc]
type = ADTimeDerivative
variable = u
[]
[]
[Materials]
[parsed]
type = ADParsedMaterial
function = '-2*(1-3*u+2*u*u)'
args = 'u'
f_name = 'L'
[]
[]
[BCs]
[both]
type = ADDirichletBC
variable = u
value = 0.51
boundary = 'left right'
[]
[]
[Postprocessors]
[mid_u]
type = PointValue
variable = u
point = '0.5 0 0'
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
num_steps = 1000
dt = 0.1
nl_abs_tol = 1e-9
[]
[Outputs]
print_linear_residuals = false
[csv]
type = CSV
file_base = 'diffusion_rate'
[]
[]
(modules/stochastic_tools/examples/parameter_study/diffusion_time.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables/T]
initial_condition = 300
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADMatDiffusion
variable = T
diffusivity = diffusivity
[]
[source]
type = ADBodyForce
variable = T
value = 100
function = 1
[]
[]
[BCs]
[left]
type = ADDirichletBC
variable = T
boundary = left
value = 300
[]
[right]
type = ADNeumannBC
variable = T
boundary = right
value = -100
[]
[]
[Materials/constant]
type = ADGenericConstantMaterial
prop_names = 'diffusivity'
prop_values = 1
[]
[Executioner]
type = Transient
solve_type = NEWTON
num_steps = 4
dt = 0.25
[]
[Postprocessors]
[T_avg]
type = ElementAverageValue
variable = T
execute_on = 'initial timestep_end'
[]
[q_left]
type = ADSideDiffusiveFluxAverage
variable = T
boundary = left
diffusivity = diffusivity
execute_on = 'initial timestep_end'
[]
[]
[VectorPostprocessors]
[T_vec]
type = LineValueSampler
variable = T
start_point = '0 0.5 0'
end_point = '1 0.5 0'
num_points = 11
sort_by = x
execute_on = 'initial timestep_end'
[]
[]
[Controls/stochastic]
type = SamplerReceiver
[]
[Outputs]
[]
(test/tests/kernels/ad_jacobians/adfunction.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
displacements = 'disp_x disp_y'
[]
[Problem]
kernel_coverage_check = false
[]
[Variables]
[c]
[]
[disp_x]
[]
[disp_y]
[]
[]
[Kernels]
[source]
type = ADBodyForce
variable = c
function = source_func
use_displaced_mesh = true
displacements = ''
#displacements = 'disp_x disp_y'
[]
[dt]
type = ADTimeDerivative
variable = c
[]
[]
[Functions]
[source_func]
type = ParsedFunction
value = 'x + y^2'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
num_steps = 1
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/phase_field_kernels/ADAllenCahn.i)
#
# Test the forward automatic differentiation Allen-Cahn Bulk kernel
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
xmax = 12
ymax = 12
elem_type = QUAD4
[]
[Variables]
[./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 = ADTimeDerivative
variable = eta
[../]
[./ACBulk]
type = ADAllenCahn
variable = eta
f_name = F
[../]
[./ACInterface]
type = ADACInterface
variable = eta
kappa_name = 1
variable_L = false
[../]
[]
[Materials]
[./consts]
type = ADGenericConstantMaterial
prop_names = 'L'
prop_values = '1'
[../]
[./free_energy]
type = ADTestDerivativeFunction
function = F1
f_name = F
op = 'eta'
[../]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
solve_type = 'NEWTON'
num_steps = 2
dt = 0.5
[]
[Outputs]
exodus = true
[]
(modules/phase_field/test/tests/ADCHSoretDiffusion/simple_transient_diffusion_with_soret.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
[./c]
[../]
[./mu]
[../]
[]
[AuxVariables]
[./T]
[./InitialCondition]
type = RampIC
value_left = 900
value_right = 1000
[../]
[../]
[]
[Kernels]
[./conc]
type = ADCHSplitConcentration
variable = c
chemical_potential_var = mu
mobility = chemical_mobility_prop
[../]
[./chempot]
type = ADCHSplitChemicalPotential
variable = mu
chemical_potential = mu_prop
[../]
[./soret]
type = ADCHSoretMobility
variable = c
T = T
mobility = thermal_mobility_prop
[../]
[./time]
type = ADTimeDerivative
variable = c
[../]
[]
[Materials]
[./chemical_potential]
type = ADPiecewiseLinearInterpolationMaterial
property = mu_prop
variable = c
x = '0 1'
y = '0 1'
[../]
[./chemical_mobility_prop]
type = ADGenericConstantMaterial
prop_names = chemical_mobility_prop
prop_values = 0.1
[../]
[./thermal_mobility_prop]
type = ADGenericConstantMaterial
prop_names = thermal_mobility_prop
prop_values = -20
[../]
[]
[BCs]
[./leftc]
type = DirichletBC
variable = c
boundary = left
value = 0
[../]
[./rightc]
type = DirichletBC
variable = c
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -ksp_grmres_restart -sub_ksp_type -sub_pc_type -pc_asm_overlap'
petsc_options_value = 'asm 31 preonly lu 2'
dt = 0.1
num_steps = 20
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/nodalkernels/constraint_enforcement/ad-upper-and-lower-bound.i)
l=10
nx=100
num_steps=10
[Mesh]
type = GeneratedMesh
dim = 1
xmax = ${l}
nx = ${nx}
[]
[Variables]
[u]
[]
[lm_upper]
[]
[lm_lower]
[]
[]
[ICs]
[u]
type = FunctionIC
variable = u
function = 'x'
[]
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = u
[]
[diff]
type = ADDiffusion
variable = u
[]
[ffn]
type = ADBodyForce
variable = u
function = 'if(x<5,-1,1)'
[]
[]
[NodalKernels]
[upper_bound]
type = ADUpperBoundNodalKernel
variable = lm_upper
v = u
exclude_boundaries = 'left right'
upper_bound = 10
[]
[forces_from_upper]
type = ADCoupledForceNodalKernel
variable = u
v = lm_upper
coef = -1
[]
[lower_bound]
type = ADLowerBoundNodalKernel
variable = lm_lower
v = u
exclude_boundaries = 'left right'
lower_bound = 0
[]
[forces_from_lower]
type = ADCoupledForceNodalKernel
variable = u
v = lm_lower
coef = 1
[]
[]
[BCs]
[left]
type = ADDirichletBC
boundary = left
value = 0
variable = u
[]
[right]
type = ADDirichletBC
boundary = right
value = ${l}
variable = u
[]
[]
[Preconditioning]
[smp]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
num_steps = ${num_steps}
solve_type = NEWTON
dtmin = 1
petsc_options_iname = '-snes_max_linear_solve_fail -ksp_max_it -pc_type -sub_pc_factor_levels -snes_linesearch_type'
petsc_options_value = '0 30 asm 16 basic'
[]
[Outputs]
exodus = true
[csv]
type = CSV
execute_on = 'nonlinear timestep_end'
[]
[dof]
type = DOFMap
execute_on = 'initial'
[]
[]
[Debug]
show_var_residual_norms = true
[]
[Postprocessors]
[active_upper_lm]
type = GreaterThanLessThanPostprocessor
variable = lm_upper
execute_on = 'nonlinear timestep_end'
value = 1e-8
comparator = 'greater'
[]
[upper_violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = ${fparse 10+1e-8}
comparator = 'greater'
[]
[active_lower_lm]
type = GreaterThanLessThanPostprocessor
variable = lm_lower
execute_on = 'nonlinear timestep_end'
value = 1e-8
comparator = 'greater'
[]
[lower_violations]
type = GreaterThanLessThanPostprocessor
variable = u
execute_on = 'nonlinear timestep_end'
value = -1e-8
comparator = 'less'
[]
[nls]
type = NumNonlinearIterations
[]
[cum_nls]
type = CumulativeValuePostprocessor
postprocessor = nls
[]
[]
(modules/tensor_mechanics/test/tests/ad_viscoplasticity_stress_update/lps_dual.i)
# This test provides an example of combining two LPS viscoplasticity models with different stress
# exponents.
[GlobalParams]
displacements = 'disp_x disp_y'
[]
[Mesh]
type = GeneratedMesh
dim = 2
nx = 1
ny = 1
xmax = 0.002
ymax = 0.002
[]
[Variables]
[./temp]
initial_condition = 1000
[../]
[]
[Kernels]
[./dt]
type = ADTimeDerivative
variable = temp
[../]
[./diff]
type = ADDiffusion
variable = temp
[../]
[]
[Modules/TensorMechanics/Master/All]
strain = FINITE
add_variables = true
generate_output = 'strain_xx strain_yy strain_xy hydrostatic_stress vonmises_stress'
use_automatic_differentiation = true
[]
[Functions]
[./pull]
type = PiecewiseLinear
x = '0 0.1'
y = '0 1e-5'
[../]
[./tot_effective_viscoplasticity]
type = ParsedFunction
vals = 'lps_1_eff_creep_strain lps_3_eff_creep_strain'
vars = 'lps_1_eff_creep_strain lps_3_eff_creep_strain'
value = 'lps_1_eff_creep_strain+lps_3_eff_creep_strain'
[../]
[]
[Materials]
[./elasticity_tensor]
type = ADComputeIsotropicElasticityTensor
youngs_modulus = 1e10
poissons_ratio = 0.3
[../]
[./stress]
type = ADComputeMultipleInelasticStress
inelastic_models = 'one two'
outputs = all
[../]
[./porosity]
type = ADPorosityFromStrain
initial_porosity = 0.1
inelastic_strain = 'combined_inelastic_strain'
outputs = 'all'
[../]
[./one]
type = ADViscoplasticityStressUpdate
coefficient = 'coef_3'
power = 3
base_name = 'lps_1'
outputs = all
relative_tolerance = 1e-11
[../]
[./two]
type = ADViscoplasticityStressUpdate
coefficient = 1e-10
power = 1
base_name = 'lps_3'
outputs = all
relative_tolerance = 1e-11
[../]
[./coef]
type = ADParsedMaterial
f_name = coef_3
# Example of creep power law
args = temp
function = '0.5e-18 * exp(-4e4 / 1.987 / temp)'
[../]
[]
[BCs]
[./no_disp_x]
type = ADDirichletBC
variable = disp_x
boundary = left
value = 0.0
[../]
[./no_disp_y]
type = ADDirichletBC
variable = disp_y
boundary = bottom
value = 0.0
[../]
[./pull_disp_y]
type = ADFunctionDirichletBC
variable = disp_y
boundary = top
function = pull
[../]
[./temp_ramp]
type = ADFunctionDirichletBC
boundary = right
function = '1000 + 400 * t / 0.12'
variable = temp
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 0.01
end_time = 0.12
[]
[Postprocessors]
[./disp_x]
type = SideAverageValue
variable = disp_x
boundary = right
[../]
[./disp_y]
type = SideAverageValue
variable = disp_y
boundary = top
[../]
[./avg_hydro]
type = ElementAverageValue
variable = hydrostatic_stress
[../]
[./avg_vonmises]
type = ElementAverageValue
variable = vonmises_stress
[../]
[./dt]
type = TimestepSize
[../]
[./num_lin]
type = NumLinearIterations
outputs = console
[../]
[./num_nonlin]
type = NumNonlinearIterations
outputs = console
[../]
[./lps_1_eff_creep_strain]
type = ElementAverageValue
variable = lps_1_effective_viscoplasticity
[../]
[./lps_3_eff_creep_strain]
type = ElementAverageValue
variable = lps_3_effective_viscoplasticity
[../]
[./lps_1_gauge_stress]
type = ElementAverageValue
variable = lps_1_gauge_stress
[../]
[./lps_3_gauge_stress]
type = ElementAverageValue
variable = lps_3_gauge_stress
[../]
[./eff_creep_strain_tot]
type = FunctionValuePostprocessor
function = tot_effective_viscoplasticity
[../]
[./porosity]
type = ElementAverageValue
variable = porosity
[../]
[]
[Outputs]
csv = true
[]
(test/tests/materials/derivative_material_interface/ad_construction_order.i)
#
# Test the the getDefaultMaterialProperty in DerivativeMaterialInterface.
# This test should only pass, if the construction order of the Materials
# using this interface does not influence the outcome.
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 5
ny = 1
xmin = 0
xmax = 1
ymin = 0
ymax = 0.1
elem_type = QUAD4
[]
[GlobalParams]
derivative_order = 2
[]
[Variables]
[./c]
[./InitialCondition]
type = FunctionIC
function = x
[../]
[../]
[]
[Kernels]
[./dummy1]
type = ADDiffusion
variable = c
[../]
[./dummy2]
type = ADTimeDerivative
variable = c
[../]
[]
[Materials]
# derivatives used both before and after being declared
[./sum_a_1]
type = ADDerivativeSumMaterial
f_name = Fa1
sum_materials = 'Fa'
args = 'c'
outputs = exodus
[../]
[./free_energy_a]
type = ADDerivativeParsedMaterial
f_name = Fa
args = 'c'
function = 'c^4'
[../]
[./sum_a_2]
type = ADDerivativeSumMaterial
f_name = Fa2
sum_materials = 'Fa'
args = 'c'
outputs = exodus
[../]
# derivatives declared after being used
[./sum_b_1]
type = ADDerivativeSumMaterial
f_name = Fb1
sum_materials = 'Fb'
args = 'c'
outputs = exodus
[../]
[./free_energy_b]
type = ADDerivativeParsedMaterial
f_name = Fb
args = 'c'
function = 'c^4'
[../]
# derivatives declared before being used
[./free_energy_c]
type = ADDerivativeParsedMaterial
f_name = Fc
args = 'c'
function = 'c^4'
[../]
[./sum_c_2]
type = ADDerivativeSumMaterial
f_name = Fc2
sum_materials = 'Fc'
args = 'c'
outputs = exodus
[../]
# non-existing derivatives
[./free_energy_d]
type = ADParsedMaterial
f_name = Fd
args = 'c'
function = 'c^4'
[../]
[./sum_d_1]
type = ADDerivativeSumMaterial
f_name = Fd1
sum_materials = 'Fd'
args = 'c'
outputs = exodus
[../]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = 'NEWTON'
num_steps = 1
dt = 1e-5
[]
[Outputs]
execute_on = 'timestep_end'
exodus = true
[]
(test/tests/materials/derivative_sum_material/ad_random_ic.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
xmax = 250
ymax = 250
elem_type = QUAD4
[]
[Variables]
[./c]
[./InitialCondition]
type = RandomIC
[../]
[../]
[]
[Kernels]
[./w_res]
type = ADDiffusion
variable = c
[../]
[./time]
type = ADTimeDerivative
variable = c
[../]
[]
[Materials]
[./free_energy1]
type = ADDerivativeParsedMaterial
f_name = Fa
args = 'c'
function = (c-0.1)^4*(1-0.1-c)^4
[../]
[./free_energy2]
type = ADDerivativeParsedMaterial
f_name = Fb
args = 'c'
function = -0.25*(c-0.1)^4*(1-0.1-c)^4
[../]
# Fa+Fb+Fb == Fc
[./free_energy3]
type = ADDerivativeParsedMaterial
f_name = Fc
args = 'c'
function = 0.5*(c-0.1)^4*(1-0.1-c)^4
outputs = all
[../]
[./dfree_energy3]
type = ADDerivativeParsedMaterial
f_name = dFc
args = 'c'
material_property_names = 'F:=D[Fc,c]'
function = F
outputs = all
[../]
[./d2free_energy3]
type = ADDerivativeParsedMaterial
f_name = d2Fc
args = 'c'
material_property_names = 'F:=D[Fc,c,c]'
function = F
outputs = all
[../]
[./free_energy]
type = ADDerivativeSumMaterial
f_name = F_sum
sum_materials = 'Fa Fb Fb'
args = 'c'
outputs = all
[../]
[./dfree_energy]
type = ADDerivativeParsedMaterial
f_name = dF_sum
material_property_names = 'F:=D[F_sum,c]'
function = F
args = 'c'
outputs = all
[../]
[./d2free_energy]
type = ADDerivativeParsedMaterial
f_name = d2F_sum
material_property_names = 'F:=D[F_sum,c,c]'
function = F
args = 'c'
outputs = all
[../]
[]
[Executioner]
type = Transient
num_steps = 1
[]
[Postprocessors]
[./F_sum]
type = ElementAverageValue
variable = F_sum
[../]
[./F_check]
type = ElementAverageValue
variable = Fc
[../]
[./dF_sum]
type = ElementAverageValue
variable = dF_sum
[../]
[./dF_check]
type = ElementAverageValue
variable = dFc
[../]
[./d2F_sum]
type = ElementAverageValue
variable = d2F_sum
[../]
[./d2F_check]
type = ElementAverageValue
variable = d2Fc
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/misc/ad_robustness/ad_two_nl_var_transient_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[./u]
[../]
[v][]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[coupled]
type = ADCoupledValueTest
variable = u
v = v
[]
[v_diff]
type = Diffusion
variable = v
[]
[]
[DGKernels]
[dummy]
type = ADDGCoupledTest
variable = u
v = v
[]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[./left_v]
type = DirichletBC
variable = v
boundary = left
value = 0
[../]
[./right_v]
type = DirichletBC
variable = v
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 2
dt = 0.1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[dof_map]
type = DOFMap
execute_on = 'initial'
[]
[]
(modules/heat_conduction/test/tests/radiative_bcs/ad_function_radiative_bc.i)
#
# If we assume that epsilon*sigma*(T_inf^4-T_s^4) is approximately equal to
# epsilon*sigma*4*T_inf^3*(T_inf-T_s), that form is equivalent to
# h*(T_inf-T_s), the convective flux bc. So, the radiative and convective
# flux bcs should give nearly the same answer if the leading terms are equal.
#
[Mesh]
[top]
type = GeneratedMeshGenerator
dim = 3
nx = 10
bias_x = 0.8
ymin = 1.2
ymax = 2.2
boundary_name_prefix = top
[]
[bottom]
type = GeneratedMeshGenerator
dim = 3
nx = 10
bias_x = 0.8
boundary_name_prefix = bot
boundary_id_offset = 6
[]
[two_blocks]
type = MeshCollectionGenerator
inputs = 'top bottom'
[]
[]
[Variables]
[temp]
initial_condition = 600.0
[]
[]
[Kernels]
[heat_dt]
type = ADTimeDerivative
variable = temp
[]
[heat_conduction]
type = ADHeatConduction
variable = temp
[]
[]
[BCs]
[top_right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = top_right
T_infinity = 300.0
heat_transfer_coefficient = 3.0
[]
[bot_right]
type = ADFunctionRadiativeBC
variable = temp
boundary = bot_right
# htc/(stefan-boltzmann*4*T_inf^3)
emissivity_function = '3/(5.670367e-8*4*300*300*300)'
[]
[]
[Materials]
[thermal]
type = ADGenericConstantMaterial
prop_names = 'density thermal_conductivity specific_heat'
prop_values = '1 10 100'
[]
[]
[Postprocessors]
[top_left_temp]
type = SideAverageValue
variable = temp
boundary = top_left
execute_on = 'TIMESTEP_END initial'
[]
[bot_left_temp]
type = SideAverageValue
variable = temp
boundary = bot_left
execute_on = 'TIMESTEP_END initial'
[]
[top_right_temp]
type = SideAverageValue
variable = temp
boundary = top_right
[]
[bot_right_temp]
type = SideAverageValue
variable = temp
boundary = bot_right
[]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1e1
nl_abs_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(modules/stochastic_tools/examples/parameter_study/diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables/T]
initial_condition = 300
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADMatDiffusion
variable = T
diffusivity = diffusivity
[]
[source]
type = ADBodyForce
variable = T
value = 100
function = 1
[]
[]
[BCs]
[left]
type = ADDirichletBC
variable = T
boundary = left
value = 300
[]
[right]
type = ADNeumannBC
variable = T
boundary = right
value = -100
[]
[]
[Materials/constant]
type = ADGenericConstantMaterial
prop_names = 'diffusivity'
prop_values = 1
[]
[Executioner]
type = Transient
solve_type = NEWTON
num_steps = 4
dt = 0.25
[]
[Postprocessors]
[T_avg]
type = ElementAverageValue
variable = T
execute_on = 'initial timestep_end'
[]
[q_left]
type = ADSideDiffusiveFluxAverage
variable = T
boundary = left
diffusivity = diffusivity
execute_on = 'initial timestep_end'
[]
[]
[Controls/stochastic]
type = SamplerReceiver
[]
[Outputs]
[]
(modules/stochastic_tools/test/tests/transfers/batch_sampler_transfer/sub.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Postprocessors]
[average]
type = AverageNodalVariableValue
variable = u
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 0.25
solve_type = NEWTON
[]
[Controls]
[stochastic]
type = SamplerReceiver
[]
[]
[Outputs]
[]
(test/tests/materials/coupled_value_function/adjac.i)
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 5
ny = 5
[]
[]
[Variables]
[u]
initial_condition = 0.1
[]
[v]
initial_condition = 0.1
[]
[]
[Materials]
[Du]
type = ADCoupledValueFunctionMaterial
function = x
v = v
prop_name = Du
[]
[Dv]
type = ADCoupledValueFunctionMaterial
function = x^2
v = u
prop_name = Dv
[]
[]
[Kernels]
[diff_u]
type = ADMatDiffusion
diffusivity = Du
variable = u
[]
[dudt]
type = ADTimeDerivative
variable = u
[]
[diff_v]
type = ADMatDiffusion
diffusivity = Dv
variable = v
[]
[dvdt]
type = ADTimeDerivative
variable = v
[]
[]
[BCs]
[u_left]
type = DirichletBC
boundary = left
variable = u
value = 1
[]
[u_right]
type = DirichletBC
boundary = right
variable = u
value = 0.1
[]
[v_top]
type = DirichletBC
boundary = top
variable = v
value = 1
[]
[v_bottom]
type = DirichletBC
boundary = bottom
variable = v
value = 0.1
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 0.1
num_steps = 4
[]
[Outputs]
exodus = true
[]
(modules/stochastic_tools/examples/parameter_study/diffusion_vector.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables/T]
initial_condition = 300
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADMatDiffusion
variable = T
diffusivity = diffusivity
[]
[source]
type = ADBodyForce
variable = T
value = 100
function = 1
[]
[]
[BCs]
[left]
type = ADDirichletBC
variable = T
boundary = left
value = 300
[]
[right]
type = ADNeumannBC
variable = T
boundary = right
value = -100
[]
[]
[Materials/constant]
type = ADGenericConstantMaterial
prop_names = 'diffusivity'
prop_values = 1
[]
[Executioner]
type = Transient
solve_type = NEWTON
num_steps = 4
dt = 0.25
[]
[Postprocessors]
[T_avg]
type = ElementAverageValue
variable = T
execute_on = 'initial timestep_end'
[]
[q_left]
type = ADSideDiffusiveFluxAverage
variable = T
boundary = left
diffusivity = diffusivity
execute_on = 'initial timestep_end'
[]
[]
[Reporters]
[acc]
type = AccumulateReporter
reporters = 'T_avg/value q_left/value'
[]
[]
[Controls/stochastic]
type = SamplerReceiver
[]
[Outputs]
[]
(modules/stochastic_tools/test/tests/multiapps/batch_full_solve_multiapp/sub.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Postprocessors]
[average]
type = AverageNodalVariableValue
variable = u
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 0.25
solve_type = NEWTON
[]
[Outputs]
exodus = true
[]
(python/mms/test/mms_temporal.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 8
ny = 8
[]
[Variables]
[u][]
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = u
[]
[diff]
type = ADDiffusion
variable = u
[]
[force]
type = BodyForce
variable = u
function = force
[]
[]
[Functions]
[exact]
type = ParsedFunction
value = 't^3*x*y'
[]
[force]
type = ParsedFunction
value = '3*x*y*t^2'
[]
[]
[BCs]
[all]
type = FunctionDirichletBC
variable = u
function = exact
boundary = 'left right top bottom'
[]
[]
[Postprocessors]
[error]
type = ElementL2Error
function = exact
variable = u
[]
[h]
type = AverageElementSize
[]
[]
[Executioner]
type = Transient
dt = 1
end_time = 3
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
csv = true
[]
(test/tests/time_integrators/newmark-beta/ad_newmark_beta_dotdot.i)
###########################################################
# This is a simple test with a time-dependent problem
# demonstrating the use of the TimeIntegrator system.
#
# Testing that the second time derivative is calculated
# correctly using the Newmark-Beta method for an AD variable
#
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 1
ny = 1
[]
[Variables]
[./u]
[../]
[]
[Functions]
[./forcing_fn]
type = PiecewiseLinear
x = '0.0 0.1 0.2 0.3 0.4 0.5 0.6'
y = '0.0 0.0 0.0025 0.01 0.0175 0.02 0.02'
[../]
[]
[Kernels]
[./ie]
type = ADTimeDerivative
variable = u
[../]
[./diff]
type = ADDiffusion
variable = u
[../]
[]
[BCs]
[./left]
type = ADFunctionDirichletBC
variable = u
preset = false
boundary = 'left'
function = forcing_fn
[../]
[./right]
type = ADFunctionDirichletBC
variable = u
preset = false
boundary = 'right'
function = forcing_fn
[../]
[]
[Executioner]
type = Transient
# Time integrator scheme
scheme = "newmark-beta"
start_time = 0.0
num_steps = 6
dt = 0.1
[]
[Postprocessors]
[./udotdot]
type = ADElementAverageSecondTimeDerivative
variable = u
[../]
[]
[Outputs]
csv = true
[]
(test/tests/bcs/mat_neumann_bc/ad_mat_neumann.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmax = 10
ymax = 10
nx = 10
ny = 10
[]
[Variables]
[./u]
[../]
[]
[AuxVariables]
[./phi]
[../]
[]
[ICs]
[./phi_IC]
type = FunctionIC
variable = phi
function = ic_func_phi
[../]
[]
[Functions]
[./ic_func_phi]
type = ParsedFunction
value = '0.5 * (1 - tanh((x - 5) / 0.8))'
[../]
[]
[BCs]
[./top]
type = ADMatNeumannBC
variable = u
boundary = top
value = 2
boundary_material = hm
[../]
[]
[Kernels]
[./dudt]
type = ADTimeDerivative
variable = u
[../]
[./diff]
type = ADDiffusion
variable = u
[../]
[]
[Materials]
[./hm]
type = ADParsedMaterial
f_name = hm
args = 'phi'
function = '3*phi^2 - 2*phi^3'
outputs = exodus
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
end_time = 10
[]
[Outputs]
exodus = true
[]
(test/tests/reporters/iteration_info/iteration_info.i)
[Mesh]
[generate]
type = GeneratedMeshGenerator
dim = 1
nx = 10
[]
[]
[Variables/u][]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 10
[]
[]
[Executioner]
type = Transient
num_steps = 3
dt = 1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Reporters/iteration_info]
type = IterationInfo
[]
[Outputs]
[out]
type = JSON
execute_system_information_on = NONE
[]
[]
(modules/phase_field/test/tests/phase_field_kernels/ADAllenCahnVariableL.i)
#
# Test the forward automatic differentiation Allen-Cahn Bulk kernel with a
# spatially varying mobility
#
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
xmax = 12
ymax = 12
elem_type = QUAD4
[]
[AuxVariables]
[./chi]
[./InitialCondition]
type = FunctionIC
function = 'x/24+0.5'
[../]
[../]
[]
[Variables]
[./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 = ADTimeDerivative
variable = eta
[../]
[./ACBulk]
type = ADAllenCahn
variable = eta
f_name = F
[../]
[./ACInterface]
type = ADACInterface
variable = eta
kappa_name = 1
variable_L = true
args = chi
[../]
[]
[Materials]
[./L]
type = ADTestDerivativeFunction
function = F2
f_name = L
op = 'eta chi'
[../]
[./free_energy]
type = ADTestDerivativeFunction
function = F1
f_name = F
op = 'eta'
[../]
[]
[Executioner]
type = Transient
scheme = 'bdf2'
solve_type = 'NEWTON'
num_steps = 2
dt = 1
[]
[Outputs]
exodus = true
[]
(modules/stochastic_tools/examples/sobol/diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables/T]
initial_condition = 300
[]
[Kernels]
[time]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADMatDiffusion
variable = T
diffusivity = diffusivity
[]
[source]
type = ADBodyForce
variable = T
value = 100
function = 1
[]
[]
[BCs]
[left]
type = ADDirichletBC
variable = T
boundary = left
value = 300
[]
[right]
type = ADNeumannBC
variable = T
boundary = right
value = -100
[]
[]
[Materials/constant]
type = ADGenericConstantMaterial
prop_names = 'diffusivity'
prop_values = 1
[]
[Executioner]
type = Transient
solve_type = NEWTON
num_steps = 4
dt = 0.25
[]
[Postprocessors]
[T_avg]
type = ElementAverageValue
variable = T
execute_on = 'initial timestep_end'
[]
[q_left]
type = ADSideDiffusiveFluxAverage
variable = T
boundary = left
diffusivity = diffusivity
execute_on = 'initial timestep_end'
[]
[]
[Controls/stochastic]
type = SamplerReceiver
[]
[Outputs]
[]
(test/tests/kernels/ad_transient_diffusion/ad_transient_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 20
dt = 0.1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(test/tests/time_integrators/central-difference/ad_central_difference_dotdot.i)
###########################################################
# This is a simple test with a time-dependent problem
# demonstrating the use of the TimeIntegrator system.
#
# Testing that the second time derivative is calculated
# correctly using the Central Difference method for an AD
# variable.
#
###########################################################
[Mesh]
type = GeneratedMesh
dim = 2
xmin = -1
xmax = 1
ymin = -1
ymax = 1
nx = 1
ny = 1
[]
[Variables]
[./u]
[../]
[]
[Functions]
[./forcing_fn]
type = PiecewiseLinear
x = '0.0 0.1 0.2 0.3 0.4 0.5 0.6'
y = '0.0 0.0 0.0025 0.01 0.0175 0.02 0.02'
[../]
[]
[Kernels]
[./ie]
type = ADTimeDerivative
variable = u
[../]
[./diff]
type = ADDiffusion
variable = u
[../]
[]
[BCs]
[./left]
type = ADFunctionDirichletBC
variable = u
boundary = 'left'
function = forcing_fn
preset = false
[../]
[./right]
type = ADFunctionDirichletBC
variable = u
boundary = 'right'
function = forcing_fn
preset = false
[../]
[]
[Executioner]
type = Transient
[./TimeIntegrator]
type = CentralDifference
[]
start_time = 0.0
num_steps = 6
dt = 0.1
[]
[Postprocessors]
[./udotdot]
type = ElementAverageSecondTimeDerivative
variable = u
[../]
[]
[Outputs]
csv = true
[]
(modules/phase_field/test/tests/ADCHSplitChemicalPotential/simple_transient_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
[./c]
[../]
[./mu]
[../]
[]
[Kernels]
[./conc]
type = ADCHSplitConcentration
variable = c
chemical_potential_var = mu
mobility = mobility_prop
[../]
[./chempot]
type = ADCHSplitChemicalPotential
variable = mu
chemical_potential = mu_prop
[../]
[./time]
type = ADTimeDerivative
variable = c
[../]
[]
[Materials]
[./chemical_potential]
type = ADPiecewiseLinearInterpolationMaterial
property = mu_prop
variable = c
x = '0 1'
y = '0 1'
[../]
[./mobility_prop]
type = ADGenericConstantMaterial
prop_names = mobility_prop
prop_values = 0.1
[../]
[]
[BCs]
[./leftc]
type = DirichletBC
variable = c
boundary = left
value = 0
[../]
[./rightc]
type = DirichletBC
variable = c
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-pc_type -ksp_grmres_restart -sub_ksp_type -sub_pc_type -pc_asm_overlap'
petsc_options_value = 'asm 31 preonly lu 2'
dt = 0.1
num_steps = 20
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/ad_jacobians/test.i)
[Mesh]
type = GeneratedMesh
dim = 2
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./damage_dt]
type = ADTimeDerivative
variable = u
[../]
[./damage]
type = ADBodyForce
value = 1
variable = u
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
num_steps = 1
[]
(modules/phase_field/test/tests/anisotropic_mobility/ad_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 15
ny = 15
xmax = 15.0
ymax = 15.0
[]
[Variables]
[./c]
[./InitialCondition]
type = CrossIC
x1 = 0.0
x2 = 30.0
y1 = 0.0
y2 = 30.0
[../]
[../]
[]
[Kernels]
[./cres]
type = ADMatAnisoDiffusion
diffusivity = D
variable = c
[../]
[./time]
type = ADTimeDerivative
variable = c
[../]
[]
[Materials]
[./D]
type = ADConstantAnisotropicMobility
tensor = '0.1 0 0
0 1 0
0 0 0'
M_name = D
[../]
[]
[Executioner]
type = Transient
scheme = 'BDF2'
solve_type = 'NEWTON'
l_max_its = 30
l_tol = 1.0e-4
nl_max_its = 50
nl_rel_tol = 1.0e-10
dt = 10.0
num_steps = 2
[]
[Outputs]
exodus = true
print_linear_residuals = true
perf_graph = true
[]
(modules/phase_field/examples/anisotropic_interfaces/ad_snow.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 14
ny = 14
xmax = 9
ymax = 9
uniform_refine = 3
[]
[Variables]
[./w]
[../]
[./T]
[../]
[]
[ICs]
[./wIC]
type = SmoothCircleIC
variable = w
int_width = 0.1
x1 = 4.5
y1 = 4.5
radius = 0.07
outvalue = 0
invalue = 1
[../]
[]
[Kernels]
[./w_dot]
type = ADTimeDerivative
variable = w
[../]
[./anisoACinterface1]
type = ADACInterfaceKobayashi1
variable = w
mob_name = adM
[../]
[./anisoACinterface2]
type = ADACInterfaceKobayashi2
variable = w
mob_name = adM
[../]
[./AllenCahn]
type = AllenCahn
variable = w
mob_name = M
f_name = fbulk
args = T
[../]
[./T_dot]
type = ADTimeDerivative
variable = T
[../]
[./CoefDiffusion]
type = ADDiffusion
variable = T
[../]
[./w_dot_T]
type = ADCoefCoupledTimeDerivative
variable = T
v = w
coef = -1.8
[../]
[]
[Materials]
[./free_energy]
type = DerivativeParsedMaterial
f_name = fbulk
args = 'w T'
constant_names = pi
constant_expressions = 4*atan(1)
function = 'm:=0.9 * atan(10 * (1 - T)) / pi; 1/4*w^4 - (1/2 - m/3) * w^3 + (1/4 - m/2) * w^2'
derivative_order = 2
outputs = exodus
[../]
[./material]
type = ADInterfaceOrientationMaterial
op = w
[../]
[./consts1]
type = ADGenericConstantMaterial
prop_names = 'adM'
prop_values = '3333.333'
[../]
[./consts2]
type = GenericConstantMaterial
prop_names = 'M'
prop_values = '3333.333'
[../]
[]
[Preconditioning]
[./SMP]
type = SMP
full = true
[../]
[]
[Executioner]
type = Transient
scheme = bdf2
solve_type = NEWTON
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu '
nl_abs_tol = 1e-10
nl_rel_tol = 1e-08
l_max_its = 30
end_time = 1
[./TimeStepper]
type = IterationAdaptiveDT
optimal_iterations = 6
iteration_window = 2
dt = 0.0005
growth_factor = 1.1
cutback_factor = 0.75
[../]
[./Adaptivity]
initial_adaptivity = 3 # Number of times mesh is adapted to initial condition
refine_fraction = 0.7 # Fraction of high error that will be refined
coarsen_fraction = 0.1 # Fraction of low error that will coarsened
max_h_level = 5 # Max number of refinements used, starting from initial mesh (before uniform refinement)
weight_names = 'w T'
weight_values = '1 0.5'
[../]
[]
[Outputs]
interval = 5
exodus = true
[]
(test/tests/outputs/xml/xml_iterations.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables/u]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[Functions/function]
type = ParsedFunction
value = 2*x
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 2
[]
[]
[Executioner]
type = Transient
num_steps = 2
solve_type = NEWTON
[]
[VectorPostprocessors]
[line]
type = LineFunctionSampler
functions = function
start_point = '0 0 0'
end_point = '1 0 0'
num_points = 5
sort_by = x
execute_on = 'LINEAR'
[]
[]
[Outputs]
[out]
type = XMLOutput
execute_on = 'LINEAR NONLINEAR'
[]
[]
(test/tests/multiapps/full_solve_multiapp_reset/sub.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[../]
[td]
type = ADTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Executioner]
type = Transient
num_steps = 4
dt = 0.25
solve_type = 'NEWTON'
[]
[Outputs]
exodus = true
[]
(modules/xfem/test/tests/moving_interface/ad_phase_transition_2d.i)
[GlobalParams]
order = FIRST
family = LAGRANGE
[]
[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
nx = 11
ny = 1
xmin = 0.0
xmax = 20.0
ymin = 0.0
ymax = 5.0
elem_type = QUAD4
[]
[]
[XFEM]
qrule = volfrac
output_cut_plane = true
[]
[UserObjects]
[velocity]
type = XFEMPhaseTransitionMovingInterfaceVelocity
diffusivity_at_positive_level_set = 5
diffusivity_at_negative_level_set = 1
equilibrium_concentration_jump = 1
value_at_interface_uo = value_uo
[]
[value_uo]
type = NodeValueAtXFEMInterface
variable = 'u'
interface_mesh_cut_userobject = 'cut_mesh'
execute_on = TIMESTEP_END
level_set_var = ls
[]
[cut_mesh]
type = InterfaceMeshCut2DUserObject
mesh_file = flat_interface_1d.e
interface_velocity_uo = velocity
heal_always = true
[]
[]
[Variables]
[u]
[]
[]
[ICs]
[ic_u]
type = FunctionIC
variable = u
function = 'if(x<5.01, 2, 1)'
[]
[]
[AuxVariables]
[ls]
order = FIRST
family = LAGRANGE
[]
[]
[Constraints]
[u_constraint]
type = XFEMEqualValueAtInterface
geometric_cut_userobject = 'cut_mesh'
use_displaced_mesh = false
variable = u
value = 2
alpha = 1e6
[]
[]
[Kernels]
[diff]
type = ADMatDiffusion
variable = u
diffusivity = diffusion_coefficient
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[AuxKernels]
[ls]
type = MeshCutLevelSetAux
mesh_cut_user_object = cut_mesh
variable = ls
execute_on = 'TIMESTEP_BEGIN'
[]
[]
[Materials]
[diffusivity_A]
type = ADGenericConstantMaterial
prop_names = A_diffusion_coefficient
prop_values = 5
[]
[diffusivity_B]
type = ADGenericConstantMaterial
prop_names = B_diffusion_coefficient
prop_values = 1
[]
[diff_combined]
type = ADLevelSetBiMaterialReal
levelset_positive_base = 'A'
levelset_negative_base = 'B'
level_set_var = ls
prop_name = diffusion_coefficient
[]
[]
[BCs]
# Define boundary conditions
[left_u]
type = ADDirichletBC
variable = u
value = 2
boundary = left
[]
[right_u]
type = ADNeumannBC
variable = u
boundary = right
value = 0
[]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu'
line_search = 'none'
l_tol = 1e-3
nl_max_its = 15
nl_rel_tol = 1e-8
nl_abs_tol = 1e-8
start_time = 0.0
dt = 1
num_steps = 5
max_xfem_update = 1
[]
[Outputs]
file_base = phase_transition_2d_out
execute_on = timestep_end
exodus = true
perf_graph = true
csv = true
[]
(modules/heat_conduction/test/tests/ad_convective_heat_flux/coupled.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
[]
[Variables]
[./temp]
initial_condition = 200.0
[../]
[]
[Kernels]
[./heat_dt]
type = ADTimeDerivative
variable = temp
[../]
[./heat_conduction]
type = Diffusion
variable = temp
[../]
[./heat]
type = ADBodyForce
variable = temp
value = 0
[../]
[]
[BCs]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = T_inf
heat_transfer_coefficient = htc
[../]
[]
[Materials]
[chf_mat]
type = ADConvectiveHeatFluxTest
temperature = temp
boundary = 'right'
[]
[]
[Postprocessors]
[./left_temp]
type = SideAverageValue
variable = temp
boundary = left
execute_on = 'TIMESTEP_END initial'
[../]
[./right_temp]
type = SideAverageValue
variable = temp
boundary = right
[../]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1
nl_abs_tol = 1e-12
[]
[Outputs]
[./out]
type = CSV
interval = 10
[../]
[]
(modules/phase_field/test/tests/anisotropic_interfaces/adkobayashi.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 32
ny = 32
xmax = 0.7
ymax = 0.7
[]
[Variables]
[w]
[]
[T]
[]
[]
[ICs]
[wIC]
type = SmoothCircleIC
variable = w
int_width = 0.1
x1 = 0.35
y1 = 0.35
radius = 0.08
outvalue = 0
invalue = 1
[]
[]
[Kernels]
[w_dot]
type = TimeDerivative
variable = w
[]
[anisoACinterface1]
type = ADACInterfaceKobayashi1
variable = w
mob_name = M
[]
[anisoACinterface2]
type = ADACInterfaceKobayashi2
variable = w
mob_name = M
[]
[AllenCahn]
type = ADAllenCahn
variable = w
mob_name = M
f_name = fbulk
[]
[T_dot]
type = ADTimeDerivative
variable = T
[]
[CoefDiffusion]
type = ADDiffusion
variable = T
[]
[w_dot_T]
type = ADCoefCoupledTimeDerivative
variable = T
v = w
coef = -1.8 #This is -K from kobayashi's paper
[]
[]
[Materials]
[free_energy]
type = ADDerivativeParsedMaterial
f_name = fbulk
args = 'w T'
constant_names = 'alpha gamma T_e pi'
constant_expressions = '0.9 10 1 4*atan(1)'
function = 'm:=alpha/pi * atan(gamma * (T_e - T)); 1/4*w^4 - (1/2 - m/3) * w^3 + (1/4 - m/2) * '
'w^2'
derivative_order = 1
outputs = exodus
[]
[material]
type = ADInterfaceOrientationMaterial
op = w
[]
[consts]
type = ADGenericConstantMaterial
prop_names = 'M'
prop_values = '3333.333'
[]
[]
[Executioner]
type = Transient
solve_type = NEWTON
scheme = bdf2
petsc_options_iname = '-pc_type'
petsc_options_value = 'lu '
nl_rel_tol = 1e-08
l_tol = 1e-4
l_max_its = 30
dt = 0.001
num_steps = 6
[]
[Outputs]
exodus = true
perf_graph = true
execute_on = 'INITIAL FINAL'
[]
(modules/level_set/test/tests/functions/olsson_bubble/olsson_bubble_adjac.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
displacements = 'disp_x disp_y'
[]
[Problem]
kernel_coverage_check = false
[]
[Variables]
[bubble]
[]
[disp_x]
[]
[disp_y]
[]
[]
[Kernels]
[bubble]
type = ADBodyForce
variable = bubble
function = bubble_func
use_displaced_mesh = true
[]
[dt]
type = ADTimeDerivative
variable = bubble
[]
[]
[Functions]
[bubble_func]
type = LevelSetOlssonBubble
center = '0.5 0.5 0'
radius = 0.4
epsilon = 0.05
[]
[]
[Executioner]
solve_type = NEWTON
type = Transient
num_steps = 1
[]
[Outputs]
exodus = true
[]
(test/tests/misc/ad_robustness/ad_two_var_transient_diffusion.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
[]
[AuxVariables]
[v][]
[]
[Variables]
[./u]
[../]
[]
[Kernels]
[./diff]
type = CoefDiffusion
variable = u
coef = 0.1
[../]
[./time]
type = ADTimeDerivative
variable = u
[../]
[coupled]
type = ADCoupledValueTest
variable = u
v = v
[]
[]
[DGKernels]
[dummy]
type = ADDGCoupledTest
variable = u
v = v
[]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 2
dt = 0.1
solve_type = NEWTON
petsc_options_iname = '-pc_type -pc_hypre_type'
petsc_options_value = 'hypre boomeramg'
[]
[Outputs]
exodus = true
[]
(modules/thermal_hydraulics/test/tests/components/hs_boundary_external_app_temperature/phy.master.i)
# This tests a transfer of temperature values computed by master app and used by a slave app
# as a heat structure boundary condition
[Mesh]
type = GeneratedMesh
dim = 1
xmax = 1
nx = 10
[]
[Functions]
[T_bc_fn]
type = ParsedFunction
value = '300+t*x*10'
[]
[T_ffn]
type = ParsedFunction
value = 'x*10'
[]
[]
[Variables]
[T]
[]
[]
[ICs]
[T_ic]
type = ConstantIC
variable = T
value = 300
[]
[]
[Kernels]
[td]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADDiffusion
variable = T
[]
[ffn]
type = BodyForce
variable = T
function = T_ffn
[]
[]
[BCs]
[left]
type = FunctionDirichletBC
variable = T
boundary = 'left right'
function = T_bc_fn
[]
[]
[Executioner]
type = Transient
dt = 0.1
num_steps = 2
nl_abs_tol = 1e-10
abort_on_solve_fail = true
solve_type = NEWTON
[]
[MultiApps]
[thm]
type = TransientMultiApp
app_type = ThermalHydraulicsApp
input_files = phy.slave.i
execute_on = 'initial timestep_end'
[]
[]
[Transfers]
[T_to_thm]
type = MultiAppNearestNodeTransfer
to_multi_app = thm
source_variable = T
variable = T_ext
target_boundary = 'hs:outer'
[]
[]
[Outputs]
exodus = true
show = 'T'
[]
(modules/phase_field/test/tests/automatic_differentiation/admatreaction.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
[]
[Variables]
[./u]
[./InitialCondition]
type = FunctionIC
function = cos(x*2*pi)
[../]
[../]
[./v]
[./InitialCondition]
type = FunctionIC
function = sin(x*2*pi)
[../]
[../]
[]
[Kernels]
[./dudt]
type = ADTimeDerivative
variable = u
[../]
[./dvdt]
type = ADTimeDerivative
variable = v
[../]
[./u]
type = ADMatReaction
variable = u
#v = v
[../]
[./v]
type = ADMatReaction
variable = v
[../]
[]
[Materials]
[./L]
type = ADTestDerivativeFunction
function = F3
f_name = L
op = 'u v'
[../]
[]
[Preconditioning]
[./smp]
type = SMP
full = true
[../]
[]
[Debug]
show_var_residual_norms = true
[]
[Executioner]
type = Transient
solve_type = NEWTON
dt = 0.1
num_steps = 5
[]
[Outputs]
exodus = true
[]
(modules/thermal_hydraulics/test/tests/components/hs_boundary_external_app_convection/plate.master.i)
# This tests a temperature and heat transfer coefficient using the MultiApp system.
# Simple heat conduction problem with heat source is solved,
# to obtain an analytical solution:
# T(x,t) = 300 + 20t(x-1)^2
# The temperature is picked up by the slave
# side of the solve, to use as ambiant temperature in a convective BC.
htc = 100
[Mesh]
type = GeneratedMesh
dim = 1
xmax = 1
nx = 10
[]
[Functions]
[left_bc_fn]
type = PiecewiseLinear
x = '0 10'
y = '300 500'
[]
[]
[Variables]
[T]
[]
[]
[AuxVariables]
[htc_ext]
initial_condition = ${htc}
[]
[]
[Functions]
[source_term]
type = ParsedFunction
value = '20 * x * x - 40 * x + 20 - 40 * t'
[]
[]
[ICs]
[T_ic]
type = ConstantIC
variable = T
value = 300
[]
[]
[Kernels]
[td]
type = ADTimeDerivative
variable = T
[]
[diff]
type = ADDiffusion
variable = T
[]
[source]
type = BodyForce
function = 'source_term'
variable = T
[]
[]
[BCs]
[left]
type = FunctionDirichletBC
variable = T
boundary = left
function = left_bc_fn
[]
[right]
type = NeumannBC
variable = T
boundary = right
value = 0
[]
[]
[Executioner]
type = Transient
dt = 0.5
end_time = 10
nl_abs_tol = 1e-10
abort_on_solve_fail = true
[]
[MultiApps]
[thm]
type = TransientMultiApp
app_type = ThermalHydraulicsApp
input_files = plate.i
execute_on = TIMESTEP_END
[]
[]
[Transfers]
[T_to_slave]
type = MultiAppNearestNodeTransfer
to_multi_app = thm
source_variable = T
variable = T_ext
[]
[htc_to_slave]
type = MultiAppNearestNodeTransfer
to_multi_app = thm
source_variable = htc_ext
variable = htc_ext
[]
[]
[Outputs]
exodus = true
[]
(modules/xfem/test/tests/moving_interface/moving_ad_diffusion.i)
[GlobalParams]
order = FIRST
family = LAGRANGE
[]
[XFEM]
qrule = volfrac
output_cut_plane = true
[]
[UserObjects]
[./level_set_cut_uo]
type = LevelSetCutUserObject
level_set_var = ls
heal_always = true
[../]
[]
[Mesh]
type = GeneratedMesh
dim = 2
nx = 5
ny = 3
xmin = 0.0
xmax = 1
ymin = 0.0
ymax = 1
elem_type = QUAD4
[]
[AuxVariables]
[./ls]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxKernels]
[./ls_function]
type = FunctionAux
variable = ls
function = ls_func
[../]
[]
[Variables]
[./u]
[../]
[]
[Functions]
[./ls_func]
type = ParsedFunction
value = 'x-0.76+0.21*t'
[../]
[]
[Kernels]
[./diff]
type = ADMatDiffusion
variable = u
diffusivity = diffusion_coefficient
[../]
[./time_deriv]
type = ADTimeDerivative
variable = u
[../]
[]
[Constraints]
[./u_constraint]
type = XFEMSingleVariableConstraint
use_displaced_mesh = false
variable = u
jump = 0
use_penalty = true
alpha = 1e5
geometric_cut_userobject = 'level_set_cut_uo'
[../]
[]
[BCs]
[./right_u]
type = ADDirichletBC
variable = u
boundary = left
value = 0
[../]
[./left_u]
type = ADDirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[Materials]
[./diffusivity_A]
type = ADGenericConstantMaterial
prop_names = A_diffusion_coefficient
prop_values = 5
[../]
[./diffusivity_B]
type = ADGenericConstantMaterial
prop_names = B_diffusion_coefficient
prop_values = 0.5
[../]
[./diff_combined]
type = ADLevelSetBiMaterialReal
levelset_positive_base = 'A'
levelset_negative_base = 'B'
level_set_var = ls
prop_name = diffusion_coefficient
[../]
[]
[Executioner]
type = Transient
solve_type = NEWTON
petsc_options_iname = '-ksp_gmres_restart -pc_type -pc_hypre_type -pc_hypre_boomeramg_max_iter'
petsc_options_value = '201 hypre boomeramg 8'
l_max_its = 20
l_tol = 1e-8
nl_max_its = 15
nl_rel_tol = 2e-12
nl_abs_tol = 1e-50
start_time = 0.0
dt = 1
end_time = 2
max_xfem_update = 1
[]
[Outputs]
exodus = true
execute_on = timestep_end
csv = true
file_base = moving_diffusion_out
perf_graph = true
[./console]
type = Console
output_linear = true
[../]
[]
(modules/stochastic_tools/examples/batch/sub.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
ny = 10
nz = 10
[]
[Variables]
[u]
[]
[]
[Kernels]
[diff]
type = ADDiffusion
variable = u
[]
[time]
type = ADTimeDerivative
variable = u
[]
[]
[BCs]
[left]
type = DirichletBC
variable = u
boundary = left
value = 0
[]
[right]
type = DirichletBC
variable = u
boundary = right
value = 1
[]
[]
[Postprocessors]
[average]
type = AverageNodalVariableValue
variable = u
[]
[]
[Executioner]
type = Transient
num_steps = 1
dt = 0.25
solve_type = NEWTON
[]
[Controls]
[receiver]
type = SamplerReceiver
[]
[]
[Outputs]
[]
(modules/heat_conduction/test/tests/ad_convective_heat_flux/equilibrium.i)
[Mesh]
type = GeneratedMesh
dim = 3
nx = 10
[]
[Variables]
[./temp]
initial_condition = 200.0
[../]
[]
[Kernels]
[./heat_dt]
type = ADTimeDerivative
variable = temp
[../]
[./heat_conduction]
type = ADDiffusion
variable = temp
[../]
[]
[BCs]
[./right]
type = ADConvectiveHeatFluxBC
variable = temp
boundary = 'right'
T_infinity = 100.0
heat_transfer_coefficient = 1
[../]
[]
[Postprocessors]
[./left_temp]
type = SideAverageValue
variable = temp
boundary = left
execute_on = 'TIMESTEP_END initial'
[../]
[./right_temp]
type = SideAverageValue
variable = temp
boundary = right
[../]
[./right_flux]
type = SideDiffusiveFluxAverage
variable = temp
boundary = right
diffusivity = 1
[../]
[]
[Executioner]
type = Transient
num_steps = 10
dt = 1e1
nl_abs_tol = 1e-12
[]
[Outputs]
[./out]
type = CSV
interval = 10
[../]
[]
(modules/navier_stokes/include/kernels/INSADHeatConductionTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "ADTimeDerivative.h"
class INSADHeatConductionTimeDerivative : public ADTimeDerivative
{
public:
static InputParameters validParams();
INSADHeatConductionTimeDerivative(const InputParameters & parameters);
protected:
ADReal precomputeQpResidual() override;
const ADMaterialProperty<Real> & _temperature_td_strong_residual;
};
(modules/heat_conduction/include/kernels/ADHeatConductionTimeDerivative.h)
// This file is part of the MOOSE framework
// https://www.mooseframework.org
//
// All rights reserved, see COPYRIGHT for full restrictions
// https://github.com/idaholab/moose/blob/master/COPYRIGHT
//
// Licensed under LGPL 2.1, please see LICENSE for details
// https://www.gnu.org/licenses/lgpl-2.1.html
#pragma once
#include "ADTimeDerivative.h"
class ADHeatConductionTimeDerivative : public ADTimeDerivative
{
public:
static InputParameters validParams();
ADHeatConductionTimeDerivative(const InputParameters & parameters);
protected:
virtual ADReal precomputeQpResidual() override;
/// Specific heat material property
const ADMaterialProperty<Real> & _specific_heat;
/// Density material property
const ADMaterialProperty<Real> & _density;
};