API Reference

Functions

class Solverz.sym_algebra.functions.sin(*args)

The sine function.

class Solverz.sym_algebra.functions.cos(*args)

The cosine function.

class Solverz.sym_algebra.functions.exp(*args)

The exponential function, \(e^x\).

class Solverz.sym_algebra.functions.Abs(*args)

The element-wise absolute value function:

\[\operatorname{Abs}(x)=|x|\]

with derivative

\[\frac{\mathrm{d}}{\mathrm{d}x}{\operatorname{Abs}}(x)=\operatorname{Sign}(x).\]
class Solverz.sym_algebra.functions.Sign(*args)

The element-wise indication of the sign of a number

\[\begin{split}\operatorname{Sign}(x)= \begin{cases} 1&x> 0\\ 0& x==0\\ -1&x< 0 \end{cases}\end{split}\]
class Solverz.sym_algebra.functions.AntiWindUp(u, umin, umax, e)

For PI controller

\[\begin{split}\begin{aligned} \dot{z}(t) & =e(t) \\ u_{d e s}(t) & =K_p e(t)+K_i z(t) \end{aligned},\end{split}\]

we limit the integrator output by setting, in the \(K_i>0\) case,

\[\begin{split}\dot{z}(t)= \begin{cases} 0 & \text { if } u_{\text {des }}(t) \geq u_{\max } \text { and } e(t) \geq 0 \\ 0 & \text { if } u_{\text {des }}(t) \leq u_{\min } \text { and } e(t) \leq 0 \\ e(t) & \text { otherwise } \end{cases}.\end{split}\]

So the AntiWindUp function reads

\[\dot{z}(t)= \operatorname{AntiWindUp}(u_{\text {des }}(t), u_{\min }, u_{\max }, e(t)).\]
class Solverz.sym_algebra.functions.Min(x, y)

The element-wise minimum of two numbers

\[\begin{split}\begin{split}\min(x,y)= \begin{cases} x&x\leq y\\ y& x>y \end{cases}\end{split}\end{split}\]
class Solverz.sym_algebra.functions.MatVecMul(*args)

Legacy sparse matrix-vector multiply. Use Mat_Mul instead.

MatVecMul decomposes the sparse matrix into CSC components and calls a Numba-compatible csc_matvec. New code should use Mat_Mul(A, x) which uses scipy.sparse directly and supports full matrix calculus.

class Solverz.sym_algebra.functions.Saturation(*args)

The element-wise saturation a number

\[\begin{split}\operatorname{Saturation}(v, v_\min, v_\max)= \begin{cases} v_\max&v> v_\max\\ v& v_\min\leq v\leq v_\max\\ v_\min&v< v_\min \end{cases}\end{split}\]
class Solverz.sym_algebra.functions.heaviside(*args)

The heaviside step function

\[\begin{split}\operatorname{Heaviside}(x)= \begin{cases} 1 & x >= 0\\ 0 & x < 0\\ \end{cases}\end{split}\]

which should be distinguished from sympy.Heaviside

class Solverz.sym_algebra.functions.ln(*args)

The ln function, \(ln(x)\).

Utilities

Solverz.utilities.miscellaneous.derive_incidence_matrix(indices, m, n) coo_array

Convert indices to m*n incidence matrix

Solvers

AE solver

Solverz.solvers.nlaesolver.nr.nr_method(eqn: nAE, y: ndarray, opt: Opt = None)

The Newton-Raphson method.

Parameters:

eqn : nAE

Numerical AE object.

y : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • ite_tol: 1e-5(default)|float

    The iteration error tolerance.

Returns:

sol : aesol

The aesol object.

References

Solverz.solvers.nlaesolver.cnr.continuous_nr(eqn: nAE, y: ndarray, opt: Opt = None)

The Continuous Newton-Raphson method.

Parameters:

eqn : nAE

Numerical AE object.

y : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • ite_tol: 1e-5(default)|float

    The iteration error tolerance.

Returns:

sol : aesol

The aesol object.

References

[R2]

F. Milano, ‘Continuous Newton’s Method for Power Flow Analysis’, IEEE Transactions on Power Systems, vol. 24, no. 1, pp. 50–57, Feb. 2009, doi: 10.1109/TPWRS.2008.2004820.

Solverz.solvers.nlaesolver.lm.lm(eqn: nAE, y: ndarray, opt: Opt = None)

The Levenberg-Marquardt method by minpack, wrapped in scipy.

Warning

Note that this function uses only dense Jacobian.

Parameters:

eqn : nAE

Numerical AE object.

y : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • ite_tol: 1e-5(default)|float

    The iteration error tolerance.

Returns:

sol : aesol

The aesol object.

References

[R3]

More, Jorge J., Burton S. Garbow, and Kenneth E. Hillstrom. 1980. User Guide for MINPACK-1.

Solverz.solvers.nlaesolver.sicnm.sicnm(ae: nAE, y0: ndarray, opt: Opt = None)

The semi-implicit continuous Newton method [R5]. The philosophy is to rewrite the algebraic equation

\[0=g(y)\]

as the differential algebraic equations

\[\begin{split}\begin{aligned} \dot{y}&=z \\ 0&=J(y)z+g(y) \end{aligned}\end{split}\]

with \(y_0\) being the initial value guess and \(z_0=-J(y_0)^{-1}g(y_0)\), where \(z\) is an intermediate variable introduced. Then the DAEs are solved by Rodas. SICNM is found to be more robust than the Newton’s method, for which the theoretical proof can be found in my paper [R5]. In addition, the non-iterative nature of Rodas guarantees the efficiency.

