Aggregator

The Aggregator class defined below is the beating heart of STREAM.

It defines a methodology where several (separate) calculations can be simulated in a coupled manner. Such coupling is defined through a graph, and is simulated through several backends (see Backends). The main scheme is neatly presented in Documentation.

Aggregator.compute

The main function of the Aggregator - bridging between the solver and the calculations.

Aggregator.load

Given a description of the system state either at one time or at many, returns a vector or a solution that fits this information.

Aggregator.save

Given either a vector solution of the system or a Solution object, creates a human-readable state description of the solution.

Aggregator.var_index

Return the Place at which a given variable lies

Aggregator.at_times

Given a transient solution and a variable, returns the variable at the calculated times.

Aggregator.solve

For a Differential Algebraic set of eqs.

Aggregator.solve_steady

Solving an Algebraic Equation \(0=F(y)\) using algebraic

Example

Consider the attached image. This graph describes the dependencies of each calculation: A, B, C. Let’s consider the Calculation A. A owns a subset of the entire vector \(\vec{y}_A\in\vec{y}\), which it receives as input. However, A depends on values owned by other Calculations: \(\left(\vec{y}_i, \vec{y}_j\right)\), namely B and C. Thus, B (e.g.) must be able to show the Aggregator the place in which it stores \(\vec{y}_i\). In this manner, there exists no communication between Calculations, and B remains undisturbed regardless of other coupling to \(\vec{y}_i\).

The above analysis is done at construction time. At evaluation time, the Aggregator divides the vector according to the references it was given, and collects the results.

Required variables which are not present on the graph must be provided as time-only functions or constants.

This can easily become a rather complex system:

class Aggregator(graph, funcs=None)[source]

Bases: object

Collects and calls calculations in order to solve a large coupled generalized equation (presented in the calculation module). This object connects the Solver to the Calculations. It receives input from the Solver, which is then distributed to the different Calculations. The Calculations then compute their allotted piece of the functional. The Aggregator then passes the results back to the Solver.

sections

A mapping of each Calculation to its given slice of \(\vec{y}\).

Type:

dict[Calculation, slice]

mass

The gathered mass_vector.

Type:

Array1D

external

For each calculation, maps the variable names which are passed from other calculations, and for each variable name maps the calculations which pass such a variable and the places in \(\vec{y}\) from which they are to be taken

Type:

dict[Calculation, dict[str, dict[Calculation, Place]]]

See also

Calculation

Create an instance of Aggregator. The following are required:

Parameters:
  • graph (DiGraph) – containing Calculations at nodes and variable coupling at the edges.

  • funcs (ExternalFunctions | None) – time-only-dependent functions, which are user controlled.

classmethod connect(a, b, *edges)[source]

Connect two Aggregator objects. In case of a clash, the second object prevails. If edges contains an edge already in either a.graph or b.graph, it is updated, not overridden.

Tip

The two inputs may share nodes. This is very useful!

Parameters:
Returns:

new – A new Aggregator whose graph and functions are composed out of a,b.

Return type:

Aggregator

classmethod from_CalculationGraph(a)[source]

Creates an Aggregator using a CalculationGraph object.

Parameters:

a (CalculationGraph) – Graph to use for the creation.

Return type:

Aggregator

classmethod from_decoupled(*nodes, funcs=None)[source]

Instantiate an Aggregator from calculations, which are not connected.

Parameters:
  • nodes (Calculation) – Calculations to ba added

  • funcs (ExternalFunctions or None) – to be passed to the aggregator

Returns:

agr – Which contains a disconnected-nodes-graph

Return type:

Aggregator

at_times(solution, node, var_name)[source]

Given a transient solution and a variable, returns the variable at the calculated times.

Parameters:
  • solution (Array2D) – The solution matrix, i.e. the variable vector (columns) at different times (rows)

  • node (Calculation) – The inquired calculation

  • var_name (str) – the inquired variable of node.

Returns:

output – A slice of the solution at the correct var_index.

Return type:

Array2D

compute(y, t=0)[source]

The main function of the Aggregator - bridging between the solver and the calculations.

The graph attribute of Aggregator contains a set of Calculation objects, each representing different equations which come together to form the full Functional \(\vec{F}(\vec{y},t)\), returned by this method. Therefore, this method orchestrates the inputs and outputs to each calculation and aggregates the result.

Parameters:
  • y (Sequence[float]) – guess/ result from solver

  • t (Second) – time

