![]() |
MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
|
Functions | |
void | magma_print_environment () |
Print MAGMA version, CUDA version, LAPACK/BLAS library version, available GPU devices, number of threads, date, etc. | |
template<typename FloatT > | |
void | magma_generate_matrix (magma_opts &opts, Matrix< FloatT > &A, Vector< typename blas::traits< FloatT >::real_t > &sigma) |
Generate an m-by-n test matrix A. | |
void magma_print_environment | ( | ) |
Print MAGMA version, CUDA version, LAPACK/BLAS library version, available GPU devices, number of threads, date, etc.
Used in testing.
void magma_generate_matrix | ( | magma_opts & | opts, |
Matrix< FloatT > & | A, | ||
Vector< typename blas::traits< FloatT >::real_t > & | sigma ) |
Generate an m-by-n test matrix A.
Similar to but does not use LAPACK's libtmg.
[in] | opts | MAGMA options. Uses matrix, cond, condD; see further details. |
[out] | A | Complex array, dimension (lda, n). On output, the m-by-n test matrix A in an lda-by-n array. |
[in,out] | sigma | Real array, dimension (min(m,n)) For matrix with "_specified", on input contains user-specified singular or eigenvalues. On output, contains singular or eigenvalues, if known, else set to NaN. sigma is not necesarily sorted. |
The --matrix
command line option specifies the matrix name according to the tables below. Where indicated, names take an optional distribution suffix (%) and an optional scaling suffix (^). The default distribution is rand
.
The --cond
and --condD
command line options specify condition numbers as described below. By default, cond = sqrt( 1/eps ) = 6.7e7 for double. condD is for specialized eigenvalue and SVD tests; by default, condD = 1.
Examples:
./testing_zgemm --matrix rand_small ./testing_zgemm --matrix svd_arith --cond 1e6 ./testing_zgemm --matrix svd_logrand_ufl --cond 1e6
In descriptions below:
See LAPACK Working Note (LAWN) 41:
Matrix | Description |
---|---|
zero | all entries are 0 |
ones | all entries are 1 |
identity | diagonal entries are 1 |
jordan | diagonal and first subdiagonal entries are 1 |
kronecker | \(A_{ij} = 1 + (m/cond) \delta_{ij} \) |
– – – | – – – |
rand^ | matrix entries random uniform on (0, 1) [note 1] |
rands^ | matrix entries random uniform on (-1, 1) [note 1] |
randn^ | matrix entries random normal with mean 0, std 1 [note 1,2] |
– – – | – – – |
diag%^ | \(A = \Sigma \) |
svd%^ | \(A = U \Sigma V^H \) |
poev%^ | \(A = V \Sigma V^H \) (eigenvalues positive [note 3], i.e., matrix SPD) |
spd%^ | alias for poev |
heev%^ | \(A = V \Lambda V^H\) (eigenvalues mixed signs) |
syev%^ | alias for heev |
geev%^ | \(A = V T V^H, \) with Schur-form \(T\), \(V\) orthogonal [not yet implemented] |
geevx%^ | \(A = X T X^{-1}, \) with Schur-form \(T\), \(X\) ill-conditioned [not yet implemented] |
[1] rand, rands, randn do not use cond, so the condition number is arbitrary.
[2] For randn, \(Expected( \log( cond ) ) = \log( 4.65 n )\) [Edelman, 1988].
Suffix | Description |
---|---|
_rand | \(\sigma_i\) random uniform on (0, 1) [default] |
_rands | \(\sigma_i\) random uniform on (-1, 1); [note 3] |
_randn | \(\sigma_i\) random normal with mean 0, std 1; [note 3] |
– – – | – – – |
_logrand | \(\log( \sigma_i )\) uniform on \((\log(1/cond), \log(1))\) |
_arith | \(\sigma_i = 1 - (i - 1)/(n - 1)(1 - 1/cond); \sigma_{i+1} - \sigma_i \) is constant |
_geo | \(\sigma_i = (cond)^{ -(i-1)/(n-1) }; \sigma_{i+1} / \sigma_i \) is constant |
_cluster0 | \(\sigma = [ 1, 1/cond, ..., 1/cond ]; \) has 1 unit value, \(n-1\) small values |
_cluster1 | \(\sigma = [ 1, ..., 1, 1/cond ]; \) has \(n-1\) unit values, 1 small value |
_rarith | _arith , reversed order |
_rgeo | _geo , reversed order |
_rcluster0 | _cluster0 , reversed order |
_rcluster1 | _cluster1 , reversed order |
_specified | user specified \(\Sigma\) on input |
[3] For _rands
, _randn
, \(\Sigma\) contains negative values.
Suffix | Description |
---|---|
_ufl | scale near underflow = 1e-308 for double |
_ofl | scale near overflow = 2e+308 for double |
_small | scale near sqrt( underflow ) = 1e-154 for double |
_large | scale near sqrt( overflow ) = 6e+153 for double |
_dominant | diagonally dominant |
For _dominant
, set diagonal entries to row or column sum:
\[ A_{ii} = \max( \sum_j | A_{ij} |, \sum_k | A_{ki} | ). \]
Note _dominant
changes the singular or eigenvalues, and cond.
Here, scaling by \(D\) is implemented, scaling by \(K\) is not yet implemented. If condD != 1, then: For SVD,
\[ A_0 = U \Sigma V^H, \\ A = A_0 K D, \]
where \(K\) is diagonal such that columns of \(U \Sigma V^H K\) have unit norm, hence \(A^T A\) has unit diagonal, and \(D\) has log-random entries in \(( \log(1/condD), \log(1) )\).
For heev,
\[ A_0 = U \Lambda U^H, \\ A = D K A_0 K D, \]
where \(K\) is diagonal such that \(K A_0 K\) has unit diagonal, and \(D\) is as above.
Note using condD changes the singular or eigenvalues; on output, \(\Sigma\) contains the singular or eigenvalues of \(A_0\), not of \(A\).
See: Demmel and Veselic, Jacobi's method is more accurate than QR, 1992.