SSJ
3.3.1
Stochastic Simulation in Java
|
This class implements a few methods for matrix calculations with double
numbers.
More...
Public Member Functions | |
DMatrix (int r, int c) | |
Creates a new DMatrix with r rows and c columns. More... | |
DMatrix (double[][] data, int r, int c) | |
Creates a new DMatrix with r rows and c columns using the data in data . More... | |
DMatrix (DMatrix that) | |
Copy constructor. More... | |
String | toString () |
Creates a String containing all the data of the DMatrix . More... | |
int | numRows () |
Returns the number of rows of the DMatrix . More... | |
int | numColumns () |
Returns the number of columns of the DMatrix . More... | |
double | get (int row, int column) |
Returns the matrix element in the specified row and column. More... | |
void | set (int row, int column, double value) |
Sets the value of the element in the specified row and column. More... | |
DMatrix | transpose () |
Returns the transposed matrix. More... | |
Static Public Member Functions | |
static void | CholeskyDecompose (double[][] M, double[][] L) |
Given a symmetric positive-definite matrix \(M\), performs the Cholesky decomposition of \(M\) and returns the result as a lower triangular matrix \(L\), such that \(M = L L^T\). More... | |
static DoubleMatrix2D | CholeskyDecompose (DoubleMatrix2D M) |
Given a symmetric positive-definite matrix \(M\), performs the Cholesky decomposition of \(M\) and returns the result as a lower triangular matrix \(L\), such that \(M = L L^T\). More... | |
static void | PCADecompose (double[][] M, double[][] A, double[] lambda) |
Computes the principal components decomposition \(M\) = \(U\Lambda U^{\mathsf{t}}\) by using the singular value decomposition of matrix \(M\). More... | |
static DoubleMatrix2D | PCADecompose (DoubleMatrix2D M, double[] lambda) |
Computes the principal components decomposition \(M\) = \(U\Lambda U^{\mathsf{t}}\) by using the singular value decomposition of matrix \(M\). More... | |
static double [] | solveLU (double[][] A, double[] b) |
Solves the matrix equation \(Ax = b\) using LU decomposition. More... | |
static void | solveTriangular (DoubleMatrix2D U, DoubleMatrix2D B, DoubleMatrix2D X) |
Solve the triangular matrix equation \(UX = B\) for \(X\). More... | |
static double [][] | exp (double[][] A) |
Similar to #exp(final DoubleMatrix2D). More... | |
static DoubleMatrix2D | exp (final DoubleMatrix2D A) |
Returns \(e^A\), the exponential of the square matrix \(A\). More... | |
static DoubleMatrix2D | expBidiagonal (final DoubleMatrix2D A) |
Returns \(e^A\), the exponential of the bidiagonal square matrix \(A\). More... | |
static DoubleMatrix1D | expBidiagonal (final DoubleMatrix2D A, final DoubleMatrix1D b) |
Computes \(c = e^A b\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\). More... | |
static DoubleMatrix2D | expmiBidiagonal (final DoubleMatrix2D A) |
Computes \(e^A - I\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\). More... | |
static DoubleMatrix1D | expmiBidiagonal (final DoubleMatrix2D A, final DoubleMatrix1D b) |
Computes \(c = (e^A - I)b\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\). More... | |
static void | copy (double[][] M, double[][] R) |
Copies the matrix \(M\) into \(R\). More... | |
static String | toString (double[][] M) |
Returns matrix \(M\) as a string. More... | |
Static Package Functions | |
static void | addMultTriang (final DoubleMatrix2D A, DoubleMatrix1D b, double h) |
static void | multBand (final DoubleMatrix2D A, int sa, final DoubleMatrix2D B, int sb, DoubleMatrix2D C) |
Matrix multiplication \(C = AB\). More... | |
static void | multBand (DoubleMatrix2D A, int sa, double h) |
Multiplication of the matrix \(A\) by the scalar \(h\). More... | |
static void | addMultBand (double g, DoubleMatrix2D A, int sa, double h, final DoubleMatrix2D B, int sb) |
Addition of the matrices \(gA + hB\) after multiplication with the scalars \(g\) and \(h\). More... | |
Package Attributes | |
int | c |
This class implements a few methods for matrix calculations with double
numbers.
DMatrix | ( | int | r, |
int | c | ||
) |
Creates a new DMatrix
with r
rows and c
columns.
r | the number of rows |
c | the number of columns |
DMatrix | ( | double | data[][], |
int | r, | ||
int | c | ||
) |
|
staticpackage |
Addition of the matrices \(gA + hB\) after multiplication with the scalars \(g\) and \(h\).
\(A\) is a square banded upper triangular matrix. It has a non-zero diagonal, sa
superdiagonals, and thus a bandwidth of sa + 1
. Similarly for \(B\) which has a bandwidth of sb + 1
. The result is put back in \(A\).
g | coefficient multiplying A |
A | input and output matrix |
sa | number of superdiagonals of A |
h | coefficient multiplying B |
B | input matrix |
sb | number of superdiagonals of B |
|
static |
Given a symmetric positive-definite matrix \(M\), performs the Cholesky decomposition of \(M\) and returns the result as a lower triangular matrix \(L\), such that \(M = L L^T\).
M | the input matrix |
L | the Cholesky lower triangular matrix |
|
static |
Given a symmetric positive-definite matrix \(M\), performs the Cholesky decomposition of \(M\) and returns the result as a lower triangular matrix \(L\), such that \(M = L L^T\).
M | the input matrix |
|
static |
Copies the matrix \(M\) into \(R\).
M | original matrix |
R | output matrix |
|
static |
Similar to #exp(final DoubleMatrix2D).
A | input matrix |
|
static |
Returns \(e^A\), the exponential of the square matrix \(A\).
The scaling and squaring method [85] is used with Padé approximants to compute the exponential.
A | input matrix |
|
static |
Returns \(e^A\), the exponential of the bidiagonal square matrix \(A\).
The only non-zero elements of \(A\) are on the diagonal and on the first superdiagonal. This method is faster than #exp(final DoubleMatrix2D) because of the special form of \(A\).
A | bidiagonal matrix |
|
static |
Computes \(c = e^A b\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\).
The only non-zero elements of \(A\) are on the diagonal and on the first superdiagonal. Uses the scaling and squaring method [85] with Padé approximants. Returns \(c\).
A | bidiagonal matrix |
b | vector |
|
static |
Computes \(e^A - I\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\).
The only non-zero elements of \(A\) are on the diagonal and on the first superdiagonal. Uses the scaling and squaring method [217], [85] with Padé approximants. Returns \(e^A - I\).
A | bidiagonal matrix |
|
static |
Computes \(c = (e^A - I)b\), where \(e^A\) is the exponential of the bidiagonal square matrix \(A\).
The only non-zero elements of \(A\) are on the diagonal and on the first superdiagonal. Uses the scaling and squaring method [217], [85] with a Taylor expansion. Returns \(c\).
A | bidiagonal matrix |
b | vector |
double get | ( | int | row, |
int | column | ||
) |
Returns the matrix element in the specified row and column.
row | the row of the selected element |
column | the column of the selected element |
IndexOutOfBoundsException | if the selected element would be outside the DMatrix |
|
staticpackage |
Matrix multiplication \(C = AB\).
All three matrices are square, banded, and upper triangular. \(A\) has a non-zero diagonal, sa
non-zero superdiagonals, and thus a bandwidth of sa + 1
. The non-zero elements of \(A_{ij}\) are those for which \(j - s_a \le i \le j\). Similarly for \(B\) which has a bandwidth of sb + 1
. The resulting matrix \(C\) has sa + sb
non-zero superdiagonals, and a bandwidth of sa + sb + 1
.
A | input matrix |
sa | number of superdiagonals of A |
B | input matrix |
sb | number of superdiagonals of B |
C | result |
|
staticpackage |
Multiplication of the matrix \(A\) by the scalar \(h\).
\(A\) is a square banded upper triangular matrix. It has a non-zero diagonal, sa
superdiagonals, and thus a bandwidth of sa + 1
. The result of the multiplication is put back in \(A\).
A | input and output matrix |
sa | number of superdiagonals of A |
h | scalar |
int numColumns | ( | ) |
Returns the number of columns of the DMatrix
.
int numRows | ( | ) |
Returns the number of rows of the DMatrix
.
|
static |
Computes the principal components decomposition \(M\) = \(U\Lambda U^{\mathsf{t}}\) by using the singular value decomposition of matrix \(M\).
Puts the eigenvalues, which are the diagonal elements of matrix \(\Lambda\), sorted by decreasing size, in vector lambda
, and puts matrix \(A = U\sqrt{\Lambda}\) in A
.
M | input matrix |
A | matrix square root of M |
lambda | the eigenvalues |
|
static |
Computes the principal components decomposition \(M\) = \(U\Lambda U^{\mathsf{t}}\) by using the singular value decomposition of matrix \(M\).
Puts the eigenvalues, which are the diagonal elements of matrix \(\Lambda\), sorted by decreasing size, in vector lambda
. Returns matrix \(A = U\sqrt{\Lambda}\).
M | input matrix |
lambda | the eigenvalues |
void set | ( | int | row, |
int | column, | ||
double | value | ||
) |
Sets the value of the element in the specified row and column.
row | the row of the selected element |
column | the column of the selected element |
value | the new value of the element |
IndexOutOfBoundsException | if the selected element would be outside the DMatrix |
|
static |
Solves the matrix equation \(Ax = b\) using LU decomposition.
\(A\) is a square matrix, \(b\) and \(x\) are vectors. Returns the solution \(x\).
A | square matrix |
b | right side vector |
|
static |
Solve the triangular matrix equation \(UX = B\) for \(X\).
\(U\) is a square upper triangular matrix. \(B\) and \(X\) must have the same number of columns.
U | input matrix |
B | right-hand side matrix |
X | output matrix |
|
static |
Returns matrix \(M\) as a string.
It is displayed in matrix form, with each row on a line.
String toString | ( | ) |
DMatrix transpose | ( | ) |
Returns the transposed matrix.
The rows and columns are interchanged.