Simulation

class pathsim.simulation.Simulation(blocks=None, connections=None, events=None, dt=0.01, dt_min=1e-16, dt_max=None, Solver=<class 'pathsim.solvers.ssprk22.SSPRK22'>, tolerance_fpi=1e-10, iterations_max=200, log=True, **solver_kwargs)[source]

Bases: object

Class that performs transient analysis of the dynamical system, defined by the blocks and connecions. It manages all the blocks and connections and the timestep update.

The global system equation is evaluated by fixed point iteration, so the information from each timestep gets distributed within the entire system and is available for all blocks at all times.

The minimum number of fixed-point iterations ‘iterations_min’ is set to ‘None’ by default and then the length of the longest internal signal path (with passthrough) is used as the estimate for minimum number of iterations needed for the information to reach all instant time blocks in each timestep. Dont change this unless you know that the actual path is shorter or something similar that prohibits instant time information flow.

Convergence check for the fixed-point iteration loop with ‘tolerance_fpi’ is based on max absolute error (max-norm) to previous iteration and should not be touched.

Multiple numerical integrators are implemented in the ‘pathsim.solvers’ module. The default solver is a fixed timestep 2nd order Strong Stability Preserving Runge Kutta (SSPRK22) method which is quite fast and has ok accuracy, especially if you are forced to take small steps to cover the behaviour of forcing functions. Adaptive timestepping and implicit integrators are also available.

Manages an event handling system based on zero crossing detection. Uses ‘Event’ objects to monitor solver states of stateful blocks and applys transformations on the state in case an event is detected.

Example

This is how to setup a simple system simulation using the ‘Simulation’ class:

import numpy as np

from pathsim import Simulation, Connection
from pathsim.blocks import Source, Integrator, Scope

src = Source(lambda t: np.cos(2*np.pi*t))
itg = Integrator()
sco = Scope(labels=["source", "integrator"])

sim = Simulation(
    blocks=[src, itg, sco],
    connections=[
        Connection(src[0], itg[0], sco[0]),
        Connection(itg[0], sco[1])
        ],
    dt=0.01
    )

sim.run(4)
sim.plot()
Parameters:
  • blocks (list[Block]) – blocks that define the system

  • connections (list[Connection]) – connections that connect the blocks

  • events (list[Event]) – list of event trackers (zero crossing detection, schedule, etc.)

  • dt (float) – transient simulation timestep in time units, default see ´SIM_TIMESTEP´ in ´_constants.py´

  • dt_min (float) – lower bound for transient simulation timestep, default see ´SIM_TIMESTEP_MIN´ in ´_constants.py´

  • dt_max (float) – upper bound for transient simulation timestep, default see ´SIM_TIMESTEP_MAX´ in ´_constants.py´

  • Solver (Solver) – ODE solver class for numerical integration from ´pathsim.solvers´, default is ´pathsim.solvers.ssprk22.SSPRK22´ (2nd order expl. Runge Kutta)

  • tolerance_fpi (float) – absolute tolerance for convergence of algebraic loops and internal optimizers of implicit ODE solvers, default see ´SIM_TOLERANCE_FPI´ in ´_constants.py´

  • iterations_max (int) – maximum allowed number of iterations for implicit ODE solver optimizers and algebraic loop solver, default see ´SIM_ITERATIONS_MAX´ in ´_constants.py´

  • log (bool | string) – flag to enable logging, default see ´LOG_ENABLE´ in ´_constants.py´ (alternatively a path to a log file can be specified)

  • solver_kwargs (dict) – additional parameters for numerical solvers such as absolute (´tolerance_lte_abs´) and relative (´tolerance_lte_rel´) tolerance, defaults are defined in ´_constants.py´

time

global simulation time, starting at ´0.0´

Type:

float

graph

internal graph representation for fast system funcion evluations using DAG with algebraic depths

Type:

Graph

boosters

list of boosters (fixed point accelerators) that wrap algebraic loop closing connections assembled from the system graph

Type:

None | list[ConnectionBooster]

engine

global integrator (ODE solver) instance serving as a dummy to get attributes and access to intermediate evaluation stages

Type:

Solver

logger

global simulation logger

Type:

logging.Logger

_blocks_dyn

blocks with internal ´Solver´ instances (stateful)

