Solver Base

class pathsim.solvers._solver.Solver(initial_value=0, parent=None, tolerance_lte_abs=1e-08, tolerance_lte_rel=0.0001)[source]

Bases: object

Base skeleton class for solver definition. Defines the basic solver methods and the metadata.

Specific solvers need to implement (some of) the base class methods defined here. This depends on the type of solver (implicit/explicit, multistage, adaptive).

Parameters:
  • initial_value (float, array) – initial condition / integration constant

  • tolerance_lte_abs (float) – absolute tolerance for local truncation error (for solvers with error estimate)

  • tolerance_lte_rel (float) – relative tolerance for local truncation error (for solvers with error estimate)

  • parent (None | Solver) – parent solver instance that manages the intermediate stages, stage counter, etc.

x

internal ‘working’ state

Type:

numeric, array[numeric]

history

internal history of past results

Type:

deque[numeric]

n

order of integration scheme

Type:

int

s

number of internal intermediate stages

Type:

int

_stage

counter for current intermediate stage

Type:

int

eval_stages

rations for evaluation times of intermediate stages

Type:

list[float]

property stage

stage property management to interface with parent solver

Returns:

stage – current intermediate evaluation stage of solver

Return type:

int

is_first_stage()[source]
is_last_stage()[source]
stages(t, dt)[source]

Generator that yields the intermediate evaluation time during the timestep ‘t + ratio * dt’ and also updates the current stage number for internal use.

Parameters:
  • t (float) – evaluation time

  • dt (float) – integration timestep

get()[source]

Returns current internal state of the solver.

Returns:

x – current internal state of the solver

Return type:

numeric, array[numeric]

set(x)[source]

Sets the internal state of the integration engine.

This method is required for event based simulations, and to handle discontinuities in state variables.

Parameters:

x (numeric, array[numeric]) – new internal state of the solver

reset()[source]

“Resets integration engine to initial value

buffer(dt)[source]

Saves the current state to an internal state buffer which is especially relevant for multistage and implicit solvers.

Multistep solver implement rolling buffers for the states and timesteps.

Resets the stage counter.

Parameters:

dt (float) – integration timestep

classmethod cast(other, parent, **solver_kwargs)[source]

Cast the integration engine to the new type and initialize with previous solver arguments so it can continue from where the ‘old’ solver stopped.

Parameters:
  • other (Solver) – solver instance to cast to new solver type

  • parent (None | Solver) – solver instance to use as parent

  • solver_kwargs (dict) – additional args for the new solver

Returns:

engine – new solver instance cast from other

Return type:

Solver

error_controller()[source]

Returns the estimated local truncation error (abs and rel) and scaling factor for the timestep, only relevant for adaptive timestepping methods.

Returns:

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

  • error (float) – estimated error of the internal error controller

  • scale (float) – estimated timestep rescale factor for error control

revert()[source]

Revert integration engine to previous timestep.

This is only relevant for adaptive methods where the simulation timestep ‘dt’ is rescaled and the engine step is recomputed with the smaller timestep.

step(f, dt)[source]

Performs the explicit timestep for (t+dt) based on the state and input at (t).

Returns the local truncation error estimate and the rescale factor for the timestep if the solver is adaptive.

Parameters:
  • f (numeric, array[numeric]) – evaluation of rhs function

  • dt (float) – integration timestep

Returns:

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

  • error (float) – estimated error of the internal error controller

  • scale (float) – estimated timestep rescale factor for error control

interpolate(r, dt)[source]

Interpolate solution after successful timestep as a ratio in the interval [t, t+dt].

This is especially relevant for Runge-Kutta solvers that have a higher order interpolant. Otherwise this is just linear interpolation using the buffered state.

Parameters:
  • r (float) – ration for interpolation within timestep

  • dt (float) – integration timestep

Returns:

x – interpolated state

Return type:

numeric, array[numeric]

class pathsim.solvers._solver.ExplicitSolver(*solver_args, **solver_kwargs)[source]

Bases: Solver

Base class for explicit solver definition.

x_0

internal ‘working’ initial value

Type:

numeric, array[numeric]

x

internal ‘working’ state

Type:

numeric, array[numeric]

n

order of integration scheme

Type:

int

s

number of internal intermediate stages

Type:

int

stage

counter for current intermediate stage

Type:

int

eval_stages

rations for evaluation times of intermediate stages

Type:

list[float]

integrate_singlestep(func, time=0.0, dt=0.01)[source]

Directly integrate the function for a single timestep ‘dt’ with explicit solvers. This method is primarily intended for testing purposes.

Parameters:
  • func (callable) – function to integrate f(x, t)

  • time (float) – starting time for timestep

  • dt (float) – integration timestep

Returns:

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

  • error_norm (float) – estimated error of the internal error controller

  • scale (float) – estimated timestep rescale factor for error control

