tdcpy.common package#

The tdcpy.common package contains common utilities for the tdcpy package, including functions for working with quasipolynomials, discretization, and composition of systems.

Submodules#

tdcpy.common.closed_loop module#

tdcpy.common.closed_loop.controller_reprezentation(order, n_inputs, n_outputs, hA=None, hB=None, hC=None, hD=None, **kwargs)#

Creates an empty controller reprezentation

Parameters:
  • order (int) – Controller order (0 for static controller)

  • n_inputs (int) – Number of controller inputs

  • n_outputs (int) – Number of controller outputs

  • hA (npt.ndarray, optional) – Vector of delays for A matrix

  • hB (npt.ndarray, optional) – Vector of delays for B matrix

  • hC (npt.ndarray, optional) – Vector of delays for C matrix

  • hD (npt.ndarray, optional) – Vector of delays for D matrix

  • **kwargs (dict, optional) – Additional arguments (not used)

Returns:

  • E (npt.NDArray) – Descriptor matrix of the controller

  • K (npt.NDArray) – System matrix of the controller

  • hK (npt.NDArray) – Vector of delays for K matrix

Notes

Given the order of the controller, number of inputs and outputs, and optional delay vectors for matrices A, B, C, D, this function creates an empty controller representation of a dynamic controller (order > 0) or static controller (order = 0).

