Builtin Matrices

TypedMatrices.BinomialType

Binomial Matrix

The binomial matrix is a multiple of an involutory matrix.

Input Options

  • dim: the dimension of the matrix.

References

G. Boyd, C. A. Micchelli, G. Strang and D. X. Zhou, Binomial matrices, Adv. Comput. Math., 14 (2001), pp. 379-391, https://doi.org/10.1023/A:1012207124894.

source
TypedMatrices.CauchyType

Cauchy Matrix

Given two vectors x and y, the (i,j) entry of the Cauchy matrix is 1/(x[i]+y[j]).

Input Options

  • x, y: two vectors.
  • x, y: two integers, as vectors 1:x and 1:y`.
  • x: an integer, as vectors 1:xand1:x``.
  • x: a vector. y defaults to x.

References

N. J. Higham, Accuracy and Stability of Numerical Algorithms, second edition, Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2002, https://doi.org/10.1137/1.9780898718027. See sect. 28.1.

source
TypedMatrices.ChebSpecType

Chebyshev Spectral Differentiation Matrix

The Chebyshev Spectral Differentiation Matrix is used to approximate numerically the derivatives of a function evaluated at Chebyshev nodes.

If k = 0,the generated matrix is nilpotent and a vector with all one entries is a null vector. If k = 1, the generated matrix is nonsingular and well-conditioned. Its eigenvalues have negative real parts.

Input Options

  • dim, k: dim is the dimension of the matrix and k is either 0 or 1.
  • dim: k=0.

References

L. N. Trefethen and M. R. Trummer, An instability phenomenon in spectral methods, SIAM J. Numer. Anal., 24 (1987), pp. 1008-1023, https://doi.org/10.1137/0724066.

source
TypedMatrices.ChowType

Chow Matrix

The Chow matrix is a singular Toeplitz lower-Hessenberg matrix.

Input Options

  • dim, alpha, delta: dim is dimension of the matrix. alpha, delta are scalars such that A[i,i] = alpha + delta and A[i,j] = alpha^(i - j + 1) for j + 1 <= i.
  • dim: alpha = 1, delta = 0.

References

T. S. Chow, A class of Hessenberg matrices with known eigenvalues and inverses, SIAM Rev., 11 (1969), pp. 391-395, https://doi.org/10.1137/1011065.

source
TypedMatrices.CirculantType

Circulant Matrix

A circulant matrix has the property that each row is obtained by cyclically permuting the entries of the previous row one step forward.

Input Options

  • vec: a vector.
  • dim: an integer, as vector 1:dim.

References

P. J. Davis, Circulant Matrices, John Wiley, 1977.

source
TypedMatrices.ClementType

Clement Matrix

The Clement matrix is a tridiagonal matrix with zero diagonal entries. If k = 1, the matrix is symmetric.

Input Options

  • dim, k: dim is the dimension of the matrix. If k = 0, the matrix is of type Tridiagonal. If k = 1, the matrix is of type SymTridiagonal.
  • dim: k = 0.

References

P. A. Clement, A class of triple-diagonal matrices for test purposes, SIAM Rev., 1 (1959), pp. 50-52, https://doi.org/10.1137/1001006.

source
TypedMatrices.CompanionType

Companion Matrix

The companion matrix to the monic polynomial

a(x) = a_0 + a_1x + ... + a_{n-1}x^{n-1} + x^n

is the n-by-n matrix with ones on the first subdiagonal and the last column given by the coefficients of a(x).

Input Options

  • vec: vec is a vector of coefficients.
  • dim: vec = [1:dim;]. dim is the dimension of the matrix.
  • polynomial: polynomial is a polynomial. Last column will contain

its coefficients.

References

N. J. Higham, What Is the Companion Matrix?, https://nhigham.com/2021/03/23/what-is-a-companion-matrix/

source
TypedMatrices.ComparisonType

Comparison Matrix

The comparison matrix of a given matrix.

Input Options

  • B, k: B is a matrix. If k = 0, then C(i,j) = abs(B(i,j)) for i ≠ j and C(i,i) = -abs(B(i,i)). If k = 1, then C(i,i) = abs(B(i,i)) and C(i,j) = -max(abs(B(i,:))) for i ≠ j.
  • B: B is a matrix and k = 1.

N. J. Higham, Efficient algorithms for computing the condition number of a tridiagonal matrix, SIAM J. Sci. Stat. Comput., 7 (1986), pp. 150-165, https://doi.org/10.1137/0907011.

source
TypedMatrices.CycolType

Cycol Matrix

This matrix has columns that repeat cyclically.

Input Options

  • m, n, k: m and n are size of the matrix. The repeated columns are randn(m, k).
  • n, k: n is size of the matrix. The repetition is randn(n, k).
  • n: n is size of the matrix. k = round(n/4)
source
TypedMatrices.DingDongType

Dingdong Matrix

The Dingdong matrix is a Hankel matrix due to F. N. Ris of IBM Thomas J. Watson Research Centre. The eigenvalues cluster around π/2 and -π/2.

Input Options

  • dim: the dimension of the matrix.

References

J. C. Nash, Compact Numerical Methods for Computers: Linear Algebra and Function Minimisation, second edition, Adam Hilger, Bristol, 1990 (Appendix 1).

source
TypedMatrices.DorrType

Dorr Matrix

The Dorr Matrix is a diagonally dominant, ill-conditioned, tridiagonal sparse matrix.

Input Options

  • dim, theta: dim is the dimension of the matrix and theta is the parameter of the matrix.
  • dim: theta = 0.01.

References

F. W. Dorr, An example of ill-conditioning in the numerical solution of singular perturbation problems, Math. Comp., 25 (1971), pp. 271-283, https://doi.org/10.1090/S0025-5718-1971-0297142-0.

source
TypedMatrices.DramadahType

Dramadah Matrix

The Dramadah matrix is a matrix of 0s and 1s whose inverse has comparatively large Frobenius norm.

Input Options

  • dim, k: the dimension of the matrix and k. k = 1 abs(det(A)) = 1, the inverse has integer entries. k = 2 the inverse has integer entries. k = 3 det(A) is equal to nth Fibonacci number.
  • dim: k = 1.

References

R.L. Graham and N.J.A. Sloane, Anti-Hadamard matrices, Linear Algebra and Appl., 62 (1984), pp. 113-137. https://doi.org/10.1016/0024-3795(84)90090-9

source
TypedMatrices.FiedlerType

Fiedler Matrix

The Fiedler matrix has exactly one positive eigenvalue, the dominant one. All the other eigenvalues are negative.

Input Options

  • vec: a vector.
  • dim: dim is the dimension of the matrix. vec=[1:dim;].

References

A. C. Schaeffer and G. Szegö, Solution to problem 3705, Amer. Math. Monthly, 43 (1936), pp. 246-259, https://doi.org/10.1090/S0002-9947-1941-0005164-7.

J. Todd, Basic Numerical Mathematics, Vol. 2: Numerical Algebra, Birkhauser, Basel, and Academic Press, New York, 1977, p. 159.

source
TypedMatrices.ForsytheType

Forsythe Matrix

The Forsythe matrix is an n-by-n perturbed Jordan block. This generator is adapted from N. J. Higham's Test Matrix Toolbox.

Input Options

  • dim, alpha, lambda: dim is the dimension of the matrix. alpha and lambda are scalars.
  • dim: alpha = sqrt(eps(type)) and lambda = 0.
source
TypedMatrices.FrankType

Frank Matrix

The Frank matrix is an upper Hessenberg matrix with determinant 1. The eigenvalues are real, positive, and very ill conditioned.

Input Options

  • dim, k: dim is the dimension of the matrix, k = 0 or 1. If k = 1 the matrix reflect about the anti-diagonal.
  • dim: the dimension of the matrix.

References

W. L. Frank, Computing eigenvalues of complex matrices by determinant evaluation and by methods of Danilewski and Wielandt, J. Soc. Indust. Appl. Math., 6 (1958), pp. 378-392, https://doi.org/10.1137/0106026. See pp. 385 and 388.

source
TypedMatrices.GCDMatType

GCDMat Matrix

A matrix whose (i,j) entry is gcd(i,j). It is a symmetric positive definite matrix.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.GearMatType

Gear Matrix

The Gear matrix has ones on the first subdiagonal and superdiagonal, and has two additional entries of value ±1. Given the two integers -n ≤ i ≤ n and -n ≤ j ≤ n, the matrix has the elements sign(i) in position (1, abs(i)) and sign(j) in position (n, n+1-abs(j)). The other elements are zeros.

Input Options

  • dim, i, j: the dimension of the matrix and the position of the 1s.
  • dim: the dimension of the matrix. i = n and j = -n by default.

References

C. W. Gear, A simple set of test matrices for eigenvalue programs, Math. Comp., 23 (1969), pp. 119-125, https://doi.org/10.1090/S0025-5718-1969-0238477-8.

source
TypedMatrices.GolubType

Golub Matrix

The Golub matrix is the product of two random matrices, the first is unit lower triangular and the second is upper triangular. The LU factorization without pivoting fails to reveal that such matrices are badly conditioned.

Input Options

  • dim: the dimension of the matrix.

References

D. Viswanath and N. Trefethen. Condition numbers of random triangular matrices, SIAM J. Matrix Anal. Appl., 19 (1998), 564-581, https://doi.org/10.1137/S0895479896312869.

source
TypedMatrices.GrcarType

Grcar Matrix

The Grcar matrix is a Toeplitz matrix with sensitive eigenvalues.

Input Options

  • dim, k: dim is the dimension of the matrix and k is the number of superdiagonals.
  • dim: the dimension of the matrix.

References

J. F. Grcar, Operator coefficient methods for linear equations, Report SAND89-8691, Sandia National Laboratories, Albuquerque, New Mexico, 1989 (Appendix 2).

source
TypedMatrices.HadamardType

Hadamard Matrix

The Hadamard matrix is a square matrix of order a power of 2, whose entries are 1 or –1. It was named after Jacques Hadamard. The rows of a Hadamard matrix are orthogonal.

Input Options

  • dim: the dimension of the matrix, dim is a power of 2.

References

S. W. Golomb and L. D. Baumert, The search for Hadamard matrices, Amer. Math. Monthly, 70 (1963) pp. 12-17, https://doi.org/10.1080/00029890.1963.11990035.

source
TypedMatrices.HankelType

Hankel Matrix

A Hankel matrix is constant across the anti-diagonals. It is symmetric.

Input Options

  • vc, vr: vc and vc are the first column and last row of the matrix. If the last element of vc differs from the first element of vr, the last element of rc prevails.
  • v: a vector, as vc = v and vr will be zeros.
  • dim: dim is the dimension of the matrix. v = [1:dim;].
source
TypedMatrices.HanowaType

Hanowa Matrix

The Hanowa matrix is a matrix whose eigenvalues lie on a vertical line in the complex plane.

Input Options

  • dim: the dimension of the matrix and alpha = -1.
  • dim, alpha: the dimension and alpha.
source
TypedMatrices.HilbertType

Hilbert Matrix

The Hilbert matrix has (i,j) element 1/(i+j-1). It is notorious for being ill conditioned. It is symmetric positive definite and totally positive.

See also InverseHilbert.

Input Options

  • dim: the dimension of the matrix.
  • row_dim, col_dim: the row and column dimensions.

References

M. D. Choi, Tricks or treats with the Hilbert matrix, Amer. Math. Monthly, 90 (1983), pp. 301-312, https://doi.org/10.1080/00029890.1983.11971218.

N. J. Higham, Accuracy and Stability of Numerical Algorithms, second edition, Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2002. See sect. 28.1.

source
TypedMatrices.InverseHilbertType

Inverse of the Hilbert Matrix

This is the inverse of the Hilbert matrix.

See also Hilbert.

Input Options

  • dim: the dimension of the matrix.

References

M. D. Choi, Tricks or treats with the Hilbert matrix, Amer. Math. Monthly, 90 (1983), pp. 301-312, https://doi.org/10.1080/00029890.1983.11971218.

N. J. Higham, Accuracy and Stability of Numerical Algorithms, second edition, Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2002. See sect. 28.1.

source
TypedMatrices.InvhessType

Invhess Matrix

This matrix is thenverse of an upper Hessenberg matrix.

Input Options

  • dim: the dimension of the matrix. x = [1:dim;].
  • x: x vector. y = -x[1:end-1].
  • x, y: x and y vectors.
source
TypedMatrices.InvolutoryType

Involutory Matrix

An involutory matrix is a matrix that is its own inverse.

Input Options

  • dim: dim is the dimension of the matrix.

References

A. S. Householder and J. A. Carpenter, The singular values of involutory and of idempotent matrices, Numer. Math. 5 (1963), pp. 234-237, https://doi.org/10.1007/BF01385894.

source
TypedMatrices.IpjfactType

Ipjfact Matrix

Hankel matrix with factorial elements.

Input Options

  • dim: the dimension of the matrix.
  • dim, k: k = 0 element (i, j) is factorial(i + j). k = 1 element (i, j) is 1 / factorial(i + j).

References

**K. Habermann*, An explicit formula for the inverse of a factorial Hankel matrix, Australasian J. Comb., 79 (2021), pp. 250-255. https://ajc.maths.uq.edu.au/pdf/79/ajcv79p250.pdf

source
TypedMatrices.JordBlocType

Jordan Block Matrix

Jordan block corresponding to the eigenvalue λ.

Input Options

  • dim: dimension of the matrix. lambda = 1.
  • dim, lambda: dimension of the matrix and the lambda.
source
TypedMatrices.KahanType

Kahan Matrix

The Kahan matrix is an upper trapezoidal matrix, i.e., the (i,j) element is equal to 0 if i > j. The useful range of θ is 0 < θ < π.

The diagonal is perturbed by pert*eps()*diagm([n:-1:1;]).

Input Options

  • rowdim, coldim, θ, pert: rowdim and coldim are the row and column dimensions of the matrix. θ and pert are scalars.
  • dim, θ, pert: dim is the dimension of the matrix.
  • dim: θ = 1.2, pert = 25.

References

W. Kahan, Numerical linear algebra, Canadian Math. Bulletin, 9 (1966), pp. 757-801, https://doi.org/10.4153/CMB-1966-083-2.

source
TypedMatrices.KMSType

Kac-Murdock-Szego Toeplitz matrix

Input Options

  • dim, rho: dim is the dimension of the matrix, rho is a scalar such that A[i,j] = rho^(abs(i-j)).
  • dim: rho = 0.5.

References

W. F. Trench, Numerical solution of the eigenvalue problem for Hermitian Toeplitz matrices, SIAM J. Matrix Anal. Appl., 10 (1989), pp. 135-146, https://doi.org/10.1137/0610010.

N. J. Higham, What Is the Kac-Murdock-Szegö Matrix?, https://nhigham.com/2021/07/06/what-is-the-kac-murdock-szego-matrix/

source
TypedMatrices.KrylovType

Krylov Matrix

The basis of a Krylow subspace. The matrix has columns [x, A*x, A^2*x, ..., A^(k-1)*x].

Input Options

  • dim: dimension of the matrix. A = randn(dim, dim). x = ones(dim). k = dim.
  • dim, x: dimension of the matrix and x.
  • dim, x, k: dimension of the matrix, x and k.
  • A: matrix. x = ones(size(A, 1)). k = size(A, 1).
  • A, x: matrix and x. k = size(A, 1).
  • A, x, k: matrix, x, and k.
source
TypedMatrices.LauchliType

Lauchli Matrix

A matrix with ones on the first row, mu on the subdiagonal, and zeros elsewhere.

Input Options

  • dim: the dimension of the matrix. mu = sqrt(eps()) by default.
  • dim, mu: the dimension and subdiagonal value of the matrix.

References

P. Lauchli, Jordan-Elimination und Ausgleichung nach kleinsten Quadraten, Numer. Math, 3 (1961), pp. 226-240. https://doi.org/10.1007/BF01386022

source
TypedMatrices.LehmerType

Lehmer Matrix

The Lehmer matrix is a symmetric positive definite matrix. It is totally nonnegative. The inverse is tridiagonal and explicitly known.

Input Options

  • dim: the dimension of the matrix.

References

M. Newman and J. Todd, The evaluation of matrix inversion programs, J. Soc. Indust. Appl. Math., 6 (1958), pp. 466-476, https://doi.org/10.1137/0106030.

D. H. Lehmer, Problem E710: The inverse of a matrix, Amer. Math. Monthly, 53 (1946), p. 97, https://doi.org/10.2307/2305463.

Solutions by D. M. Smiley and M. F. Smiley, and by John Williamson. Amer. Math. Monthly, 53 (1946), pp. 534-535, https://doi.org/10.2307/2305078.

source
TypedMatrices.LeslieType

Leslie Matrix

Matrix for birth numbers and survival rates in the Leslie population model.

Input Options

  • dim: the dimension of the matrix. x = ones(n) and y = ones(n - 1) by default.
  • x, y: x and y.
source
TypedMatrices.LespType

Lesp Matrix

A matrix with eigenvalues smoothly distributed in the interval [-2*n-3.5,-4.5].

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.LotkinType

Lotkin Matrix

The Lotkin matrix is the Hilbert matrix with its first row altered to all ones. It is ill conditioned and has many negative eigenvalues of small magnitude.

Input Options

  • dim: dim is the dimension of the matrix.

References

M. Lotkin, A set of test matrices, Math. Tables Aid Comput., 9 (1955), pp. 153-161, https://doi.org/10.2307/2002051.

source
TypedMatrices.MagicType

Magic Square Matrix

The magic matrix is a matrix with integer entries such that the row elements, column elements, diagonal elements and anti-diagonal elements all add up to the same number.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.MinijType

MIN[I,J] Matrix

A matrix with (i,j) entry min(i,j). It is a symmetric positive definite matrix. The eigenvalues and eigenvectors are known explicitly. Its inverse is tridiagonal.

Input Options

  • dim: the dimension of the matrix.

References

J. Fortiana and C. M. Cuadras, A family of matrices, the discretized Brownian bridge, and distance-based regression, Linear Algebra Appl., 264 (1997), pp. 173-188, https://doi.org/10.1016/S0024-3795(97)00051-7.

source
TypedMatrices.MolerType

Moler Matrix

The Moler matrix is a symmetric positive definite matrix. It has one small eigenvalue.

Input Options

  • dim, alpha: dim is the dimension of the matrix, alpha is a scalar.
  • dim: alpha = -1.

References

J.C. Nash, Compact Numerical Methods for Computers: Linear Algebra and Function Minimisation, second edition, Adam Hilger, Bristol, 1990 (Appendix 1).

source
TypedMatrices.NeumannType

Neumann Matrix

A singular matrix from the discrete Neumann problem. The matrix is sparse and the null space is formed by a vector of ones

Input Options

  • dim: the dimension of the matrix, must be a perfect square integer.

References

R. J. Plemmons, Regular splittings and the discrete Neumann problem, Numer. Math., 25 (1976), pp. 153-161, https://doi.org/10.1007/BF01462269.

source
TypedMatrices.OrthogType

Orthogonal Matrix

Orthogonal and nearly orthogonal matrices.

Input Options

  • dim: the dimension of the matrix. k = 1 by default.
  • dim, k: the dimension and type of the matrix.
source
TypedMatrices.OscillateType

Oscillating Matrix

A matrix A is called oscillating if A is totally nonnegative and if there exists an integer q > 0 such that A^q is totally positive.

Input Options

  • Σ: the singular value spectrum of the matrix.
  • dim, mode: dim is the dimension of the matrix. mode = 1: geometrically distributed singular values. mode = 2: arithmetrically distributed singular values.
  • dim: mode = 1.

References

P. C. Hansen, Test matrices for regularization methods, SIAM J. Sci. Comput., 16 (1995), pp. 506-512, https://doi.org/10.1137/0916032. .

source
TypedMatrices.ParterType

Parter Matrix

The Parter matrix is a Toeplitz and Cauchy matrix with singular values near π.

Input Options

  • dim: the dimension of the matrix.

References

The MathWorks Newsletter, Volume 1, Issue 1, March 1986, page 2.

S. V. Parter, On the distribution of the singular values of Toeplitz matrices, Linear Algebra Appl., 80 (1986), pp. 115-130, https://doi.org/10.1016/0024-3795(86)90280-6.

source
TypedMatrices.PascalType

Pascal Matrix

The Pascal matrix’s anti-diagonals form the Pascal’s triangle.

Input Options

  • dim: the dimension of the matrix.

References

R. Brawer and M. Pirovino, The linear algebra of the Pascal matrix, Linear Algebra and Appl., 174 (1992), pp. 13-23 (this paper gives a factorization of L = PASCAL(N,1) and a formula for the elements of L^k).

N. J. Higham, Accuracy and Stability of Numerical Algorithms, second edition, Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2002. See sect. 28.4.

source
TypedMatrices.PeiType

Pei Matrix

The Pei matrix is a symmetric matrix with known inversion.

Input Options

  • dim, alpha: dim is the dimension of the matrix. alpha is a scalar.
  • dim: the dimension of the matrix.

References

M. L. Pei, A test matrix for inversion procedures, Comm. ACM, 5 (1962), p. 508, https://doi.org/10.1145/368959.368975.

source
TypedMatrices.PoissonType

Poisson Matrix

A block tridiagonal matrix from Poisson's equation. This matrix is sparse, symmetric positive definite, and has known eigenvalues.

Input Options

  • dim: the dimension of the matrix.

References

G. H. Golub and C. F. Van Loan, Matrix Computations, second edition, Johns Hopkins University Press, Baltimore, Maryland, 1989. See sect. 4.5.4.

source
TypedMatrices.ProlateType

Prolate Matrix

A prolate matrix is a symmetirc, ill-conditioned Toeplitz matrix.

Input Options

  • dim, alpha: dim is the dimension of the matrix. w is a real scalar.
  • dim: the case when w = 0.25.

References

J. M. Varah. The Prolate matrix. Linear Algebra Appl., 187 (1993), pp. 267-278, https://doi.org/10.1016/0024-3795(93)90142-B.

source
TypedMatrices.RandcoluType

Randcolu Matrix

Random matrix with normalized columns and given singular values.

Input Options

  • dim: the dimension of the matrix, x will be generated randomly.
  • n, m: the size of the matrix.
  • n, m, k: size and k flag. Enable initial transformation if k = 0.
  • x: the x vector.
  • x, m: the x vector and m.
  • x, m, k: the x vector, m, and k flag.
source
TypedMatrices.RandcorrType

Random Correlation Matrix

A random correlation matrix is a symmetric positive semidefinite matrix with 1s on the diagonal.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.RandoType

Random Matrix with Element -1, 0, 1

Input Options

  • rowdim, coldim, k: row_dim and col_dim are row and column dimensions, k = 1: entries are 0 or 1. k = 2: entries are -1 or 1. k = 3: entries are -1, 0 or 1.
  • dim, k: row_dim = col_dim = dim.
  • dim: k = 1.
source
TypedMatrices.RandSVDType

Random Matrix with Pre-assigned Singular Values

Input Options

  • rowdim, coldim, kappa, mode: row_dim and col_dim are the row and column dimensions. kappa is the condition number of the matrix. mode = 1: one large singular value. mode = 2: one small singular value. mode = 3: geometrically distributed singular values. mode = 4: arithmetrically distributed singular values. mode = 5: random singular values with unif. dist. logarithm.
  • dim, kappa, mode: row_dim = col_dim = dim.
  • dim, kappa: mode = 3.
  • dim: kappa = sqrt(1/eps()), mode = 3.

References

N. J. Higham, Accuracy and Stability of Numerical Algorithms, second edition, Society for Industrial and Applied Mathematics, Philadelphia, PA, USA, 2002. See sect. 28.3.

source
TypedMatrices.RedheffType

Redheffer Matrix

The Redheffer matrix contains only 1s and 0s.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.RiemannType

Riemann Matrix

A matrix associated with the Riemann hypothesis.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.RISType

RIS Matrix

The RIS matrix has (i,j) element 0.5/(n-i-j+1.5). It is symmetric.

Input Options

  • dim: the dimension of the matrix.
source
TypedMatrices.RohessType

Random Orthogonal Upper Hessenberg Matrix

The matrix is constructed via a product of Givens rotations.

Input Options

  • dim: the dimension of the matrix.

References

W. B. Gragg, The QR algorithm for unitary Hessenberg matrices, J. Comp. Appl. Math., 16 (1986), pp. 1-8, https://doi.org/10.1016/0377-0427(86)90169-X.

source
TypedMatrices.RosserType

Rosser Matrix

The Rosser matrix’s eigenvalues are very close together so it is a challenging matrix for many eigenvalue algorithms.

Input Options

  • dim, a, b: dim is the dimension of the matrix. dim must be a power of 2. a and b are scalars. For dim = 8, a = 2 and b = 1, the generated matrix is the test matrix used by Rosser.
  • dim: a = rand(1:5), b = rand(1:5).

References

J. B. Rosser, C. Lanczos, M. R. Hestenes, and W. Karush, Separation of close eigenvalues of a real symmetric matrix, J. Research National Bureau Standards, 47 (1951), pp. 291-297.

source
TypedMatrices.SamplingType

Matrix with Application in Sampling Theory

A nonsymmetric matrix with eigenvalues 0, 1, 2, ... n-1.

Input Options

  • vec: vec is a vector with no repeated elements.
  • dim: dim is the dimension of the matrix. vec = [1:dim;]/dim.

References

L. Bondesson and I. Traat, A nonsymmetric matrix with integer eigenvalues, Linear Multilinear Algebra, 55 (2007), pp. 239-247, https://doi.org/10.1080/03081080600906455.

source
TypedMatrices.SmokeType

Smoke Matrix

Complex matrix with a "smoke ring" pseudospectrum. The matrix has ones on the superdiagonal, and cos(w) + sin(w) * im on the diagonal. The A(n, 1)` entry is 1 if k = 0, 0 if k = 1.

