Class CommonOps_FSCC
FMatrixSparseCSC
. For example, addition, matrix multiplication, ... etc-
Method Summary
Modifier and TypeMethodDescriptionstatic FMatrixSparseCSC
add
(float alpha, FMatrixSparseCSC A, float beta, FMatrixSparseCSC B, @Nullable FMatrixSparseCSC outputC, @Nullable IGrowArray gw, @Nullable FGrowArray gx) Performs matrix addition:
C = αA + βBstatic FMatrixSparseCSC
apply
(FMatrixSparseCSC input, FOperatorUnary func) static FMatrixSparseCSC
apply
(FMatrixSparseCSC input, FOperatorUnary func, @Nullable FMatrixSparseCSC output) This applies a given unary function on every value stored in the matrix B = f(A).static FMatrixSparseCSC
applyColumnIdx
(FMatrixSparseCSC input, FOperatorBinaryIdx func, @Nullable FMatrixSparseCSC output) This applies a given unary function on every non-zero column, value stored in the matrix B = f(A).static FMatrixSparseCSC
applyRowIdx
(FMatrixSparseCSC input, FOperatorBinaryIdx func, @Nullable FMatrixSparseCSC output) This applies a given unary function on every nz row,value stored in the matrix B = f(A).static void
changeSign
(FMatrixSparseCSC A, FMatrixSparseCSC outputB) B = -A.static boolean
Checks for duplicate elements.static boolean
Checks to see if row indicies are sorted into ascending order.static boolean
static boolean
static FMatrixSparseCSC
concatColumns
(FMatrixSparseCSC left, FMatrixSparseCSC right, @Nullable FMatrixSparseCSC out) Concats two matrices along their columns (horizontal).static FMatrixSparseCSC
concatRows
(FMatrixSparseCSC top, FMatrixSparseCSC bottom, @Nullable FMatrixSparseCSC out) Concats two matrices along their rows (vertical).static float
Returns the determinant of the matrix.static FMatrixSparseCSC
diag
(float... values) Returns a diagonal matrix with the specified diagonal elements.static FMatrixSparseCSC
diag
(@Nullable FMatrixSparseCSC A, float[] values, int offset, int length) Creates a diagonal matrix from an array.static void
divide
(float scalar, FMatrixSparseCSC A, FMatrixSparseCSC outputB) B = scalar/A.static void
divide
(FMatrixSparseCSC A, float scalar, FMatrixSparseCSC outputB) B = A/scalar.static void
divideColumns
(FMatrixSparseCSC A, float[] values, int offset) Divides all elements of column 'i' by values[i].static void
divideRows
(float[] diag, int offset, FMatrixSparseCSC A) Divides all elements of row 'i' by value[i].static void
divideRowsCols
(float[] diagA, int offsetA, FMatrixSparseCSC B, float[] diagC, int offsetC) Equivalent to multiplying a matrix B by the inverse of two diagonal matrices.static float
dotInnerColumns
(FMatrixSparseCSC A, int colA, FMatrixSparseCSC B, int colB, @Nullable IGrowArray gw, @Nullable FGrowArray gx) Computes the inner product of two column vectors taken from the input matrices.static void
duplicatesAdd
(FMatrixSparseCSC A, @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.static float
Returns the value of the element with the largest valuestatic float
Returns the value of the element with the largest abs()static float
Returns the value of the element with the minimum valuestatic float
Returns the value of the element with the smallest abs()static FMatrixSparseCSC
elementMult
(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable FMatrixSparseCSC output, @Nullable IGrowArray gw, @Nullable FGrowArray gx) Performs an element-wise multiplication.
output[i,j] = A[i,j]*B[i,j]
All matrices must have the same shape.static float
Sum of all elementsstatic void
extract
(FMatrixSparseCSC src, int srcY0, int srcY1, int srcX0, int srcX1, FMatrixSparseCSC dst, int dstY0, int dstX0) Extracts a submatrix from 'src' and inserts it in a submatrix in 'dst'.static FMatrixSparseCSC
extractColumn
(FMatrixSparseCSC A, int column, @Nullable FMatrixSparseCSC out) Extracts a column from A and stores it into out.static void
extractDiag
(FMatrixSparseCSC A, FMatrixRMaj outputB) Extracts the diagonal elements 'src' write it to the 'dst' vector.static void
extractDiag
(FMatrixSparseCSC A, FMatrixSparseCSC outputB) Extracts the diagonal elements 'src' write it to the 'dst' vector.static FMatrixSparseCSC
extractRows
(FMatrixSparseCSC A, int row0, int row1, @Nullable FMatrixSparseCSC out) Creates a submatrix by extracting the specified rows from A.static void
fill
(FMatrixSparseCSC A, float value) Sets every element in the matrix to the specified value.static FMatrixSparseCSC
identity
(int length) static FMatrixSparseCSC
identity
(int numRows, int numCols) static boolean
invert
(FMatrixSparseCSC A, FMatrixRMaj inverse) Performs a matrix inversion operation that does not modify the original and stores the results in another matrix.static void
maxAbsCols
(FMatrixSparseCSC A, @Nullable FMatrixRMaj outputB) Finds the maximum abs in each column of A and stores it into valuesstatic FMatrixRMaj
maxCols
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output) Computes the maximums of each column in the input matrix and returns the results in a vector:
bj = max(i=1:m ; aij)static FMatrixRMaj
maxRows
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output, @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)static FMatrixRMaj
minCols
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output) Computes the minimum of each column in the input matrix and returns the results in a vector:
bj = min(i=1:m ; aij)static FMatrixRMaj
minRows
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output, @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)static FMatrixRMaj
mult
(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable FMatrixRMaj outputC) Performs matrix multiplication.static FMatrixSparseCSC
mult
(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable FMatrixSparseCSC outputC) static FMatrixSparseCSC
mult
(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable FMatrixSparseCSC outputC, @Nullable IGrowArray gw, @Nullable FGrowArray gx) Performs matrix multiplication.static void
multAdd
(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC) C = C + A*Bstatic void
multAddTransA
(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC, @Nullable FGrowArray work) C = C + AT*Bstatic void
multAddTransAB
(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC) C = C + AT*BTstatic void
multAddTransB
(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC, @Nullable FGrowArray work) C = C + A*BTstatic void
multColumns
(FMatrixSparseCSC A, float[] values, int offset) Multiply all elements of column 'i' by value[i].static void
multRows
(float[] diag, int offset, FMatrixSparseCSC A) Multiply all elements of row 'i' by value[i].static void
multRowsCols
(float[] diagA, int offsetA, FMatrixSparseCSC B, float[] diagC, int offsetC) Equivalent to multiplying a matrix B by two diagonal matrices.static FMatrixRMaj
multTransA
(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable FMatrixRMaj outputC, @Nullable FGrowArray work) Performs matrix multiplication.static FMatrixRMaj
multTransAB
(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC) Performs matrix multiplication.static FMatrixRMaj
multTransB
(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable FMatrixRMaj outputC, @Nullable FGrowArray work) Performs matrix multiplication.static int[]
permutationInverse
(int[] original, int length) static void
permutationInverse
(int[] original, int[] inverse, int length) Computes the inverse permutation vectorstatic FMatrixSparseCSC
permutationMatrix
(int[] p, boolean inverse, int N, @Nullable FMatrixSparseCSC P) Converts the permutation vector into a matrix.static void
permutationVector
(FMatrixSparseCSC P, int[] vector) Converts the permutation matrix into a vectorstatic void
permute
(@org.jetbrains.annotations.Nullable int[] permRowInv, FMatrixSparseCSC input, @org.jetbrains.annotations.Nullable int[] permCol, FMatrixSparseCSC 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.static void
permute
(int[] perm, float[] input, float[] output, int N) Permutes a vector.static void
permuteInv
(int[] perm, float[] input, float[] output, int N) Permutes a vector in the inverse.static void
permuteRowInv
(int[] permInv, FMatrixSparseCSC input, FMatrixSparseCSC output) Applies the row permutation specified by the vector to the input matrix and save the results in the output matrix.static void
permuteSymmetric
(FMatrixSparseCSC input, int[] permInv, FMatrixSparseCSC output, @Nullable IGrowArray gw) Applies the permutation to upper triangular symmetric matrices.static FMatrixRMaj
reduceColumnWise
(FMatrixSparseCSC input, float initValue, FOperatorBinary func, @Nullable FMatrixRMaj output, @Nullable Mask mask) This accumulates the values per column to a scalar valuestatic FMatrixRMaj
reduceRowWise
(FMatrixSparseCSC input, float initValue, FOperatorBinary func, @Nullable FMatrixRMaj output) This accumulates the values per row to a scalar valuestatic float
reduceScalar
(FMatrixSparseCSC input, float initValue, FOperatorBinary func) This accumulates the matrix values to a scalar value.static float
reduceScalar
(FMatrixSparseCSC input, FOperatorBinary func) static void
removeZeros
(FMatrixSparseCSC A, float tol) Removes all elements from the matrix that are > tol.static void
removeZeros
(FMatrixSparseCSC input, FMatrixSparseCSC output, float tol) Copies all elements from input into output which are > tol.static void
scale
(float scalar, FMatrixSparseCSC A, FMatrixSparseCSC outputB) B = scalar*A.static FMatrixSparseCSC
select
(FMatrixSparseCSC A, IPredicateBinary selector, @Nullable FMatrixSparseCSC output) Select entries from A and save them in C.static void
static boolean
solve
(FMatrixSparseCSC a, FMatrixRMaj b, FMatrixRMaj x) Solves for x in the following equation:
A*x = bstatic boolean
Solves for x in the following equation:
A*x = bstatic FMatrixRMaj
sumCols
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output) Computes the sum of each column in the input matrix and returns the results in a vector:
bj = sum(i=1:m ; aij)static FMatrixRMaj
sumRows
(FMatrixSparseCSC input, @Nullable FMatrixRMaj output) Computes the sum of each row in the input matrix and returns the results in a vector:
bj = sum(i=1:n ; aji)static void
symmLowerToFull
(FMatrixSparseCSC A, FMatrixSparseCSC outputB, @Nullable IGrowArray gw) Given a symmetric matrix, which is represented by a lower triangular matrix, convert it back into a full symmetric matrixstatic float
This computes the trace of the matrix:
trace = ∑i=1:n { aii }
where n = min(numRows,numCols)static FMatrixSparseCSC
transpose
(FMatrixSparseCSC A, @Nullable FMatrixSparseCSC A_t, @Nullable IGrowArray gw) Perform matrix transposestatic void
zero
(FMatrixSparseCSC A, int row0, int row1, int col0, int col1) Zeros an inner rectangle inside the matrix.
-
Method Details
-
checkIndicesSorted
Checks to see if row indicies are sorted into ascending order. O(N)- Returns:
- true if sorted and false if not
-
checkStructure
-
checkSortedFlag
-
checkDuplicateElements
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 FMatrixSparseCSC transpose(FMatrixSparseCSC A, @Nullable @Nullable FMatrixSparseCSC A_t, @Nullable @Nullable IGrowArray gw) Perform matrix transpose- Parameters:
A
- Input matrix. Not modifiedA_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 FMatrixSparseCSC mult(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable @Nullable FMatrixSparseCSC outputC) -
mult
public static FMatrixSparseCSC mult(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable @Nullable FMatrixSparseCSC outputC, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable FGrowArray 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 FMatrixRMaj mult(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable @Nullable FMatrixRMaj outputC) Performs matrix multiplication. C = A*B- Parameters:
A
- MatrixB
- Dense MatrixoutputC
- Dense Matrix
-
multAdd
C = C + A*B
-
multTransA
public static FMatrixRMaj multTransA(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable @Nullable FMatrixRMaj outputC, @Nullable @Nullable FGrowArray work) Performs matrix multiplication. C = AT*B- Parameters:
A
- MatrixB
- Dense MatrixoutputC
- Dense Matrix
-
multAddTransA
public static void multAddTransA(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC, @Nullable @Nullable FGrowArray work) C = C + AT*B
-
multTransB
public static FMatrixRMaj multTransB(FMatrixSparseCSC A, FMatrixRMaj B, @Nullable @Nullable FMatrixRMaj outputC, @Nullable @Nullable FGrowArray work) Performs matrix multiplication. C = A*BT- Parameters:
A
- MatrixB
- Dense MatrixoutputC
- Dense Matrix
-
multAddTransB
public static void multAddTransB(FMatrixSparseCSC A, FMatrixRMaj B, FMatrixRMaj outputC, @Nullable @Nullable FGrowArray work) C = C + A*BT
-
multTransAB
Performs matrix multiplication. C = AT*BT- Parameters:
A
- MatrixB
- Dense MatrixoutputC
- Dense Matrix
-
multAddTransAB
C = C + AT*BT
-
symmLowerToFull
public static void symmLowerToFull(FMatrixSparseCSC A, FMatrixSparseCSC 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 matrixoutputB
- (Output) Symmetric matrix.gw
- (Optional) Workspace. Can be null.
-
add
public static FMatrixSparseCSC add(float alpha, FMatrixSparseCSC A, float beta, FMatrixSparseCSC B, @Nullable @Nullable FMatrixSparseCSC outputC, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable FGrowArray gx) Performs matrix addition:
C = αA + βB- Parameters:
alpha
- scalar value multiplied against AA
- Matrixbeta
- scalar value multiplied against BB
- MatrixoutputC
- Output matrix.gw
- (Optional) Storage for internal workspace. Can be null.gx
- (Optional) Storage for internal workspace. Can be null.
-
identity
-
identity
-
setIdentity
-
scale
B = scalar*A. A and B can be the same instance.- Parameters:
scalar
- (Input) Scalar valueA
- (Input) Matrix. Not modified.outputB
- (Output) Matrix. Modified.
-
divide
B = A/scalar. A and B can be the same instance.- Parameters:
scalar
- (Input) Scalar valueA
- (Input) Matrix. Not modified.outputB
- (Output) Matrix. Modified.
-
divide
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 valueoutputB
- (Output) Matrix. Modified.
-
changeSign
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
Returns the value of the element with the smallest abs()- Parameters:
A
- (Input) Matrix. Not modified.- Returns:
- scalar
-
elementMaxAbs
Returns the value of the element with the largest abs()- Parameters:
A
- (Input) Matrix. Not modified.- Returns:
- scalar
-
elementMin
Returns the value of the element with the minimum value- Parameters:
A
- (Input) Matrix. Not modified.- Returns:
- scalar
-
elementMax
Returns the value of the element with the largest value- Parameters:
A
- (Input) Matrix. Not modified.- Returns:
- scalar
-
elementSum
Sum of all elements- Parameters:
A
- (Input) Matrix. Not modified.- Returns:
- scalar
-
elementMult
public static FMatrixSparseCSC elementMult(FMatrixSparseCSC A, FMatrixSparseCSC B, @Nullable @Nullable FMatrixSparseCSC output, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable FGrowArray 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) Matrixoutput
- (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
Finds the maximum abs in each column of A and stores it into values- Parameters:
A
- (Input) MatrixoutputB
- (Output) storage for column max abs
-
multColumns
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 columnoffset
- (Input) first index in values to start at
-
divideColumns
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 columnoffset
- (Input) first index in values to start at
-
multRowsCols
public static void multRowsCols(float[] diagA, int offsetA, FMatrixSparseCSC B, float[] 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.numRowsoffsetA
- First index in AB
- Rectangular matrixdiagC
- Array of length indexC + B.numColsoffsetC
- First index in C
-
divideRowsCols
public static void divideRowsCols(float[] diagA, int offsetA, FMatrixSparseCSC B, float[] 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.numRowsoffsetA
- First index in AB
- Rectangular matrixdiagC
- Array of length indexC + B.numColsoffsetC
- First index in C
-
diag
Returns a diagonal matrix with the specified diagonal elements.- Parameters:
values
- values of diagonal elements- Returns:
- A diagonal matrix
-
diag
public static FMatrixSparseCSC diag(@Nullable @Nullable FMatrixSparseCSC A, float[] 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 matirxlength
- Length of the diagonal matrixoffset
- First index in values- Returns:
- The diagonal matrix
-
extractDiag
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
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 FMatrixSparseCSC permutationMatrix(int[] p, boolean inverse, int N, @Nullable @Nullable FMatrixSparseCSC P) Converts the permutation vector into a matrix. B = P*A. B[p[i],:] = A[i,:]- Parameters:
p
- (Input) Permutation vectorinverse
- (Input) If it is the inverse. B[i,:] = A[p[i],:)P
- (Output) Permutation matrix
-
permutationVector
Converts the permutation matrix into a vector- Parameters:
P
- (Input) Permutation matrixvector
- (Output) Permutation vector
-
permutationInverse
public static void permutationInverse(int[] original, int[] inverse, int length) Computes the inverse permutation vector- Parameters:
original
- Original permutation vectorinverse
- It's inverse
-
permutationInverse
public static int[] permutationInverse(int[] original, int length) -
permuteRowInv
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 permutedoutput
- (Output) Matrix which has the permutation stored in it. Is reshaped.
-
permute
public static void permute(@Nullable @org.jetbrains.annotations.Nullable int[] permRowInv, FMatrixSparseCSC input, @Nullable @org.jetbrains.annotations.Nullable int[] permCol, FMatrixSparseCSC 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 permutedpermCol
- (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, float[] input, float[] output, int N) Permutes a vector. output[i] = input[perm[i]]- Parameters:
perm
- (Input) permutation vectorinput
- (Input) Vector which is to be permutedoutput
- (Output) Where the permuted vector is stored.N
- Number of elements in the vector.
-
permuteInv
public static void permuteInv(int[] perm, float[] input, float[] output, int N) Permutes a vector in the inverse. output[perm[k]] = input[k]- Parameters:
perm
- (Input) permutation vectorinput
- (Input) Vector which is to be permutedoutput
- (Output) Where the permuted vector is stored.N
- Number of elements in the vector.
-
permuteSymmetric
public static void permuteSymmetric(FMatrixSparseCSC input, int[] permInv, FMatrixSparseCSC 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 FMatrixSparseCSC concatRows(FMatrixSparseCSC top, FMatrixSparseCSC bottom, @Nullable @Nullable FMatrixSparseCSC out) Concats two matrices along their rows (vertical).- Parameters:
top
- Matrix on the topbottom
- Matrix on the bototmout
- (Output) (Optional) Storage for combined matrix. Resized.- Returns:
- Combination of the two matrices
-
concatColumns
public static FMatrixSparseCSC concatColumns(FMatrixSparseCSC left, FMatrixSparseCSC right, @Nullable @Nullable FMatrixSparseCSC out) Concats two matrices along their columns (horizontal).- Parameters:
left
- Matrix on the leftright
- Matrix on the rightout
- (Output) (Optional) Storage for combined matrix. Resized.- Returns:
- Combination of the two matrices
-
extractColumn
public static FMatrixSparseCSC extractColumn(FMatrixSparseCSC A, int column, @Nullable @Nullable FMatrixSparseCSC out) Extracts a column from A and stores it into out.- Parameters:
A
- (Input) Source matrix. not modified.column
- The column in Aout
- (Output, Optional) Storage for column vector- Returns:
- The column of A.
-
extractRows
public static FMatrixSparseCSC extractRows(FMatrixSparseCSC A, int row0, int row1, @Nullable @Nullable FMatrixSparseCSC out) Creates a submatrix by extracting the specified rows from A. rows = {row0 %le; i %le; row1}.- Parameters:
A
- (Input) matrixrow0
- First row. Inclusiverow1
- Last row+1.out
- (Output, Option) Storage for output matrix- Returns:
- The submatrix
-
extract
public static void extract(FMatrixSparseCSC src, int srcY0, int srcY1, int srcX0, int srcX1, FMatrixSparseCSC 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 FMatrixSparseCSC select(FMatrixSparseCSC A, IPredicateBinary selector, @Nullable @Nullable FMatrixSparseCSC 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 selectedoutput
- (Optional/Output) Matrix to use for the output. Can be the same as A- Returns:
- Matrix storing the selected entries of A
-
fill
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
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 matrixoutput
- Optional storage for output. Reshaped into a row vector. Modified.- Returns:
- Vector containing the sum of each column
-
minCols
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 matrixoutput
- Optional storage for output. Reshaped into a row vector. Modified.- Returns:
- Vector containing the minimums of each column
-
maxCols
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 matrixoutput
- Optional storage for output. Reshaped into a row vector. Modified.- Returns:
- Vector containing the maximums of each column
-
sumRows
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 matrixoutput
- Optional storage for output. Reshaped into a column vector. Modified.- Returns:
- Vector containing the sum of each row
-
minRows
public static FMatrixRMaj minRows(FMatrixSparseCSC input, @Nullable @Nullable FMatrixRMaj 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 matrixoutput
- Optional storage for output. Reshaped into a column vector. Modified.gw
- work space- Returns:
- Vector containing the minimum of each row
-
maxRows
public static FMatrixRMaj maxRows(FMatrixSparseCSC input, @Nullable @Nullable FMatrixRMaj 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 matrixoutput
- Optional storage for output. Reshaped into a column vector. Modified.gw
- work space- Returns:
- Vector containing the maximum of each row
-
zero
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 float dotInnerColumns(FMatrixSparseCSC A, int colA, FMatrixSparseCSC B, int colB, @Nullable @Nullable IGrowArray gw, @Nullable @Nullable FGrowArray gx) Computes the inner product of two column vectors taken from the input matrices.dot = A(:,colA)'*B(:,colB)
- Parameters:
A
- MatrixcolA
- Column in AB
- MatrixcolB
- Column in B- Returns:
- Dot product
-
solve
Solves for x in the following equation:
A*x = bIf 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_FSCC
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
Solves for x in the following equation:
A*x = bIf 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_FSCC
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
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-1If 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_FSCC
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
Returns the determinant of the matrix. If the inverse of the matrix is also needed, then usingLUDecomposition_F32
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
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
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
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
Multiply all elements of row 'i' by value[i]. A[i,:] *= values[i]- Parameters:
diag
- (Input) multiplication factorsoffset
- (Input) First index in valuesA
- (Input/Output) Matrix. Modified.
-
divideRows
Divides all elements of row 'i' by value[i]. A[i,:] /= values[i]- Parameters:
diag
- (Input) division factorsoffset
- (Input) First index in valuesA
- (Input/Output) Matrix. Modified.
-
trace
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 FMatrixSparseCSC apply(FMatrixSparseCSC input, FOperatorUnary func, @Nullable @Nullable FMatrixSparseCSC 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 modifiedfunc
- Unary function accepting a floatoutput
- (Output) Matrix. Modified.- Returns:
- The output matrix
-
apply
-
applyRowIdx
public static FMatrixSparseCSC applyRowIdx(FMatrixSparseCSC input, FOperatorBinaryIdx func, @Nullable @Nullable FMatrixSparseCSC 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 modifiedfunc
- Binary function accepting (row-index, value)output
- (Output) Matrix. Modified.- Returns:
- The output matrix
-
applyColumnIdx
public static FMatrixSparseCSC applyColumnIdx(FMatrixSparseCSC input, FOperatorBinaryIdx func, @Nullable @Nullable FMatrixSparseCSC 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 modifiedfunc
- Binary function accepting (column-index, value)output
- (Output) Matrix. Modified.- Returns:
- The output matrix
-
reduceScalar
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 modifiedinitValue
- initial value for accumulatorfunc
- Accumulator function defining "+" for accumulator += cellValue- Returns:
- accumulated value
-
reduceScalar
-
reduceColumnWise
public static FMatrixRMaj reduceColumnWise(FMatrixSparseCSC input, float initValue, FOperatorBinary func, @Nullable @Nullable FMatrixRMaj output, @Nullable @Nullable Mask mask) This accumulates the values per column to a scalar valuefor each column `j`, result = initialValue for each row 'i' result = func( result, A[i,j] ) output[j] = result
- Parameters:
input
- (Input) input matrix. Not modifiedinitValue
- initial value for accumulatorfunc
- Accumulator function defining "+" for accumulator += cellValueoutput
- output (Output) Vector, where result can be stored inmask
- (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 FMatrixRMaj reduceRowWise(FMatrixSparseCSC input, float initValue, FOperatorBinary func, @Nullable @Nullable FMatrixRMaj output) This accumulates the values per row to a scalar valuefor each column `j`, result = initialValue for each row 'i' result = func( result, A[i,j] ) output[j] = result
- Parameters:
input
- (Input) input matrix. Not modifiedinitValue
- initial value for accumulatorfunc
- Accumulator function defining "+" for accumulator += cellValueoutput
- 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`
-