API

The following functions are designed for public use.

Network Construction API

NetworkDynamics.NetworkType
Network([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 defined graphelement. See vidx and src/dst keywors for VertexModel and EdgeModel constructors respectively.

  • vertexm: A single VertexModel or a vector of VertexModel objects. The order of the vertex models must mirror the order of the vertices(g) iterator.

  • edgem: A single EdgeModel or a vector of EdgeModel objects. The order of the edge models must mirror the order of the edges(g) iterator.

Optional keyword arguments:

  • execution=SequentialExecution{true}(): Execution model of the network. E.g. SequentialExecution, KAExecution, PolyesterExecution or ThreadedExecution.
  • aggregator=execution isa SequentialExecution ? SequentialAggregator(+) : PolyesterAggregator(+): Aggregation function applied to the edge models. E.g. SequentialAggregator, PolyesterAggregator, ThreadedAggregator, SparseAggregator.
  • check_graphelement=true: Check if the graphelement 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.
source
Network(nw::Network; g, vertexm, edgem, kwargs...)

Rebuild the Network with same graph and vertex/edge models but possibly different kwargs.

source
NetworkDynamics.dimMethod
dim(nw::Network)

Returns the number of dynamic states in the network, corresponts to the length of the flat state vector.

source
NetworkDynamics.pdimMethod
pdim(nw::Network)

Returns the number of parameters in the network, corresponts to the length of the flat parameter vector.

source

Component Models

NetworkDynamics.VertexModelMethod
VertexModel(; kwargs...)

Build a VertexModel according to the keyword arguments.

Main Arguments:

  • f=nothing: Dynamic function of the component. Can be nothing if dim is 0.
  • g: Output function of the component. Usefull helpers: StateMask
  • sym/dim: Symbolic names of the states. If dim is provided, sym is set automaticially.
  • outsym/outdim: Symbolic names of the outputs. If outdim is provided, outsym is set automaticially. Can be infered automaticially if g isa StateMask.
  • psym/pdim=0: Symbolic names of the parameters. Ifpdimis provided,psym` is set automaticially.
  • mass_matrix=I: Mass matrix of component. Can be a vector v and is then interpreted as Diagonal(v).
  • name=dim>0 ? :VertexM : :StaticVertexM: Name of the component.

Optional Arguments:

  • insym/indim: Symbolic names of the inputs. If indim 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 from g automaticially.
  • obssym/obsf: Define additional "observable" states.
  • symmetadata/metadata: Provide prefilled metadata dictionaries.
  • extin=nothing: Define "external" inputs for the model with Network indices, i.e. extin=[VIndex(7,:x), ..]. Those inputs will be provided as another input vector f(x, in, extin, p, t) and g(y, x, in, extin, p, t).

All Symbol arguments can be used to set default values, i.e. psym=[:K=>1, :p].

source
NetworkDynamics.EdgeModelMethod
EdgeModel(; kwargs...)

Build a EdgeModel according to the keyword arguments.

Main Arguments:

  • f=nothing: Dynamic function of the component. Can be nothing if dim is 0.
  • g: Output function of the component. Usefull helpers: AntiSymmetric, Symmetric, Fiducial, Directed and StateMask.
  • sym/dim: Symbolic names of the states. If dim is provided, sym is set automaticially.
  • outsym/outdim: Symbolic names of the outputs. If outdim is provided, outsym is set automaticially. In general, outsym for edges isa named tuple (; src, dst). However, depending on the g function, it might be enough to provide a single vector or even nothing (e.g. AntiSymmetric(StateMask(1:2))). See Building EdgeModels for examples.
  • psym/pdim=0: Symbolic names of the parameters. Ifpdimis provided,psym` is set automaticially.
  • mass_matrix=I: Mass matrix of component. Can be a vector v and is then interpreted as Diagonal(v).
  • name=dim>0 ? :EdgeM : :StaticEdgeM: Name of the component.

Optional Arguments:

  • insym/indim: Symbolic names of the inputs. If indim 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 from g automaticially.
  • obssym/obsf: Define additional "observable" states.
  • symmetadata/metadata: Provide prefilled metadata dictionaries.
  • extin=nothing: Define "external" inputs for the model with Network indices, i.e. extin=[VIndex(7,:x), ..]. Those inputs will be provided as another input vector f(x, insrc, indst, extin, p, t) and g(ysrc, ydst, x, insrc, indst, extin, p, t).

All Symbol arguments can be used to set default values, i.e. psym=[:K=>1, :p].

source

Component Models with MTK

NetworkDynamics.VertexModelMethod
VertexModel(sys::ODESystem, inputs, outputs;
            verbose=false, name=getname(sys), extin=nothing, 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 states
  • outputs: names of variables in you equation representing the node output

Additional kw arguments:

  • name: Set name of the component model. Will be lifted from the ODESystem name.
  • extin=nothing: Provide external inputs as pairs, i.e. extin=[:extvar => VIndex(1, :a)] will bound the variable extvar(t) in the equations to the state a of the first vertex.
  • ff_to_constraint=true: Controlls, whether output transformations g which depend on inputs should be transformed into constraints. Defaults to true since ND.jl does not handle vertices with FF yet.
source
NetworkDynamics.EdgeModelMethod
EdgeModel(sys::ODESystem, srcin, srcout, dstin, dstout;
          verbose=false, name=getname(sys), extin=nothing, 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 source
  • dstin: names of variables in you equation representing the node state at the destination
  • srcout: names of variables in you equation representing the output at the source
  • dstout: names of variables in you equation representing the output at the destination

Additional kw arguments:

  • name: Set name of the component model. Will be lifted from the ODESystem name.
  • extin=nothing: Provide external inputs as pairs, i.e. extin=[:extvar => VIndex(1, :a)] will bound the variable extvar(t) in the equations to the state a of the first vertex.
  • ff_to_constraint=false: Controlls, whether output transformations g which depend on inputs should be transformed into constraints.
source
NetworkDynamics.EdgeModelMethod
EdgeModel(sys::ODESystem, srcin, dstin, AntiSymmetric(dstout); kwargs...)

Create a EdgeModel object from a given ODESystem created with ModelingToolkit for single sided models.

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), or
  • Directed(dstout).

Additional kwargs are the same as for the double-sided EdgeModel MTK constructor.

source

Output Function Helpers/Wrappers

NetworkDynamics.StateMaskType
StateMask(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, StateMasks 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.

source

Accessors for Component Properties

NetworkDynamics.outdimFunction
outdim(c::VertexModel)::Int
outdim(c::EdgeModel)::@NamedTuple(src::Int, dst::Int)

Retrieve the output dimension of the component

source
NetworkDynamics.outsymFunction

outsym(c::VertexModel)::Vector{Symbol} outsym(c::EdgeModel)::@NamedTuple{src::Vector{Symbol}, dst::Vector{Symbol}}

Retrieve the output symbols of the component.

source
NetworkDynamics.insymFunction
insym(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.

source

FeedForwardType-Traits

Symbolic Indexing API

Network Parameter Object

NetworkDynamics.NWParameterType
NWParameter(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).

source
NetworkDynamics.NWParameterMethod
NWParameter(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.

source

Network State Object

NetworkDynamics.NWStateType
NWState(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).

source
NetworkDynamics.NWStateMethod
NWState(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.

source
NetworkDynamics.NWStateMethod
NWState(p::NWState; utype=typeof(uflat(s)), ptype=typeof(pflat(s)))

Create NWState based on other state object, just convert types.

source
NetworkDynamics.NWStateMethod
NWState(p::NWParameter; utype=Vector{Float64}, ufill=filltype(utype), default=true)

Create NWState based on existing NWParameter object.

source
NetworkDynamics.pflatFunction
pflat(p::NWParameter)
pflat(s::NWState)

Retrieve the wrapped flat array representation of the parameters.

source

Symbolic Indices

NetworkDynamics.VIndexType
VIndex{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 ints
  • sub: 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
VIndex(2)          # references the second vertex model

Can be used to index into objects supporting the SymbolicIndexingInterface, e.g. NWState, NWParameter or ODESolution.

See also: EIndex, VPIndex, EPIndex

source
NetworkDynamics.EIndexType
EIndex{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 ints
  • sub: 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
EIndex(2)          # references the second edge model

Can be used to index into objects supporting the SymbolicIndexingInterface, e.g. NWState, NWParameter or ODESolution.

See also: VIndex, VPIndex, EPIndex

source
NetworkDynamics.VPIndexType
VPIndex{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 ints
  • sub: 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.

See also: EPIndex, VIndex, EIndex

source
NetworkDynamics.EPIndexType
EPIndex{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 ints
  • sub: 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.

See also: VPIndex, VIndex, EIndex

source
NetworkDynamics.@obsexMacro
@obsex([name =] expression)

Define observable expressions, which are simple combinations of knonw states/parameters/observables. @obsex(...) returns an ObservableExpression which can be used as an symbolic index. This is mainly intended for quick plotting or export of common "derived" variables, such as the argument of a 2-component complex state. For example:

sol(t; idxs=@obsex(arg = atan(VIndex(1,:u_i), VIndex(1,:u_r))]
sol(t; idxs=@obsex(δrel = VIndex(1,:δ) - VIndex(2,:δ)))
source

Index generators

NetworkDynamics.vidxsFunction
vidxs([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"
source
NetworkDynamics.eidxsFunction
vidxs([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"
source
NetworkDynamics.vpidxsFunction
vpidxs([inpr], components=:, variables=:) :: Vector{VPIndex}

Generate vector of symbolic indexes for parameters. See vidxs for more information.

source
NetworkDynamics.epidxsFunction
epidxs([inpr], components=:, variables=:) :: Vector{EPIndex}

Generate vector of symbolic indexes for parameters. See eidxs for more information.

source

Metadata API

Component Metadata API

NetworkDynamics.get_graphelementFunction
get_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).

source
NetworkDynamics.set_graphelement!Function
set_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.

source

Per-Symbol Metadata API

NetworkDynamics.set_metadata!Method
set_metadata!(c::ComponentModel, sym::Symbol, key::Symbol, value)
set_metadata!(c::ComponentModel, sym::Symbol, pair)

Sets the metadata key for symbol sym to value.

source

Callbacks API

Define Callbacks

NetworkDynamics.ComponentCallbackType
abstract type ComponentCallback end

Abstract type for a component based callback. A component callback bundles a ComponentCondition as well as a ComponentAffect which can be then tied to a component model using add_callback! or set_callback!.

On a Network level, you can automaticially create network wide CallbackSets using get_callbacks.

See ContinousComponentCallback and VectorContinousComponentCallback for concrete implemenations of this abstract type.

source
NetworkDynamics.VectorContinousComponentCallbackType
VectorContinousComponentCallback(condition, affect, len; kwargs...)

Connect a ComponentCondition and a [ComponentAffect)[@ref] to a continous callback which can be attached to a component model using add_callback! or set_callback!. This vector version allows for condions which have len output dimensions. The affect will be triggered with the additional event_idx argument to know in which dimension the zerocrossing was detected.

The kwargs will be forwarded to the VectorContinuousCallback when the component based callbacks are collected for the whole network using get_callbacks(::Network). DiffEq.jl docs for available options.

source
NetworkDynamics.DiscreteComponentCallbackType
DiscreteComponentCallback(condition, affect; kwargs...)

Connect a ComponentCondition and a [ComponentAffect)[@ref] to a discrete callback which can be attached to a component model using add_callback! or set_callback!.

Note that the condition function returns a boolean value, as the discrete callback perform no rootfinding.

The kwargs will be forwarded to the DiscreteCallback when the component based callbacks are collected for the whole network using get_callbacks(::Network). DiffEq.jl docs for available options.

source
NetworkDynamics.PresetTimeComponentCallbackType
PresetTimeComponentCallback(ts, affect; kwargs...)

Tirgger a ComponentAffect at given timesteps ts in discrete callback, which can be attached to a component model using add_callback! or set_callback!.

The kwargs will be forwarded to the PresetTimeCallback when the component based callbacks are collected for the whole network using get_callbacks(::Network).

The PresetTimeCallback will take care of adding the timesteps to the solver, ensuring to exactly trigger at the correct times.

source
NetworkDynamics.ComponentConditionType
ComponentCondition(f::Function, sym, psym)

Creates a callback condition for a [ComponentCallback].

  • f: The condition function. Must be a function of the form out=f(u, p, t) when used for ContinousComponentCallback or DiscreteComponentCallback and f!(out, u, p, t) when used for VectorContinousComponentCallback.
    • Arguments of f
      • u: The current value of the selecte sym states, provided as a SymbolicView object.
      • p: The current value of the selected psym parameters.
      • t: The current simulation time.
  • sym: A vector or tuple of symbols, which represent states (including inputs, outputs, observed) of the component model. Determines, which states will be available thorugh parameter u in the callback condition function f.
  • psym: A vector or tuple of symbols, which represetn parameters of the component mode. Determines, which parameters will be available in the condition function f

Example

Consider a component model with states [:u1, :u2], inputs [:i], outputs [:o] and parameters [:p1, :p2].

ComponentCondition([:u1, :o], [:p1]) do (u, p, t)
    # access states symbolicially or via int index
    u[:u1] == u[1]
    u[:o] == u[2]
    p[:p1] == p[1]
    # the states/prameters `:u2`, `:i` and `:p2` are not available as
    # they are not listed in the `sym` and `psym` arguments.
end
source
NetworkDynamics.ComponentAffectType
ComponentAffect(f::Function, sym, psym)

Creates a callback condition for a [ComponentCallback].

  • f: The affect function. Must be a function of the form f(u, p, [event_idx], ctx) where event_idx is only available in VectorContinousComponentCallback.
    • Arguments of f
      • u: The current (mutable) value of the selected sym states, provided as a SymbolicView object.
      • p: The current (mutalbe) value of the selected psym parameters.
      • event_idx: The current event index, i.e. which out element triggerd in case of VectorContinousComponentCallback.
      • ctx::NamedTuple a named tuple with context variables.
        • ctx.model: a referenc to the ocmponent model
        • ctx.vidx/ctx.eidx: The index of the vertex/edge model.
        • ctx.src/ctx.dst: src and dst indices (only for edge models).
        • ctx.integrator: The integrator object.
        • ctx.t=ctx.integrator.t: The current simulation time.
  • sym: A vector or tuple of symbols, which represent states (excluding inputs, outputs, observed) of the component model. Determines, which states will be available thorugh parameter u in the callback condition function f.
  • psym: A vector or tuple of symbols, which represetn parameters of the component mode. Determines, which parameters will be available in the condition function f

Example

Consider a component model with states [:u1, :u2], inputs [:i], outputs [:o] and parameters [:p1, :p2].

ComponentAffect([:u1, :o], [:p1]) do (u, p, ctx)
    u[:u1] = 0 # change the state
    p[:p1] = 1 # change the parameter
    @info "Changed :u1 and :p1 on vertex $(ctx.vidx)" # access context
end
source
NetworkDynamics.SymbolicViewType
SymbolicView{N,VT} <: AbstractVetor{VT}

Is a (smallish) fixed size vector type with named dimensions. Its main purpose is to allow named acces to variables in ComponentCondition and ComponentAffect functions.

I.e. when the ComponentAffect declared sym=[:x, :y], you can acces u[:x] and u[:y] inside the condition function.

source
NetworkDynamics.get_callbacksMethod
get_callbacks(nw::Network)::CallbackSet

Returns a CallbackSet composed of all the "component-based" callbacks in the metadata of the Network components.

source

Attach Callbacks to Edge/VertexModels

NetworkDynamics.get_callbacksMethod
get_callback(c::ComponentModel)

Gets all callback functions for the component. Wraps in tuple, even if there is only a single one.

source

Initialization

NetworkDynamics.find_fixpointFunction
find_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.

source
NetworkDynamics.initialize_component!Function
initialize_component!(cf::ComponentModel; verbose=true, apply_bound_transformation=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.

Bounds of free variables

When encountering any bounds in the free variables, NetworkDynamics will try to conserve them by applying a coordinate transforamtion. This behavior can be supressed by setting apply_bound_transformation. The following transformations are used:

  • (a, b) intervals where both a and b are positive are transformed to u^2/sqrt(u)
  • (a, b) intervals where both a and b are negative are transformed to -u^2/sqrt(-u)
source
NetworkDynamics.init_residualFunction
init_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!.

source
NetworkDynamics.get_initial_stateFunction
get_initial_state(c::ComponentModel, syms; missing_val=nothing)

Returns the initial state for symbol sym (single symbol of vector) of the component model c. Returns missing_val if the symbol is not initialized. Also works for observed symbols.

See also: dump_initial_state.

source
NetworkDynamics.set_defaults!Function
set_defaults!(nw::Network, s::NWState)

Set the default values of the network to the values of the given state. Can be used to "store" the found fixpoint in the network metadata.

source

Execution Types

NetworkDynamics.ExecutionStyleType
abstract 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.
source

Aggregators

NetworkDynamics.AggregatorType
abstract 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(+)
source

Utils

NetworkDynamics.save_parameters!Function
save_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.

source
NetworkDynamics.ff_to_constraintFunction
ff_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.

source
Base.copyMethod
copy(c::NetworkDynamics.ComponentModel)

Shallow copy of the component model. Creates a deepcopy of metadata and symmetadata but references the same objects everywhere else.

source