One can change the rodas scheme according to the ones implemented in the DAE version of Rodas.

SICNM used the Hessian-vector product (HVP) interface of the algebraic equation models, the make_hvp flag should be set to True when printing numerical modules.

An illustrative example of SICNM can be found in the power flow section of Solverz’ cookbook.

Parameters:

ae : nAE

Numerical AE object.

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • ite_tol: 1e-5(default)|float

    The iteration error tolerance.

Returns:

sol : aesol

The aesol object.

References

[R5] (1,2,3)
  1. Yu, W. Gu, S. Lu, and Y. Xu, “Semi-implicit continuous newton method for power flow analysis,” 2023, arXiv:2312.02809.

FDAE solver

Solverz.solvers.fdesolver.fdae_solver(fdae: nFDAE, tspan: List | ndarray, y0: ndarray, opt: Opt = None, **kwargs)

The general solver of FDAE.

Parameters:

fdae : nFDAE

Numerical FDAE object.

tspan : List | np.ndarray

An array specifying t0 and tend

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • step_size: 1e-3(default)|float

    The step size

  • ite_tol: 1e-8(default)|float

    The error tolerance of inner Newton iterations.

kwargs : dict

The value of y at previous time nodes. For example, if an FDAE uses \(y_0\), \(y_{-1}\) and \(y_{-2}\), then we can call FDAE solver with

sol = fdae_solver(mdl,
                  [0, 120],
                  y0,
                  Opt(step_size=60),
                  y1=y1,
                  y2=y2)

where y1 denotes \(y_{-1}\) and y2 denotes \(y_{-2}\).

Returns:

sol : daesol

The daesol object.

DAE solver

Solverz.solvers.daesolver.beuler.backward_euler(dae: nDAE, tspan: List | ndarray, y0: ndarray, opt: Opt = None)

The fixed step backward Euler (implicit Euler) method.

Parameters:

dae : nDAE

Numerical DAE object.

tspan : List | np.ndarray

An array specifying t0 and tend

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • step_size: 1e-3(default)|float

    The step size

Returns:

sol : daesol

The daesol object.

References

Solverz.solvers.daesolver.trapezoidal.implicit_trapezoid(dae: nDAE, tspan: List | ndarray, y0: ndarray, opt: Opt = None)

The fixed step implicit trapezoidal method.

Parameters:

dae : nDAE

Numerical DAE object.

tspan : List | np.ndarray

An array specifying t0 and tend

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • step_size: 1e-3(default)|float

    The step size

Returns:

sol : daesol

The daesol object.

References

Solverz.solvers.daesolver.rodas.Rodas(dae: nDAE, tspan: List | ndarray, y0: ndarray, opt: Opt = None)

The stiffly accurate Rosenbrock methods including Rodas4 [R8], Rodasp [R9], Rodas5p [R10].

Parameters:

dae : nDAE

Numerical DAE object.

tspan : List | np.ndarray

Either - a list specifying [t0, tend], or - a np.ndarray specifying the time nodes that you are concerned about

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • scheme: ‘rodas4’(default)|’rodasp’|’rodas5p’

    The rodas scheme

  • rtol: 1e-3(default)|float

    The relative error tolerance

  • atol: 1e-6(default)|float

    The absolute error tolerance

  • event: Callable

    The simulation events, with \(t\) and \(y\) being args, and value, is_terminal and direction being outputs

  • fixed_h: False(default)|bool

    To use fixed step size

  • hinit: None(default)|float

    Initial step size

  • hmax: None(default)|float

    Maximum step size.

  • pbar: False(default)|bool

    To display progress bar

Returns:

sol : daesol

The daesol object.

References

[R8] (1,2)

Hairer and Wanner, Solving Ordinary Differential Equations II , 2nd ed. Berlin, Heidelberg, Germany: Springer-Verlag, 1996.

[R9] (1,2)

Steinebach, Order-reduction of ROW-methods for DAEs and method of lines applications. Preprint-Nr. 1741, FB Mathematik, TH Darmstadt, 1995

[R10] (1,2)

Steinebach, “Construction of rosenbrock-wanner method rodas5p and numerical benchmarks within the julia differential equations package,” BIT, vol. 63, no. 27, Jun 2023.

Solverz.solvers.daesolver.ode15s.ode15s(dae: nDAE, tspan: List | ndarray, y0: ndarray, opt: Opt = None)

The python implementation of MATLAB Ode15s without event detection [R11].

Parameters:

dae : nDAE

Numerical DAE object.

tspan : List | np.ndarray

Either - a list specifying [t0, tend], or - a np.ndarray specifying the time nodes that you are concerned about

y0 : np.ndarray

The initial values of variables

opt : Opt

The solver options, including:

  • rtol: 1e-3(default)|float

    The relative error tolerance

  • atol: 1e-6(default)|float

    The absolute error tolerance

  • hinit: None(default)|float

    Initial step size

  • hmax: None(default)|float

    Maximum step size.

  • pbar: False(default)|bool

    To display progress bar

Returns:

sol : daesol

The daesol object.

References

[R11] (1,2)
    1. Shampine and M. W. Reichelt, “The MATLAB ODE Suite,” SIAM J. Sci. Comput., vol. 18, no. 1, pp. 1–22, Jan. 1997, doi: 10.1137/S1064827594276424.