Skip to contents

Compute the poles and zeroes of a rational matrix. For polynomial matrices and rational matrices in left matrix fraction form the poles (and zeroes) are computed via the (reciprocals of the) eigenvalues of the associated companion matrices, see also companion_matrix. For statespace realizations the poles are computed via (the reciprocals of) the eigenvalues of the state transition matrix \(A\) and for the zeroes the eigenvalues of the state transition matrix \(A-BD^{-1}C\) of the inverse are used.

Usage

poles(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for polm
zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for lpolm
zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for lmfd
zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for rmfd
zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for stsp
zeroes(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for polm
poles(x, ...)

# S3 method for lmfd
poles(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for rmfd
poles(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

# S3 method for stsp
poles(x, tol = sqrt(.Machine$double.eps), print_message = TRUE, ...)

Arguments

x

an object which represents a rational matrix (i.e. a polm, lmfd or stsp object).

tol

Double. Default set to sqrt(.Machine$double.eps). Required to decide on when a root is to be considered "at infinity".

print_message

Boolean. Default set to TRUE. Prints a message if roots "at infinity " are discarded.

...

not used.

Value

Vector of poles, respectively zeroes.

Details

The methods do not return numerically reliable and correct results in all cases. For some more details see the vignette Rational Matrices.

  • Zeroes are only computed for square, non singular matrices which have no zero at \(z=0\). If the matrix evaluated at \(z=0\) is close to singular, the results may be unreliable.

  • The procedures use a threshold tol in order to decide whether a small eigenvalue returned by eigen corresponds to a "true zero" eigenvalue or not.

  • If the pair \(a,b\) of polynomials of the LMFD is not left coprime then a pole/zero cancellation occurs. This is not taken into account by the procedures. Hence, in this case, the results also contain some spurious poles/zeroes. This happens also for non minimal state space realizations.

Examples

 
# zeroes of polynomial matrices #############################################

# scalar polynomial ###
(a = polm(c(1, 0, 0, 0.5, 0)))
#> ( 1 x 1 ) matrix polynomial with degree <= 4 
#>      z^0 [,1] z^1 [,1] z^2 [,1] z^3 [,1] z^4 [,1]
#> [1,]        1        0        0      0.5        0
(z = zeroes(a))
#> [1] -1.259921+0.000000i  0.629961-1.091124i  0.629961+1.091124i

# compare with the result of "polyroot"
all.equal(sort(z), sort(polyroot(as.vector(a))))
#> [1] TRUE

# zero degree polynomial (have no zeroes) ###
zeroes(polm(diag(3)))
#> numeric(0)

# (2 x 2) polynomial of degree 2 ### 
a = polm(dbind(d = 3, diag(2), test_array(dim = c(2,2,2))))
(z = zeroes(a))
#> [1] -0.002994225  0.265114779 -0.991220211 -1.270900343

# check the rank of a(z) at the computed zeroes 
az = zvalues(a, z)
apply(az, MARGIN = 3, FUN = function(x) {d = svd(x)$d; min(d)/max(d)})
#> [1] 1.403359e-16 8.297615e-17 1.374392e-14 0.000000e+00

if (FALSE) {
# the following examples throw an error
zeroes(polm(c(0, 0, 0, 0.5))) # constant term is zero
zeroes(polm(test_array(dim = c(2, 1, 3)))) # non-square polynomial
zeroes(polm(test_array(dim = c(2, 2, 0)))) # zero polynomial
}
 
# zeroes of a Laurent polynomial #################################
(lp = lpolm(1:5, min_deg = -7))
#> ( 1 x 1 ) Laurent polynomial matrix with degree <= -3, and minimal degree >= -7
#>      z^-7 [,1] z^-6 [,1] z^-5 [,1] z^-4 [,1] z^-3 [,1]
#> [1,]         1         2         3         4         5
(p = polm(1:5))
#> ( 1 x 1 ) matrix polynomial with degree <= 4 
#>      z^0 [,1] z^1 [,1] z^2 [,1] z^3 [,1] z^4 [,1]
#> [1,]        1        2        3        4        5
zeroes(p)
#> [1] -0.5378323-0.3582847i -0.5378323+0.3582847i  0.1378323-0.6781544i
#> [4]  0.1378323+0.6781544i
zeroes(lp)
#> [1] -0.5378323-0.3582847i -0.5378323+0.3582847i  0.1378323-0.6781544i
#> [4]  0.1378323+0.6781544i
 
# zeroes of a rational matrix in LMFD form #################################

c = lmfd(test_polm(dim = c(2,2), degree = 3, random = TRUE),
         test_polm(dim = c(2,2), degree = 1, random = TRUE))
(z = zeroes(c))
#> [1] -0.2049177  1.2448546
all.equal(z, zeroes(c$b))
#> [1] TRUE
 
# zeroes of a rational matrix in RMFD form #################################

k = rmfd(c = test_polm(dim = c(2,2), degree = 3, random = TRUE),
         d = test_polm(dim = c(2,2), degree = 1, random = TRUE))
(z = zeroes(k))
#> [1]   0.6372353 -11.4740943
all.equal(z, zeroes(k$d))
#> [1] TRUE
 
# zeroes of a rational matrix in statespace form ###########################

k = stsp(A = matrix(rnorm(3*3), nrow = 3, ncol = 3),
         B = matrix(rnorm(3*2), nrow = 3, ncol = 2),
         C = matrix(rnorm(3*2), nrow = 2, ncol = 3),
         D = matrix(rnorm(2*2), nrow = 2, ncol = 2))
(z = zeroes(k, tol = 0))
#> [1] 0.1837212+0.0000000i 0.4998418-0.3803373i 0.4998418+0.3803373i
all.equal(z, 1/(eigen(k$A - k$B %*% solve(k$D, k$C), only.values = TRUE)$values))
#> [1] TRUE

if (FALSE) {
k = stsp(k$A, k$B, k$C, 
         D = matrix(rnorm(2*1), nrow = 2, ncol = 1)[,c(1,1)])  # D is singular
zeroes(k)                                                      # zeroes() throws an error

k = stsp(k$A, k$B[,1,drop = FALSE], k$C, k$D[,1,drop = FALSE]) # (2 x 1) rational matrix
zeroes(k)                          # throws an error, since k is not square
}

# poles of polynomial matrices #############################################

# polynomials have no poles ###
poles(test_polm(dim = c(2,1), degree = 2, random = TRUE)) 
#> numeric(0)

# poles of a rational matrix in LMFD form ##################################

(z = poles(c))
#> [1]  0.2813151+0.0000000i -0.6982848-0.7111637i -0.6982848+0.7111637i
#> [4]  1.0177287+0.0000000i  0.3638302-1.6370207i  0.3638302+1.6370207i
all.equal(z, zeroes(c$a))
#> [1] TRUE

# poles of a rational matrix in RMFD form ##################################

(z = poles(c))
#> [1]  0.2813151+0.0000000i -0.6982848-0.7111637i -0.6982848+0.7111637i
#> [4]  1.0177287+0.0000000i  0.3638302-1.6370207i  0.3638302+1.6370207i
all.equal(z, zeroes(c$a))
#> [1] TRUE

# poles of a rational matrix in statespace form ###########################

(z = poles(k, tol = 0))
#> [1] 0.4729606+0.0000000i 0.7541721-0.9958197i 0.7541721+0.9958197i
all.equal(z, 1/(eigen(k$A, only.values = TRUE)$values))
#> [1] TRUE