integrate(func, time_start=0.0, time_end=1.0, dt=0.01, dt_min=1e-16, dt_max=None, adaptive=True)[source]

Directly integrate the function ‘func’ from ‘time_start’ to ‘time_end’ with timestep ‘dt’ for explicit solvers.

This method is primarily intended for testing purposes or for use as a standalone numerical integrator.

Example

This is how to directly use the solver to integrate an ODE:

#1st order linear ODE
def f(x, u, t):
    return -x

#initial condition
x0 = 1

#initialize ODE solver
sol = Solver(x0)

#integrate from 0 to 5 with timestep 0.1
t, x = sol.integrate(f, time_end=5, dt=0.1)
Parameters:
  • func (callable) – function to integrate f(x, t)

  • time_start (float) – starting time for integration

  • time_end (float) – end time for integration

  • dt (float) – timestep or initial timestep for adaptive solvers

  • dt_min (float) – lower bound for timestep, default ‘0.0’

  • dt_max (float) – upper bound for timestep, default ‘None’

  • adaptive (bool) – use adaptive timestepping if available

Returns:

  • outout_times (array[float]) – time points of the solution

  • output_states (array[numeric], array[array[numeric]]) – state values at solution time points

class pathsim.solvers._solver.ImplicitSolver(*solver_args, **solver_kwargs)[source]

Bases: Solver

Base class for implicit solver definition.

x_0

internal ‘working’ initial value

Type:

numeric, array[numeric]

x

internal ‘working’ state

Type:

numeric, array[numeric]

n

order of integration scheme

Type:

int

s

number of internal intermediate stages

Type:

int

stage

counter for current intermediate stage

Type:

int

eval_stages

rations for evaluation times of intermediate stages

Type:

list[float]

opt

optimizer instance to solve the implicit update equation

Type:

NewtonAnderson, Anderson, etc.

buffer(dt)[source]

Saves the current state to an internal state buffer which is especially relevant for multistage and implicit solvers.

Resets the stage counter and the optimizer of implicit methods.

Parameters:

dt (float) – integration timestep

solve(j, J, dt)[source]

Advances the solution of the implicit update equation of the solver with the optimizer of the engine and tracks the evolution of the solution by providing the residual norm of the fixed-point solution.

Parameters:
  • f (numeric, array[numeric]) – evaluation of rhs function

  • J (array[numeric]) – evaluation of jacobian of rhs function

  • dt (float) – integration timestep

Returns:

err – residual error of the fixed point update equation

Return type:

float

integrate_singlestep(func, jac, time=0.0, dt=0.01, tolerance_fpi=1e-09, max_iterations=200)[source]

Directly integrate the function ‘func’ for a single timestep ‘dt’ with implicit solvers. This method is primarily intended for testing purposes.

Parameters:
  • func (callable) – function to integrate f(x, t)

  • jac (callable) – jacobian of f w.r.t. x

  • time_start (float) – starting time for timestep

  • dt (float) – integration timestep

  • tolerance_fpi (float) – convergence criterion for implicit update equation

  • max_iterations (int) – maximum numer of iterations for optimizer to solve implicit update equation

Returns:

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

  • error_norm (float) – estimated error of the internal error controller or solver when not converged

  • scale (float) – estimated timestep rescale factor for error control

integrate(func, jac, time_start=0.0, time_end=1.0, dt=0.01, dt_min=1e-16, dt_max=None, adaptive=True, tolerance_fpi=1e-09, max_iterations=200)[source]

Directly integrate the function ‘func’ from ‘time_start’ to ‘time_end’ with timestep ‘dt’ for implicit solvers.

This method is primarily intended for testing purposes or for use as a standalone numerical integrator.

Example

This is how to directly use the solver to integrate an ODE:

#1st order linear ODE
def f(x, t):
    return -x

#initial condition
x0 = 1

#initialize ODE solver
sol = Solver(x0)

#integrate from 0 to 5 with timestep 0.1
t, x = sol.integrate(f, time_end=5, dt=0.1)
Parameters:
  • func (callable) – function to integrate f(x, t)

  • jac (callable) – jacobian of f w.r.t. x

  • time_start (float) – starting time for integration

  • time_end (float) – end time for integration

  • dt (float) – timestep or initial timestep for adaptive solvers

  • dt_min (float) – lower bound for timestep, default ‘0.0’

  • dt_max (float) – upper bound for timestep, default ‘None’

  • adaptive (bool) – use adaptive timestepping if available

  • tolerance_fpi (float) – convergence criterion for implicit update equation

  • max_iterations (int) – maximum numer of iterations for optimizer to solve implicit update equation

Returns:

  • outout_times (array[float]) – time points of the solution

  • output_states (array[numeric], array[array[numeric]]) – state values at solution time points