tdcpy.ddae module#

DDAE implementation

class tdcpy.ddae.DDAE(A, hA, E=None, B=None, hB=None, C=None, hC=None, D=None, hD=None, **kwargs)#

Bases: TDSBase

Delay Differential-Algebraic Equation (DDAE)

\[E \dot{x}(t) = \sum_{k=1}^{m_A} A_k x(t - h_{A,k}) + \sum_{k=1}^{m_B} B_k u(t - h_{B,k})\]
\[y(t) = \sum_{k=1}^{m_C} C_k x(t - h_{C,k}) + \sum_{k=1}^{m_D} D_k u(t - h_{D,k})\]

where \(x(t) \in \mathbb{R}^n\) is the state vector, \(u(t) \in \mathbb{R}^p\) the input vector, and \(y(t) \in \mathbb{R}^q\) the output vector.

Parameters:
  • A (npt.NDArray) – 3D array of shape (n, n, mA) containing system matrices A_k

  • hA (npt.NDArray) – 1D array of shape (mA,) containing system delays h_{A,k}

  • E (npt.NDArray, optional) – 2D array of shape (n, n) containing descriptor matrix E, by default None, which is equivalent to identity matrix

  • B (npt.NDArray, optional) – 3D array of shape (n, p, mB) containing input matrices B_k, by default None

  • hB (npt.NDArray, optional) – 1D array of shape (mB,) containing input delays h_{B,k}, by default None

  • C (npt.NDArray, optional) – 3D array of shape (q, n, mC) containing output matrices C_k, by default None

  • hC (npt.NDArray, optional) – 1D array of shape (mC,) containing output delays h_{C,k}, by default None

  • D (npt.NDArray, optional) – 3D array of shape (q, p, mD) containing feed-through matrices D_k, by default None

  • hD (npt.NDArray, optional) – 1D array of shape (mD,) containing feed-through delays h_{D,k}, by default None

  • **kwargs (additional arguments) –

    tol_singularfloat, optional

    tolerance for considering matrix singular in null space computations, by default 1e-12

Returns:

DDAE object

Return type:

DDAE

Notes

The DDAE class represents a linear time-invariant delay differential-algebraic equation (DDAE) of the form described above.

Examples

>>> import numpy as np
>>> from tdcpy.ddae import DDAE
>>> A = np.array([[[1, 0], [0, 1]], [[0, 0], [0, 0]]], dtype=float).transpose((1,2,0))
>>> hA = np.array([0.0, 1.0], dtype=float)
>>> E = np.array([[1, 0], [0, 1]], dtype=float)
>>> B = np.array([[[1], [0]]], dtype=float).transpose((1,2,0))
>>> hB = np.array([0.0], dtype=float)
>>> C = np.array([[[1, 0], [0, 1]]], dtype=float).transpose((1,2,0))
>>> hC = np.array([0.0], dtype=float)
>>> D = np.array([[[0], [0]]], dtype=float).transpose((1,2,0))
>>> hD = np.array([0.0], dtype=float)
>>> ddae = DDAE(A=A, hA=hA, E=E, B=B, hB=hB, C=C, hC=hC, D=D, hD=hD)
>>> ddae.n
2
>>> ddae.n_inputs
1
>>> ddae.n_outputs
2
property n: int#

number of variables

property n_inputs: int#

number of inputs

property n_outputs: int#

number of outputs

property E: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

E matrix from TDS representation, shape (n, n)

property A: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

dynamics matrices

property hA: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

dynamic delays

property mA: int#

number of delays

property B: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

input matrices

property hB: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

input delays

property mB: int#

number of input delays

property C: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

output matrices

property hC: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

output delays

property mC: int#

number of output delays

property D: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

feed-through matrices

property hD: ndarray[tuple[Any, ...], dtype[_ScalarT]] | None#

feed-through delays

property mD: int#

number of feed-through delays

property uE: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

orthonormal basis for left null space of E

property vE: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

orthonormal basis for right null space of E

property is_logical: bool#

property from original MATLAB package, unused as of now

property is_compressed: bool#

Cheks if DDAE is in compressed form (no duplicates in hA)

property is_sorted: bool#

Checks if DDAe is in sorted form (ascending hA)

property is_lti: bool#

Checks if DDAE is Linear Time-invariant

property is_real: bool#

Checks if DDAE uses complex storage for any of defining arrays

property is_delay_difference_equation: bool#

Checks if DDAE is a delay difference equation, i.e. E==0

property is_normalized_delay_difference_equation: bool#

Check if DDAE is a normalized delay difference equation.

A DDAE is normalized if it satisfies all of the following conditions:

  • DDAE is not logical

  • E == 0

  • A[0] == I

Tolerances used for comparisons are absolute tolerance of 1e-12.

property is_essentially_retarded: bool#

Checks if DDAE is essentially retarded

DDAE is essentially retarded IFF characteristic equation of underlying delay-difference equation

det | SUM A[i] * exp(-s*hA[i]) |

does not depend on complex argument s, i.e. delay difference equation has to take form:

0 = A[0] * x(t)

property is_essentially_neutral#

Checks if DDAE is essentialy netural

get_delay_difference_equation(**kwargs)#

Creates Delay-Difference Equation from DDAE

For a DDAE, the associated delay difference equation is given by: add equation where U and V are orthogonal matrices whose columns form a basis for the null add equation

Parameters:

**kwargs

Additional keyword arguments.

tolfloat, optional

Norm tolerance for considering matrix vanish. Default is 1e-14.

rcondfloat, optional

Relative condition number. Singular values s smaller than rcond * max(s) are considered zero in null space construction. Default is 1e-12.

Returns:

DDAE representing delay difference equation or None if E is non-singular (there is no delay difference equation).

Return type:

DDAE or None

Raises:

ValueError – If the DDAE is logical. (not supported as of now, but is ready for future)

to_asymptotic_transfer_function(**kwargs)#

TODO - implement this function

Return type:

DDAE

sort(inplace=False)#

Sorts arrays containing matrices and delays (ascending order)

Parameters:

inplace (bool) – wheter to modify existing DDEA or create a new one, default False

Returns:

sorted representation None: if inplace=True (current object is updated)

Return type:

DDAE

compress(inplace=False, rtol=1e-05, atol=1e-08)#

Removes delay duplicates, sorts delays into ascending order

Parameters:
  • inplace (bool) – wheter to modify existing DDEA or create a new one, default False

  • rtol (float) – relative tolerance for determining matrix element is zero, default 1e-5

  • atol (float) – absolute tolerance for determining matrix element is zero, default 1e-8

Returns:

compressed representation None: if inplace=True (current object is updated)

Return type:

DDAE

eval_char_matrix(s)#

Evaluate characteristic matrix at s

M(s) = E*s - A[0]*exp(-s*hA[0]) - … - A[mA]*exp(-s*hA[mA])

Parameters:

s (complex, float, int) – s from complex plane

Returns:

characteristic matrix M evaluated at s

Return type:

M (array)

eval_char_matrix_derivative(s)#

Derivative of characteristic matrix with respect to s evaluated at s

dM(s) = E + hA[0]*A[0]*exp(-s*hA[0]) + … + hA[mA]*A[mA]*exp(-s*hA[mA])

Parameters:

s (complex, float, int) – s from complex plane

Returns:

derivative of characteristic matrix M evaluated at s

Return type:

dM (array)