Type:

set[Block]

_blocks_evt

blocks with internal events (discrete time, eventful)

Type:

set[Block]

_active

flag for setting the simulation as active, used for interrupts

Type:

bool

property size

Get size information of the simulation, such as total number of blocks and dynamic states, with recursive retrieval from subsystems

Returns:

sizes – size of simulation (number of blocks) and number of internal states (from internal engines)

Return type:

tuple[int]

plot(*args, **kwargs)[source]

Plot the simulation results by calling all the blocks that have visualization capabilities such as the ‘Scope’ and ‘Spectrum’.

This is a quality of life method. Blocks can be visualized individually due to the object oriented nature, but it might be nice to just call the plot metho globally and look at all the results at once. Also works for models loaded from an external file.

Parameters:
  • args (tuple) – args for the plot methods

  • kwargs (dict) – kwargs for the plot method

save(path='', **metadata)[source]

Save the dictionary representation of the simulation instance to an external file

Parameters:
  • path (str) – filepath to save data to

  • metadata (dict) – metadata for the simulation model

classmethod load(path='', **kwargs)[source]

Load and instantiate a Simulation from an external file in json format

Parameters:
  • path (str) – filepath to load data from

  • kwargs (dict) – additional args for the simulation, overwriting metadata

Returns:

out – reconstructed object from dict representation

Return type:

Simulation

to_dict(**metadata)[source]

Convert simulation to a complete model representation as a dict with additional metadata.

Parameters:

metadata (dict) – metadata for the simulation model

Returns:

data – dict that describes the simulation model

Return type:

dict

classmethod from_dict(data, **kwargs)[source]

Create simulation from model data dict

Parameters:
  • data (dict) – model definition in json format

  • kwargs (dict) – additional args for the simulation, overwriting metadata

Returns:

simulation – instance of the Simulation class with mode definition

Return type:

Simulation

add_block(block, _defer_graph=False)[source]

Adds a new block to the simulation, initializes its local solver instance and collects internal events of the new block.

This works dynamically for running simulations.

Parameters:
  • block (Block) – block to add to the simulation

  • _defer_graph (bool) – flag for defering graph construction to a later stage

add_connection(connection, _defer_graph=False)[source]

Adds a new connection to the simulaiton and checks if the new connection overwrites any existing connections.

This works dynamically for running simulations.

Parameters:
  • connection (Connection) – connection to add to the simulation

  • _defer_graph (bool) – flag for defering graph construction to a later stage

add_event(event)[source]

Checks and adds a new event to the simulation.

This works dynamically for running simulations.

Parameters:

event (Event) – event to add to the simulation

reset(time=0.0)[source]

Reset the blocks to their initial state and the global time of the simulation.

For recording blocks such as ‘Scope’, their recorded data is also reset.

Resets linearization automatically, since resetting the blocks resets their internal operators.

Afterwards the system function is evaluated with ‘_update’ to update the block inputs and outputs.

Parameters:

time (float) – simulation time for reset

linearize()[source]

Linearize the full system in the current simulation state at the current simulation time.

This is achieved by linearizing algebraic and dynamic operators of the internal blocks. See definition of the ‘Block’ class.

Before linearization, the global system function is evaluated to get the blocks into the current simulation state. This is only really relevant if no solving attempt has been happened before.

delinearize()[source]

Revert the linearization of the full system.

steadystate(reset=False)[source]

Find steady state solution (DC operating point) of the system by switching all blocks to steady state solver, solving the fixed point equations, then switching back.

The steady state solver forces all the temporal derivatives, i.e. the right hand side equation (including external inputs) of the engines of dynamic blocks to zero.

Note

This is really a sort of pseudo-steady-state solve. It does NOT compute the limit \(t\rightarrow\infty\) but rather forces all time derivatives to zero at a given moment in time.

This means, for a given t it computes the block states x such that:

\[0 = f(x, t)\]

instead of the real steady state:

\[\lim_{t \rightarrow \infty} x(t)\]
Parameters:

reset (bool) – reset the simulation before solving for steady state (default False)

timestep_fixed_explicit(dt=None)[source]

Advances the simulation by one timestep ‘dt’ for explicit fixed step solvers.

If discrete events are detected, they are resolved immediately within the timestep.

