API
The following functions are designed for public use.
Network Construction API
NetworkDynamics.Network
— TypeNetwork([g,] vertexf, edgef; kwarg...)
Construct a Network
object from a graph g
and edge and component models vertexf
and edgef
.
Arguments:
g::AbstractGraph
: The graph on which the network is defined. Optional, can be ommittet if all component models have a definedgraphelement
. Seevidx
andsrc
/dst
keywors forVertexModel
andEdgeModel
constructors respectively.vertexm
: A singleVertexModel
or a vector ofVertexModel
objects. The order of the vertex models must mirror the order of thevertices(g)
iterator.edgem
: A singleEdgeModel
or a vector ofEdgeModel
objects. The order of the edge models must mirror the order of theedges(g)
iterator.
Optional keyword arguments:
execution=SequentialExecution{true}()
: Execution model of the network. E.g.SequentialExecution
,KAExecution
,PolyesterExecution
orThreadedExecution
.aggregator=execution isa SequentialExecution ? SequentialAggregator(+) : PolyesterAggregator(+)
: Aggregation function applied to the edge models. E.g.SequentialAggregator
,PolyesterAggregator
,ThreadedAggregator
,SparseAggregator
.check_graphelement=true
: Check if thegraphelement
metadata is consistent with the graph.dealias=false
Check if the components alias eachother and create copies if necessary. This is necessary if the same component model is referenced in multiple places in the Network but you want to dynamicially asign metadata, such as initialization information to specific instances.verbose=false
: Show additional information during construction.
Network(nw::Network; g, vertexm, edgem, kwargs...)
Rebuild the Network with same graph and vertex/edge models but possibly different kwargs.
NetworkDynamics.dim
— Methoddim(nw::Network)
Returns the number of dynamic states in the network, corresponts to the length of the flat state vector.
NetworkDynamics.pdim
— Methodpdim(nw::Network)
Returns the number of parameters in the network, corresponts to the length of the flat parameter vector.
Component Models
NetworkDynamics.VertexModel
— MethodVertexModel(; kwargs...)
Build a VertexModel
according to the keyword arguments.
Main Arguments:
f=nothing
: Dynamic function of the component. Can be nothing ifdim
is 0.g
: Output function of the component. Usefull helpers:StateMask
sym
/dim
: Symbolic names of the states. Ifdim
is provided,sym
is set automaticially.outsym
/outdim
: Symbolic names of the outputs. Ifoutdim
is provided,outsym
is set automaticially. Can be infered automaticially ifg
isaStateMask
.- psym
/
pdim=0: Symbolic names of the parameters. If
pdimis provided,
psym` is set automaticially. mass_matrix=I
: Mass matrix of component. Can be a vectorv
and is then interpreted asDiagonal(v)
.name=dim>0 ? :VertexM : :StaticVertexM
: Name of the component.
Optional Arguments:
insym
/indim
: Symbolic names of the inputs. Ifindim
is provided,insym
is set automaticially.vidx
: Index of the vertex in the graph, enables graphless constructor.ff
:FeedForwardType
of component. Will be typically infered fromg
automaticially.obssym
/obsf
: Define additional "observable" states.symmetadata
/metadata
: Provide prefilled metadata dictionaries.
All Symbol arguments can be used to set default values, i.e. psym=[:K=>1, :p]
.
NetworkDynamics.EdgeModel
— MethodEdgeModel(; kwargs...)
Build a EdgeModel
according to the keyword arguments.
Main Arguments:
f=nothing
: Dynamic function of the component. Can be nothing ifdim
is 0.g
: Output function of the component. Usefull helpers:AntiSymmetric
,Symmetric
,Fiducial
,Directed
andStateMask
.sym
/dim
: Symbolic names of the states. Ifdim
is provided,sym
is set automaticially.outsym
/outdim
: Symbolic names of the outputs. Ifoutdim
is provided,outsym
is set automaticially. In general, outsym for edges isa named tuple(; src, dst)
. However, depending on theg
function, it might be enough to provide a single vector or even nothing (e.g.AntiSymmetric(StateMask(1:2))
). See BuildingEdgeModel
s for examples.- psym
/
pdim=0: Symbolic names of the parameters. If
pdimis provided,
psym` is set automaticially. mass_matrix=I
: Mass matrix of component. Can be a vectorv
and is then interpreted asDiagonal(v)
.name=dim>0 ? :EdgeM : :StaticEdgeM
: Name of the component.
Optional Arguments:
insym
/indim
: Symbolic names of the inputs. Ifindim
is provided,insym
is set automaticially. For edges,insym
is a named tuple(; src, dst)
. If give as vector tuple is created automaticially.src
/dst
: Index or name of the vertices at src and dst end. Enables graphless constructor.ff
:FeedForwardType
of component. Will be typically infered fromg
automaticially.obssym
/obsf
: Define additional "observable" states.symmetadata
/metadata
: Provide prefilled metadata dictionaries.
All Symbol arguments can be used to set default values, i.e. psym=[:K=>1, :p]
.
Component Models with MTK
NetworkDynamics.VertexModel
— MethodVertexModel(sys::ODESystem, inputs, outputs; ff_to_constraint=true, kwargs...)
Create a VertexModel
object from a given ODESystem
created with ModelingToolkit. You need to provide 2 lists of symbolic names (Symbol
or Vector{Symbols}
):
inputs
: names of variables in you equation representing the aggregated edge statesoutputs
: names of variables in you equation representing the node output
ff_to_constraint
controlls, whether output transformations g
which depend on inputs should be
NetworkDynamics.EdgeModel
— MethodEdgeModel(sys::ODESystem, srcin, srcout, dstin, dstout; ff_to_constraint=false, kwargs...)
Create a EdgeModel
object from a given ODESystem
created with ModelingToolkit. You need to provide 4 lists of symbolic names (Symbol
or Vector{Symbols}
):
srcin
: names of variables in you equation representing the node state at the sourcedstin
: names of variables in you equation representing the node state at the destinationsrcout
: names of variables in you equation representing the output at the sourcedstout
: names of variables in you equation representing the output at the destination
ff_to_constraint
controlls, whether output transformations g
which depend on inputs should be transformed into constraints.
NetworkDynamics.EdgeModel
— MethodEdgeModel(sys::ODESystem, srcin, dstin, AntiSymmetric(dstout); ff_to_constraint=false, kwargs...)
Create a EdgeModel
object from a given ODESystem
created with ModelingToolkit.
Here you only need to provide one list of output symbols: dstout
. To make it clear how to handle the single-sided output definiton, you musst wrap the symbol vector in
AntiSymmetric(dstout)
,Symmetric(dstout)
, orDirected(dstout)
.
ff_to_constraint
controlls, whether output transformations g
which depend on inputs should be transformed into constraints.
Output Function Helpers/Wrappers
NetworkDynamics.StateMask
— TypeStateMask(i::AbstractArray)
StateMaks(i::Number)
A StateMask
is a predefined output function. It can be used to define the output of a component model by picking from the internal state.
I.e. g=StateMask(2:3)
in a vertex function will output the internal states 2 and 3. In many contexts, StateMask
s can be constructed implicitly by just providing the indices, e.g. g=1:2
.
For EdgeModel
this needs to be combined with a Directed
, Symmetric
, AntiSymmetric
or Fiducial
coupling, e.g. g=Fiducial(1:2, 3:4)
forwards states 1:2 to dst and states 3:4 to src.
NetworkDynamics.Symmetric
— TypeSymmetric(g)
Wraps a single-sided output function g
turns it into a double sided output function which applies
y_dst = g(...)
y_src = y_dst
g
can be a Number
/AbstractArray
to impicitly wrap the corresponding StateMask
.
See also AntiSymmetric
, Directed
, Fiducial
and StateMask
.
NetworkDynamics.AntiSymmetric
— TypeAntiSymmetric(g_dst)
Wraps a single-sided output function g_dst
turns it into a double sided output function which applies
y_dst = g_dst(...)
y_src = -y_dst
g_dst
can be a Number
/AbstractArray
to impicitly wrap the corresponding StateMask
.
NetworkDynamics.Directed
— TypeDirected(g_dst)
Wraps a single-sided output function g_dst
turns it into a double sided output function which applies
y_dst = g_dst(...)
With Directed
there is no output for the src
side. g_dst
can be a Number
/AbstractArray
to impicitly wrap the corresponding StateMask
.
See also AntiSymmetric
, Symmetric
, Fiducial
and StateMask
.
NetworkDynamics.Fiducial
— TypeFiducial(g_src, g_dst)
Wraps two single-sided output function g_src
and g_dst
and turns them into a double sided output function which applies
y_dst = g_src(...)
y_src = g_dst(...)
g
can be a Number
/AbstractArray
to impicitly wrap the corresponding StateMask
.
See also AntiSymmetric
, Directed
, Fiducial
and StateMask
.
Accessors for Component Properties
NetworkDynamics.fftype
— Functionfftype(x)
Retrieve the feed forward trait of x
.
NetworkDynamics.dim
— Methoddim(c::ComponentModel)::Int
Retrieve the dimension of the component.
NetworkDynamics.sym
— Functionsym(c::ComponentModel)::Vector{Symbol}
Retrieve the symbols of the component.
NetworkDynamics.outdim
— Functionoutdim(c::VertexModel)::Int
outdim(c::EdgeModel)::@NamedTuple(src::Int, dst::Int)
Retrieve the output dimension of the component
NetworkDynamics.outsym
— Functionoutsym(c::VertexModel)::Vector{Symbol} outsym(c::EdgeModel)::@NamedTuple{src::Vector{Symbol}, dst::Vector{Symbol}}
Retrieve the output symbols of the component.
NetworkDynamics.pdim
— Methodpdim(c::ComponentModel)::Int
Retrieve the parameter dimension of the component.
NetworkDynamics.psym
— Functionpsym(c::ComponentModel)::Vector{Symbol}
Retrieve the parameter symbols of the component.
NetworkDynamics.obssym
— Functionobssym(c::ComponentModel)::Vector{Symbol}
Retrieve the observation symbols of the component.
NetworkDynamics.hasinsym
— Functionhasinsym(c::ComponentModel)
Checks if the optioan field insym
is present in the component model.
NetworkDynamics.insym
— Functioninsym(c::VertexModel)::Vector{Symbol}
insym(c::EdgeModel)::@NamedTuple{src::Vector{Symbol}, dst::Vector{Symbol}}
Musst be called after hasinsym
/hasindim
returned true. Gives the insym
vector(s). For vertex model just a single vector, for edges it returns a named tuple (; src, dst)
with two symbol vectors.
NetworkDynamics.hasindim
— Functionhasindim(c::ComponentModel)
Checks if the optioan field insym
is present in the component model.
NetworkDynamics.indim
— Functionindim(c::VertexModel)::Int
indim(c::EdgeModel)::@NamedTuple{src::Int,dst::Int}
Musst be called after hasinsym
/hasindim
returned true. Gives the input dimension(s).
FeedForwardType
-Traits
NetworkDynamics.FeedForwardType
— Typeabstract type FeedForwardType end
Abstract supertype for the FeedForwardType traits.
NetworkDynamics.PureFeedForward
— TypePureFeedForward <: FeedForwardType
Trait for component output functions g
that have pure feed forward behavior (do not depend on x):
g!(outs..., ins..., p, t)
See also FeedForward
, NoFeedForward
and PureStateMap
.
NetworkDynamics.FeedForward
— TypeFeedForward <: FeedForwardType
Trait for component output functions g
that have feed forward behavior. May depend on everything:
g!(outs..., x, ins..., p, t)
See also PureFeedForward
, NoFeedForward
and PureStateMap
.
NetworkDynamics.NoFeedForward
— TypeNoFeedForward <: FeedForwardType
Trait for component output functions g
that have no feed forward behavior (do not depend on inputs):
g!(outs..., x, p, t)
See also PureFeedForward
, FeedForward
and PureStateMap
.
NetworkDynamics.PureStateMap
— TypePureStateMap <: FeedForwardType
Trait for component output functions g
that only depends on state:
g!(outs..., x)
See also PureFeedForward
, FeedForward
and NoFeedForward
.
Symbolic Indexing API
Network Parameter Object
NetworkDynamics.NWParameter
— TypeNWParameter(nw_or_nw_wraper, pflat)
Indexable wrapper for flat parameter array pflat
. Needs Network or wrapper of Network, e.g. ODEProblem
.
p = NWParameter(nw)
p.v[idx, :sym] # get parameter :sym of vertex idx
p.e[idx, :sym] # get parameter :sym of edge idx
p[s::Union{VPIndex, EPIndex}] # get parameter for specific index
Get flat array representation using pflat(p)
.
NetworkDynamics.NWParameter
— MethodNWParameter(nw_or_nw_wraper;
ptype=Vector{Float64}, pfill=filltype(ptype), default=true)
Creates "empty" NWParameter
object for the Network/Wrapper nw
with flat type ptype
. The array will be prefilled with pfill
(defaults to NaN).
If default=true
the default parameter values attached to the network components will be loaded.
NetworkDynamics.NWParameter
— MethodNWParameter(p::NWParameter; ptype=typeof(p.pflat))
Create NWParameter
based on other parameter object, just convert type.
NetworkDynamics.NWParameter
— MethodNWParameter(int::SciMLBase.DEIntegrator)
Create NWParameter
object from integrator
.
Network State Object
NetworkDynamics.NWState
— TypeNWState(nw_or_nw_wrapper, uflat, [pflat], [t])
Indexable wrapper for flat state & parameter array. Needs Network or wrapper of Network, e.g. ODEProblem
.
s = NWState(nw)
s.v[idx, :sym] # get state :sym of vertex idx
s.e[idx, :sym] # get state :sym of edge idx
s.p.v[idx, :sym] # get parameter :sym of vertex idx
s.p.e[idx, :sym] # get parameter :sym of edge idx
s[s::Union{VIndex, EIndex, EPIndex, VPIndex}] # get parameter for specific index
Get flat array representation using uflat(s)
and pflat(s)
.
NetworkDynamics.NWState
— MethodNWState(nw_or_nw_wrapper;
utype=Vector{Float64}, ufill=filltype(utype),
ptype=Vector{Float64}, pfill=filltype(ptype), default=true)
Creates "empty" NWState
object for the Network/Wrapper nw
with flat types utype
& ptype
. The arrays will be prefilled with ufill
and pfill
respectively (defaults to NaN).
If default=true
the default state & parameter values attached to the network components will be loaded.
NetworkDynamics.NWState
— MethodNWState(p::NWState; utype=typeof(uflat(s)), ptype=typeof(pflat(s)))
Create NWState
based on other state object, just convert types.
NetworkDynamics.NWState
— MethodNWState(p::NWParameter; utype=Vector{Float64}, ufill=filltype(utype), default=true)
Create NWState
based on existing NWParameter
object.
NetworkDynamics.NWState
— MethodNWState(int::SciMLBase.DEIntegrator)
Create NWState
object from integrator
.
NetworkDynamics.uflat
— Functionuflat(s::NWState)
Retrieve the wrapped flat array representation of the state.
NetworkDynamics.pflat
— Functionpflat(p::NWParameter)
pflat(s::NWState)
Retrieve the wrapped flat array representation of the parameters.
Symbolic Indices
NetworkDynamics.VIndex
— TypeVIndex{C,S} <: SymbolicStateIndex{C,S}
idx = VIndex(comp, sub)
A symbolic index for a vertex state variable.
comp
: the component index, either int or a collection of intssub
: the subindex, either int, symbol or a collection of those.
VIndex(1, :P) # vertex 1, variable :P
VIndex(1:5, 1) # first state of vertices 1 to 5
VIndex(7, (:x,:y)) # states :x and :y of vertex 7
Can be used to index into objects supporting the SymbolicIndexingInterface
, e.g. NWState
, NWParameter
or ODESolution
.
NetworkDynamics.EIndex
— TypeEIndex{C,S} <: SymbolicStateIndex{C,S}
idx = EIndex(comp, sub)
A symbolic index for an edge state variable.
comp
: the component index, either int or a collection of intssub
: the subindex, either int, symbol or a collection of those.
EIndex(1, :P) # edge 1, variable :P
EIndex(1:5, 1) # first state of edges 1 to 5
EIndex(7, (:x,:y)) # states :x and :y of edge 7
Can be used to index into objects supporting the SymbolicIndexingInterface
, e.g. NWState
, NWParameter
or ODESolution
.
NetworkDynamics.VPIndex
— TypeVPIndex{C,S} <: SymbolicStateIndex{C,S}
idx = VPIndex(comp, sub)
A symbolic index into the parameter a vertex:
comp
: the component index, either int or a collection of intssub
: the subindex, either int, symbol or a collection of those.
Can be used to index into objects supporting the SymbolicIndexingInterface
, e.g. NWParameter
or ODEProblem
.
NetworkDynamics.EPIndex
— TypeEPIndex{C,S} <: SymbolicStateIndex{C,S}
idx = VEIndex(comp, sub)
A symbolic index into the parameter a vertex:
comp
: the component index, either int or a collection of intssub
: the subindex, either int, symbol or a collection of those.
Can be used to index into objects supporting the SymbolicIndexingInterface
, e.g. NWParameter
or ODEProblem
.
Index generators
NetworkDynamics.vidxs
— Functionvidxs([inpr], components=:, variables=:) :: Vector{VIndex}
Generate vector of symbolic indexes for vertices.
inpr
: Only needed for name matching or:
access. Can be Network, sol, prob, ...components
: Number/Vector,:
,Symbol
(name matches),String
/Regex
(name contains)variables
: Symbol/Number/Vector,:
,String
/Regex
(all sym containing)
Examples:
vidxs(nw) # all vertex state indices
vidxs(1:2, :u) # [VIndex(1, :u), VIndex(2, :u)]
vidxs(nw, :, [:u, :v]) # [VIndex(i, :u), VIndex(i, :v) for i in 1:nv(nw)]
vidxs(nw, "ODEVertex", :) # all symbols of all vertices with name containing "ODEVertex"
NetworkDynamics.eidxs
— Functionvidxs([inpr], components=:, variables=:) :: Vector{EIndex}
Generate vector of symbolic indexes for edges.
inpr
: Only needed for name matching or:
access. Can be Network, sol, prob, ...components
: Number/Vector,:
,Symbol
(name matches),String
/Regex
(name contains)variables
: Symbol/Number/Vector,:
,String
/Regex
(all sym containing)
Examples:
eidxs(nw) # all edge state indices
eidxs(1:2, :u) # [EIndex(1, :u), EIndex(2, :u)]
eidxs(nw, :, [:u, :v]) # [EIndex(i, :u), EIndex(i, :v) for i in 1:ne(nw)]
eidxs(nw, "FlowEdge", :) # all symbols of all edges with name containing "FlowEdge"
NetworkDynamics.vpidxs
— Functionvpidxs([inpr], components=:, variables=:) :: Vector{VPIndex}
Generate vector of symbolic indexes for parameters. See vidxs
for more information.
NetworkDynamics.epidxs
— Functionepidxs([inpr], components=:, variables=:) :: Vector{EPIndex}
Generate vector of symbolic indexes for parameters. See eidxs
for more information.
Metadata API
Component Metadata API
NetworkDynamics.metadata
— FunctionNetworkDynamics.has_metadata
— Methodhas_metadata(c::ComponentModel, key::Symbol)
Checks if metadata key
is present for the component.
NetworkDynamics.get_metadata
— Methodget_metadata(c::ComponentModel, key::Symbol)
Retrieves the metadata key
for the component.
NetworkDynamics.set_metadata!
— Methodset_metadata!(c::ComponentModel, key::Symbol, value)
Sets the metadata key
for the component to value
.
NetworkDynamics.has_graphelement
— Functionhas_graphelement(c)
Checks if the edge or vetex function function has the graphelement
metadata.
NetworkDynamics.get_graphelement
— Functionget_graphelement(c::EdgeModel)::@NamedTuple{src::T, dst::T}
get_graphelement(c::VertexModel)::Int
Retrieves the graphelement
metadata for the component model. For edges this returns a named tupe (;src, dst)
where both are either integers (vertex index) or symbols (vertex name).
NetworkDynamics.set_graphelement!
— Functionset_graphelement!(c::EdgeModel, src, dst)
set_graphelement!(c::VertexModel, vidx)
Sets the graphelement
metadata for the edge model. For edges this takes two arguments src
and dst
which are either integer (vertex index) or symbol (vertex name). For vertices it takes a single integer vidx
.
Per-Symbol Metadata API
NetworkDynamics.symmetadata
— Functionsymmetadata(c::ComponentModel)::Dict{Symbol,Dict{Symbol,Any}}
Retrieve the metadata dictionary for the symbols. Keys are the names of the symbols as they appear in sym
, psym
, obssym
and insym
.
See also symmetadata
NetworkDynamics.get_metadata
— Methodget_metadata(c::ComponentModel, sym::Symbol, key::Symbol)
Retrievs the metadata key
for symbol sym
.
NetworkDynamics.has_metadata
— Methodhas_metadata(c::ComponentModel, sym::Symbol, key::Symbol)
Checks if symbol metadata key
is present for symbol sym
.
NetworkDynamics.set_metadata!
— Methodset_metadata!(c::ComponentModel, sym::Symbol, key::Symbol, value)
set_metadata!(c::ComponentModel, sym::Symbol, pair)
Sets the metadata key
for symbol sym
to value
.
NetworkDynamics.has_default
— Functionhas_default(c::ComponentModel, sym::Symbol)
Checks if a default
value is present for symbol sym
.
See also get_default
, set_default!
.
NetworkDynamics.get_default
— Functionget_default(c::ComponentModel, sym::Symbol)
Returns the default
value for symbol sym
.
See also has_default
, set_default!
.
NetworkDynamics.set_default!
— Functionset_default(c::ComponentModel, sym::Symbol, value)
Sets the default
value for symbol sym
to value
.
See also has_default
, get_default
.
NetworkDynamics.has_guess
— Functionhas_guess(c::ComponentModel, sym::Symbol)
Checks if a guess
value is present for symbol sym
.
See also get_guess
, set_guess!
.
NetworkDynamics.get_guess
— Functionget_guess(c::ComponentModel, sym::Symbol)
Returns the guess
value for symbol sym
.
See also has_guess
, set_guess!
.
NetworkDynamics.set_guess!
— Functionset_guess(c::ComponentModel, sym::Symbol, value)
Sets the guess
value for symbol sym
to value
.
NetworkDynamics.has_init
— Functionhas_init(c::ComponentModel, sym::Symbol)
Checks if a init
value is present for symbol sym
.
NetworkDynamics.get_init
— Functionget_init(c::ComponentModel, sym::Symbol)
Returns the init
value for symbol sym
.
NetworkDynamics.set_init!
— Functionset_init(c::ComponentModel, sym::Symbol, value)
Sets the init
value for symbol sym
to value
.
NetworkDynamics.has_bounds
— Functionhas_bounds(c::ComponentModel, sym::Symbol)
Checks if a bounds
value is present for symbol sym
.
See also get_bounds
, set_bounds!
.
NetworkDynamics.get_bounds
— Functionget_bounds(c::ComponentModel, sym::Symbol)
Returns the bounds
value for symbol sym
.
See also has_bounds
, set_bounds!
.
NetworkDynamics.set_bounds!
— Functionset_bounds(c::ComponentModel, sym::Symbol, value)
Sets the bounds
value for symbol sym
to value
.
See also has_bounds
, get_bounds
.
Initialization
NetworkDynamics.find_fixpoint
— Functionfind_fixpoint(nw::Network, [x0::NWState=NWState(nw)], [p::NWParameter=x0.p]; kwargs...)
find_fixpoint(nw::Network, x0::AbstractVector, p::AbstractVector; kwargs...)
find_fixpoint(nw::Network, x0::AbstractVector; kwargs...)
Convenience wrapper around SteadyStateProblem
from SciML-ecosystem. Constructs and solves the steady state problem, returns found value wrapped as NWState
.
NetworkDynamics.initialize_component!
— Functioninitialize_component!(cf::ComponentModel; verbose=true, kwargs...)
Initialize a ComponentModel
by solving the corresponding NonlinearLeastSquaresProblem
. During initialization, everyting which has a default
value (see Metadata) is considered "fixed". All other variables are considered "free" and are solved for. The initial guess for each variable depends on the guess
value in the Metadata.
The result is stored in the ComponentModel
itself. The values of the free variables are stored in the metadata field init
.
The kwargs
are passed to the nonlinear solver.
NetworkDynamics.init_residual
— Functioninit_residual(cf::T; t=NaN, recalc=false)
Calculates the residual |du| for the given component model for the values provided via default
and init
Metadata.
If recalc=false just return the residual determined in the actual initialization process.
See also initialize_component!
.
Execution Types
NetworkDynamics.ExecutionStyle
— Typeabstract type ExecutionStyle{buffered::Bool} end
Abstract type for execution style. The coreloop dispatches based on the Execution style stored in the network object.
buffered=true
means that the edge input es explicitly gathered, i.e. the vertex outputs in the output buffer will be copied into a dedicated input buffer for the edges.buffered=false
means, that the edge inputs are not explicitly gathered, but the corloop will perform a redirected lookup into the output buffer.
NetworkDynamics.SequentialExecution
— Typestruct SequentialExecution{buffered::Bool}
Sequential execution, no parallelism. For buffered
see ExecutionStyle
.
NetworkDynamics.PolyesterExecution
— Typestruct PolyesterExecution{buffered}
Parallel execution using Polyester.jl
. For buffered
see ExecutionStyle
.
NetworkDynamics.ThreadedExecution
— Typestruct ThreadedExecution{buffered}
Parallel execution using Julia threads. For buffered
see ExecutionStyle
.
NetworkDynamics.KAExecution
— Typestruct KAExecution{buffered}
Parallel execution using KernelAbstractions.jl
. Works with GPU and CPU arrays. For buffered
see ExecutionStyle
.
Aggregators
NetworkDynamics.Aggregator
— Typeabstract type Aggregator end
Abstract sypertype for aggregators. Aggregators operate on the output buffer of all components and fill the aggregation buffer with the aggregatated edge values per vertex.
All aggregators have the constructor
Aggegator(aggfun)
for example
SequentialAggreator(+)
NetworkDynamics.SequentialAggregator
— TypeSequentialAggregator(aggfun)
Sequential aggregation.
NetworkDynamics.SparseAggregator
— TypeSparseAggregator(+)
Only works with additive aggregation +
. Aggregates via sparse inplace matrix multiplication. Works with GPU Arrays.
NetworkDynamics.ThreadedAggregator
— TypeThreadedAggregator(aggfun)
Parallel aggregation using Julia threads.
NetworkDynamics.PolyesterAggregator
— TypePolyesterAggregator(aggfun)
Parallel aggregation using Polyester.jl
.
NetworkDynamics.KAAggregator
— TypeKAAggregator(aggfun)
Parallel aggregation using KernelAbstractions.jl
. Works with both GPU and CPU arrays.
Utils
NetworkDynamics.save_parameters!
— Functionsave_parameters!(integrator::SciMLBase.DEIntegrator)
Save the current parameter values in the integrator. Call this function inside callbacks if the parameter values have changed. This will store a timeseries of said parameters in the solution object, thus alowing us to recosntruct observables which depend on time-dependet variables.
NetworkDynamics.ff_to_constraint
— Functionff_to_constraint(v::VertexModel)
Takes VertexModel
v
with feed forward and turns all algebraic output states into internal states by defining algebraic constraints contraints 0 = out - g(...)
. The new output function is just a StateMask
into the extended internal state vector.
Returns the transformed VertexModel
.
Base.copy
— Methodcopy(c::NetworkDynamics.ComponentModel)
Shallow copy of the component model. Creates a deepcopy of metadata
and symmetadata
but references the same objects everywhere else.