- Notifications
You must be signed in to change notification settings - Fork38
Arrays with arbitrarily nested named components.
License
SciML/ComponentArrays.jl
Folders and files
Name | Name | Last commit message | Last commit date | |
---|---|---|---|---|
Repository files navigation
Documentation | Build Status |
---|---|
The main export of this package is theComponentArray
type. "Components" ofComponentArray
sare really just array blocks that can be accessed through a named index. This will create a newComponentArray
whose data is a view into the original,allowing for standalone models to be composed together by simple function composition. Inessence,ComponentArray
s allow you to do the things you would usually need a modelinglanguage for, but without actually needing a modeling language. The main targets are for useinDifferentialEquations.jl andOptim.jl, but anything that requiresflat vectors is fair game.
Check out theNEWS for new features by minor release version.
The easiest way to construct 1-dimensionalComponentArray
s (aliased asComponentVector
) is as if they wereNamedTuple
s. In fact, a good way to think about them is as arbitrarily nested, mutableNamedTuple
s that can be passed through a solver.
julia> c= (a=2, b=[1,2]);julia> x=ComponentArray(a=5, b=[(a=20., b=0), (a=33., b=0), (a=44., b=3)], c=c)ComponentVector{Float64}(a=5.0, b= [(a=20.0, b=0.0), (a=33.0, b=0.0), (a=44.0, b=3.0)], c= (a=2.0, b= [1.0,2.0]))julia> x.c.a=400; xComponentVector{Float64}(a=5.0, b= [(a=20.0, b=0.0), (a=33.0, b=0.0), (a=44.0, b=3.0)], c= (a=400.0, b= [1.0,2.0]))julia> x[8]400.0julia>collect(x)10-element Vector{Float64}:5.020.00.033.00.044.03.0400.01.02.0julia>typeof(similar(x, Int32))===typeof(ComponentVector{Int32}(a=5, b=[(a=20., b=0), (a=33., b=0), (a=44., b=3)], c=c))true
ComponentArray
s can be constructed from existingComponentArray
s (currently nested fields cannot be changed this way):
julia> x=ComponentVector(a=1, b=2, c=3);julia>ComponentVector(x; a=11, new=42)ComponentVector{Int64}(a=11, b=2, c=3, new=42)
Higher dimensionalComponentArray
s can be created too, but it's a little messy at the moment. The nice thing for modeling is that dimension expansion through broadcasted operations can create higher-dimensionalComponentArray
s automatically, so Jacobian cache arrays that are created internally withfalse .* x .* x'
will be two-dimensionalComponentArray
s (aliased asComponentMatrix
) with proper axes. Check out theODE with Jacobian example in the examples folder to see how this looks in practice.
julia> x=ComponentArray(a=1, b=[2,1,4.0], c=c)ComponentVector{Float64}(a=1.0, b= [2.0,1.0,4.0], c= (a=2.0, b= [1.0,2.0]))julia> x2= x.* x'7×7 ComponentMatrix{Float64} with axesAxis(a=1, b=2:4, c=ViewAxis(5:7,Axis(a=1, b=2:3)))×Axis(a=1, b=2:4, c=ViewAxis(5:7,Axis(a=1, b=2:3)))1.02.01.04.02.01.02.02.04.02.08.04.02.04.01.02.01.04.02.01.02.04.08.04.016.08.04.08.02.04.02.08.04.02.04.01.02.01.04.02.01.02.02.04.02.08.04.02.04.0julia> x2[:c,:c]3×3 ComponentMatrix{Float64} with axesAxis(a=1, b=2:3)×Axis(a=1, b=2:3)4.02.04.02.01.02.04.02.04.0julia> x2[:a,:a]1.0julia>@view x2[:a,:c]ComponentVector{Float64,SubArray...}(a=2.0, b= [1.0,2.0])julia> x2[:b,:c]3×3 ComponentMatrix{Float64} with axesFlatAxis()×Axis(a=1, b=2:3)4.02.04.02.01.02.08.04.08.0
This example uses@unpack
fromParameters.jlfor nice syntax. Example taken from:SciML/ModelingToolkit.jl#36 (comment)
using ComponentArraysusing DifferentialEquationsusing Parameters:@unpacktspan= (0.0,20.0)## Lorenz systemfunctionlorenz!(D, u, p, t; f=0.0)@unpack σ, ρ, β= p@unpack x, y, z= u D.x= σ*(y- x) D.y= x*(ρ- z)- y- f D.z= x*y- β*zreturnnothingendlorenz_p= (σ=10.0, ρ=28.0, β=8/3)lorenz_ic=ComponentArray(x=0.0, y=0.0, z=0.0)lorenz_prob=ODEProblem(lorenz!, lorenz_ic, tspan, lorenz_p)## Lotka-Volterra systemfunctionlotka!(D, u, p, t; f=0.0)@unpack α, β, γ, δ= p@unpack x, y= u D.x= α*x- β*x*y+ f D.y=-γ*y+ δ*x*yreturnnothingendlotka_p= (α=2/3, β=4/3, γ=1.0, δ=1.0)lotka_ic=ComponentArray(x=1.0, y=1.0)lotka_prob=ODEProblem(lotka!, lotka_ic, tspan, lotka_p)## Composed Lorenz and Lotka-Volterra systemfunctioncomposed!(D, u, p, t) c= p.c#coupling parameter@unpack lorenz, lotka= ulorenz!(D.lorenz, lorenz, p.lorenz, t, f=c*lotka.x)lotka!(D.lotka, lotka, p.lotka, t, f=c*lorenz.x)returnnothingendcomp_p= (lorenz=lorenz_p, lotka=lotka_p, c=0.01)comp_ic=ComponentArray(lorenz=lorenz_ic, lotka=lotka_ic)comp_prob=ODEProblem(composed!, comp_ic, tspan, comp_p)## Solve problem# We can solve the composed system...comp_sol=solve(comp_prob)# ...or we can unit test one of the component systemslotka_sol=solve(lotka_prob)
Notice how cleanly thecomposed!
function can pass variables from one function to another with no array index juggling in sight. This is especially useful for large models as it becomes harder to keep track top-level model array position when adding new or deleting old components from the model. We could go further and composecomposed!
with other components ad (practically) infinitum with no mental bookkeeping.
The main benefit, however, is now our differential equations are unit testable. Bothlorenz
andlotka
can be run as their ownODEProblem
withf
set to zero to see the unforced response.
About
Arrays with arbitrarily nested named components.