Parameters:

dt (float) – timestep

Returns:

  • success (bool) – indicator if the timestep was successful

  • max_error (float) – maximum local truncation error from integration

  • scale (float) – rescale factor for timestep

  • total_evals (int) – total number of system evaluations

  • total_solver_its (int) – total number of implicit solver iterations

timestep_fixed_implicit(dt=None)[source]

Advances the simulation by one timestep ‘dt’ for implicit fixed step solvers.

If discrete events are detected, they are resolved immediately within the timestep.

Parameters:

dt (float) – timestep

Returns:

  • success (bool) – indicator if the timestep was successful

  • max_error (float) – maximum local truncation error from integration

  • scale (float) – rescale factor for timestep

  • total_evals (int) – total number of system evaluations

  • total_solver_its (int) – total number of implicit solver iterations

timestep_adaptive_explicit(dt=None)[source]

Advances the simulation by one timestep ‘dt’ for explicit adaptive solvers.

If the local truncation error of the solver exceeds the tolerances set in the ‘solver_kwargs’, the simulation state is reverted to the state that was buffered (_buffer(time, dt)) at the beginning of the timestep.

If discrete events are detected, the chronologically first event is handled only. The event location (in time) is approached adaptively by reverting the step and adjusting the stepsize (this is equivalent to the secant method for finding zeros of the event function) until the tolerance of the event is satisfied (close==True).

Parameters:

dt (float) – timestep

Returns:

  • success (bool) – indicator if the timestep was successful

  • max_error (float) – maximum local truncation error from integration

  • scale (float) – rescale factor for timestep

  • total_evals (int) – total number of system evaluations

  • total_solver_its (int) – total number of implicit solver iterations

timestep_adaptive_implicit(dt=None)[source]

Advances the simulation by one timestep ‘dt’ for implicit adaptive solvers.

If the local truncation error of the solver exceeds the tolerances set in the ‘solver_kwargs’, the simulation state is reverted to the state that was buffered (_buffer(time, dt)) at the beginning of the timestep.

If the solution of the implicit update equation in ‘solve’ doesnt converge, the timestep is also considered unsuccessful. Then it is reverted and the timestep is halfed.

If discrete events are detected, the chronologically first event is handled only. The event location (in time) is approached adaptively by reverting the step and adjusting the stepsize (this is equivalent to the secant method for finding zeros of the event function) until the tolerance of the event is satisfied (close==True).

Parameters:

dt (float) – timestep

Returns:

  • success (bool) – indicator if the timestep was successful

  • max_error (float) – maximum local truncation error from integration

  • scale (float) – rescale factor for timestep

  • total_evals (int) – total number of system evaluations

  • total_solver_its (int) – total number of implicit solver iterations

timestep(dt=None, adaptive=True)[source]

Advances the transient simulation by one timestep ‘dt’.

Automatic stepping method selection based on selected Solver.

Parameters:
  • dt (float) – timestep size for transient simulation

  • adaptive (bool) – explicitly select the addaptive timestepping branch

Returns:

  • success (bool) – indicator if the timestep was successful

  • max_error (float) – maximum local truncation error from integration

  • scale (float) – rescale factor for timestep

  • total_evals (int) – total number of system evaluations

  • total_solver_its (int) – total number of implicit solver iterations

step(dt=None, adaptive=True)[source]

Wraps ‘Simulation.timestep’ for backward compatibility

stop()[source]

Set the flag for active simulation to ‘False’, intended to be called from the outside (for example by events) to interrupt the timestepping loop in ‘run’.

run(duration=10, reset=False, adaptive=True)[source]

Perform multiple simulation timesteps for a given ‘duration’.

Tracks the total number of block evaluations (proxy for function calls, although larger, since one function call of the system equation consists of many block evaluations) and the total number of solver iterations for implicit solvers.

Additionally the progress of the simulation is tracked by a custom ‘ProgressTracker’ class that is a dynamic generator and interfaces the logging system.

Parameters:
  • duration (float) – simulation time (in time units)

  • reset (bool) – reset the simulation before running (default False)

  • adaptive (bool) – use adaptive timesteps if solver is adaptive (default True)

Returns:

stats – stats of simulation run tracked by the ´ProgressTracker´

Return type:

dict