Class CommonOps_DSCC

java.lang.Object
org.ejml.sparse.csc.CommonOps_DSCC

public class CommonOps_DSCC extends Object
Most common operations on DMatrixSparseCSC. For example, addition, matrix multiplication, ... etc
  • Method Details

    • checkIndicesSorted

      public static boolean checkIndicesSorted(DMatrixSparseCSC A)
      Checks to see if row indicies are sorted into ascending order. O(N)
      Returns:
      true if sorted and false if not
    • checkStructure

      public static boolean checkStructure(DMatrixSparseCSC A)
    • checkSortedFlag

      public static boolean checkSortedFlag(DMatrixSparseCSC A)
    • checkDuplicateElements

      public static boolean checkDuplicateElements(DMatrixSparseCSC A)
      Checks for duplicate elements. A is sorted
      Parameters:
      A - Matrix to be tested.
      Returns:
      true if duplicates or false if false duplicates
    • transpose

      public static DMatrixSparseCSC transpose(DMatrixSparseCSC A, @Nullable @Nullable DMatrixSparseCSC A_t, @Nullable @Nullable IGrowArray gw)
      Perform matrix transpose
      Parameters:
      A - Input matrix. Not modified
      A_t - Storage for transpose of 'a'. Must be correct shape. data length might be adjusted.
      gw - (Optional) Storage for internal workspace. Can be null.
      Returns:
      The transposed matrix
    • mult

      public static DMatrixSparseCSC mult(DMatrixSparseCSC A, DMatrixSparseCSC B, @Nullable @Nullable DMatrixSparseCSC outputC)
    • mult

      public static DMatrixSparseCSC mult(DMatrixSparseCSC A, DMatrixSparseCSC B, @Nullable @Nullable DMatrixSparseCSC outputC, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable DGrowArray gx)
      Performs matrix multiplication. C = A*B
      Parameters:
      A - (Input) Matrix. Not modified.
      B - (Input) Matrix. Not modified.
      outputC - (Output) Storage for results. Data length is increased if insufficient.
      gw - (Optional) Storage for internal workspace. Can be null.
      gx - (Optional) Storage for internal workspace. Can be null.
    • mult

      public static DMatrixRMaj mult(DMatrixSparseCSC A, DMatrixRMaj B, @Nullable @Nullable DMatrixRMaj outputC)
      Performs matrix multiplication. C = A*B
      Parameters:
      A - Matrix
      B - Dense Matrix
      outputC - Dense Matrix
    • multAdd

      public static void multAdd(DMatrixSparseCSC A, DMatrixRMaj B, DMatrixRMaj outputC)

      C = C + A*B

    • multTransA

      public static DMatrixRMaj multTransA(DMatrixSparseCSC A, DMatrixRMaj B, @Nullable @Nullable DMatrixRMaj outputC, @Nullable @Nullable DGrowArray work)
      Performs matrix multiplication. C = AT*B
      Parameters:
      A - Matrix
      B - Dense Matrix
      outputC - Dense Matrix
    • multAddTransA

      public static void multAddTransA(DMatrixSparseCSC A, DMatrixRMaj B, DMatrixRMaj outputC, @Nullable @Nullable DGrowArray work)

      C = C + AT*B

    • multTransB

      public static DMatrixRMaj multTransB(DMatrixSparseCSC A, DMatrixRMaj B, @Nullable @Nullable DMatrixRMaj outputC, @Nullable @Nullable DGrowArray work)
      Performs matrix multiplication. C = A*BT
      Parameters:
      A - Matrix
      B - Dense Matrix
      outputC - Dense Matrix
    • multAddTransB

      public static void multAddTransB(DMatrixSparseCSC A, DMatrixRMaj B, DMatrixRMaj outputC, @Nullable @Nullable DGrowArray work)

      C = C + A*BT

    • multTransAB

      public static DMatrixRMaj multTransAB(DMatrixSparseCSC A, DMatrixRMaj B, DMatrixRMaj outputC)
      Performs matrix multiplication. C = AT*BT
      Parameters:
      A - Matrix
      B - Dense Matrix
      outputC - Dense Matrix
    • multAddTransAB

      public static void multAddTransAB(DMatrixSparseCSC A, DMatrixRMaj B, DMatrixRMaj outputC)

      C = C + AT*BT

    • symmLowerToFull

      public static void symmLowerToFull(DMatrixSparseCSC A, DMatrixSparseCSC outputB, @Nullable @Nullable IGrowArray gw)
      Given a symmetric matrix, which is represented by a lower triangular matrix, convert it back into a full symmetric matrix
      Parameters:
      A - (Input) Lower triangular matrix
      outputB - (Output) Symmetric matrix.
      gw - (Optional) Workspace. Can be null.
    • add

      public static DMatrixSparseCSC add(double alpha, DMatrixSparseCSC A, double beta, DMatrixSparseCSC B, @Nullable @Nullable DMatrixSparseCSC outputC, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable DGrowArray gx)
      Performs matrix addition:
      C = αA + βB
      Parameters:
      alpha - scalar value multiplied against A
      A - Matrix
      beta - scalar value multiplied against B
      B - Matrix
      outputC - Output matrix.
      gw - (Optional) Storage for internal workspace. Can be null.
      gx - (Optional) Storage for internal workspace. Can be null.
    • identity

      public static DMatrixSparseCSC identity(int length)
    • identity

      public static DMatrixSparseCSC identity(int numRows, int numCols)
    • setIdentity

      public static void setIdentity(DMatrixSparseCSC A)
    • scale

      public static void scale(double scalar, DMatrixSparseCSC A, DMatrixSparseCSC outputB)
      B = scalar*A. A and B can be the same instance.
      Parameters:
      scalar - (Input) Scalar value
      A - (Input) Matrix. Not modified.
      outputB - (Output) Matrix. Modified.
    • divide

      public static void divide(DMatrixSparseCSC A, double scalar, DMatrixSparseCSC outputB)
      B = A/scalar. A and B can be the same instance.
      Parameters:
      scalar - (Input) Scalar value
      A - (Input) Matrix. Not modified.
      outputB - (Output) Matrix. Modified.
    • divide

      public static void divide(double scalar, DMatrixSparseCSC A, DMatrixSparseCSC outputB)
      B = scalar/A. A and B can be the same instance. Only non-zero values are affected
      Parameters:
      A - (Input) Matrix. Not modified.
      scalar - (Input) Scalar value
      outputB - (Output) Matrix. Modified.
    • changeSign

      public static void changeSign(DMatrixSparseCSC A, DMatrixSparseCSC outputB)
      B = -A. Changes the sign of elements in A and stores it in B. A and B can be the same instance.
      Parameters:
      A - (Input) Matrix. Not modified.
      outputB - (Output) Matrix. Modified.
    • elementMinAbs

      public static double elementMinAbs(DMatrixSparseCSC A)
      Returns the value of the element with the smallest abs()
      Parameters:
      A - (Input) Matrix. Not modified.
      Returns:
      scalar
    • elementMaxAbs

      public static double elementMaxAbs(DMatrixSparseCSC A)
      Returns the value of the element with the largest abs()
      Parameters:
      A - (Input) Matrix. Not modified.
      Returns:
      scalar
    • elementMin

      public static double elementMin(DMatrixSparseCSC A)
      Returns the value of the element with the minimum value
      Parameters:
      A - (Input) Matrix. Not modified.
      Returns:
      scalar
    • elementMax

      public static double elementMax(DMatrixSparseCSC A)
      Returns the value of the element with the largest value
      Parameters:
      A - (Input) Matrix. Not modified.
      Returns:
      scalar
    • elementSum

      public static double elementSum(DMatrixSparseCSC A)
      Sum of all elements
      Parameters:
      A - (Input) Matrix. Not modified.
      Returns:
      scalar
    • elementMult

      public static DMatrixSparseCSC elementMult(DMatrixSparseCSC A, DMatrixSparseCSC B, @Nullable @Nullable DMatrixSparseCSC output, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable DGrowArray gx)
      Performs an element-wise multiplication.
      output[i,j] = A[i,j]*B[i,j]
      All matrices must have the same shape.
      Parameters:
      A - (Input) Matrix.
      B - (Input) Matrix
      output - (Output) Matrix. data array is grown to min(A.nz_length,B.nz_length), resulting a in a large speed boost.
      gw - (Optional) Storage for internal workspace. Can be null.
      gx - (Optional) Storage for internal workspace. Can be null.
    • maxAbsCols

      public static void maxAbsCols(DMatrixSparseCSC A, @Nullable @Nullable DMatrixRMaj outputB)
      Finds the maximum abs in each column of A and stores it into values
      Parameters:
      A - (Input) Matrix
      outputB - (Output) storage for column max abs
    • multColumns

      public static void multColumns(DMatrixSparseCSC A, double[] values, int offset)
      Multiply all elements of column 'i' by value[i]. A[:,i] *= values[i].
      Equivalent to A = A*diag(values)
      Parameters:
      A - (Input/Output) Matrix. Modified.
      values - (Input) multiplication factor for each column
      offset - (Input) first index in values to start at
    • divideColumns

      public static void divideColumns(DMatrixSparseCSC A, double[] values, int offset)
      Divides all elements of column 'i' by values[i]. A[:,i] /= values[i].
      Equivalent to A = A*inv(diag(values))
      Parameters:
      A - (Input/Output) Matrix. Modified.
      values - (Input) multiplication factor for each column
      offset - (Input) first index in values to start at
    • multRowsCols

      public static void multRowsCols(double[] diagA, int offsetA, DMatrixSparseCSC B, double[] diagC, int offsetC)
      Equivalent to multiplying a matrix B by two diagonal matrices. B = A*B*C, where A=diag(a) and C=diag(c).
      Parameters:
      diagA - Array of length offsteA + B.numRows
      offsetA - First index in A
      B - Rectangular matrix
      diagC - Array of length indexC + B.numCols
      offsetC - First index in C
    • divideRowsCols

      public static void divideRowsCols(double[] diagA, int offsetA, DMatrixSparseCSC B, double[] diagC, int offsetC)
      Equivalent to multiplying a matrix B by the inverse of two diagonal matrices. B = inv(A)*B*inv(C), where A=diag(a) and C=diag(c).
      Parameters:
      diagA - Array of length offsteA + B.numRows
      offsetA - First index in A
      B - Rectangular matrix
      diagC - Array of length indexC + B.numCols
      offsetC - First index in C
    • diag

      public static DMatrixSparseCSC diag(double... values)
      Returns a diagonal matrix with the specified diagonal elements.
      Parameters:
      values - values of diagonal elements
      Returns:
      A diagonal matrix
    • diag

      public static DMatrixSparseCSC diag(@Nullable @Nullable DMatrixSparseCSC A, double[] values, int offset, int length)
      Creates a diagonal matrix from an array. Elements in the array can be offset.
      Parameters:
      A - (Optional) Storage for diagonal matrix. If null a new one will be declared.
      values - First index in the diagonal matirx
      length - Length of the diagonal matrix
      offset - First index in values
      Returns:
      The diagonal matrix
    • extractDiag

      public static void extractDiag(DMatrixSparseCSC A, DMatrixSparseCSC outputB)

      Extracts the diagonal elements 'src' write it to the 'dst' vector. 'dst' can either be a row or column vector.

      Parameters:
      A - Matrix whose diagonal elements are being extracted. Not modified.
      outputB - A vector the results will be written into. Modified.
    • extractDiag

      public static void extractDiag(DMatrixSparseCSC A, DMatrixRMaj outputB)

      Extracts the diagonal elements 'src' write it to the 'dst' vector. 'dst' can either be a row or column vector.

      Parameters:
      A - Matrix whose diagonal elements are being extracted. Not modified.
      outputB - A vector the results will be written into. Modified.
    • permutationMatrix

      public static DMatrixSparseCSC permutationMatrix(int[] p, boolean inverse, int N, @Nullable @Nullable DMatrixSparseCSC P)
      Converts the permutation vector into a matrix. B = P*A. B[p[i],:] = A[i,:]
      Parameters:
      p - (Input) Permutation vector
      inverse - (Input) If it is the inverse. B[i,:] = A[p[i],:)
      P - (Output) Permutation matrix
    • permutationVector

      public static void permutationVector(DMatrixSparseCSC P, int[] vector)
      Converts the permutation matrix into a vector
      Parameters:
      P - (Input) Permutation matrix
      vector - (Output) Permutation vector
    • permutationInverse

      public static void permutationInverse(int[] original, int[] inverse, int length)
      Computes the inverse permutation vector
      Parameters:
      original - Original permutation vector
      inverse - It's inverse
    • permutationInverse

      public static int[] permutationInverse(int[] original, int length)
    • permuteRowInv

      public static void permuteRowInv(int[] permInv, DMatrixSparseCSC input, DMatrixSparseCSC output)
      Applies the row permutation specified by the vector to the input matrix and save the results in the output matrix. output[perm[j],:] = input[j,:]
      Parameters:
      permInv - (Input) Inverse permutation vector. Specifies new order of the rows.
      input - (Input) Matrix which is to be permuted
      output - (Output) Matrix which has the permutation stored in it. Is reshaped.
    • permute

      public static void permute(@Nullable @org.jetbrains.annotations.Nullable int[] permRowInv, DMatrixSparseCSC input, @Nullable @org.jetbrains.annotations.Nullable int[] permCol, DMatrixSparseCSC output)
      Applies the forward column and inverse row permutation specified by the two vector to the input matrix and save the results in the output matrix. output[permRow[j],permCol[i]] = input[j,i]
      Parameters:
      permRowInv - (Input) Inverse row permutation vector. Null is the same as passing in identity.
      input - (Input) Matrix which is to be permuted
      permCol - (Input) Column permutation vector. Null is the same as passing in identity.
      output - (Output) Matrix which has the permutation stored in it. Is reshaped.
    • permute

      public static void permute(int[] perm, double[] input, double[] output, int N)
      Permutes a vector. output[i] = input[perm[i]]
      Parameters:
      perm - (Input) permutation vector
      input - (Input) Vector which is to be permuted
      output - (Output) Where the permuted vector is stored.
      N - Number of elements in the vector.
    • permuteInv

      public static void permuteInv(int[] perm, double[] input, double[] output, int N)
      Permutes a vector in the inverse. output[perm[k]] = input[k]
      Parameters:
      perm - (Input) permutation vector
      input - (Input) Vector which is to be permuted
      output - (Output) Where the permuted vector is stored.
      N - Number of elements in the vector.
    • permuteSymmetric

      public static void permuteSymmetric(DMatrixSparseCSC input, int[] permInv, DMatrixSparseCSC output, @Nullable @Nullable IGrowArray gw)
      Applies the permutation to upper triangular symmetric matrices. Typically a symmetric matrix only stores the upper triangular part, so normal permutation will have undesirable results, e.g. the zeros will get mixed in and will no longer be symmetric. This algorithm will handle the implicit lower triangular and construct new upper triangular matrix.

      See page cs_symperm() on Page 22 of "Direct Methods for Sparse Linear Systems"

      Parameters:
      input - (Input) Upper triangular symmetric matrix which is to be permuted. Entries below the diagonal are ignored.
      permInv - (Input) Inverse permutation vector. Specifies new order of the rows and columns.
      output - (Output) Upper triangular symmetric matrix which has the permutation stored in it. Reshaped.
      gw - (Optional) Storage for internal workspace. Can be null.
    • concatRows

      public static DMatrixSparseCSC concatRows(DMatrixSparseCSC top, DMatrixSparseCSC bottom, @Nullable @Nullable DMatrixSparseCSC out)
      Concats two matrices along their rows (vertical).
      Parameters:
      top - Matrix on the top
      bottom - Matrix on the bototm
      out - (Output) (Optional) Storage for combined matrix. Resized.
      Returns:
      Combination of the two matrices
    • concatColumns

      public static DMatrixSparseCSC concatColumns(DMatrixSparseCSC left, DMatrixSparseCSC right, @Nullable @Nullable DMatrixSparseCSC out)
      Concats two matrices along their columns (horizontal).
      Parameters:
      left - Matrix on the left
      right - Matrix on the right
      out - (Output) (Optional) Storage for combined matrix. Resized.
      Returns:
      Combination of the two matrices
    • extractColumn

      public static DMatrixSparseCSC extractColumn(DMatrixSparseCSC A, int column, @Nullable @Nullable DMatrixSparseCSC out)
      Extracts a column from A and stores it into out.
      Parameters:
      A - (Input) Source matrix. not modified.
      column - The column in A
      out - (Output, Optional) Storage for column vector
      Returns:
      The column of A.
    • extractRows

      public static DMatrixSparseCSC extractRows(DMatrixSparseCSC A, int row0, int row1, @Nullable @Nullable DMatrixSparseCSC out)
      Creates a submatrix by extracting the specified rows from A. rows = {row0 %le; i %le; row1}.
      Parameters:
      A - (Input) matrix
      row0 - First row. Inclusive
      row1 - Last row+1.
      out - (Output, Option) Storage for output matrix
      Returns:
      The submatrix
    • extract

      public static void extract(DMatrixSparseCSC src, int srcY0, int srcY1, int srcX0, int srcX1, DMatrixSparseCSC dst, int dstY0, int dstX0)

      Extracts a submatrix from 'src' and inserts it in a submatrix in 'dst'.

      si-y0 , j-x0 = oij for all y0 ≤ i < y1 and x0 ≤ j < x1

      where 'sij' is an element in the submatrix and 'oij' is an element in the original matrix.

      WARNING: This is a very slow operation for sparse matrices. The current implementation is simple but involves excessive memory copying.

      Parameters:
      src - The original matrix which is to be copied. Not modified.
      srcX0 - Start column.
      srcX1 - Stop column+1.
      srcY0 - Start row.
      srcY1 - Stop row+1.
      dst - Where the submatrix are stored. Modified.
      dstY0 - Start row in dst.
      dstX0 - start column in dst.
    • select

      public static DMatrixSparseCSC select(DMatrixSparseCSC A, IPredicateBinary selector, @Nullable @Nullable DMatrixSparseCSC output)
      Select entries from A and save them in C. A more generic but probably also slower version of `extract`

      Can be used to select f.i. the lower or upper triangle of a matrix

      Simplified version of: GxB_Select

      Parameters:
      A - (Input) Matrix. Not modified.
      selector - Function to decide whether an entry gets selected
      output - (Optional/Output) Matrix to use for the output. Can be the same as A
      Returns:
      Matrix storing the selected entries of A
    • fill

      public static void fill(DMatrixSparseCSC A, double value)

      Sets every element in the matrix to the specified value. This can require a very large amount of memory and might exceed the maximum array size

      Aij = value

      Parameters:
      A - A matrix whose elements are about to be set. Modified.
      value - The value each element will have.
    • sumCols

      public static DMatrixRMaj sumCols(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output)

      Computes the sum of each column in the input matrix and returns the results in a vector:

      bj = sum(i=1:m ; aij)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a row vector. Modified.
      Returns:
      Vector containing the sum of each column
    • minCols

      public static DMatrixRMaj minCols(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output)

      Computes the minimum of each column in the input matrix and returns the results in a vector:

      bj = min(i=1:m ; aij)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a row vector. Modified.
      Returns:
      Vector containing the minimums of each column
    • maxCols

      public static DMatrixRMaj maxCols(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output)

      Computes the maximums of each column in the input matrix and returns the results in a vector:

      bj = max(i=1:m ; aij)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a row vector. Modified.
      Returns:
      Vector containing the maximums of each column
    • sumRows

      public static DMatrixRMaj sumRows(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output)

      Computes the sum of each row in the input matrix and returns the results in a vector:

      bj = sum(i=1:n ; aji)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a column vector. Modified.
      Returns:
      Vector containing the sum of each row
    • minRows

      public static DMatrixRMaj minRows(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output, @Nullable @Nullable IGrowArray gw)

      Computes the minimum of each row in the input matrix and returns the results in a vector:

      bj = min(i=1:n ; aji)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a column vector. Modified.
      gw - work space
      Returns:
      Vector containing the minimum of each row
    • maxRows

      public static DMatrixRMaj maxRows(DMatrixSparseCSC input, @Nullable @Nullable DMatrixRMaj output, @Nullable @Nullable IGrowArray gw)

      Computes the maximum of each row in the input matrix and returns the results in a vector:

      bj = max(i=1:n ; aji)

      Parameters:
      input - Input matrix
      output - Optional storage for output. Reshaped into a column vector. Modified.
      gw - work space
      Returns:
      Vector containing the maximum of each row
    • zero

      public static void zero(DMatrixSparseCSC A, int row0, int row1, int col0, int col1)
      Zeros an inner rectangle inside the matrix.
      Parameters:
      A - Matrix that is to be modified.
      row0 - Start row.
      row1 - Stop row+1.
      col0 - Start column.
      col1 - Stop column+1.
    • dotInnerColumns

      public static double dotInnerColumns(DMatrixSparseCSC A, int colA, DMatrixSparseCSC B, int colB, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable DGrowArray gx)
      Computes the inner product of two column vectors taken from the input matrices.

      dot = A(:,colA)'*B(:,colB)

      Parameters:
      A - Matrix
      colA - Column in A
      B - Matrix
      colB - Column in B
      Returns:
      Dot product
    • solve

      public static boolean solve(DMatrixSparseCSC a, DMatrixRMaj b, DMatrixRMaj x)

      Solves for x in the following equation:

      A*x = b

      If the system could not be solved then false is returned. If it returns true that just means the algorithm finished operating, but the results could still be bad because 'A' is singular or nearly singular.

      If repeat calls to solve are being made then one should consider using LinearSolverFactory_DSCC instead.

      It is ok for 'b' and 'x' to be the same matrix.

      Parameters:
      a - (Input) A matrix that is m by n. Not modified.
      b - (Input) A matrix that is n by k. Not modified.
      x - (Output) A matrix that is m by k. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • solve

      public static boolean solve(DMatrixSparseCSC a, DMatrixSparseCSC b, DMatrixSparseCSC x)

      Solves for x in the following equation:

      A*x = b

      If the system could not be solved then false is returned. If it returns true that just means the algorithm finished operating, but the results could still be bad because 'A' is singular or nearly singular.

      If repeat calls to solve are being made then one should consider using LinearSolverFactory_DSCC instead.

      It is ok for 'b' and 'x' to be the same matrix.

      Parameters:
      a - (Input) A matrix that is m by n. Not modified.
      b - (Input) A matrix that is n by k. Not modified.
      x - (Output) A matrix that is m by k. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • invert

      public static boolean invert(DMatrixSparseCSC A, DMatrixRMaj inverse)

      Performs a matrix inversion operation that does not modify the original and stores the results in another matrix. The two matrices must have the same dimension.

      B = A-1

      If the algorithm could not invert the matrix then false is returned. If it returns true that just means the algorithm finished. The results could still be bad because the matrix is singular or nearly singular.

      For medium to large matrices there might be a slight performance boost to using LinearSolverFactory_DSCC instead.

      Parameters:
      A - (Input) The matrix that is to be inverted. Not modified.
      inverse - (Output) Where the inverse matrix is stored. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • det

      public static double det(DMatrixSparseCSC A)
      Returns the determinant of the matrix. If the inverse of the matrix is also needed, then using LUDecomposition_F64 directly (or any similar algorithm) can be more efficient.
      Parameters:
      A - The matrix whose determinant is to be computed. Not modified.
      Returns:
      The determinant.
    • removeZeros

      public static void removeZeros(DMatrixSparseCSC input, DMatrixSparseCSC output, double tol)
      Copies all elements from input into output which are > tol.
      Parameters:
      input - (Input) input matrix. Not modified.
      output - (Output) Output matrix. Modified and shaped to match input.
      tol - Tolerance for defining zero
    • removeZeros

      public static void removeZeros(DMatrixSparseCSC A, double tol)
      Removes all elements from the matrix that are > tol. The modification is done in place and no temporary storage is declared.
      Parameters:
      A - (Input/Output) input matrix. Modified.
      tol - Tolerance for defining zero
    • duplicatesAdd

      public static void duplicatesAdd(DMatrixSparseCSC A, @Nullable @Nullable IGrowArray work)
      For each duplicate element in the matrix it will remove the duplicates and replace them with a single element that is the sum of all the duplicates. The result will be a valid matrix. This is done inplace and does not require the matrix to be initially sorted.
      Parameters:
      A - (Input/Output) input matrix. Modified.
      work - Nullable. Internal workspace array.
    • multRows

      public static void multRows(double[] diag, int offset, DMatrixSparseCSC A)
      Multiply all elements of row 'i' by value[i]. A[i,:] *= values[i]
      Parameters:
      diag - (Input) multiplication factors
      offset - (Input) First index in values
      A - (Input/Output) Matrix. Modified.
    • divideRows

      public static void divideRows(double[] diag, int offset, DMatrixSparseCSC A)
      Divides all elements of row 'i' by value[i]. A[i,:] /= values[i]
      Parameters:
      diag - (Input) division factors
      offset - (Input) First index in values
      A - (Input/Output) Matrix. Modified.
    • trace

      public static double trace(DMatrixSparseCSC A)

      This computes the trace of the matrix:

      trace = ∑i=1:n { aii }
      where n = min(numRows,numCols)

      Parameters:
      A - (Input) Matrix. Not modified.
    • apply

      public static DMatrixSparseCSC apply(DMatrixSparseCSC input, DOperatorUnary func, @Nullable @Nullable DMatrixSparseCSC output)
      This applies a given unary function on every value stored in the matrix B = f(A). A and B can be the same instance.
      Parameters:
      input - (Input) input matrix. Not modified
      func - Unary function accepting a double
      output - (Output) Matrix. Modified.
      Returns:
      The output matrix
    • apply

      public static DMatrixSparseCSC apply(DMatrixSparseCSC input, DOperatorUnary func)
    • applyRowIdx

      public static DMatrixSparseCSC applyRowIdx(DMatrixSparseCSC input, DOperatorBinaryIdx func, @Nullable @Nullable DMatrixSparseCSC output)
      This applies a given unary function on every nz row,value stored in the matrix B = f(A). A and B can be the same instance.
      Parameters:
      input - (Input) input matrix. Not modified
      func - Binary function accepting (row-index, value)
      output - (Output) Matrix. Modified.
      Returns:
      The output matrix
    • applyColumnIdx

      public static DMatrixSparseCSC applyColumnIdx(DMatrixSparseCSC input, DOperatorBinaryIdx func, @Nullable @Nullable DMatrixSparseCSC output)
      This applies a given unary function on every non-zero column, value stored in the matrix B = f(A). A and B can be the same instance.
      Parameters:
      input - (Input) input matrix. Not modified
      func - Binary function accepting (column-index, value)
      output - (Output) Matrix. Modified.
      Returns:
      The output matrix
    • reduceScalar

      public static double reduceScalar(DMatrixSparseCSC input, double initValue, DOperatorBinary func)
      This accumulates the matrix values to a scalar value.
       result = initialValue
       for all (i,j) result = func( result, A[i,j] )
       
      Parameters:
      input - (Input) input matrix. Not modified
      initValue - initial value for accumulator
      func - Accumulator function defining "+" for accumulator += cellValue
      Returns:
      accumulated value
    • reduceScalar

      public static double reduceScalar(DMatrixSparseCSC input, DOperatorBinary func)
    • reduceColumnWise

      public static DMatrixRMaj reduceColumnWise(DMatrixSparseCSC input, double initValue, DOperatorBinary func, @Nullable @Nullable DMatrixRMaj output, @Nullable @Nullable Mask mask)
      This accumulates the values per column to a scalar value
       for each column `j`,
         result = initialValue
         for each row 'i' result = func( result, A[i,j] )
         output[j] = result
       
      Parameters:
      input - (Input) input matrix. Not modified
      initValue - initial value for accumulator
      func - Accumulator function defining "+" for accumulator += cellValue
      output - output (Output) Vector, where result can be stored in
      mask - (Optional) Mask for specifying which entries should be overwritten
      Returns:
      a column-vector, where v[i] == values of column i reduced to scalar based on `func`
    • reduceRowWise

      public static DMatrixRMaj reduceRowWise(DMatrixSparseCSC input, double initValue, DOperatorBinary func, @Nullable @Nullable DMatrixRMaj output)
      This accumulates the values per row to a scalar value
       for each column `j`,
         result = initialValue
         for each row 'i' result = func( result, A[i,j] )
         output[j] = result
       
      Parameters:
      input - (Input) input matrix. Not modified
      initValue - initial value for accumulator
      func - Accumulator function defining "+" for accumulator += cellValue
      output - output (Output) Vector, where result can be stored in
      Returns:
      a row-vector, where v[i] == values of row i reduced to scalar based on `func`