Input Options

  • dim: dimension of the matrix. k = 0.
  • dim, k: dimension of the matrix and the k.
source
TypedMatrices.ToeplitzType

Toeplitz Matrix

A Toeplitz matrix is a matrix in which each descending diagonal from left to right is constant.

Input Options

  • vc, vr: vc and vr are the first column and row of the matrix.
  • v: symmatric case, i.e., vc = vr = v.
  • dim: dim is the dimension of the matrix. v = [1:dim;] is the first row and column vector.
source
TypedMatrices.TriwType

Triw Matrix

Upper triangular matrices discussed by Wilkinson and others.

Input Options

  • rowdim, coldim, α, k: row_dim and col_dim are row and column dimension of the matrix. α is a scalar representing the entries on the superdiagonals. k is the number of superdiagonals.
  • dim: the dimension of the matrix.

References

G. H. Golub and J. H. Wilkinson, Ill-conditioned eigensystems and the computation of the Jordan canonical form, SIAM Rev., 18 (1976), pp. 578-619, https://doi.org/10.1137/1018113.

source
TypedMatrices.WathenType

Wathen Matrix

The Wathen Matrix is the consistent mass matrix of a regular nx-by-ny` grid of 8 nodes in the finite element method. The matrix is a sparse, symmetric positive definite, and has random entries.

Input Options

  • [type,] nx, ny: the dimension of the matrix is equal to 3 * nx * ny + 2 * nx * ny + 1.
  • [type,] n: nx = ny = n.

References

A. J. Wathen, Realistic eigenvalue bounds for the Galerkin mass matrix, IMA J. Numer. Anal., 7 (1987), pp. 449-457, https://doi.org/10.1093/imanum/7.4.449.

source
TypedMatrices.WilkinsonType

Wilkinson Matrix

The Wilkinson matrix is a symmetric tridiagonal matrix with pairs of nearly equal eigenvalues. The most frequently used ordre is 21.

Input Options

  • dim: the dimension of the matrix.

References

J. H. Wilkinson, Error analysis of direct methods of matrix inversion, J. Assoc. Comput. Mach., 8 (1961), pp. 281-330, https://doi.org/10.1145/321075.321076.

source