Returns:

Functional – the differential-algebraic functional f(y,t) which is provided in parts from the different calculations.

Return type:

Array1D

draw(node_options=None, edge_options=None)[source]

Method equivalent of draw_aggregator

load(s: dict[str, dict[str, float | ndarray]]) ndarray[source]
load(s: dict[float, State]) Solution

Given a description of the system state either at one time or at many, returns a vector or a solution that fits this information.

Parameters:

s – The system description to parse.

save(solution: Solution) dict[float, State][source]
save(solution: Sequence[float], t: float | ndarray = 0, strict: bool = False) State

Given either a vector solution of the system or a Solution object, creates a human-readable state description of the solution.

Parameters:
  • solution – The solution to make human-readable.

  • t – The time at which the solution is given (if it is a single vector)

  • strict – Whether information beyond the vector state variables should be added.

solve(y0, time, yp0=None, eq_type=None, *, progressbar=False, **options)[source]

For a Differential Algebraic set of eqs. (DAE), the chosen solver is IDA from the LLNL SUNDIALS suite, which is kindly wrapped by Scikits.Odes, originally written in C with DASPK (Fortran) usages. This solver performs (among many other capabilities) integration by variable-order, variable-coefficient BDF. Newton iteration is used to find a solution.

Parameters:
  • y0 (Array1D or DictState) – Initial values or guess. Can either be an array or a State, in the latter case load will be used to obtain the desired array.

  • time (Sequence[float]) – Return results at these time points.

  • yp0 (Array1D or None) – Initial derivatives. It helps if they’re known (in the DAE case), but by default the consistent yp0 is found from y0.

  • eq_type ('ODE', 'DAE', 'ALG' or None) – A solver may be chosen deliberately from [ODE, DAE, ALG]. If None, the method is set by looking at the mass matrix and whether time is none.

  • progressbar (ProgressBarLike or bool) – Whether to use a progressbar, and if so, which one. If True, use use progressbar.ProgressBar

  • options – Other options

Returns:

solution – Calculated vector at requested times: [time, variable].

Return type:

Solution

References

Scikits.Odes documentation

solve_steady(guess, **options)[source]

Solving an Algebraic Equation \(0=F(y)\) using algebraic

Parameters:
  • guess (Array1D or DictState) – Initial guess. Can either be an array or a State, in the latter case load will be used to obtain the desired array.

  • options – Solver options

Returns:

solution – Calculated vector.

Return type:

Array1D

var_index(node, var_name)[source]

Return the Place at which a given variable lies

Parameters:
  • node (Calculation) – The calculation whose variable is requested

  • var_name (str) – Variable name

Returns:

The Place in aggregator vector where this variable resides

Return type:

Place

Solvers

The Aggregator defines a functional \(\vec{F}(\vec{y},t)\) which can be utilized to solve several problems. Below are several backends used for different classes of problems.

Backends

Each eq_type in solve relies on a different backend solver, for which different options are relevant. For the sake of brevity, these options are not specified here, and the user is encouraged to look at each backend, as well as the defaults set in aggregator to achieve a higher level of proficiency.

eq_type

Solver

Function Name

ODE

Scipy.integrate.solve_ivp

differential

DAE

Scikits.odes.ida

differential_algebraic

ALG

Scipy.optimize.root

algebraic

exception TransientRuntimeError(t, y, ydot, message, *args)[source]

RuntimeError which occurred during a transient simulation, mostly because of solver convergence problems.

Data for debugging the error may be extracted by:

try:
    # error raising simulation
except TransientRuntimeError as e:
    t, y, ydot = e.t, e.y, e.ydot
    results = Solution(t, y)
    raise e
Parameters:

message (str)

algebraic(F, y0, time=None, R=None, **options)[source]

Solving an Algebraic Equation \(0=F(y, t)\)

Parameters:
  • F (Functional) – The main right-hand side function \(F(y,t)\).

  • y0 (Array1D) – Initial Guess.

  • time (Sequence[float] | None) – If None, the root of the functional is found. Else, the root is found at the specified time points, given sequential initial guesses. It is a quasi-static simulation, if one wills it. The first vector is then the initial guess.

  • R (Functional | None) – A function controlling transient simulation stop events.

  • options – Other options to be passed to the Scipy.optimize.root solver.

Returns:

solution – The solution matrix at requested times: [time, variable].

Return type:

Array