\[E x' = \sum_{k=0}^{hK} K_k x(t - hK_k)\]

where K is constructed for controller matrices A, B, C, D as:

\[\begin{split}K = \begin{bmatrix} A & B \\ C & D \end{bmatrix}\end{split}\]

and hK is the concatenated vector of delays i.e. \(hK = [hA_c \quad hB_c \quad hC_c \quad hD_c]\).

  • For static controller (order=0) only hD is used, other delay vectors are ignored (and a warning is issued if they are provided)

  • If order > 0 and user does not provide delay vectors, they are set to [0.0] by default

  • The resulting controller is fully connected (all entries in A, B, C, D matrices are True)

  • The resulting controller is in the form suitable for creating ClosedLoop object

Examples

>>> import numpy as np
>>> from tdcpy.common.closed_loop import controller_reprezentation
>>> E, K, hK = controller_reprezentation(order=1, n_inputs=2, n_outputs=1)
>>> E.shape
(2, 3)
>>> K.shape
(2, 3, 1)
>>> hK
array([0.], dtype=float32)
>>> E, K, hK = controller_reprezentation(order=0, n_inputs=2, n_outputs=1, hD=np.array([0.0, 1.0]))
>>> E.shape
(1, 2)
>>> K.shape
(1, 2, 2)
>>> hK
array([0., 1.])

tdcpy.common.composition module#

Set of functions for TDS composition#

The composition module contains functions required for creating a packed representation of a controller or a system itself. The packed representation can later be used in the closed-loop interconnection with the system.

Implemented functions:

  1. concatenate_2x2_by_delays: creates a concatenated DDAE (E*, A*, hA*) given a DDAE representation (E,A,B,C,D) and the repective delay matrices (hA, hB, hC, hD).

tdcpy.common.composition.concatenate_2x2_by_delays(E, A, B, C, D, hA, hB, hC, hD, EE=None, AA=None, hAA=None)#

Concatenates system into compact form respecting delay vectors.

Parameters:
  • E (array) – LHS matrix

  • A (array) – left hand side matrices of DDAE, assumed non-empty

  • B (array) – 3D array of representing input matrices

  • C (array) – 3D array of representing output matrices

  • D (array) – left hand side matrices of DDAE, assumed non-empty

  • hA (array) – vector of delays associated with array A

  • hB (array) – vector of delays associated with array B

  • hC (array) – vector of delays associated with array C

  • hD (array) – vector of delays associated with array D

  • EE (array, optional) – if provided, used as LHS matrix of the concatenated system

  • AA (array, optional) – if provided, used as RHS 3D array of the concatenated system

  • hAA (array, optional) – if provided, used as delay vector of the concatenated system

Returns:

A tuple containing

EEarray

2d array of concatenated LHS

AAarray

3d array of concatenated RHS

hAAarray

1d vector of concatenated delays associated with RHS

Return type:

tuple

Notes

Assumes the system is defined as

\[ \begin{align}\begin{aligned}E \dot{x}(t) &= A_0 x(t - h_{A,0}) + ... + A_{n} x(t - h_{A,n}) + + B_0 u(t - h_{B,0}) + ... + B_{m} u(t - h_{B,m})\\ y(t) &= C_0 x(t - h_{C,0}) + ... + C_{p} x(t - h_{C,p}) + + D_0 u(t - h_{D,0}) + ... + D_{q} u(t - h_{D,q})\end{aligned}\end{align} \]

Concatenates the system into:

\[E^* \dot{x}_1(t) = A_0^* x_2(t - h A_0^*) + \dotsb + A_{n^*}^* x_2(t - h A_{n^*}^*)\]

with

\(x_1 := \begin{bmatrix} x \\ y \end{bmatrix}\), \(x_2 := \begin{bmatrix} x \\ u \end{bmatrix}\), \(hA^* = [hA \quad hB \quad hC \quad hD]\), and \(n^* = n+m+p+q\).

\[\begin{split}x_1 := \begin{bmatrix} x \\ y \end{bmatrix}, \quad x_2 := \begin{bmatrix} x \\ u \end{bmatrix}\end{split}\]

and therefore:

\[h_A^* = [hA \quad hB \quad hC \quad hD], \quad n^* = n+m+p+q\]

with the new left-hand side matrix:

\[\begin{split}E^* := \begin{bmatrix} E & 0 \\ 0 & 0\end{bmatrix}\end{split}\]

right-hand side array:

\(A^*\{:,:, :n\} := \begin{bmatrix} A & 0 \\ 0 & 0\end{bmatrix}\), \(A^*\{:,:, n:n+m\} := \begin{bmatrix} 0 & B \\ 0 & 0\end{bmatrix}\), \(A^*\{:,:, n+m:n+m+p\} := \begin{bmatrix} 0 & 0 \\ C & 0\end{bmatrix}\), \(A^*\{:,:, n+m+p:\} = \begin{bmatrix} 0 & 0 \\ 0 & D \end{bmatrix}\).

  • Assumes all input arrays are non-empty

  • If EE, AA, hAA are provided, they are updated in place and returned

  • If EE, AA, hAA are not provided, they are created and returned

  • The resulting system is in the form suitable for creating a ClosedLoop object

Examples

>>> import numpy as np
>>> from tdcpy.common.composition import concatenate_2x2_by_delays
>>> E = np.array([[1, 0], [0, 0]])
>>> A = np.array([[[0, -1], [1, 0]], [[0, 0], [0, 0]]])
>>> B = np.array([[[0], [1]], [[0], [0]]])
>>> C = np.array([[[1, 0]], [[0, 0]]])
>>> D = np.array([[[0]], [[1]]])
>>> hA = np.array([0., 1.])
>>> hB = np.array([0.])
>>> hC = np.array([0.,0.])
>>> hD = np.array([0.])
>>> EE, AA, hAA = concatenate_2x2_by_delays(E, A, B, C, D, hA, hB, hC, hD)
>>> EE
array([[1, 0, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0]])
>>> AA[:,:,0]
array([[0, 1, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0]])
>>> AA[:,:,1]
array([[-1,  0,  0,  0],
       [ 0,  0,  0,  0],
       [ 0,  0,  0,  0],
       [ 0,  0,  0,  0]])
>>> AA[:,:,2]
array([[0, 0, 0, 1],
       [0, 0, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0]])
>>> hAA
array([0., 1., 0., 0., 0., 0.])
tdcpy.common.composition.interconnect()#

tdcpy.common.compress module#

Set of function for representation compressions#

The compress module consists of functions for obtaining a minimal representaation of a given set of matrices and delays.

compression := obtaining minimal sorted representation of “something”

Implemented functions:

  1. compress_matrices_delays: creates a compressed matrix representation (A*, hA*) by removing duplicate delays

  2. compress_bool_matrices_delays: removes delay duplicates and sorts matrices in ascending of delays

  3. sort_matrices_delays: sorts the matrices in ascending order of the delays

Notes

tdcpy.common.compress.compress_matrices_delays(A, hA, rtol=1e-05, atol=1e-08)#

Compresses matrices - delays representation

Removes delay duplicates, sorts delays into ascending order and removes matrices close to zero, i.e. converts the representation (A, hA):

\[A_0 x(t - h_{A,0}) + ... + A_{mA} x(t - h_{A,mA})\]

into representation (A*, hA*), where:

  1. matrices A*[i] are NOT close to zero

  2. hA* does not contain duplicates

Parameters:
  • A (array) – 3D array of stacked matrices (axis 2)

  • hA (array) – 1D array (vector) of delays

  • 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:

A tuple containing

compressed_Aarray

compressed representation of A

compressed_hAarray

compressed vector of delays

Return type:

tuple

Notes

  1. the sort is “stable” (see numpy.argsort implementation)

  2. if all matrices are close to zero, returns empty array with

    shape (A.shape[0], A.shape[1], 0) and empty hA array

  3. if A is empty, returns A and hA unchanged

  4. if hA is empty, returns A and hA unchanged

  5. if A and hA have inconsistent shapes, raises ValueError

  6. if rtol or atol are negative, raises ValueError

  7. if A is not 3D array, raises ValueError

  8. if hA is not 1D array, raises ValueError

  9. if A.shape[2] != hA.shape[0], raises ValueError

Examples

>>> import numpy as np
>>> from tdcpy.common.compress import compress_matrices_delays
>>> A0 = np.array([[0, 1],[1, 0]])
>>> A1 = np.array([[1, 0],[0, 0]])
>>> A2 = np.array([[0, 0],[0, 1]])
>>> A = np.stack((A0, A1, A2), axis=2)
>>> hA = np.array([0., 2., 0.])
>>> compress_matrices_delays(A, hA)
(array([[[0., 1.],
    [1., 0.]],
   [[1., 0.],
    [1., 0.]]]), array([0., 2.]))
tdcpy.common.compress.compress_bool_matrices_delays(A, hA, keep_zeros=False)#

Compresses boolean matrices - delays representation

Removes delay duplicates, and sorts delays into ascending order, i.e. converts the representation (A, hA):

\[A_0 x(t - h_{A,0}) + ... + A_{mA} x(t - h_{A,mA})\]

into the representation (A*, hA*), where:

  1. matrices A*[i] are NOT zero

  2. hA* does not contain duplicates

Parameters:
  • A (array) – 3D array of stacked boolean matrices (axis 2)

  • hA (array) – 1D array (vector) of delays

  • keep_zeros (bool) – whether to keep zero matrices in the compressed representation, default False

Returns:

A tuple containing

compressed_Aarray

compressed representation of A

compressed_hAarray

compressed vector of delays

Return type:

tuple

Notes

  1. the sort is “stable” (see numpy.argsort implementation)

  2. if all matrices are zero, returns empty array with shape (A.shape[0], A.shape[1], 0) and empty hA array

  3. if A is empty, returns A and hA unchanged

  4. if hA is empty, returns A and hA unchanged

  5. if A and hA have inconsistent shapes, raises ValueError

  6. if A is not 3D array, raises ValueError

  7. if hA is not 1D array, raises ValueError

  8. if A.shape[2] != hA.shape[0], raises ValueError

Examples

>>> import numpy as np
>>> from tdcpy.common.compress import compress_bool_matrices_delays
>>> A0 = np.array([[0, 1],[1, 0]])
>>> A1 = np.array([[1, 0],[0, 0]])
>>> A2 = np.array([[0, 0],[0, 1]])
>>> A = np.stack((A0, A1, A2), axis=2)
>>> hA = np.array([0., 2., 0.])
>>> compress_bool_matrices_delays(A, hA)
(array([[[0., 0.],
    [1., 0.]],
   [[1., 0.],
    [1., 0.]]]), array([0., 2.]))
tdcpy.common.compress.sort_matrices_delays(A, hA)#

Sorts matrices - delays representation (ascending order by delays)

Sorts delays into ascending order, i.e. changes the representation (A, hA):

\[A_0 x(t - h_{A,0}) + ... + A_{mA} x(t - h_{A,mA})\]

into representation (A*, hA*), where:

  1. hA* is now in ascending order

  2. shapes are preserved

Parameters:
  • A (array) – 3D array of stacked matrices (axis 2)

  • hA (array) – 1D array (vector) of delays

Returns:

A tuple containing

compressed_Aarray

compressed representation of A

compressed_hAarray

compressed vector of delays

Return type:

tuple

Notes

  1. the sort is “stable” (see numpy.argsort implementation)

  2. if A is empty, returns A and hA unchanged

  3. if hA is empty, returns A and hA unchanged

  4. if A and hA have inconsistent shapes, raises ValueError

  5. if A is not 3D array, raises ValueError

  6. if hA is not 1D array, raises ValueError

  7. if A.shape[2] != hA.shape[0], raises ValueError

  8. if hA is already sorted, returns A and hA unchanged

Examples

>>> import numpy as np
>>> from tdcpy.common.compress import sort_matrices_delays
>>> A0 = np.array([[0, 1],[1, 0]])
>>> A1 = np.array([[1, 0],[0, 0]])
>>> A2 = np.array([[0, 0],[0, 1]])
>>> A = np.stack((A0, A1, A2), axis=2)
>>> hA = np.array([0., 2., 0.])
>>> sort_matrices_delays(A, hA)
(array([[[0, 0, 1],
        [1, 0, 0]],

        [[1, 0, 0],
        [0, 1, 0]]]), array([0., 0., 2.]))
tdcpy.common.compress.compress_ddae(E, A, hA, B, hB, C, hC, D, hD, rtol=1e-05, atol=1e-08)#

Compresses DDAE representation

Parameters:

tdcpy.common.delay_difference_equation module#

Set of functions for obtaining and manipulation of delay difference equations#

Implemented functions:

  1. ddae_to_diff: extracts the delay-difference equation represented by (D,hD) from the delay-differential algebraic equation represented by (E,A,hA)

  2. ndde_to_diff: extracts the associated delay-difference equation represented by (D,hD) given an NDDE represented by matrices (H,hH)

  3. normalize_diff: obtained the normalized ADDE with the leading matrix corresponding to the zero-delay term equal to the identity matrix

Notes:

  1. these functions are internal, they do not operate via high level API

  2. these functions assume correct inputs, input types, etc. (that is to save computation time), use these function only if you know what you are doing

tdcpy.common.delay_difference_equation.ddae_to_diff(E, A, hA, uE=None, vE=None, **kwargs)#

Converts delay differential algebraic equation (DDAE) to delay difference equation

DDAE dynamics represented by

\[E \dot{x}(t) = A_0 x(t - h_{A,0}) + ... + A_{m_A} x(t - h_{A,m_A})\]

into delay difference equation represented by

\[D_0 x(t - h_{D,0}) + ... + D_{m_D} x(t - h_{D,m_D}) = 0\]
Parameters:
  • E (array) – left-side matrix of shape (n,n)

  • A (array) – right-side matrices of shape (n,n,m)

  • hA (array) – vector of delays of shape (m,)

  • uE (array, optional) – orthonormal basis for left null space of E, optional, default None means uE will be calculated via SVD

  • vE (array, optional) – orthonormal basis for right null space of E, optional, default None means uE will be calculated via SVD

  • **kwargs – tol: norm tolerance for considering matrix vanish, default 1e-14 rcond (float): relative condition number. Singular values s smaller than rcond * max(s) are considered zero in null space construction, default 1e-12 is_compressed: if True, hA is assumed to be in compressed form, i.e. hA[0] == 0, default True

Returns:

A tuple containing

Darray

right-side delay difference equation matrices of shape (p,p,q)

hDarray

right-side delay difference equation delays of shape (q,)

Return type:

tuple

Notes

  1. n, m are assumed to be > 1

  2. A, hA is assumed to be in compressed form, i.e. hA[0] == 0

  3. if E is non-singular, D, hD are returned as empty arrays

  4. D, hD are returned in compressed form, i.e. no zero matrices in D (and associated delays in hD)

  5. if E is singular, the size of D, hD depends on the rank of E and the number of non-vanishing matrices uE.T @ A[i] @ vE

Examples

>>> import numpy as np
>>> from tdcpy.common.delay_difference_equation import ddae_to_diff
>>> E = np.array([[1,0,0],[0,0,0],[0,0,1]])
>>> A = np.zeros(shape=(3,3,3))
>>> A[:,:,0] = np.array([[0,1,0],[0,0,0],[0,0,0]])
>>> A[:,:,1] = np.array([[0,0,0],[1,0,0],[0,0,1]])
>>> A[:,:,2] = np.array([[0,0,1],[0,0,0],[0,1,0]])
>>> hA = np.array([0,1,2])
>>> D,hD = ddae_to_diff(E,A,hA)
>>> D
    array([], shape=(1, 1, 0), dtype=float64)
>>> hD
array([], dtype=int64)
tdcpy.common.delay_difference_equation.ndde_to_diff(H, hH, **kwargs)#

Converts NDDE to delay difference equation

Parameters:
  • H (array) – right-side matrices of shape (n,n,mH)

  • hH (array) – vector of delays of shape (mH,)

  • **kwargs – tol: norm tolerance for considering matrix vanish, default 1e-14

Returns:

A tuple containing

Darray

right-side delay difference equation matrices of shape (n,n,mH+1)

hDarray

right-side delay difference equation matrices and delays of shape (n,n,mH+1) and (mH+1,)

Return type:

tuple

Notes

For a NDDAE, the associated delay difference equation is given by

\[I x(t) + H_0 x(t-hH_0) + ... + H_{mH} x(t-hH_{mH}) = 0 (1)\]
  1. n, mH are assumed to be > 0

  2. hH is assumed to be non-zero delays

  3. D, hD are returned in compressed form, i.e. no zero matrices in D (and associated delays in hD)

  4. D, hD are returned in normalized form, i.e. D[0] == I is omitted

Examples

>>> import numpy as np
>>> from tdcpy.common.delay_difference_equation import ndde_to_diff
>>> H = np.zeros(shape=(2,2,2))
>>> H[:,:,0] = np.array([[0,1],[0,0]])
>>> H[:,:,1] = np.array([[0,0],[1,0]])
>>> hH = np.array([1,2])
>>> D,hD = ndde_to_diff(H,hH)
>>> D
array([[[1., 0., 0.],
        [0., 1., 0.]],

        [[0., 0., 1.],
         [1., 0., 0.]]])
>>> hD
array([0, 1, 2])
tdcpy.common.delay_difference_equation.normalize_diff(D, hD)#

Normalizes delay difference equation

Transforms the delay difference equation such that the leading zero delay matrix D[0] equals identity (and can be omitted).

Parameters:
  • D (array) – right-side matrices of shape (n,n,m)

  • hD (array) – vector of delays of shape (m,)

Returns:

A tuple containing

DDarray

3D array representing matrices: [inv(D[0])*D[1], … , inv(D[0])*D[m-1]]

hDDarray

array of non-zero delays of shape (m-1,)

Return type:

tuple

Notes

  1. m > 2 is assumed

  2. D[0] is assumed to be invertible

  3. hD[0] == 0

Examples

>>> import numpy as np
>>> from tdcpy.common.delay_difference_equation import normalize_diff
>>> D = np.zeros(shape=(2,2,3))
>>> D[:,:,0] = np.array([[1,0],[0,1]])
>>> D[:,:,1] = np.array([[0,1],[1,0]])
>>> D[:,:,2] = np.array([[1,1],[0,0]])
>>> hD = np.array([0,1,2])
>>> DD,hDD = normalize_diff(D,hD)
>>> DD
array([[[0., 1.],
        [1., 1.]],

        [[1., 0.],
         [0., 0.]]])
>>> hDD
array([1, 2])

tdcpy.common.discretization module#

TODO Set of functions for discretizing a DDAE into a DAE —————————————————

Implemented functions:

  1. discretize_ddae: returns a discretized DAE represented by (E*,A*) given a ddae represented by (E,A,hA)

tdcpy.common.discretization.discretize_ddae(E, A, hA, discretization, s0=0j, method='cheb')#

Discretizes DDAE into DAE

DDAE of form:

\[E \dot{x}(t) = A_0 x(t) + A_1 x(t - h_{A,1}) + ... + A_{m} x(t - h_{A,m}),\]

is discretized into DAE of form:

\[E \dot{x}(t) = A x(t),\]

When method == ‘cheb’, the code uses the companion-type reformulation of the spectral discretisaion of the infinitesimal generator of the solution operator underlying the DDE as presented in [1] Section 2.2.

When method == ‘legendre’, the code uses a similar companion-type reformulation but now based on Legendre polynomials instead of Chebyshev polynomials.

Parameters:
  • E (array) – right hand side matrix of DDAE, assumed non-empty

  • A (array) – left hand side matrices of DDAE, assumed non-empty

  • hA (array) – vector of delays, assumed non-empty, hA[0] == 0

  • discretization (int) – degree of discretization > 0

  • s0 (complex, optional) – point discretization is done around, default 0

  • method (str, optional) – type of approximation, default ‘cheb’, allowed ‘cheb’, ‘legendre’

Returns:

A tuple containing

Earray

left hand-side matrix of DAE

Aarray

right hand-side matrix of DAE

Return type:

tuple

Notes

  1. assumes hA is in compressed form, i.e. hA[0] == 0

  2. if s0 != 0, the matrices A are shifted accordingly before and after discretization

  3. if method is not recognized, raises ValueError

References

Examples

>>> import numpy as np
>>> from tdcpy.common.discretization import discretize_ddae
>>> E = np.array([[1,0],[0,1]])
>>> A = np.zeros(shape=(2,2,2))
>>> A[:,:,0] = np.array([[0,1],[0,0]])
>>> A[:,:,1] = np.array([[0,0],[1,0]])
>>> hA =  np.array([0,1])
>>> E_dae, A_dae = discretize_ddae(E,A,hA,discretization=2, method="cheb")
>>> E_dae
array([[ 0.5  ,  0.   ,  0.   ,  0.   , -0.25 , -0.   ],
       [ 0.   ,  0.5  ,  0.   ,  0.   , -0.   , -0.25 ],
       [ 0.   ,  0.   ,  0.125,  0.   ,  0.   ,  0.   ],
       [ 0.   ,  0.   ,  0.   ,  0.125,  0.   ,  0.   ],
       [ 1.   ,  0.   ,  1.   ,  0.   ,  1.   ,  0.   ],
       [ 0.   ,  1.   ,  0.   ,  1.   ,  0.   ,  1.   ]])
>>> A_dae
array([[ 0.+0.j,  0.+0.j,  1.+0.j,  0.+0.j,  0.+0.j,  0.+0.j],
       [ 0.+0.j,  0.+0.j,  0.+0.j,  1.+0.j,  0.+0.j,  0.+0.j],
       [ 0.+0.j,  0.+0.j,  0.+0.j,  0.+0.j,  1.+0.j,  0.+0.j],
       [ 0.+0.j,  0.+0.j,  0.+0.j,  0.+0.j,  0.+0.j,  1.+0.j],
       [ 0.+0.j,  1.+0.j,  0.+0.j,  1.+0.j,  0.+0.j,  1.+0.j],
       [ 1.+0.j,  0.+0.j, -1.+0.j,  0.+0.j,  1.+0.j,  0.+0.j]])

tdcpy.common.quasipoly module#

Set of functions for necessary quasipolynomial manipulation#

Implemented functions:

  1. compress_qp: converts a quasipolynomial to minimal form

  2. qp_to_ndde: returns an ndde represented by (H,hH,A,hA) given a quasipolynomial of the representation (coeffs,delays)

tdcpy.common.quasipoly.sort_qp()#

I do not think this function is needed, but will implement it in future if necessary.

tdcpy.common.quasipoly.compress_qp(coefs, delays, atol=None, rtol=None)#

Converts quasipolynomial to minimal form

Minimal form of QP: No zero rows in coefs, last coef column is not zero column, delays are sorted in ascending order.

Quasipolynomial is represented via matrix coefs of shape (m, n+1), and vector of delays delays of shape (m), where the resulting quasipolynomial is defined as:

\[QP(s) = \sum\limits_{i=0}^{m-1} exp(-delays[i]*s) \sum\limits_{j=0}^{n} coefs[i,j] * s^j\]
Parameters:
  • coeffs (array) – matrix definition of polynomial coefficients (each row represents polynomial coefficients corresponding to delay)

  • delays (array) – vector definition of associated delays (each delay corresponds to row in coefs)

  • atol (float, optional) – absolute tolerance for determining if coefficient is sufficiently close to zero, default None, see numpy.isclose

  • rtol (float, optional) – relative tolerance for determining if coefficient is sufficiently close to zero, default None, see numpy.isclose

  • coefs (ndarray[tuple[Any, ...], dtype[_ScalarT]])

Returns:

A tuple containing

new_coefsarray

matrix definition of polynomial coefficients

new_delaysarray

vector definition of associated delays

Return type:

tuple

Notes

  1. if all coefficients are close to zero, returns empty arrays

  2. if coefs is empty, returns coefs and delays unchanged

  3. if delays is empty, returns coefs and delays unchanged

  4. if coefs and delays have inconsistent shapes, raises ValueError

  5. if rtol or atol are negative, raises ValueError

  6. if coefs is not 2D array, raises ValueError

  7. if delays is not 1D array, raises ValueError

  8. if coefs.shape[0] != delays.shape[0], raises ValueError

Examples

>>> import numpy as np
>>> from tdcpy.common.quasipoly import compress_qp
>>> coefs = np.array([[0.0, 1.0], [0.0, 2.0], [0.0, 0.0]])
>>> delays = np.array([0.0, 1.0, 2.0])
>>> new_coefs, new_delays = compress_qp(coefs, delays)
>>> new_coefs
array([[0., 1.],
       [0., 2.]])
>>> new_delays
array([0., 1.])
>>> coefs = np.array([[0.0, 0.0], [0.0, 0.0]])
>>> delays = np.array([0.0, 1.0])
>>> new_coefs, new_delays = compress_qp(coefs, delays)
>>> new_coefs
array([], shape=(0, 2), dtype=float64)
>>> new_delays
array([], dtype=float64)
tdcpy.common.quasipoly.qp_to_ndde(coefs, delays, ascending=True)#

Converts quasipolynomial into neutral delay differential equation

Converts quasipolynomial defined via coefs and delays

\[QP(s) = \sum\limits_{i=0}^{m-1} exp(-delays[i]*s) \sum\limits_{j=0}^{n} coefs[i,j] * s^j\]

into NDDE represented via arrays A, hA, H, hH

\[\dot{x}(t) = A_0*x(t - hA_0) + ... + A_{mA}*x(t - hA_{mA}) - H_0* \dot{x}(t - hH_0) - ... - H_{mH}* \dot{x}(t - hH_{mH})\]

with mA number of delays associated with A and mH number of delays associated with H.

Parameters:
  • coeffs (array) – matrix definition of polynomial coefficients (each row represents polynomial coefficients corresponding to delay)

  • delays (array) – vector definition of associated delays (each delay corresponds to row in coefs)

  • ascending (bool, optional) – ordering of powers of s in each row, default ascending meaning that coefs[i,j] is associated to ith polynomial and jth power of s, setting this to False will default to original MATLAB behaviour, where coefs[i,j] is associated to ith polynomial and (n-j)th power of s

Returns:

A tuple containing

Aarray

matrices defining delay differential equation, with shape (n,n,mA)

hAarray

vector of delays associated with A

Harray

matrices defining delay difference equation, with shape (n,n,mH)

hHarray

vector of delays associated with H

Return type:

tuple

Notes

  1. if all coefficients are close to zero, returns empty arrays

  2. if coefs is empty, returns empty arrays

  3. if delays is empty, returns empty arrays

  4. if coefs and delays have inconsistent shapes, raises ValueError

  5. if system is of advanced type (delay[0] != 0.0 or coefs[0,-1] == 0.0), raises ValueError

Examples

>>> import numpy as np
>>> from tdcpy.common.quasipoly import qp_to_ndde
>>> coefs = np.array([[1.0, 0.0], [0.0, 2.0]])
>>> delays = np.array([0.0, 1.0])
>>> A, hA, H, hH = qp_to_ndde(coefs, delays, ascending=False)
>>> A
array([[[-2.]]])
>>> hA
array([1.])
>>> H
array([], shape=(1, 1, 0), dtype=float64)
>>> hH
array([], dtype=float64)
>>> coefs = np.array([[0.0, 1.0], [0.0, 2.0], [0.0, 3.0]])
>>> delays = np.array([0.0, 1.0, 2.0])
>>> A, hA, H, hH = qp_to_ndde(coefs, delays)
>>> A
array([[[-0., -0.]]])
>>> A.shape
(1, 1, 2)
>>> hA
array([1., 2.])
>>> H
array([[[2., 3.]]])
>>> H.shape
(1, 1, 2)
>>> hH
array([1., 2.])

Module contents#