Class IMatrix

java.lang.Object
smile.math.matrix.IMatrix
All Implemented Interfaces:
Serializable, Cloneable
Direct Known Subclasses:
BandMatrix, BigMatrix, Matrix, SparseMatrix, SymmMatrix

public abstract class IMatrix extends Object implements Cloneable, Serializable
Matrix base class. The most important method is the matrix vector multiplication, which is the only operation needed in many iterative matrix algorithms, e.g. biconjugate gradient method for solving linear equations and power iteration and Lanczos algorithm for eigen decomposition, which are usually very efficient for very large and sparse matrices.

A matrix is a rectangular array of numbers. An item in a matrix is called an entry or an element. Entries are often denoted by a variable with two subscripts. Matrices of the same size can be added and subtracted entrywise and matrices of compatible size can be multiplied. These operations have many of the properties of ordinary arithmetic, except that matrix multiplication is not commutative, that is, AB and BA are not equal in general.

Matrices are a key tool in linear algebra. One use of matrices is to represent linear transformations and matrix multiplication corresponds to composition of linear transformations. Matrices can also keep track of the coefficients in a system of linear equations. For a square matrix, the determinant and inverse matrix (when it exists) govern the behavior of solutions to the corresponding system of linear equations, and eigenvalues and eigenvectors provide insight into the geometry of the associated linear transformation.

There are several methods to render matrices into a more easily accessible form. They are generally referred to as matrix transformation or matrix decomposition techniques. The interest of all these decomposition techniques is that they preserve certain properties of the matrices in question, such as determinant, rank or inverse, so that these quantities can be calculated after applying the transformation, or that certain matrix operations are algorithmically easier to carry out for some types of matrices.

The LU decomposition factors matrices as a product of lower (L) and an upper triangular matrices (U). Once this decomposition is calculated, linear systems can be solved more efficiently, by a simple technique called forward and back substitution. Likewise, inverses of triangular matrices are algorithmically easier to calculate. The QR decomposition factors matrices as a product of an orthogonal (Q) and a right triangular matrix (R). QR decomposition is often used to solve the linear least squares problem, and is the basis for a particular eigenvalue algorithm, the QR algorithm. Singular value decomposition expresses any matrix A as a product UDV', where U and V are unitary matrices and D is a diagonal matrix. The eigendecomposition or diagonalization expresses A as a product VDV-1, where D is a diagonal matrix and V is a suitable invertible matrix. If A can be written in this form, it is called diagonalizable.

