- variableThe name of the variable that this residual object operates on
C++ Type:NonlinearVariableName
Unit:(no unit assumed)
Controllable:No
Description:The name of the variable that this residual object operates on
ODETimeDerivative
Returns the time derivative contribution to the residual for a scalar variable.
A scalar variable can be set to the solution of an ordinary differential equation (ODE), as specified in the Scalar Kernels syntax page. This kernel adds a time derivative term. The time integration scheme will be shared with the other non-linear variables. To use a different time integrating scheme, the ODETimeDerivative
scalar kernel should be replaced with a custom implementation.
Example input syntax
In this example, the scalar variables x
and y
are the solutions to the coupled ODE problem:
The time derivative terms are added for each variable using two ODETimeDerivative
scalar kernels.
[ScalarKernels<<<{"href": "../../syntax/ScalarKernels/index.html"}>>>]
[./td1]
type = ODETimeDerivative<<<{"description": "Returns the time derivative contribution to the residual for a scalar variable.", "href": "ODETimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = x
[../]
[./ode1]
type = ParsedODEKernel<<<{"description": "Parsed expression ODE kernel.", "href": "ParsedODEKernel.html"}>>>
expression<<<{"description": "function expression"}>>> = '-3*x - 2*y'
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = x
coupled_variables<<<{"description": "Scalar variables coupled in the parsed expression."}>>> = y
[../]
[./td2]
type = ODETimeDerivative<<<{"description": "Returns the time derivative contribution to the residual for a scalar variable.", "href": "ODETimeDerivative.html"}>>>
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = y
[../]
[./ode2]
type = ParsedODEKernel<<<{"description": "Parsed expression ODE kernel.", "href": "ParsedODEKernel.html"}>>>
expression<<<{"description": "function expression"}>>> = '-4*x - y'
variable<<<{"description": "The name of the variable that this residual object operates on"}>>> = y
coupled_variables<<<{"description": "Scalar variables coupled in the parsed expression."}>>> = x
[../]
[]
(test/tests/kernels/ode/parsedode_sys_impl_test.i)Input Parameters
- matrix_onlyFalseWhether this object is only doing assembly to matrices (no vectors)
Default:False
C++ Type:bool
Controllable:No
Description:Whether this object is only doing assembly to matrices (no vectors)
Optional Parameters
- absolute_value_vector_tagsThe tags for the vectors this residual object should fill with the absolute value of the residual contribution
C++ Type:std::vector<TagName>
Controllable:No
Description:The tags for the vectors this residual object should fill with the absolute value of the residual contribution
- 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
Contribution To Tagged Field Data 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.
- 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
- 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
Input Files
- (test/tests/controls/time_periods/scalarkernels/scalarkernels.i)
- (test/tests/kernels/bad_scaling_scalar_kernels/ill_conditioned_field_scalar_system.i)
- (test/tests/time_integrators/scalar/stiff.i)
- (test/tests/time_integrators/explicit_ssp_runge_kutta/explicit_ssp_runge_kutta.i)
- (test/tests/outputs/variables/output_vars_test.i)
- (examples/ex18_scalar_kernel/ex18_parsed.i)
- (test/tests/controls/conditional_functional_enable/conditional_function_enable.i)
- (test/tests/adaptivity/scalar/scalar_adaptivity.i)
- (test/tests/kernels/ode/ode_sys_impl_test.i)
- (test/tests/outputs/variables/output_vars_nonexistent.i)
- (test/tests/misc/check_error/scalar_kernel_with_var.i)
- (test/tests/postprocessors/scalar_variable/scalar_variable_pps.i)
- (test/tests/time_integrators/scalar/scalar.i)
- (test/tests/outputs/variables/output_vars_hidden_shown_check.i)
- (test/tests/scaling/ignore-variables/ignore.i)
- (examples/ex18_scalar_kernel/ex18.i)
- (modules/thermal_hydraulics/test/tests/scalarkernels/postprocessor_source/postprocessor_source.i)
- (test/tests/bcs/periodic/no_add_scalar.i)
- (test/tests/time_integrators/actually_explicit_euler_verification/ee-ode.i)
- (test/tests/ics/function_scalar_ic/function_scalar_ic.i)
- (test/tests/scalar_kernels/ad_coupled_scalar/ad_coupled_scalar.i)
- (test/tests/postprocessors/scalar_coupled_postprocessor/scalar_coupled_postprocessor_test.i)
- (test/tests/tag/scalar_tag_vector.i)
- (test/tests/kernels/ode/parsedode_sys_impl_test.i)
- (test/tests/kernels/ode/parsedode_pp_test.i)
Child Objects
(test/tests/kernels/ode/parsedode_sys_impl_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
[]
[Functions]
[./f_fn]
type = ParsedFunction
expression = -4
[../]
[./bc_all_fn]
type = ParsedFunction
expression = x*x+y*y
[../]
# ODEs
[./exact_x_fn]
type = ParsedFunction
expression = (-1/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[]
# NL
[Variables]
[./u]
family = LAGRANGE
order = FIRST
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff]
type = Diffusion
variable = u
[../]
[./uff]
type = BodyForce
variable = u
function = f_fn
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ParsedODEKernel
expression = '-3*x - 2*y'
variable = x
coupled_variables = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ParsedODEKernel
expression = '-4*x - y'
variable = y
coupled_variables = x
[../]
[]
[BCs]
[./all]
type = FunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = bc_all_fn
[../]
[]
[Postprocessors]
active = 'exact_x l2err_x'
[./exact_x]
type = FunctionValuePostprocessor
function = exact_x_fn
execute_on = 'initial timestep_end'
point = '0 0 0'
[../]
[./l2err_x]
type = ScalarL2Error
variable = x
function = exact_x_fn
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Transient
start_time = 0
dt = 0.01
num_steps = 100
solve_type = 'PJFNK'
[]
[Outputs]
file_base = ode_sys_impl_test_out
exodus = true
[]
(test/tests/controls/time_periods/scalarkernels/scalarkernels.i)
# This tests controllability of the enable parameter of scalar kernels.
#
# There are 2 scalar variables, {u, v}, with the ODEs:
# du/dt = 1 u(0) = 0
# v = u v(0) = -10
# A control switches the ODE 'v = u' to the following ODE when t >= 2:
# dv/dt = 2
#
# 5 time steps (of size dt = 1) will be taken, and the predicted values are as follows:
# t u v
# ------------------
# 0 0 -10
# 1 1 1
# 2 2 2
# 3 3 4
# 4 4 6
# 5 5 8
u_initial = 0
u_growth = 1
v_initial = -10
v_growth = 2
t_transition = 2
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
[./u]
family = SCALAR
order = FIRST
[../]
[./v]
family = SCALAR
order = FIRST
[../]
[]
[ICs]
[./u_ic]
type = ScalarConstantIC
variable = u
value = ${u_initial}
[../]
[./v_ic]
type = ScalarConstantIC
variable = v
value = ${v_initial}
[../]
[]
[ScalarKernels]
[./u_time]
type = ODETimeDerivative
variable = u
[../]
[./u_src]
type = ParsedODEKernel
variable = u
expression = '-${u_growth}'
[../]
[./v_time]
type = ODETimeDerivative
variable = v
enable = false
[../]
[./v_src]
type = ParsedODEKernel
variable = v
expression = '-${v_growth}'
enable = false
[../]
[./v_constraint]
type = ParsedODEKernel
variable = v
coupled_variables = 'u'
expression = 'v - u'
[../]
[]
[Controls]
[./time_period_control]
type = TimePeriod
end_time = ${t_transition}
enable_objects = 'ScalarKernel::v_constraint'
disable_objects = 'ScalarKernel::v_time ScalarKernel::v_src'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[]
[Executioner]
type = Transient
scheme = implicit-euler
dt = 1
num_steps = 5
abort_on_solve_fail = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-8
[]
[Outputs]
csv = true
[]
(test/tests/kernels/bad_scaling_scalar_kernels/ill_conditioned_field_scalar_system.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 2
[]
[Variables]
[./u]
[../]
[v]
family = SCALAR
initial_condition = 1
[]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[scalar]
type = ScalarLagrangeMultiplier
variable = u
lambda = v
[]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = left
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = right
value = 1
[../]
[]
[ScalarKernels]
[reaction]
type = ParsedODEKernel
expression = '10^20 * v'
variable = v
[]
[time]
type = ODETimeDerivative
variable = v
[]
[]
[Executioner]
type = Transient
num_steps = 1
dtmin = 1
solve_type = NEWTON
petsc_options = '-pc_svd_monitor -ksp_view_pmat -snes_converged_reason -ksp_converged_reason'
petsc_options_iname = '-pc_type -snes_stol'
petsc_options_value = 'svd 0'
[]
[Outputs]
exodus = true
[]
(test/tests/time_integrators/scalar/stiff.i)
# This is a linear model problem described in Frank et al, "Order
# results for implicit Runge-Kutta methods applied to stiff systems",
# SIAM J. Numerical Analysis, vol. 22, no. 3, 1985, pp. 515-534.
#
# Problems "PL" and "PNL" from page 527 of the paper:
# { dy1/dt = lambda*y1 + y2**p, y1(0) = -1/(lambda+p)
# { dy2/dt = -y2, y2(0) = 1
#
# The exact solution is:
# y1 = -exp(-p*t)/(lambda+p)
# y2 = exp(-t)
#
# According to the following paragraph from the reference above, the
# p=1 version of this problem should not exhibit order reductions
# regardless of stiffness, while the nonlinear version (p>=2) will
# exhibit order reductions down to the "stage order" of the method for
# lambda large, negative.
# Use Dollar Bracket Expressions (DBEs) to set the value of LAMBDA in
# a single place. You can also set this on the command line with
# e.g. LAMBDA=-4, but note that this does not seem to override the
# value set in the input file. This is a bit different from the way
# that command line values normally work...
# Note that LAMBDA == Y2_EXPONENT is not allowed!
# LAMBDA = -10
# Y2_EXPONENT = 2
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 1
ny = 1
elem_type = QUAD4
[]
[Variables]
[./y1]
family = SCALAR
order = FIRST
[../]
[./y2]
family = SCALAR
order = FIRST
[../]
[]
[ICs]
[./y1_init]
type = FunctionScalarIC
variable = y1
function = y1_exact
[../]
[./y2_init]
type = FunctionScalarIC
variable = y2
function = y2_exact
[../]
[]
[ScalarKernels]
[./y1_time]
type = ODETimeDerivative
variable = y1
[../]
[./y1_space]
type = ParsedODEKernel
variable = y1
expression = '-(${LAMBDA})*y1 - y2^${Y2_EXPONENT}'
coupled_variables = 'y2'
[../]
[./y2_time]
type = ODETimeDerivative
variable = y2
[../]
[./y2_space]
type = ParsedODEKernel
variable = y2
expression = 'y2'
[../]
[]
[Executioner]
type = Transient
[./TimeIntegrator]
type = LStableDirk2
[../]
start_time = 0
end_time = 1
dt = 0.125
solve_type = 'PJFNK'
nl_max_its = 6
nl_abs_tol = 1.e-13
nl_rel_tol = 1.e-32 # Force nl_abs_tol to be used.
line_search = 'none'
[]
[Functions]
[./y1_exact]
type = ParsedFunction
expression = '-exp(-${Y2_EXPONENT}*t)/(lambda+${Y2_EXPONENT})'
symbol_names = 'lambda'
symbol_values = ${LAMBDA}
[../]
[./y2_exact]
type = ParsedFunction
expression = exp(-t)
[../]
[]
[Postprocessors]
[./error_y1]
type = ScalarL2Error
variable = y1
function = y1_exact
execute_on = 'initial timestep_end'
[../]
[./error_y2]
type = ScalarL2Error
variable = y2
function = y2_exact
execute_on = 'initial timestep_end'
[../]
[./max_error_y1]
# Estimate ||e_1||_{\infty}
type = TimeExtremeValue
value_type = max
postprocessor = error_y1
execute_on = 'initial timestep_end'
[../]
[./max_error_y2]
# Estimate ||e_2||_{\infty}
type = TimeExtremeValue
value_type = max
postprocessor = error_y2
execute_on = 'initial timestep_end'
[../]
[./value_y1]
type = ScalarVariable
variable = y1
execute_on = 'initial timestep_end'
[../]
[./value_y2]
type = ScalarVariable
variable = y2
execute_on = 'initial timestep_end'
[../]
[./value_y1_abs_max]
type = TimeExtremeValue
value_type = abs_max
postprocessor = value_y1
execute_on = 'initial timestep_end'
[../]
[./value_y2_abs_max]
type = TimeExtremeValue
value_type = abs_max
postprocessor = value_y2
execute_on = 'initial timestep_end'
[../]
[]
[Outputs]
csv = true
[]
(test/tests/time_integrators/explicit_ssp_runge_kutta/explicit_ssp_runge_kutta.i)
# This test solves the following IVP:
# du/dt = f(u(t), t), u(0) = 1
# f(u(t), t) = -u(t) + t^3 + 3t^2
# The exact solution is the following:
# u(t) = exp(-t) + t^3
[Mesh]
[./mesh]
type = GeneratedMeshGenerator
dim = 1
nx = 1
[../]
[]
[Variables]
[./u]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[]
[ScalarKernels]
[./time_derivative]
type = ODETimeDerivative
variable = u
[../]
[./source_part1]
type = ParsedODEKernel
variable = u
expression = 'u'
[../]
[./source_part2]
type = PostprocessorSinkScalarKernel
variable = u
postprocessor = sink_pp
[../]
[]
[Functions]
[./sink_fn]
type = ParsedFunction
expression = '-t^3 - 3*t^2'
[../]
[]
[Postprocessors]
[./sink_pp]
type = FunctionValuePostprocessor
function = sink_fn
execute_on = 'LINEAR NONLINEAR'
[../]
[./l2_err]
type = ScalarL2Error
variable = u
function = ${fparse exp(-0.5) + 0.5^3}
[../]
[]
[Executioner]
type = Transient
[./TimeIntegrator]
type = ExplicitSSPRungeKutta
order = 1
[../]
end_time = 0.5
dt = 0.1
[]
[Outputs]
file_base = 'first_order'
[./csv]
type = CSV
show = 'u'
execute_on = 'FINAL'
[../]
[]
(test/tests/outputs/variables/output_vars_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
elem_type = QUAD9
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = SECOND
family = LAGRANGE
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[AuxVariables]
[./elemental]
order = CONSTANT
family = MONOMIAL
[../]
[./elemental_restricted]
order = CONSTANT
family = MONOMIAL
[../]
[./nodal]
order = FIRST
family = LAGRANGE
[../]
[./nodal_restricted]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff_u]
type = Diffusion
variable = u
[../]
[./conv_u]
type = CoupledForce
variable = u
v = v
[../]
[./diff_v]
type = Diffusion
variable = v
[../]
[]
[AuxKernels]
[./elemental]
type = ConstantAux
variable = elemental
value = 1
[../]
[./elemental_restricted]
type = ConstantAux
variable = elemental_restricted
value = 1
[../]
[./nodal]
type = ConstantAux
variable = elemental
value = 2
[../]
[./nodal_restricted]
type = ConstantAux
variable = elemental_restricted
value = 2
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ImplicitODEx
variable = x
y = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ImplicitODEy
variable = y
x = x
[../]
[]
[BCs]
active = 'left_u right_u left_v'
[./left_u]
type = DirichletBC
variable = u
boundary = 1
value = 1
[../]
[./right_u]
type = DirichletBC
variable = u
boundary = 3
value = 9
[../]
[./left_v]
type = DirichletBC
variable = v
boundary = 1
value = 5
[../]
[./right_v]
type = DirichletBC
variable = v
boundary = 2
value = 2
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
dt = 0.01
num_steps = 1
[]
[Outputs]
show = 'x u nodal elemental'
[./out]
type = Exodus
elemental_as_nodal = true
scalar_as_nodal = true
[../]
[]
(examples/ex18_scalar_kernel/ex18_parsed.i)
#
# Example 18 modified to use parsed ODE kernels.
#
# The ParsedODEKernel takes expression expressions in the input file and computes
# Jacobian entries via automatic differentiation. It allows for rapid development
# of new models without the need for code recompilation.
#
# This input file should produce the exact same result as ex18.i
#
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
elem_type = QUAD4
[]
[Functions]
# ODEs
[./exact_x_fn]
type = ParsedFunction
expression = (-1/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[./exact_y_fn]
type = ParsedFunction
expression = (2/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[]
[Variables]
[./diffused]
order = FIRST
family = LAGRANGE
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = diffused
[../]
[./diff]
type = Diffusion
variable = diffused
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
#
# This parsed expression ODE Kernel behaves exactly as the ImplicitODEx kernel
# in the main example. Checkout ImplicitODEx::computeQpResidual() in the
# source code file ImplicitODEx.C to see the matching residual function.
#
# The ParsedODEKernel automaticaly generates the On- and Off-Diagonal Jacobian
# entries.
#
[./ode1]
type = ParsedODEKernel
expression = '-3*x - 2*y'
variable = x
coupled_variables = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
#
# This parsed expression ODE Kernel behaves exactly as the ImplicitODEy Kernel
# in the main example.
#
[./ode2]
type = ParsedODEKernel
expression = '-4*x - y'
variable = y
coupled_variables = x
[../]
[]
[BCs]
[./right]
type = ScalarDirichletBC
variable = diffused
boundary = 1
scalar_var = x
[../]
[./left]
type = ScalarDirichletBC
variable = diffused
boundary = 3
scalar_var = y
[../]
[]
[Postprocessors]
# to print the values of x, y into a file so we can plot it
[./x_pp]
type = ScalarVariable
variable = x
execute_on = timestep_end
[../]
[./y_pp]
type = ScalarVariable
variable = y
execute_on = timestep_end
[../]
[./exact_x]
type = FunctionValuePostprocessor
function = exact_x_fn
execute_on = timestep_end
[../]
[./exact_y]
type = FunctionValuePostprocessor
function = exact_y_fn
execute_on = timestep_end
point = '0 0 0'
[../]
# Measure the error in ODE solution for 'x'.
[./l2err_x]
type = ScalarL2Error
variable = x
function = exact_x_fn
[../]
# Measure the error in ODE solution for 'y'.
[./l2err_y]
type = ScalarL2Error
variable = y
function = exact_y_fn
[../]
[]
[Executioner]
type = Transient
start_time = 0
dt = 0.01
num_steps = 10
solve_type = 'PJFNK'
[]
[Outputs]
file_base = 'ex18_out'
exodus = true
[]
(test/tests/controls/conditional_functional_enable/conditional_function_enable.i)
# This tests controllability of the enable parameter of a MOOSE object via a
# conditional function.
#
# There are 2 scalar variables, {u, v}, with the ODEs:
# du/dt = 1 u(0) = 0
# v = u v(0) = -10
# A control switches the ODE 'v = u' to the following ODE when u >= 1.99:
# dv/dt = 2
#
# 5 time steps (of size dt = 1) will be taken, and the predicted values are as follows:
# t u v
# ------------------
# 0 0 -10
# 1 1 1
# 2 2 2
# 3 3 4
# 4 4 6
# 5 5 8
u_initial = 0
u_growth = 1
u_threshold = 1.99
v_initial = -10
v_growth = 2
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Variables]
[./u]
family = SCALAR
order = FIRST
[../]
[./v]
family = SCALAR
order = FIRST
[../]
[]
[ICs]
[./u_ic]
type = ScalarConstantIC
variable = u
value = ${u_initial}
[../]
[./v_ic]
type = ScalarConstantIC
variable = v
value = ${v_initial}
[../]
[]
[ScalarKernels]
[./u_time]
type = ODETimeDerivative
variable = u
[../]
[./u_src]
type = ParsedODEKernel
variable = u
expression = '-${u_growth}'
[../]
[./v_time]
type = ODETimeDerivative
variable = v
enable = false
[../]
[./v_src]
type = ParsedODEKernel
variable = v
expression = '-${v_growth}'
enable = false
[../]
[./v_constraint]
type = ParsedODEKernel
variable = v
coupled_variables = 'u'
expression = 'v - u'
[../]
[]
[Functions]
[./conditional_function]
type = ParsedFunction
symbol_names = 'u_sol'
symbol_values = 'u'
expression = 'u_sol >= ${u_threshold}'
[../]
[]
[Controls]
[./u_threshold]
type = ConditionalFunctionEnableControl
conditional_function = conditional_function
enable_objects = 'ScalarKernel::v_time ScalarKernel::v_src'
disable_objects = 'ScalarKernel::v_constraint'
execute_on = 'INITIAL TIMESTEP_END'
[../]
[]
[Executioner]
type = Transient
scheme = implicit-euler
dt = 1
num_steps = 5
abort_on_solve_fail = true
nl_rel_tol = 1e-8
nl_abs_tol = 1e-8
[]
[Outputs]
csv = true
[]
(test/tests/adaptivity/scalar/scalar_adaptivity.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 10
ny = 10
xmin = -0.5
xmax = 0.5
ymin = -0.5
ymax = 0.5
[]
[Variables]
[scalar]
order = THIRD
family = SCALAR
[]
[u]
[InitialCondition]
type = FunctionIC
function = 'x*x+y*y'
[]
[]
[]
[Kernels]
[u_dot]
type = TimeDerivative
variable = u
[]
[c_res]
type = Diffusion
variable = u
[]
[]
[ScalarKernels]
[d1]
type = ODETimeDerivative
variable = scalar
[]
[]
[BCs]
[Periodic]
[all]
auto_direction = 'x y'
variable = 'u'
[]
[]
[]
[Preconditioning]
[SMP]
type = SMP
full = true
[]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
petsc_options_iname = '-pc_type -sub_pc_type'
petsc_options_value = 'asm lu '
num_steps = 2
[]
[Adaptivity]
initial_steps = 2
max_h_level = 2
marker = EFM
[Markers]
[EFM]
type = ErrorFractionMarker
coarsen = 0.2
refine = 0.8
indicator = GJI
[]
[]
[Indicators]
[GJI]
type = GradientJumpIndicator
variable = u
[]
[]
[]
[Outputs]
exodus = true
[]
(test/tests/kernels/ode/ode_sys_impl_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
[]
[Functions]
[./f_fn]
type = ParsedFunction
expression = -4
[../]
[./bc_all_fn]
type = ParsedFunction
expression = x*x+y*y
[../]
# ODEs
[./exact_x_fn]
type = ParsedFunction
expression = (-1/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[]
# NL
[Variables]
[./u]
family = LAGRANGE
order = FIRST
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff]
type = Diffusion
variable = u
[../]
[./uff]
type = BodyForce
variable = u
function = f_fn
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ImplicitODEx
variable = x
y = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ImplicitODEy
variable = y
x = x
[../]
[]
[BCs]
[./all]
type = FunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = bc_all_fn
[../]
[]
[Postprocessors]
active = 'exact_x l2err_x'
[./exact_x]
type = FunctionValuePostprocessor
function = exact_x_fn
execute_on = 'initial timestep_end'
point = '0 0 0'
[../]
[./l2err_x]
type = ScalarL2Error
variable = x
function = exact_x_fn
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Transient
start_time = 0
dt = 0.01
num_steps = 100
solve_type = 'PJFNK'
[]
[Outputs]
exodus = true
[]
(test/tests/outputs/variables/output_vars_nonexistent.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
elem_type = QUAD9
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = SECOND
family = LAGRANGE
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[AuxVariables]
[./elemental]
order = CONSTANT
family = MONOMIAL
[../]
[./elemental_restricted]
order = CONSTANT
family = MONOMIAL
[../]
[./nodal]
order = FIRST
family = LAGRANGE
[../]
[./nodal_restricted]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff_u]
type = Diffusion
variable = u
[../]
[./conv_u]
type = CoupledForce
variable = u
v = v
[../]
[./diff_v]
type = Diffusion
variable = v
[../]
[]
[AuxKernels]
[./elemental]
type = ConstantAux
variable = elemental
value = 1
[../]
[./elemental_restricted]
type = ConstantAux
variable = elemental_restricted
value = 1
[../]
[./nodal]
type = ConstantAux
variable = elemental
value = 2
[../]
[./nodal_restricted]
type = ConstantAux
variable = elemental_restricted
value = 2
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ImplicitODEx
variable = x
y = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ImplicitODEy
variable = y
x = x
[../]
[]
[BCs]
active = 'left_u right_u left_v'
[./left_u]
type = DirichletBC
variable = u
boundary = 1
value = 1
[../]
[./right_u]
type = DirichletBC
variable = u
boundary = 3
value = 9
[../]
[./left_v]
type = DirichletBC
variable = v
boundary = 1
value = 5
[../]
[./right_v]
type = DirichletBC
variable = v
boundary = 2
value = 2
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
dt = 0.01
num_steps = 10
[]
[Outputs]
file_base = out_nonexistent
exodus = true
show = 'u elemental nodal x foo1 foo2'
[]
(test/tests/misc/check_error/scalar_kernel_with_var.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 2
ny = 2
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[]
[AuxVariables]
[./v]
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[./rea]
type = Reaction
variable = u
[../]
[]
[ScalarKernels]
[./nope]
type = ODETimeDerivative
variable = u
[../]
[]
[BCs]
[./left]
type = DirichletBC
variable = u
boundary = 1
value = 0
[../]
[./right]
type = DirichletBC
variable = u
boundary = 2
value = 1
[../]
[]
[Executioner]
type = Transient
solve_type = 'NEWTON'
[]
[Outputs]
file_base = out
[]
(test/tests/postprocessors/scalar_variable/scalar_variable_pps.i)
[Problem]
use_hash_table_matrix_assembly = true
[]
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[time]
type = TimeDerivative
variable = u
[]
[]
[ScalarKernels]
[time]
type = ODETimeDerivative
variable = v
[]
[flux_sink]
type = PostprocessorSinkScalarKernel
variable = v
postprocessor = scale_flux
[]
[]
[BCs]
[right]
type = DirichletBC
value = 0
variable = u
boundary = 'right'
[]
[left]
type = ADMatchedScalarValueBC
variable = u
v = v
boundary = 'left'
[]
[]
[Variables]
[u][]
[v]
family = SCALAR
order = FIRST
initial_condition = 1
[]
[]
[Postprocessors]
[flux]
type = SideDiffusiveFluxIntegral
variable = u
diffusivity = 1
boundary = 'left'
execute_on = 'initial nonlinear linear timestep_end'
[]
[scale_flux]
type = ScalePostprocessor
scaling_factor = -1
value = flux
execute_on = 'initial nonlinear linear timestep_end'
[]
[reporter]
type = ScalarVariable
variable = v
execute_on = 'initial timestep_end'
[]
[]
[Executioner]
type = Transient
dt = .1
end_time = 1
solve_type = PJFNK
nl_rel_tol = 1e-12
[]
[Outputs]
csv = true
[]
(test/tests/time_integrators/scalar/scalar.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 1
ny = 1
elem_type = QUAD4
[]
[Variables]
[./n]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[]
[ScalarKernels]
[./dn]
type = ODETimeDerivative
variable = n
[../]
[./ode1]
type = ParsedODEKernel
expression = '-n'
variable = n
# implicit = false
[../]
[]
[Executioner]
type = Transient
[./TimeIntegrator]
# type = ImplicitEuler
# type = BDF2
type = CrankNicolson
# type = ImplicitMidpoint
# type = LStableDirk2
# type = LStableDirk3
# type = LStableDirk4
# type = AStableDirk4
#
# Explicit methods
# type = ExplicitEuler
# type = ExplicitMidpoint
# type = Heun
# type = Ralston
[../]
start_time = 0
end_time = 1
dt = 0.001
dtmin = 0.001 # Don't allow timestep cutting
solve_type = 'PJFNK'
nl_max_its = 2
nl_abs_tol = 1.e-12 # This is an ODE, so nl_abs_tol makes sense.
[]
[Functions]
[./exact_solution]
type = ParsedFunction
expression = exp(t)
[../]
[]
[Postprocessors]
[./error_n]
# Post processor that computes the difference between the computed
# and exact solutions. For the exact solution used here, the
# error at the final time should converge at O(dt^p), where p is
# the order of the method.
type = ScalarL2Error
variable = n
function = exact_solution
# final is not currently supported for Postprocessor execute_on...
# execute_on = 'final'
[../]
[]
[Outputs]
csv = true
[]
(test/tests/outputs/variables/output_vars_hidden_shown_check.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
elem_type = QUAD9
[]
[Variables]
[./u]
order = FIRST
family = LAGRANGE
[../]
[./v]
order = SECOND
family = LAGRANGE
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[AuxVariables]
[./elemental]
order = CONSTANT
family = MONOMIAL
[../]
[./elemental_restricted]
order = CONSTANT
family = MONOMIAL
[../]
[./nodal]
order = FIRST
family = LAGRANGE
[../]
[./nodal_restricted]
order = FIRST
family = LAGRANGE
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff_u]
type = Diffusion
variable = u
[../]
[./conv_u]
type = CoupledForce
variable = u
v = v
[../]
[./diff_v]
type = Diffusion
variable = v
[../]
[]
[AuxKernels]
[./elemental]
type = ConstantAux
variable = elemental
value = 1
[../]
[./elemental_restricted]
type = ConstantAux
variable = elemental_restricted
value = 1
[../]
[./nodal]
type = ConstantAux
variable = elemental
value = 2
[../]
[./nodal_restricted]
type = ConstantAux
variable = elemental_restricted
value = 2
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ImplicitODEx
variable = x
y = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ImplicitODEy
variable = y
x = x
[../]
[]
[BCs]
active = 'left_u right_u left_v'
[./left_u]
type = DirichletBC
variable = u
boundary = 1
value = 1
[../]
[./right_u]
type = DirichletBC
variable = u
boundary = 3
value = 9
[../]
[./left_v]
type = DirichletBC
variable = v
boundary = 1
value = 5
[../]
[./right_v]
type = DirichletBC
variable = v
boundary = 2
value = 2
[../]
[]
[Executioner]
type = Transient
solve_type = 'PJFNK'
dt = 0.01
num_steps = 10
[]
[Outputs]
file_base = out_hidden
exodus = true
hide = 'u elemental nodal x'
show = u
[]
(test/tests/scaling/ignore-variables/ignore.i)
[Mesh]
type = GeneratedMesh
dim = 1
nx = 3
[]
[Variables]
[u][]
[v][]
[x]
family = SCALAR
type = MooseVariableBase
[]
[y]
family = SCALAR
[]
[]
[Kernels]
[dt_u]
type = TimeDerivative
variable = u
[]
[diff_u]
type = Diffusion
variable = u
[]
[dt_v]
type = TimeDerivative
variable = v
[]
[diff_v]
type = MatDiffusion
variable = v
diffusivity = 1e-3
[]
[]
[ScalarKernels]
[dt_x]
type = ODETimeDerivative
variable = x
[]
[ode_x]
type = ParsedODEKernel
variable = x
coupled_variables = y
expression = '-3*x - 2*y'
[]
[dt_y]
type = ODETimeDerivative
variable = y
[]
[ode_y ]
type = ParsedODEKernel
variable = y
expression = '10*y'
[]
[]
[Executioner]
type = Transient
num_steps = 2
automatic_scaling = true
compute_scaling_once = false
ignore_variables_for_autoscaling = 'v y'
solve_type = NEWTON
verbose = true
[]
(examples/ex18_scalar_kernel/ex18.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 10
ny = 10
elem_type = QUAD4
[]
[Functions]
# ODEs
[./exact_x_fn]
type = ParsedFunction
expression = (-1/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[./exact_y_fn]
type = ParsedFunction
expression = (2/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[]
[Variables]
[./diffused]
order = FIRST
family = LAGRANGE
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = diffused
[../]
[./diff]
type = Diffusion
variable = diffused
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ImplicitODEx
variable = x
y = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ImplicitODEy
variable = y
x = x
[../]
[]
[BCs]
[./right]
type = ScalarDirichletBC
variable = diffused
boundary = 1
scalar_var = x
[../]
[./left]
type = ScalarDirichletBC
variable = diffused
boundary = 3
scalar_var = y
[../]
[]
[Postprocessors]
# to print the values of x, y into a file so we can plot it
[./x_pp]
type = ScalarVariable
variable = x
execute_on = timestep_end
[../]
[./y_pp]
type = ScalarVariable
variable = y
execute_on = timestep_end
[../]
[./exact_x]
type = FunctionValuePostprocessor
function = exact_x_fn
execute_on = timestep_end
point = '0 0 0'
[../]
[./exact_y]
type = FunctionValuePostprocessor
function = exact_y_fn
execute_on = timestep_end
point = '0 0 0'
[../]
# Measure the error in ODE solution for 'x'.
[./l2err_x]
type = ScalarL2Error
variable = x
function = exact_x_fn
[../]
# Measure the error in ODE solution for 'y'.
[./l2err_y]
type = ScalarL2Error
variable = y
function = exact_y_fn
[../]
[]
[Executioner]
type = Transient
start_time = 0
dt = 0.01
num_steps = 10
#Preconditioned JFNK (default)
solve_type = 'PJFNK'
[]
[Outputs]
exodus = true
[]
(modules/thermal_hydraulics/test/tests/scalarkernels/postprocessor_source/postprocessor_source.i)
# This input file tests PostprocessorSourceScalarKernel.
#
# The following initial value problem is modeled here:
# du/dt = t, u(0) = 0
# Using backward Euler time integration with dt=1, the solution values should
# be as follows:
# u(0) = 0
# u(1) = 1
# u(2) = 3
[Mesh]
type = GeneratedMesh
dim = 1
nx = 1
[]
[Variables]
[u]
family = SCALAR
order = FIRST
[]
[]
[ICs]
[ic_u]
type = ScalarConstantIC
variable = u
value = 0
[]
[]
[ScalarKernels]
[sk_time]
type = ODETimeDerivative
variable = u
[]
[sk_source]
type = PostprocessorSourceScalarKernel
variable = u
pp = pp_source
[]
[]
[Functions]
[fn_source]
type = ParsedFunction
expression = 't'
[]
[]
[Postprocessors]
[pp_source]
type = FunctionValuePostprocessor
function = fn_source
execute_on = 'LINEAR NONLINEAR'
[]
[]
[Executioner]
type = Transient
scheme = implicit-euler
dt = 1
num_steps = 2
[]
[Outputs]
csv = true
show = 'u'
execute_on = 'INITIAL TIMESTEP_END'
[]
(test/tests/bcs/periodic/no_add_scalar.i)
# Test to make sure that periodic boundaries
# are not applied to scalar variables.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 20
[]
[Variables]
[./c]
[./InitialCondition]
type = FunctionIC
function = x
[../]
[../]
[./scalar]
family = SCALAR
[../]
[]
[BCs]
[./Periodic]
[./all]
auto_direction = x
[../]
[../]
[]
[Kernels]
[./dt]
type = TimeDerivative
variable = c
[../]
[./diff]
type = Diffusion
variable = c
[../]
[]
[ScalarKernels]
[./scalar]
type = ODETimeDerivative
variable = scalar
[../]
[]
[Executioner]
type = Transient
dt = 0.1
num_steps = 3
[]
[Outputs]
exodus = true
[]
(test/tests/time_integrators/actually_explicit_euler_verification/ee-ode.i)
# Tests that ActuallyExplicitEuler works with scalar variables.
#
# The ODE and IC used are the following:
# du/dt = 2, u(0) = 0
# Thus the solution is u(t) = 2*t.
[Mesh]
type = GeneratedMesh
dim = 1
nx = 1
[]
[Variables]
[./u]
family = SCALAR
order = FIRST
initial_condition = 0
[../]
[]
[ScalarKernels]
[./time]
type = ODETimeDerivative
variable = u
[../]
[./source]
type = ParsedODEKernel
variable = u
expression = -2
[../]
[]
[Executioner]
type = Transient
[./TimeIntegrator]
type = ActuallyExplicitEuler
[../]
dt = 1
num_steps = 5
[]
[Outputs]
csv = true
[]
(test/tests/ics/function_scalar_ic/function_scalar_ic.i)
[Mesh]
# a dummy mesh
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 1
ny = 1
elem_type = QUAD4
[]
[Variables]
[./n]
family = SCALAR
order = FIRST
[../]
[]
[Functions]
[./f]
type = ParsedFunction
expression = cos(t)
[../]
[]
[ICs]
[./f]
type = FunctionScalarIC
variable = n
function = f
[../]
[]
[ScalarKernels]
[./dn]
type = ODETimeDerivative
variable = n
[../]
[./ode1]
type = ParsedODEKernel
expression = '-n'
variable = n
[../]
[]
[Executioner]
type = Transient
start_time = 0
end_time = 1
dt = 0.01
scheme = bdf2
solve_type = 'PJFNK'
timestep_tolerance = 1e-12
[]
[Outputs]
csv = true
[]
(test/tests/scalar_kernels/ad_coupled_scalar/ad_coupled_scalar.i)
[Problem]
use_hash_table_matrix_assembly = true
[]
[Mesh]
type = GeneratedMesh
dim = 1
nx = 10
[]
[Kernels]
[diff]
type = Diffusion
variable = u
[]
[time]
type = TimeDerivative
variable = u
[]
[]
[ScalarKernels]
[time]
type = ODETimeDerivative
variable = v
[]
[flux_sink]
type = PostprocessorSinkScalarKernel
variable = v
postprocessor = scale_flux
[]
[]
[BCs]
[right]
type = DirichletBC
value = 0
variable = u
boundary = 'right'
[]
[left]
type = ADMatchedScalarValueBC
variable = u
v = v
boundary = 'left'
[]
[]
[Variables]
[u][]
[v]
family = SCALAR
order = FIRST
initial_condition = 1
[]
[]
[Postprocessors]
[flux]
type = SideDiffusiveFluxIntegral
variable = u
diffusivity = 1
boundary = 'left'
execute_on = 'initial nonlinear linear timestep_end'
[]
[scale_flux]
type = ScalePostprocessor
scaling_factor = -1
value = flux
execute_on = 'initial nonlinear linear timestep_end'
[]
[]
[Executioner]
type = Transient
dt = .1
end_time = 1
solve_type = PJFNK
nl_rel_tol = 1e-12
[]
[Outputs]
exodus = true
[]
(test/tests/postprocessors/scalar_coupled_postprocessor/scalar_coupled_postprocessor_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
nx = 5
ny = 5
xmax = 1
ymax = 1
elem_type = QUAD4
[]
[Variables]
[./u]
initial_condition = 1
[../]
[./scalar_variable]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./diff]
type = Diffusion
variable = u
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = scalar_variable
[../]
[]
[BCs]
[./leftDirichlet]
type = DirichletBC
variable = u
boundary = 'left'
value = 1
[../]
[./rightDirichlet]
type = DirichletBC
variable = u
boundary = 'right'
value = 0
[../]
[]
[Postprocessors]
[./totalFlux]
type = ScalarCoupledPostprocessor
variable = u
coupled_scalar = scalar_variable
boundary = left
[../]
[]
[Executioner]
type = Transient
dt = 1
num_steps = 1
solve_type = JFNK
l_max_its = 30
l_tol = 1e-6
nl_max_its = 20
nl_rel_tol = 1e-5
[]
[Outputs]
csv = true
[]
(test/tests/tag/scalar_tag_vector.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 1
ny = 1
elem_type = QUAD4
[]
[Variables]
[./n]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[]
[AuxVariables]
[./tag_vector_var1]
family = SCALAR
order = FIRST
[../]
[./tag_vector_var2]
family = SCALAR
order = FIRST
[../]
[./tag_matrix_var2]
family = SCALAR
order = FIRST
[../]
[]
[ScalarKernels]
[./dn]
type = ODETimeDerivative
variable = n
extra_matrix_tags = 'mat_tag1 mat_tag2'
extra_vector_tags = 'vec_tag1'
[../]
[./ode1]
type = ParsedODEKernel
expression = '-n'
variable = n
extra_matrix_tags = 'mat_tag1'
extra_vector_tags = 'vec_tag1'
[../]
[./ode2]
type = ParsedODEKernel
expression = '-n'
variable = n
vector_tags = 'vec_tag2'
matrix_tags = 'mat_tag2'
[../]
[]
[AuxScalarKernels]
[./TagVectorAux]
type = ScalarTagVectorAux
variable = tag_vector_var1
v = n
vector_tag = vec_tag1
[../]
[./TagVectorAux2]
type = ScalarTagVectorAux
variable = tag_vector_var2
v = n
vector_tag = vec_tag2
[../]
[./TagMatrixAux2]
type = ScalarTagMatrixAux
variable = tag_matrix_var2
v = n
matrix_tag = mat_tag2
[../]
[]
[Problem]
type = TagTestProblem
test_tag_vectors = 'time nontime residual vec_tag1 vec_tag2'
test_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_matrices = 'mat_tag1 mat_tag2'
extra_tag_vectors = 'vec_tag1 vec_tag2'
[]
[Executioner]
type = Transient
start_time = 0
num_steps = 10
dt = 0.001
dtmin = 0.001 # Don't allow timestep cutting
solve_type = NEWTON
nl_max_its = 2
nl_abs_tol = 1.e-12 # This is an ODE, so nl_abs_tol makes sense.
[]
[Functions]
[./exact_solution]
type = ParsedFunction
expression = exp(t)
[../]
[]
[Postprocessors]
[./error_n]
# Post processor that computes the difference between the computed
# and exact solutions. For the exact solution used here, the
# error at the final time should converge at O(dt^p), where p is
# the order of the method.
type = ScalarL2Error
variable = n
function = exact_solution
[../]
[]
[Outputs]
csv = true
[]
(test/tests/kernels/ode/parsedode_sys_impl_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
[]
[Functions]
[./f_fn]
type = ParsedFunction
expression = -4
[../]
[./bc_all_fn]
type = ParsedFunction
expression = x*x+y*y
[../]
# ODEs
[./exact_x_fn]
type = ParsedFunction
expression = (-1/3)*exp(-t)+(4/3)*exp(5*t)
[../]
[]
# NL
[Variables]
[./u]
family = LAGRANGE
order = FIRST
[../]
# ODE variables
[./x]
family = SCALAR
order = FIRST
initial_condition = 1
[../]
[./y]
family = SCALAR
order = FIRST
initial_condition = 2
[../]
[]
[Kernels]
[./td]
type = TimeDerivative
variable = u
[../]
[./diff]
type = Diffusion
variable = u
[../]
[./uff]
type = BodyForce
variable = u
function = f_fn
[../]
[]
[ScalarKernels]
[./td1]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ParsedODEKernel
expression = '-3*x - 2*y'
variable = x
coupled_variables = y
[../]
[./td2]
type = ODETimeDerivative
variable = y
[../]
[./ode2]
type = ParsedODEKernel
expression = '-4*x - y'
variable = y
coupled_variables = x
[../]
[]
[BCs]
[./all]
type = FunctionDirichletBC
variable = u
boundary = '0 1 2 3'
function = bc_all_fn
[../]
[]
[Postprocessors]
active = 'exact_x l2err_x'
[./exact_x]
type = FunctionValuePostprocessor
function = exact_x_fn
execute_on = 'initial timestep_end'
point = '0 0 0'
[../]
[./l2err_x]
type = ScalarL2Error
variable = x
function = exact_x_fn
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Transient
start_time = 0
dt = 0.01
num_steps = 100
solve_type = 'PJFNK'
[]
[Outputs]
file_base = ode_sys_impl_test_out
exodus = true
[]
(test/tests/kernels/ode/parsedode_pp_test.i)
[Mesh]
type = GeneratedMesh
dim = 2
xmin = 0
xmax = 1
ymin = 0
ymax = 1
nx = 2
ny = 2
elem_type = QUAD4
[]
[Variables]
[./x]
family = SCALAR
order = FIRST
initial_condition = 0
[../]
[]
[ScalarKernels]
[./dt]
type = ODETimeDerivative
variable = x
[../]
[./ode1]
type = ParsedODEKernel
expression = '-mytime'
postprocessors = mytime
variable = x
[../]
[]
[Postprocessors]
[./computed_x]
type = ScalarVariable
variable = x
execute_on = 'initial timestep_end'
[../]
[./mytime]
type = FunctionValuePostprocessor
function = t
execute_on = 'initial timestep_begin'
[../]
[./exact_x]
type = FunctionValuePostprocessor
function = '0.5*t^2'
execute_on = 'initial timestep_end'
[../]
[./l2err_x]
type = ScalarL2Error
variable = x
function = '0.5*t^2'
execute_on = 'initial timestep_end'
[../]
[]
[Executioner]
type = Transient
scheme = bdf2
dt = 0.1
num_steps = 10
solve_type = 'NEWTON'
[]
[Outputs]
file_base = ode_pp_test_out
hide = 'x mytime'
csv = true
[]
(modules/thermal_hydraulics/include/scalarkernels/ODECoefTimeDerivative.h)
// This file is part of the MOOSE framework
// https://mooseframework.inl.gov
//
// 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 "ODETimeDerivative.h"
/**
* Time derivative multiplied by a coefficient for ODEs
*/
class ODECoefTimeDerivative : public ODETimeDerivative
{
public:
ODECoefTimeDerivative(const InputParameters & parameters);
protected:
virtual Real computeQpResidual() override;
virtual Real computeQpJacobian() override;
/// Coefficient that the time derivative terms is multiplied with
const Real & _coef;
public:
static InputParameters validParams();
};