|
Parallel Colt 0.6.1 | |||||||||
| PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
| SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD | |||||||||
java.lang.Objectcern.colt.PersistentObject
cern.colt.matrix.AbstractMatrix
cern.colt.matrix.AbstractMatrix2D
cern.colt.matrix.tdouble.DoubleMatrix2D
cern.colt.matrix.tdouble.impl.SparseDoubleMatrix2D
public class SparseDoubleMatrix2D
Sparse hashed 2-d matrix holding double elements. First see the package summary and javadoc tree view to get the broad picture.
Implementation:
Note that this implementation is not synchronized. Uses a
OpenIntDoubleHashMap, which is a compact and
performant hashing technique.
Memory requirements:
Cells that
trimToSize().
worst case: memory [bytes] = (1/minLoadFactor) * nonZeros * 13.
best case: memory [bytes] = (1/maxLoadFactor) * nonZeros * 13.
Where nonZeros = cardinality() is the number of non-zero cells.
Thus, a 1000 x 1000 matrix with minLoadFactor=0.25 and maxLoadFactor=0.5 and
1000000 non-zero cells consumes between 25 MB and 50 MB. The same 1000 x 1000
matrix with 1000 non-zero cells consumes between 25 and 50 KB.
Time complexity:
This class offers expected time complexity O(1) (i.e.
constant time) for the basic operations get, getQuick,
set, setQuick and size assuming the hash
function disperses the elements properly among the buckets. Otherwise,
pathological cases, although highly improbable, can occur, degrading
performance to O(N) in the worst case. As such this sparse class
is expected to have no worse time complexity than its dense counterpart
DenseDoubleMatrix2D. However, constant factors are considerably
larger.
Cells are internally addressed in row-major. Performance sensitive applications can exploit this fact. Setting values in a loop row-by-row is quicker than column-by-column, because fewer hash collisions occur. Thus
for (int row = 0; row < rows; row++) {
for (int column = 0; column < columns; column++) {
matrix.setQuick(row, column, someValue);
}
}
is quicker than
for (int column = 0; column < columns; column++) {
for (int row = 0; row < rows; row++) {
matrix.setQuick(row, column, someValue);
}
}
cern.colt.map,
OpenIntDoubleHashMap,
Serialized Form| Field Summary |
|---|
| Fields inherited from class cern.colt.PersistentObject |
|---|
serialVersionUID |
| Constructor Summary | |
|---|---|
SparseDoubleMatrix2D(double[][] values)
Constructs a matrix with a copy of the given values. |
|
SparseDoubleMatrix2D(int rows,
int columns)
Constructs a matrix with a given number of rows and columns and default memory usage. |
|
SparseDoubleMatrix2D(int rows,
int columns,
int initialCapacity,
double minLoadFactor,
double maxLoadFactor)
Constructs a matrix with a given number of rows and columns using memory as specified. |
|
| Method Summary | |
|---|---|
DoubleMatrix2D |
assign(double value)
Sets all cells to the state specified by value. |
DoubleMatrix2D |
assign(DoubleFunction function)
Assigns the result of a function to each cell; x[row,col] = function(x[row,col]). |
DoubleMatrix2D |
assign(DoubleMatrix2D source)
Replaces all cell values of the receiver with the values of another matrix. |
DoubleMatrix2D |
assign(DoubleMatrix2D y,
DoubleDoubleFunction function)
Assigns the result of a function to each cell; x[row,col] = function(x[row,col],y[row,col]). |
int |
cardinality()
Returns the number of cells having non-zero values. |
void |
dct2(boolean scale)
Computes the 2D discrete cosine transform (DCT-II) of this matrix. |
void |
dctColumns(boolean scale)
Computes the discrete cosine transform (DCT-II) of each column of this matrix. |
void |
dctRows(boolean scale)
Computes the discrete cosine transform (DCT-II) of each row of this matrix. |
void |
dht2()
Computes the 2D discrete Hartley transform (DHT) of this matrix. |
void |
dhtColumns()
Computes the discrete Hartley transform (DHT) of each column of this matrix. |
void |
dhtRows()
Computes the discrete Hartley transform (DHT) of each row of this matrix. |
void |
dst2(boolean scale)
Computes the 2D discrete sine transform (DST-II) of this matrix. |
void |
dstColumns(boolean scale)
Computes the discrete sine transform (DST-II) of each column of this matrix. |
void |
dstRows(boolean scale)
Computes the discrete sine transform (DST-II) of each row of this matrix. |
AbstractIntDoubleMap |
elements()
Returns the elements of this matrix. |
void |
ensureCapacity(int minCapacity)
Ensures that the receiver can hold at least the specified number of non-zero cells without needing to allocate new internal memory. |
void |
fft2()
Computes the 2D discrete Fourier transform (DFT) of this matrix. |
DoubleMatrix2D |
forEachNonZero(IntIntDoubleFunction function)
Assigns the result of a function to each non-zero cell; x[row,col] = function(x[row,col]). |
DComplexMatrix2D |
getFft2()
Returns new complex matrix which is the 2D discrete Fourier transform (DFT) of this matrix. |
DComplexMatrix2D |
getFftColumns()
Returns new complex matrix which is the discrete Fourier transform (DFT) of each column of this matrix. |
DComplexMatrix2D |
getFftRows()
Returns new complex matrix which is the discrete Fourier transform (DFT) of each row of this matrix. |
DComplexMatrix2D |
getIfft2(boolean scale)
Returns new complex matrix which is the 2D inverse of the discrete Fourier transform (IDFT) of this matrix. |
DComplexMatrix2D |
getIfftColumns(boolean scale)
Returns new complex matrix which is the inverse of the discrete Fourier transform (IDFT) of each column of this matrix. |
DComplexMatrix2D |
getIfftRows(boolean scale)
Returns new complex matrix which is the inverse of the discrete Fourier transform (IDFT) of each row of this matrix. |
double |
getQuick(int row,
int column)
Returns the matrix cell value at coordinate [row,column]. |
void |
idct2(boolean scale)
Computes the 2D inverse of the discrete cosine transform (DCT-III) of this matrix. |
void |
idctColumns(boolean scale)
Computes the inverse of the discrete cosine transform (DCT-III) of each column of this matrix. |
void |
idctRows(boolean scale)
Computes the inverse of the discrete cosine transform (DCT-III) of each row of this matrix. |
void |
idht2(boolean scale)
Computes the 2D inverse of the discrete Hartley transform (IDHT) of this matrix. |
void |
idhtColumns(boolean scale)
Computes the inverse of the discrete Hartley transform (IDHT) of each column of this matrix. |
void |
idhtRows(boolean scale)
Computes the inverse of the discrete Hartley transform (IDHT) of each row of this matrix. |
void |
idst2(boolean scale)
Computes the 2D inverse of the discrete sine transform (DST-III) of this matrix. |
void |
idstColumns(boolean scale)
Computes the inverse of the discrete sine transform (DST-III) of each column of this matrix. |
void |
idstRows(boolean scale)
Computes the inverse of the discrete sine transform (DST-III) of each row of this matrix. |
void |
ifft2(boolean scale)
Computes the 2D inverse of the discrete Fourier transform (IDFT) of this matrix. |
int |
index(int row,
int column)
Returns the position of the given coordinate within the (virtual or non-virtual) internal 1-dimensional array. |
DoubleMatrix2D |
like(int rows,
int columns)
Construct and returns a new empty matrix of the same dynamic type as the receiver, having the specified number of rows and columns. |
DoubleMatrix1D |
like1D(int size)
Construct and returns a new 1-d matrix of the corresponding dynamic type, entirelly independent of the receiver. |
void |
setQuick(int row,
int column,
double value)
Sets the matrix cell at coordinate [row,column] to the specified value. |
void |
trimToSize()
Releases any superfluous memory created by explicitly putting zero values into cells formerly having non-zero values; An application can use this operation to minimize the storage of the receiver. |
DoubleMatrix1D |
vectorize()
Returns a vector obtained by stacking the columns of the matrix on top of one another. |
DoubleMatrix1D |
zMult(DoubleMatrix1D y,
DoubleMatrix1D z,
double alpha,
double beta,
boolean transposeA)
Linear algebraic matrix-vector multiplication; z = alpha * A * y + beta*z. |
DoubleMatrix2D |
zMult(DoubleMatrix2D B,
DoubleMatrix2D C,
double alpha,
double beta,
boolean transposeA,
boolean transposeB)
Linear algebraic matrix-matrix multiplication; C = alpha * A x B + beta*C. |
| Methods inherited from class cern.colt.matrix.tdouble.DoubleMatrix2D |
|---|
aggregate, aggregate, aggregate, aggregate, assign, assign, assign, assign, assign, assign, copy, equals, equals, get, getMaxLocation, getMinLocation, getNegativeValues, getNonZeros, getPositiveValues, like, normalize, set, setNcolumns, setNrows, toArray, toString, viewColumn, viewColumnFlip, viewDice, viewPart, viewRow, viewRowFlip, viewSelection, viewSelection, viewSelection, viewSorted, viewStrides, zAssign8Neighbors, zMult, zMult, zSum |
| Methods inherited from class cern.colt.matrix.AbstractMatrix2D |
|---|
checkShape, checkShape, columns, columnStride, rows, rowStride, size, toStringShort |
| Methods inherited from class cern.colt.matrix.AbstractMatrix |
|---|
isView |
| Methods inherited from class cern.colt.PersistentObject |
|---|
clone |
| Methods inherited from class java.lang.Object |
|---|
getClass, hashCode, notify, notifyAll, wait, wait, wait |
| Constructor Detail |
|---|
public SparseDoubleMatrix2D(double[][] values)
The values are copied. So subsequent changes in values are not reflected in the matrix, and vice-versa.
values - The values to be filled into the new matrix.
IllegalArgumentException - if
for any 1 <= row < values.length: values[row].length != values[row-1].length.
public SparseDoubleMatrix2D(int rows,
int columns)
rows - the number of rows the matrix shall have.columns - the number of columns the matrix shall have.
IllegalArgumentException - if
rows<0 || columns<0 || (double)columns*rows > Integer.MAX_VALUE.
public SparseDoubleMatrix2D(int rows,
int columns,
int initialCapacity,
double minLoadFactor,
double maxLoadFactor)
OpenIntDoubleHashMap.
rows - the number of rows the matrix shall have.columns - the number of columns the matrix shall have.initialCapacity - the initial capacity of the hash map. If not known, set
initialCapacity=0 or small.minLoadFactor - the minimum load factor of the hash map.maxLoadFactor - the maximum load factor of the hash map.
IllegalArgumentException - if
initialCapacity < 0 || (minLoadFactor < 0.0 || minLoadFactor >= 1.0) || (maxLoadFactor <= 0.0 || maxLoadFactor >= 1.0) || (minLoadFactor >= maxLoadFactor).
IllegalArgumentException - if
rows<0 || columns<0 || (double)columns*rows > Integer.MAX_VALUE.| Method Detail |
|---|
public DoubleMatrix2D assign(DoubleFunction function)
Example:
matrix = 2 x 2 matrix
0.5 1.5
2.5 3.5
// change each cell to its sine
matrix.assign(cern.jet.math.Functions.sin);
-->
2 x 2 matrix
0.479426 0.997495
0.598472 -0.350783
For further examples, see the package doc.
assign in class DoubleMatrix2Dfunction - a function object taking as argument the current cell's value.
DoubleFunctionspublic DoubleMatrix2D assign(double value)
assign in class DoubleMatrix2Dvalue - the value to be filled into the cells.
public DoubleMatrix2D assign(DoubleMatrix2D source)
assign in class DoubleMatrix2Dsource - the source matrix to copy from (may be identical to the
receiver).
IllegalArgumentException - if
columns() != source.columns() || rows() != source.rows()
public DoubleMatrix2D assign(DoubleMatrix2D y,
DoubleDoubleFunction function)
DoubleMatrix2DExample:
// assign x[row,col] = x[row,col]<sup>y[row,col]</sup>
m1 = 2 x 2 matrix
0 1
2 3
m2 = 2 x 2 matrix
0 2
4 6
m1.assign(m2, cern.jet.math.Functions.pow);
-->
m1 == 2 x 2 matrix
1 1
16 729
For further examples, see the package doc.
assign in class DoubleMatrix2Dy - the secondary matrix to operate on.function - a function object taking as first argument the current cell's
value of this, and as second argument the current
cell's value of y,
DoubleFunctionspublic int cardinality()
cardinality in class DoubleMatrix2Dpublic void dct2(boolean scale)
DoubleMatrix2D
dct2 in class DoubleMatrix2Dscale - if true then scaling is performedpublic void dctColumns(boolean scale)
DoubleMatrix2D
dctColumns in class DoubleMatrix2Dscale - if true then scaling is performedpublic void dctRows(boolean scale)
DoubleMatrix2D
dctRows in class DoubleMatrix2Dscale - if true then scaling is performedpublic void dht2()
DoubleMatrix2D
dht2 in class DoubleMatrix2Dpublic void dhtColumns()
DoubleMatrix2D
dhtColumns in class DoubleMatrix2Dpublic void dhtRows()
DoubleMatrix2D
dhtRows in class DoubleMatrix2Dpublic void dst2(boolean scale)
DoubleMatrix2D
dst2 in class DoubleMatrix2Dscale - if true then scaling is performedpublic void dstColumns(boolean scale)
DoubleMatrix2D
dstColumns in class DoubleMatrix2Dscale - if true then scaling is performedpublic void dstRows(boolean scale)
DoubleMatrix2D
dstRows in class DoubleMatrix2Dscale - if true then scaling is performedpublic AbstractIntDoubleMap elements()
elements in class DoubleMatrix2Dpublic void ensureCapacity(int minCapacity)
This method never need be called; it is for performance tuning only. Calling this method before tt>set()ing a large number of non-zero values boosts performance, because the receiver will grow only once instead of potentially many times and hash collisions get less probable.
ensureCapacity in class AbstractMatrixminCapacity - the desired minimum number of non-zero cells.public void fft2()
DoubleMatrix2D
this[k1][2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
this[k1][2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
this[0][2*k2] = Re[0][k2] = Re[0][columns-k2],
this[0][2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
this[k1][0] = Re[k1][0] = Re[rows-k1][0],
this[k1][1] = Im[k1][0] = -Im[rows-k1][0],
this[rows-k1][1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
this[rows-k1][0] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
this[0][0] = Re[0][0],
this[0][1] = Re[0][columns/2],
this[rows/2][0] = Re[rows/2][0],
this[rows/2][1] = Re[rows/2][columns/2]
This method computes only half of the elements of the real transform. The
other half satisfies the symmetry condition. If you want the full real
forward transform, use getFft2. To get back the original
data, use ifft2.
fft2 in class DoubleMatrix2Dpublic DoubleMatrix2D forEachNonZero(IntIntDoubleFunction function)
DoubleMatrix2D
forEachNonZero in class DoubleMatrix2Dfunction - a function object taking as argument the current non-zero
cell's row, column and value.
public DComplexMatrix2D getFft2()
DoubleMatrix2D
getFft2 in class DoubleMatrix2Dpublic DComplexMatrix2D getFftColumns()
DoubleMatrix2D
getFftColumns in class DoubleMatrix2Dpublic DComplexMatrix2D getFftRows()
DoubleMatrix2D
getFftRows in class DoubleMatrix2Dpublic DComplexMatrix2D getIfft2(boolean scale)
DoubleMatrix2D
getIfft2 in class DoubleMatrix2Dpublic DComplexMatrix2D getIfftColumns(boolean scale)
DoubleMatrix2D
getIfftColumns in class DoubleMatrix2Dpublic DComplexMatrix2D getIfftRows(boolean scale)
DoubleMatrix2D
getIfftRows in class DoubleMatrix2D
public double getQuick(int row,
int column)
Provided with invalid parameters this method may return invalid objects without throwing any exception. You should only use this method when you are absolutely sure that the coordinate is within bounds. Precondition (unchecked): 0 <= column < columns() && 0 <= row < rows().
getQuick in class DoubleMatrix2Drow - the index of the row-coordinate.column - the index of the column-coordinate.
public void idct2(boolean scale)
DoubleMatrix2D
idct2 in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idctColumns(boolean scale)
DoubleMatrix2D
idctColumns in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idctRows(boolean scale)
DoubleMatrix2D
idctRows in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idht2(boolean scale)
DoubleMatrix2D
idht2 in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idhtColumns(boolean scale)
DoubleMatrix2D
idhtColumns in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idhtRows(boolean scale)
DoubleMatrix2D
idhtRows in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idst2(boolean scale)
DoubleMatrix2D
idst2 in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idstColumns(boolean scale)
DoubleMatrix2D
idstColumns in class DoubleMatrix2Dscale - if true then scaling is performedpublic void idstRows(boolean scale)
DoubleMatrix2D
idstRows in class DoubleMatrix2Dscale - if true then scaling is performedpublic void ifft2(boolean scale)
DoubleMatrix2D
this[k1][2*k2] = Re[k1][k2] = Re[rows-k1][columns-k2],
this[k1][2*k2+1] = Im[k1][k2] = -Im[rows-k1][columns-k2],
0<k1<rows, 0<k2<columns/2,
this[0][2*k2] = Re[0][k2] = Re[0][columns-k2],
this[0][2*k2+1] = Im[0][k2] = -Im[0][columns-k2],
0<k2<columns/2,
this[k1][0] = Re[k1][0] = Re[rows-k1][0],
this[k1][1] = Im[k1][0] = -Im[rows-k1][0],
this[rows-k1][1] = Re[k1][columns/2] = Re[rows-k1][columns/2],
this[rows-k1][0] = -Im[k1][columns/2] = Im[rows-k1][columns/2],
0<k1<rows/2,
this[0][0] = Re[0][0],
this[0][1] = Re[0][columns/2],
this[rows/2][0] = Re[rows/2][0],
this[rows/2][1] = Re[rows/2][columns/2]
This method computes only half of the elements of the real transform. The
other half satisfies the symmetry condition. If you want the full real
inverse transform, use getIfft2.
ifft2 in class DoubleMatrix2Dscale - if true then scaling is performed
public int index(int row,
int column)
index in class AbstractMatrix2Drow - the index of the row-coordinate.column - the index of the column-coordinate.
public DoubleMatrix2D like(int rows,
int columns)
like in class DoubleMatrix2Drows - the number of rows the matrix shall have.columns - the number of columns the matrix shall have.
public DoubleMatrix1D like1D(int size)
like1D in class DoubleMatrix2Dsize - the number of cells the matrix shall have.
public void setQuick(int row,
int column,
double value)
Provided with invalid parameters this method may access illegal indexes without throwing any exception. You should only use this method when you are absolutely sure that the coordinate is within bounds. Precondition (unchecked): 0 <= column < columns() && 0 <= row < rows().
setQuick in class DoubleMatrix2Drow - the index of the row-coordinate.column - the index of the column-coordinate.value - the value to be filled into the specified cell.public void trimToSize()
Background:
Cells that
trimToSize in class AbstractMatrixpublic DoubleMatrix1D vectorize()
vectorize in class DoubleMatrix2D
public DoubleMatrix1D zMult(DoubleMatrix1D y,
DoubleMatrix1D z,
double alpha,
double beta,
boolean transposeA)
DoubleMatrix2D
zMult in class DoubleMatrix2Dy - the source vector.z - the vector where results are to be stored. Set this parameter
to null to indicate that a new result vector shall be
constructed.
public DoubleMatrix2D zMult(DoubleMatrix2D B,
DoubleMatrix2D C,
double alpha,
double beta,
boolean transposeA,
boolean transposeB)
DoubleMatrix2D
zMult in class DoubleMatrix2DB - the second source matrix.C - the matrix where results are to be stored. Set this parameter
to null to indicate that a new result matrix shall be
constructed.
|
Parallel Colt 0.6.1 | |||||||||
| PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
| SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD | |||||||||