Symbolic Metadata

It is possible to add metadata to symbolic variables, the metadata will be displayed when calling help on a variable.

The following information can be added (note, it's possible to extend this to user-defined metadata as well)

Variable descriptions

Descriptive strings can be attached to variables using the [description = "descriptive string"] syntax:

using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D
@variables u [description = "This is my input"]
getdescription(u)
"This is my input"

When variables with descriptions are present in systems, they will be printed when the system is shown in the terminal:

@variables u(t) [description = "A short description of u"]
@parameters p [description = "A description of p"]
@named sys = ODESystem([u ~ p], t)
Model sys:
Equations (1):
  1 standard: see equations(sys)
Unknowns (1): see unknowns(sys)
  u(t): A short description of u
Parameters (1): see parameters(sys)
  p: A description of p

Calling help on the variable u displays the description, alongside other metadata:

help?> u

  A variable of type Symbolics.Num (Num wraps anything in a type that is a subtype of Real)

  Metadata
  ≡≡≡≡≡≡≡≡≡≡

  ModelingToolkit.VariableDescription: This is my input

  Symbolics.VariableSource: (:variables, :u)

Connect

Variables in connectors can have connect metadata which describes the type of connections.

Flow is used for variables that represent physical quantities that "flow" ex: current in a resistor. These variables sum up to zero in connections.

Stream can be specified for variables that flow bi-directionally.

using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D

@variables i(t) [connect = Flow]
@variables k(t) [connect = Stream]
hasconnect(i)
true
getconnect(k)
Stream

Input or output

Designate a variable as either an input or an output using the following

using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D

@variables u [input = true]
isinput(u)
true
@variables y [output = true]
isoutput(y)
true

Bounds

Bounds are useful when parameters are to be optimized, or to express intervals of uncertainty.

@variables u [bounds = (-1, 1)]
hasbounds(u)
true
getbounds(u)
(-1, 1)

Bounds can also be specified for array variables. A scalar array bound is applied to each element of the array. A bound may also be specified as an array, in which case the size of the array must match the size of the symbolic variable.

@variables x[1:2, 1:2] [bounds = (-1, 1)]
hasbounds(x)
true
getbounds(x)
([-1 -1; -1 -1], [1 1; 1 1])
getbounds(x[1, 1])
(-1, 1)
getbounds(x[1:2, 1])
([-1, -1], [1, 1])
@variables x[1:2] [bounds = (-Inf, [1.0, Inf])]
hasbounds(x)
true
getbounds(x)
([-Inf, -Inf], [1.0, Inf])
getbounds(x[2])
(-Inf, Inf)
hasbounds(x[2])
false

Guess

Specify an initial guess for custom initial conditions of an ODESystem.

@variables u [guess = 1]
hasguess(u)
true
getguess(u)
1

Mark input as a disturbance

Indicate that an input is not available for control, i.e., it's a disturbance input.

@variables u [input = true, disturbance = true]
isdisturbance(u)
true

Mark parameter as tunable

Indicate that a parameter can be automatically tuned by parameter optimization or automatic control tuning apps.

@parameters Kp [tunable = true]
istunable(Kp)
true

Probability distributions

A probability distribution may be associated with a parameter to indicate either uncertainty about its value, or as a prior distribution for Bayesian optimization.

using Distributions
d = Normal(10, 1)
@parameters m [dist = d]
hasdist(m)
getdist(m)

Irreducible

A variable can be marked irreducible to prevent it from being moved to an observed state. This forces the variable to be computed during solving so that it can be accessed in callbacks

@variables important_value [irreducible = true]
isirreducible(important_value)
true

State Priority

When a model is structurally simplified, the algorithm will try to ensure that the variables with higher state priority become states of the system. A variable's state priority is a number set using the state_priority metadata.

@variables important_dof [state_priority = 10] unimportant_dof [state_priority = -2]
state_priority(important_dof)
10.0

Units

Units for variables can be designated using symbolic metadata. For more information, please see the model validation and units section of the docs. Note that getunit is not equivalent to get_unit - the former is a metadata getter for individual variables (and is provided so the same interface function for unit exists like other metadata), while the latter is used to handle more general symbolic expressions.

using DynamicQuantities
@variables speed [unit = u"m/s"]
hasunit(speed)
true
getunit(speed)
1.0 m s⁻¹

Miscellaneous metadata

User-defined metadata can be added using the misc metadata. This can be queried using the hasmisc and getmisc functions.

@variables u [misc = :conserved_parameter] y [misc = [2, 4, 6]]
hasmisc(u)
true
getmisc(y)
3-element Vector{Int64}:
 2
 4
 6

Additional functions

For systems that contain parameters with metadata like described above, have some additional functions defined for convenience. In the example below, we define a system with tunable parameters and extract bounds vectors

@variables x(t)=0 u(t)=0 [input = true] y(t)=0 [output = true]
@parameters T [tunable = true, bounds = (0, Inf)]
@parameters k [tunable = true, bounds = (0, Inf)]
eqs = [D(x) ~ (-x + k * u) / T # A first-order system with time constant T and gain k
       y ~ x]
sys = ODESystem(eqs, t, name = :tunable_first_order)

\[ \begin{align} \frac{\mathrm{d} x\left( t \right)}{\mathrm{d}t} &= \frac{ - x\left( t \right) + k u\left( t \right)}{T} \\ y\left( t \right) &= x\left( t \right) \end{align} \]

p = tunable_parameters(sys) # extract all parameters marked as tunable
2-element Vector{SymbolicUtils.BasicSymbolic{Real}}:
 k
 T
lb, ub = getbounds(p) # operating on a vector, we get lower and upper bound vectors
(lb = [0, 0], ub = [Inf, Inf])
b = getbounds(sys) # Operating on the system, we get a dict
Dict{SymbolicUtils.BasicSymbolic{Real}, Tuple{Int64, Float64}} with 2 entries:
  k => (0, Inf)
  T => (0, Inf)

See also: ModelingToolkit.dump_variable_metadata, ModelingToolkit.dump_parameters, ModelingToolkit.dump_unknowns.

Index

Docstrings

ModelingToolkit.getboundsFunction
getbounds(sys::ModelingToolkit.AbstractSystem, p = parameters(sys))

Returns a dict with pairs p => (lb, ub) mapping parameters of sys to lower and upper bounds. Create parameters with bounds like this

@parameters p [bounds=(-1, 1)]

To obtain unknown variable bounds, call getbounds(sys, unknowns(sys))

source
ModelingToolkit.getboundsMethod
getbounds(x)

Get the bounds associated with symbolic variable x. Create parameters with bounds like this

@parameters p [bounds=(-1, 1)]
source
ModelingToolkit.getdistMethod
getdist(x)

Get the probability distribution associated with symbolic variable x. If no distribution is associated with x, nothing is returned. Create parameters with associated distributions like this

using Distributions
d = Normal(0, 1)
@parameters u [dist = d]
hasdist(u) # true
getdist(u) # retrieve distribution
source
ModelingToolkit.getguessMethod
getguess(x)

Get the guess for the initial value associated with symbolic variable x. Create variables with a guess like this

@variables x [guess=1]
source
ModelingToolkit.getunitMethod
getunit(x)

Fetch the unit associated with variable x. This function is a metadata getter for an individual variable, while get_unit is used for unit inference on more complicated sdymbolic expressions.

source
ModelingToolkit.istunableFunction
istunable(x, default = true)

Determine whether symbolic variable x is marked as a tunable for an automatic tuning algorithm.

default indicates whether variables without tunable metadata are to be considered tunable or not.

Create a tunable parameter by

@parameters u [tunable=true]

See also tunable_parameters, getbounds

source
ModelingToolkit.tunable_parametersFunction
tunable_parameters(sys, p = parameters(sys; initial_parameters = true); default=true)

Get all parameters of sys that are marked as tunable.

Keyword argument default indicates whether variables without tunable metadata are to be considered tunable or not.

Create a tunable parameter by

@parameters u [tunable=true]

For systems created with split = true (the default) and default = true passed to this function, the order of parameters returned is the order in which they are stored in the tunables portion of MTKParameters. Note that array variables will not be scalarized. To obtain the flattened representation of the tunables portion, call Symbolics.scalarize(tunable_parameters(sys)) and concatenate the resulting arrays.

See also getbounds, istunable, MTKParameters, complete

source
ModelingToolkit.dump_variable_metadataFunction
dump_variable_metadata(var)

Return all the metadata associated with symbolic variable var as a NamedTuple.

using ModelingToolkit

@parameters p::Int [description = "My description", bounds = (0.5, 1.5)]
ModelingToolkit.dump_variable_metadata(p)
source
ModelingToolkit.dump_parametersFunction
dump_parameters(sys::AbstractSystem)

Return an array of NamedTuples containing the metadata associated with each parameter in sys. Also includes the default value of the parameter, if provided.

using ModelingToolkit
using DynamicQuantities
using ModelingToolkit: t, D

@parameters p = 1.0, [description = "My parameter", tunable = false] q = 2.0, [description = "Other parameter"]
@variables x(t) = 3.0 [unit = u"m"]
@named sys = ODESystem(Equation[], t, [x], [p, q])
ModelingToolkit.dump_parameters(sys)

See also: ModelingToolkit.dump_variable_metadata, ModelingToolkit.dump_unknowns

source
ModelingToolkit.dump_unknownsFunction
dump_unknowns(sys::AbstractSystem)

Return an array of NamedTuples containing the metadata associated with each unknown in sys. Also includes the default value of the unknown, if provided.

using ModelingToolkit
using DynamicQuantities
using ModelingToolkit: t, D

@parameters p = 1.0, [description = "My parameter", tunable = false] q = 2.0, [description = "Other parameter"]
@variables x(t) = 3.0 [unit = u"m"]
@named sys = ODESystem(Equation[], t, [x], [p, q])
ModelingToolkit.dump_unknowns(sys)

See also: ModelingToolkit.dump_variable_metadata, ModelingToolkit.dump_parameters

source