Skip to content

Commit 65c9a04

Browse files
Merge pull request #2542 from AayushSabharwal/as/update-docs
docs: update FAQ section with relevant MTKv9 questions
2 parents 57ccce8 + 28ff9c0 commit 65c9a04

File tree

2 files changed

+106
-20
lines changed

2 files changed

+106
-20
lines changed

docs/src/basics/Events.md

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -164,7 +164,7 @@ documentation. In affect functions, we have that
164164
function affect!(integ, u, p, ctx)
165165
# integ.t is the current time
166166
# integ.u[u.v] is the value of the unknown `v` above
167-
# integ.p[p.q] is the value of the parameter `q` above
167+
# integ.ps[p.q] is the value of the parameter `q` above
168168
end
169169
```
170170

docs/src/basics/FAQ.md

Lines changed: 105 additions & 19 deletions
Original file line numberDiff line numberDiff line change
@@ -1,25 +1,97 @@
11
# Frequently Asked Questions
22

3+
## Why are my parameters some obscure object?
4+
5+
In ModelingToolkit.jl version 9, the parameter vector was replaced with a custom
6+
`MTKParameters` object, whose internals are intentionally undocumented and subject
7+
to change without a breaking release. This enables us to efficiently store and generate
8+
code for parameters of multiple types. To obtain parameter values use
9+
[SymbolicIndexingInterface.jl](https://github.com/SciML/SymbolicIndexingInterface.jl/) or
10+
[SciMLStructures.jl](https://github.com/SciML/SciMLStructures.jl/). For example:
11+
12+
```julia
13+
prob.ps[lorenz.β] # obtains the value of parameter `β`. Note the `.ps` instead of `.p`
14+
getβ = getp(prob, lorenz.β) # returns a function that can fetch the value of `β`
15+
getβ(sol) # can be used on any object that is based off of the same system
16+
getβ(prob)
17+
```
18+
19+
Indexes into the `MTKParameters` object take the form of `ParameterIndex` objects, which
20+
are similarly undocumented. Following is the list of behaviors that should be relied on for
21+
`MTKParameters`:
22+
23+
- It implements the SciMLStructures interface.
24+
- It can be queried for parameters using functions returned from
25+
`SymbolicIndexingInterface.getp`.
26+
- `getindex(::MTKParameters, ::ParameterIndex)` can be used to obtain the value of a
27+
parameter with the given index.
28+
- `setindex!(::MTKParameters, value, ::ParameterIndex)` can be used to set the value of a
29+
parameter with the given index.
30+
- `parameter_values(sys, sym)` will return a `ParameterIndex` object if `sys` has been
31+
`complete`d (through `structural_simplify`, `complete` or `@mtkbuild`).
32+
- `copy(::MTKParameters)` is defined and duplicates the parameter object, including the
33+
memory used by the underlying buffers.
34+
35+
Any other behavior of `MTKParameters` (other `getindex`/`setindex!` methods, etc.) is an
36+
undocumented internal and should not be relied upon.
37+
38+
## How do I use non-numeric/array-valued parameters?
39+
40+
In ModelingToolkit.jl version 9, parameters are required to have a `symtype` matching
41+
the type of their values. For example, this will error during problem construction:
42+
43+
```julia
44+
@parameters p = [1, 2, 3]
45+
```
46+
47+
Since by default parameters have a `symtype` of `Real` (which is interpreted as `Float64`)
48+
but the default value given to it is a `Vector{Int}`. For array-valued parameters, use the
49+
following syntax:
50+
51+
```julia
52+
@parameters p[1:n, 1:m]::T # `T` is the `eltype` of the parameter array
53+
@parameters p::T # `T` is the type of the array
54+
```
55+
56+
The former approach is preferred, since the size of the array is known. If the array is not
57+
a `Base.Array` or the size is not known during model construction, the second syntax is
58+
required.
59+
60+
The same principle applies to any parameter type that is not `Float64`.
61+
62+
```julia
63+
@parameters p1::Int # integer-valued
64+
@parameters p2::Bool # boolean-valued
65+
@parameters p3::MyCustomStructType # non-numeric
66+
@parameters p4::ComponentArray{...} # non-standard array
67+
```
68+
369
## Getting the index for a symbol
470

5-
Since **ordering of symbols is not guaranteed after symbolic transformations**,
6-
one should normally refer to values by their name. For example, `sol[lorenz.x]`
7-
from the solution. But what if you need to get the index? The following helper
8-
function will do the trick:
71+
Ordering of symbols is not guaranteed after symbolic transformations, and parameters
72+
are now stored in a custom `MTKParameters` object instead of a vector. Thus, values
73+
should be referred to by their name. For example `sol[lorenz.x]`. To obtain the index,
74+
use the following functions from
75+
[SymbolicIndexingInterface.jl](https://github.com/SciML/SymbolicIndexingInterface.jl/):
976

1077
```julia
11-
indexof(sym, syms) = findfirst(isequal(sym), syms)
12-
indexof(σ, parameters(sys))
78+
variable_index(sys, sym)
79+
parameter_index(sys, sym)
1380
```
1481

82+
Note that while the variable index will be an integer, the parameter index is a struct of
83+
type `ParameterIndex` whose internals should not be relied upon.
84+
1585
## Transforming value maps to arrays
1686

1787
ModelingToolkit.jl allows (and recommends) input maps like `[x => 2.0, y => 3.0]`
1888
because symbol ordering is not guaranteed. However, what if you want to get the
19-
lowered array? You can use the internal function `varmap_to_vars`. For example:
89+
lowered array? You can use the internal function `varmap_to_vars` for unknowns.
90+
and the `MTKParameters` constructor for parameters. For example:
2091

2192
```julia
22-
pnew = varmap_to_vars([β => 3.0, c => 10.0, γ => 2.0], parameters(sys))
93+
unew = varmap_to_vars([x => 1.0, y => 2.0, z => 3.0], unknowns(sys))
94+
pnew = ModelingToolkit.MTKParameters(sys, [β => 3.0, c => 10.0, γ => 2.0])
2395
```
2496

2597
## How do I handle `if` statements in my symbolic forms?
@@ -63,37 +135,51 @@ end
63135
Since `ODEProblem` on a MTK `sys` will have to generate code, this will be slower than
64136
caching the generated code, and will require automatic differentiation to go through the
65137
code generation process itself. All of this is unnecessary. Instead, generate the problem
66-
once outside the loss function, and remake the prob inside the loss function:
138+
once outside the loss function, and update the parameter values inside the loss function:
67139

68140
```julia
69141
prob = ODEProblem(sys, [], [p1 => p[1], p2 => p[2]])
70142
function loss(p)
71-
remake(prob, p = ...)
143+
# update parameters
72144
sol = solve(prob, Tsit5())
73145
sum(abs2, sol)
74146
end
75147
```
76148

77-
Now, one has to be careful with `remake` to ensure that the parameters are in the right
78-
order. One can use the previously mentioned indexing functionality to generate index
79-
maps for reordering `p` like:
149+
If a subset of the parameters are optimized, `setp` from SymbolicIndexingInterface.jl
150+
should be used to generate an efficient function for setting parameter values. For example:
80151

81152
```julia
82-
p = @parameters x y z
83-
idxs = ModelingToolkit.varmap_to_vars([p[1] => 1, p[2] => 2, p[3] => 3], p)
84-
p[idxs]
153+
using SymbolicIndexingInterface
154+
155+
prob = ODEProblem(sys, [], [p1 => p[1], p2 => p[2]])
156+
setter! = setp(sys, [p1, p2])
157+
function loss(p)
158+
setter!(prob, p)
159+
sol = solve(prob, Tsit5())
160+
sum(abs2, sol)
161+
end
85162
```
86163

87-
Using this, the fixed index map can be used in the loss function. This would look like:
164+
[SciMLStructures.jl](https://github.com/SciML/SciMLStructures.jl/) can be leveraged to
165+
obtain all the parameters for optimization using the `Tunable` portion. By default, all
166+
numeric or numeric array parameters are marked as tunable, unless explicitly marked as
167+
`tunable = false` in the variable metadata.
88168

89169
```julia
170+
using SciMLStructures: replace!, Tunable
171+
90172
prob = ODEProblem(sys, [], [p1 => p[1], p2 => p[2]])
91-
idxs = Int.(ModelingToolkit.varmap_to_vars([p1 => 1, p2 => 2], p))
92173
function loss(p)
93-
remake(prob, p = p[idxs])
174+
replace!(Tunable(), prob.p, p)
94175
sol = solve(prob, Tsit5())
95176
sum(abs2, sol)
96177
end
178+
179+
p, replace, alias = SciMLStructures.canonicalize(Tunable(), prob.p)
180+
# p is an `AbstractVector` which can be optimized
181+
# if `alias == true`, then `p` aliases the memory used by `prob.p`, so
182+
# changes to the array will be reflected in parameter values
97183
```
98184

99185
# ERROR: ArgumentError: SymbolicUtils.BasicSymbolic{Real}[xˍt(t)] are missing from the variable map.

0 commit comments

Comments
 (0)