See Also:
  • Nested Class Summary

    Nested Classes
    Modifier and Type
    Class
    Description
    static interface 
    The preconditioner matrix.
  • Constructor Summary

    Constructors
    Constructor
    Description
     
  • Method Summary

    Modifier and Type
    Method
    Description
    double
    apply(int i, int j)
    Returns A[i,j] for Scala users.
    colName(int i)
    Returns the name of i-th column.
    Returns the column names.
    void
    colNames(String[] names)
    Sets the column names.
    double[]
    Returns the diagonal elements.
    double
    eigen(double[] v)
    Returns the largest eigen pair of matrix with the power iteration under the assumptions A has an eigenvalue that is strictly greater in magnitude than its other eigenvalues and the starting vector has a nonzero component in the direction of an eigenvector associated with the dominant eigenvalue.
    double
    eigen(double[] v, double p, double tol, int maxIter)
    Returns the largest eigen pair of matrix with the power iteration under the assumptions A has an eigenvalue that is strictly greater in magnitude than its other eigenvalues and the starting vector has a nonzero component in the direction of an eigenvector associated with the dominant eigenvalue.
    double
    get(int i, int j)
    Returns A[i,j].
    Returns a simple Jacobi preconditioner matrix that is the trivial diagonal part of A in some cases.
    static IMatrix
    market(Path path)
    Reads a matrix from a Matrix Market File Format file.
    double[]
    mv(double[] x)
    Returns the matrix-vector multiplication A * x.
    void
    mv(double[] x, double[] y)
    Matrix-vector multiplication y = A * x.
    abstract void
    mv(double[] work, int inputOffset, int outputOffset)
    Matrix-vector multiplication A * x.
    void
    mv(double alpha, double[] x, double beta, double[] y)
    Matrix-vector multiplication.
    abstract void
    mv(Transpose trans, double alpha, double[] x, double beta, double[] y)
    Matrix-vector multiplication.
    abstract int
    Returns the number of columns.
    abstract int
    Returns the number of rows.
    rowName(int i)
    Returns the name of i-th row.
    Returns the row names.
    void
    rowNames(String[] names)
    Sets the row names.
    void
    set(int i, int j, double x)
    Sets A[i,j] = x.
    abstract long
    Returns the number of stored matrix elements.
    double
    solve(double[] b, double[] x)
    Solves A * x = b by iterative biconjugate gradient method with Jacobi preconditioner matrix.
    double
    solve(double[] b, double[] x, IMatrix.Preconditioner P, double tol, int itol, int maxIter)
    Solves A * x = b by iterative biconjugate gradient method.
    Returns the square matrix of A' * A or A * A', whichever is smaller.
     
    toString(boolean full)
    Returns the string representation of matrix.
    toString(int m, int n)
    Returns the string representation of matrix.
    double
    Returns the matrix trace.
    double[]
    tv(double[] x)
    Returns Matrix-vector multiplication A' * x.
    void
    tv(double[] x, double[] y)
    Matrix-vector multiplication y = A' * x.
    abstract void
    tv(double[] work, int inputOffset, int outputOffset)
    Matrix-vector multiplication A' * x.
    void
    tv(double alpha, double[] x, double beta, double[] y)
    Matrix-vector multiplication.
    void
    update(int i, int j, double x)
    Sets A[i,j] = x for Scala users.

    Methods inherited from class java.lang.Object

    clone, equals, finalize, getClass, hashCode, notify, notifyAll, wait, wait, wait
  • Constructor Details

    • IMatrix

      public IMatrix()
  • Method Details

    • nrow

      public abstract int nrow()
      Returns the number of rows.
      Returns:
      the number of rows.
    • ncol

      public abstract int ncol()
      Returns the number of columns.
      Returns:
      the number of columns.
    • size

      public abstract long size()
      Returns the number of stored matrix elements. For conventional matrix, it is simplify nrow * ncol. But it is usually much less for band, packed or sparse matrix.
      Returns:
      the number of stored matrix elements.
    • rowNames

      public String[] rowNames()
      Returns the row names.
      Returns:
      the row names.
    • rowNames

      public void rowNames(String[] names)
      Sets the row names.
      Parameters:
      names - the row names.
    • rowName

      public String rowName(int i)
      Returns the name of i-th row.
      Parameters:
      i - the row index.
      Returns:
      the name of i-th row.
    • colNames

      public String[] colNames()
      Returns the column names.
      Returns:
      the column names.
    • colNames

      public void colNames(String[] names)
      Sets the column names.
      Parameters:
      names - the column names.
    • colName

      public String colName(int i)
      Returns the name of i-th column.
      Parameters:
      i - the column index.
      Returns:
      the name of i-th column.
    • toString

      public String toString()
      Overrides:
      toString in class Object
    • toString

      public String toString(boolean full)
      Returns the string representation of matrix.
      Parameters:
      full - Print the full matrix if true. Otherwise, print only top left 7 x 7 submatrix.
      Returns:
      the string representation of matrix.
    • toString

      public String toString(int m, int n)
      Returns the string representation of matrix.
      Parameters:
      m - the number of rows to print.
      n - the number of columns to print.
      Returns:
      the string representation of matrix.
    • mv

      public abstract void mv(Transpose trans, double alpha, double[] x, double beta, double[] y)
      Matrix-vector multiplication.
      
           y = alpha * op(A) * x + beta * y
       
      where op is the transpose operation.
      Parameters:
      trans - normal, transpose, or conjugate transpose operation on the matrix.
      alpha - the scalar alpha.
      x - the input vector.
      beta - the scalar beta. When beta is supplied as zero then y need not be set on input.
      y - the input and output vector.
    • mv

      public double[] mv(double[] x)
      Returns the matrix-vector multiplication A * x.
      Parameters:
      x - the vector.
      Returns:
      the matrix-vector multiplication A * x.
    • mv

      public void mv(double[] x, double[] y)
      Matrix-vector multiplication y = A * x.
      Parameters:
      x - the input vector.
      y - the output vector.
    • mv

      public void mv(double alpha, double[] x, double beta, double[] y)
      Matrix-vector multiplication.
      
           y = alpha * A * x + beta * y
       
      Parameters:
      alpha - the scalar alpha.
      x - the input vector.
      beta - the scalar beta. When beta is supplied as zero then y need not be set on input.
      y - the input and output vector.
    • mv

      public abstract void mv(double[] work, int inputOffset, int outputOffset)
      Matrix-vector multiplication A * x.
      Parameters:
      work - the workspace for both input and output vector.
      inputOffset - the offset of input vector in workspace.
      outputOffset - the offset of output vector in workspace.
    • tv

      public double[] tv(double[] x)
      Returns Matrix-vector multiplication A' * x.
      Parameters:
      x - the vector.
      Returns:
      the matrix-vector multiplication A' * x.
    • tv

      public void tv(double[] x, double[] y)
      Matrix-vector multiplication y = A' * x.
      Parameters:
      x - the input vector.
      y - the output vector.
    • tv

      public void tv(double alpha, double[] x, double beta, double[] y)
      Matrix-vector multiplication.
      
           y = alpha * A' * x + beta * y
       
      Parameters:
      alpha - the scalar alpha.
      x - the input vector.
      beta - the scalar beta. When beta is supplied as zero then y need not be set on input.
      y - the input and output vector.
    • tv

      public abstract void tv(double[] work, int inputOffset, int outputOffset)
      Matrix-vector multiplication A' * x.
      Parameters:
      work - the workspace for both input and output vector.
      inputOffset - the offset of input vector in workspace.
      outputOffset - the offset of output vector in workspace.
    • set

      public void set(int i, int j, double x)
      Sets A[i,j] = x.
      Parameters:
      i - the row index.
      j - the column index.
      x - the matrix cell value.
    • update

      public void update(int i, int j, double x)
      Sets A[i,j] = x for Scala users.
      Parameters:
      i - the row index.
      j - the column index.
      x - the matrix cell value.
    • get

      public double get(int i, int j)
      Returns A[i,j].
      Parameters:
      i - the row index.
      j - the column index.
      Returns:
      the matrix cell value.
    • apply

      public double apply(int i, int j)
      Returns A[i,j] for Scala users.
      Parameters:
      i - the row index.
      j - the column index.
      Returns:
      the matrix cell value.
    • diag

      public double[] diag()
      Returns the diagonal elements.
      Returns:
      the diagonal elements.
    • trace

      public double trace()
      Returns the matrix trace. The sum of the diagonal elements.
      Returns:
      the matrix trace.
    • eigen

      public double eigen(double[] v)
      Returns the largest eigen pair of matrix with the power iteration under the assumptions A has an eigenvalue that is strictly greater in magnitude than its other eigenvalues and the starting vector has a nonzero component in the direction of an eigenvector associated with the dominant eigenvalue.
      Parameters:
      v - on input, it is the non-zero initial guess of the eigen vector. On output, it is the eigen vector corresponding largest eigen value.
      Returns:
      the largest eigen value.
    • eigen

      public double eigen(double[] v, double p, double tol, int maxIter)
      Returns the largest eigen pair of matrix with the power iteration under the assumptions A has an eigenvalue that is strictly greater in magnitude than its other eigenvalues and the starting vector has a nonzero component in the direction of an eigenvector associated with the dominant eigenvalue.
      Parameters:
      v - on input, it is the non-zero initial guess of the eigen vector. On output, it is the eigen vector corresponding largest eigen value.
      p - the origin in the shifting power method. A - pI will be used in the iteration to accelerate the method. p should be such that |(λ2 - p) / (λ1 - p)| < |λ2 / λ1|, where λ2 is the second largest eigenvalue in magnitude. If we known the eigenvalue spectrum of A, (λ2 + λn)/2 is the optimal choice of p, where λn is the smallest eigenvalue in magnitude. Good estimates of λ2 are more difficult to compute. However, if μ is an approximation to largest eigenvector, then using any x0 such that x0*μ = 0 as the initial vector for a few iterations may yield a reasonable estimate of λ2.
      tol - the desired convergence tolerance.
      maxIter - the maximum number of iterations in case that the algorithm does not converge.
      Returns:
      the largest eigen value.
    • market

      public static IMatrix market(Path path) throws IOException, ParseException
      Reads a matrix from a Matrix Market File Format file. For details, see http://people.sc.fsu.edu/~jburkardt/data/mm/mm.html. The returned matrix may be dense or sparse.
      Parameters:
      path - the input file path.
      Returns:
      a dense or sparse matrix.
      Throws:
      IOException - when fails to read the file.
      ParseException - when fails to parse the file.
    • square

      public IMatrix square()
      Returns the square matrix of A' * A or A * A', whichever is smaller. For SVD, we compute eigenvalue decomposition of A' * A when m >= n, or that of A * A' when m < n.
      Returns:
      the matrix of A' * A or A * A', whichever is smaller.
    • Jacobi

      public IMatrix.Preconditioner Jacobi()
      Returns a simple Jacobi preconditioner matrix that is the trivial diagonal part of A in some cases.
      Returns:
      the preconditioner matrix.
    • solve

      public double solve(double[] b, double[] x)
      Solves A * x = b by iterative biconjugate gradient method with Jacobi preconditioner matrix.
      Parameters:
      b - the right hand side of linear equations.
      x - on input, x should be set to an initial guess of the solution (or all zeros). On output, x is reset to the improved solution.
      Returns:
      the estimated error.
    • solve

      public double solve(double[] b, double[] x, IMatrix.Preconditioner P, double tol, int itol, int maxIter)
      Solves A * x = b by iterative biconjugate gradient method.
      Parameters:
      b - the right hand side of linear equations.
      x - on input, x should be set to an initial guess of the solution (or all zeros). On output, x is reset to the improved solution.
      P - The preconditioner matrix.
      tol - The desired convergence tolerance.
      itol - Which convergence test is applied. If itol = 1, iteration stops when |Ax - b| / |b| is less than the parameter tolerance. If itol = 2, the stop criterion is that |A-1 (Ax - b)| / |A-1b| is less than tolerance. If tol = 3, |xk+1 - xk|2 is less than tolerance. The setting of tol = 4 is same as tol = 3 except that the L norm instead of L2.
      maxIter - The maximum number of iterations.
      Returns:
      the estimated error.