differential(F, y0, time, **options)[source]

Solving an Ordinary Differential Equation (ODE) \(\dot{y}=F(y, t)\)

Parameters:
  • F (Functional) – The main right-hand side function \(F(y,t)\).

  • y0 (Array1D) – Initial values.

  • time (Sequence[float]) – Time points for which the simulation values should be returned.

  • options – Other options to be passed to the scipy.integrate.solve_ivp solver.

Returns:

solution – The solution matrix at requested times: [time, variable].

Return type:

Array2D

differential_algebraic(F, mass, y0, time, yp0=None, R=None, continuous=False, **options)[source]

Solving a Differential Algebraic Equation (DAE) \(M\dot{y}=F(y, t)\)

Parameters:
  • F (Functional) – The main right-hand side function \(F(y,t)\).

  • mass (Array1D) – Mass vector, determining which index is differential and which is algebraic.

  • y0 (Array1D) – Initial values.

  • time (Sequence[float]) – Time points for which the simulation values should be returned.

  • yp0 (Array1D | None) – Initial value derivatives. If None, considered 0. However, when compute_initcond = 'yp0' is set (as is default), this vector is deduced from y0.

  • R (Functional | None) – A function controlling simulation stop events. Simulation continues unless a False valued index is returned.

  • continuous (bool) – If True, simulation continues after R has yielded a non-True value, whereas the simulation stops in that case for False. An added behavior is that using continuous=True forces a restart of the simulation, such that initial steps are much smaller, and controlled by first_step_size

  • options – Other options to be passed to the scikits.odes solver.

Returns:

solution – The solution matrix: [time, variable], and the vector of times in which it was calculated.

Return type:

tuple[Array2D, Array1D]

Jacobians

Aggregator solvers can greatly benefit from calculated Jacobians. The following are some implementations for the DAE (Differential-Algebraic Equations) form and the ALG (Algebraic) form.

The main idea behind these implementations is that the graphical structure of an Aggregator object can be utilized to deduce the Jacobian sparsity, at least across Calculation boundaries.

DAE_jacobian(agr, step_strategy=<function _default_step_strategy>)[source]

A function for creating a one-sided approximation of the Jacobian of an Aggregator functional. The Jacobian is defined as follows:

\[J_{ij} = \frac{dF_i}{dy_j} \approx \frac{F_i(y + h_j) - F_i(y)}{h_j}\]

Where \(h_j\) is the jth step size and \(y + h_j\) simply means the jth component is \(y_j + h_j\). This is however a bit simplistic, since the actual required Jacobian is of the full functional:

\[G(t, y, \dot{y}) = F(y, t) - M\dot{y}\]

Meaning the Jacobian is actually

\[\frac{dG}{dy}=\frac{\partial G}{\partial\dot{y}} \frac{\partial\dot{y}}{\partial y} + \frac{\partial G}{\partial y}\]

Since \(\dot{y}\) is approximated linearly, \(\partial_y\dot{y}\equiv\sigma\) does not depend on the functions considered above, only on the approximation scheme and time step.

In our case, \(\partial_{\dot{y}_j} G_i=-M_{ij} = -m_i \delta_{i,j}\).

Parameters:
  • agr (Aggregator) – The Aggregator whose functional is desired

  • step_strategy (Callable or None) – A function for the step vector \(h = f(y, \dot{y})\)

Returns:

J_func – A function whose signature is \(J(t, y, \dot{y}, G(t, y, \dot{y}), \sigma, \text{result})\), matching the required SUNDIALS format.

Return type:

Callable

State

The solution of a steady state problem is represented as a State object, which is mostly a fancy dict-type.

The concept of a human-readable state of the system.

Basically a 2-deep nested dict with str keys and the inner values are Value objects. We provide additional useful methods in the same namespace.

class State[source]

A nested dictionary with str keys that connects calculations to a dictionary of the values of their variables.

classmethod from_dataframe(df)[source]

Creates a State from a DataFrame representation.

Parameters:

df (DataFrame) – Data to turn into a State

Return type:

State

classmethod from_liststate(s)[source]

Make a new State from its serializable list version.

Parameters:

s (dict[str, dict[str, float | list]])

Return type:

State

classmethod load(f)[source]

Load a State from a YAML IO.

Parameters:

f (IO) – Stream to read from

Return type:

State

classmethod merge(*st)[source]

Merge states together. Later states have precedence.

Parameters:

st (dict[str, dict[str, float | ndarray]])

Return type:

State

classmethod uniform(calculations, value, *variables)[source]

Create a state, or a partial state for given calculations and given target variables, which are assigned a uniform value.

Parameters:
  • calculations (Iterable[Calculation]) – Desired calculations to assign a partial state.

  • value (Value) – Value to assign to the variable states.

  • variables (str) – Variables to have states. If none are given, all variables in calculations are included.

Returns:

state – A partially legal state

Return type:

State

dump(f=None)[source]

Dump this State as YML.

Parameters:

f (Stream to dump to.)

Return type:

str | None

filter_calculations(f)[source]

Filter out state calculations predicated upon a function.

Parameters:

f (Callable) – Predicate function. False values are filtered out.

Returns:

state’ – A filtered partial state

Return type:

State

Examples

>>> State(A=dict(b=1), a=dict(B=2)).filter_calculations(str.isupper)
{'A': {'b': 1}}
filter_values(f)[source]

Filter out state variables predicated upon a function.

Parameters:

f (Callable) – Predicate function. False values are filtered out.

Returns:

state’ – A filtered partial state

Return type:

State

Examples

>>> State(A=dict(b=1), a=dict(B=2)).filter_values(lambda x: x == 2)
{'a': {'B': 2}}
filter_var_names(f)[source]

Filter out state variable names predicated upon a function.

Parameters:

f (Callable) – Predicate function. False values are filtered out.

Returns:

state’ – A filtered partial state

Return type:

State

Examples

>>> State(A=dict(b=1), a=dict(B=2)).filter_var_names(str.isupper)
{'a': {'B': 2}}
listify()[source]

Replace all arrays in self with serializable lists

Return type:

dict[str, dict[str, float | list]]

records()[source]

Generate records from this state.

Return type:

Iterable[dict[str, Any]]

to_dataframe()[source]

Represents a State as a DataFrame

Return type:

DataFrame

from_dataframe(df)[source]

Creates a State or a dictionary of {time: State} from a DataFrame representation.

Parameters:

df (DataFrame) – Data to turn into a State

Return type:

State | dict[float, State]

parse_value(records)[source]

Parse records in a dataframe into a float/array, depending on the shape they represent.

Parameters:

records (DataFrame) – The DataFrame that contains only records that belong to one variable.

Returns:

This is a float if the array is a scalar, or a numpy array if it is dimensional.

Return type:

Value

state_timeseries_from_dataframe(df)[source]

Reads a DataFrame into an Aggregator-readable StateTimeseries.

Parameters:

df (DataFrame) – The DataFrame to parse.

Raises:

ValueError – If the given DataFrame has no column named ‘time’.

Return type:

dict[float, State]

See also

state_timeseries_to_dataframe

to_dataframe(s)[source]

Transforms the state(s) into a DataFrame.

Parameters:

s (State | StateTimeseries) – State or a time keyed dictionary of states.

Returns:

A Pandas DataFrame.

Return type:

DataFrame

StateTimeseries

A mapping from time points to Aggregator State

alias of dict[float, State]

Solution

The solution at different times can be represented bya bunch of States and times, but a compact equivalent object is the Solution object.

The result of a time dependent solution of the initial value problem with an Aggregator

class Solution(time, data)[source]

The result of asking an Aggregator to solve a system of equations.

Parameters:
  • time (Array1D) – The times at which the solution was calculated.

  • data (Array2D) – The vector of values for each time in the time vector. Shaped as (len(time), len(state_vector))

Constraints

When solving using the IDA solver, it is possible to constrain the solution. Tools for this can be found in stream.aggregator See the sundials package documentation for more details.

class CONSTRAINT(*values)[source]

Possible values for IDA for sign constraints. See the sundials documentation for explanation of these values

create_constraints(agr, default_sign=CONSTRAINT.none, **kwargs)[source]

Create a constraint array, as expected by IDA Currently, we support sign constraints in DAE mode only. Defaults to no sign constraint for all variables.

Meant to be used as the contraints_type option for differential_algebraic.

Parameters:
  • agr (Aggregator) – The aggregator for which to create the constraints array.

  • default_sign (CONSTRAINT) – The default option to set all variables to, if not specified in kwargs.

  • kwargs (_ConstraintTypes) – The variables to set to each of the possible CONSTRAINT values.

Returns:

Array with the same shape as agr.graph, with sign contraints.

Return type:

np.ndarray