1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299
|
/**
\addtogroup arrayfire_func
@{
\defgroup lapack_factor_func_lu lu
\ingroup lapack_factor_mat
\brief Perform LU decomposition
This function decomposes input matrix **A** into a lower triangle **L**, an upper triangle **U** such that
\f$A = L * U\f$
For stability, a permutation array **P** is also used to modify the formula in the following manner.
\f$A(P, span) = L * U\f$
This operation can be performed in ArrayFire using the following code snippet.
\snippet test/lu_dense.cpp ex_lu_unpacked
The permuted version of the original matrix can be reconstructed using the following snippet.
\snippet test/lu_dense.cpp ex_lu_recon
The sample output for these operations can be seen below.
\code
a_orig [3 3 1 1]
0.0000 3.0000 6.0000
1.0000 4.0000 7.0000
2.0000 5.0000 8.0000
l [3 3 1 1]
1.0000 0.0000 0.0000
0.0000 1.0000 0.0000
0.5000 0.5000 1.0000
u [3 3 1 1]
2.0000 5.0000 8.0000
0.0000 3.0000 6.0000
0.0000 0.0000 0.0000
pivot [3 1 1 1]
2
0
1
a_recon [3 3 1 1]
2.0000 5.0000 8.0000
0.0000 3.0000 6.0000
1.0000 4.0000 7.0000
a_perm [3 3 1 1]
2.0000 5.0000 8.0000
0.0000 3.0000 6.0000
1.0000 4.0000 7.0000
\endcode
When memory is a concern, users can perform the LU decomposition in place as shown below.
\snippet test/lu_dense.cpp ex_lu_packed
The lower and upper triangle matrices can be obtained if necessary in the following manner.
\snippet test/lu_dense.cpp ex_lu_extract
LU decompositions has many applications including <a href="http://en.wikipedia.org/wiki/LU_decomposition#Solving_linear_equations">solving a system of linear equations</a>. Check \ref af::solveLU fore more information.
=======================================================================
\defgroup lapack_factor_func_qr qr
\ingroup lapack_factor_mat
\brief Perform QR decomposition
This function decomposes input matrix **A** into an orthogonal matrix **Q** and an upper triangular matrix **R** such that
\f$A = Q * R\f$
\f$Q * Q^T = I\f$
Where **I** is an identity matrix. The matrix **Q** is a square matrix of size **max(M, N)** where **M** and **N** are rows and columns of **A** respectively. The matrix **R** is the same size as **A*.
This operation can be performed in ArrayFire using the following code snippet.
\snippet test/qr_dense.cpp ex_qr_unpacked
The additional parameter **Tau** can be used to speed up solving over and under determined system of equations.
The original matrix can be reconstructed using the following code snippet.
\snippet test/qr_dense.cpp ex_qr_recon
When memory is a concern, users can perform QR decomposition in place as shown below.
\snippet test/qr_dense.cpp ex_qr_packed
=======================================================================
\defgroup lapack_factor_func_cholesky cholesky
\ingroup lapack_factor_mat
\brief Perform Cholesky decomposition
This function decomposes a <a href="http://en.wikipedia.org/wiki/Positive-definite_matrix">positive definite</a> matrix **A** into two triangular matrices such that
\f$A = L * U\f$
\f$L = U^T\f$
Only one of **L** and **U** is stored to conserve space when solving linear equations.
This operation can be performed in ArrayFire using the following code snippet.
\snippet test/cholesky_dense.cpp ex_chol_reg
When memory is a concern, users can perform Cholesky decomposition in place as shown below.
\snippet test/cholesky_dense.cpp ex_chol_inplace
=======================================================================
\defgroup lapack_factor_func_svd svd
\ingroup lapack_factor_mat
\brief Perform Singular Value Decomposition
This function factorizes a matrix **A** into two unitary matrices **U** and **Vt**, and a diagonal matrix **S** such that
\f$A = U * S * Vt\f$
If **A** has **M** rows and **N** columns, **U** is of the size **M x M** , **V** is of size **N x N**, and **S** is of size **M x N**
The arrayfire function only returns the non zero diagonal elements of **S**. To reconstruct the original matrix **A** from the individual factors, the following code snuppet can be used:
\snippet test/svd_dense.cpp ex_svd_reg
When memory is a concern, and **A** is dispensible, \ref svdInPlace() can be used
=======================================================================
\defgroup lapack_solve_func_gen solve
\ingroup lapack_solve_mat
\brief Solve a system of equations
This function takes a co-efficient matrix **A** and an output matrix **B** as inputs to solve the following equation for **X**
\f$A * X = B\f$
This operation can be done in ArrayFire using the following code snippet.
\snippet test/solve_dense.cpp ex_solve
The results can be verified by reconstructing the output matrix using \ref af::matmul in the following manner.
\snippet test/solve_dense.cpp ex_solve_recon
The sample output can be seen below
\code
A [3 3 1 1]
0.1000 3.1000 6.1000
1.1000 4.1000 7.0000
2.0000 5.0000 8.0000
B0 [3 1 1 1]
21.9000
30.7000
39.0000
X1 [3 1 1 1]
4.0000
3.0000
2.0000
B1 [3 1 1 1]
21.9000
30.7000
39.0000
\endcode
If the coefficient matrix is known to be a triangular matrix, \ref AF_MAT_LOWER or \ref AF_MAT_UPPER can be passed to make solve faster.
The sample code snippets for solving a lower triangular matrix can be seen below.
\snippet test/solve_dense.cpp ex_solve_lower
Similarily, the code snippet for solving an upper triangular matrix can be seen below.
\snippet test/solve_dense.cpp ex_solve_upper
See also: \ref af::solveLU
=======================================================================
\defgroup lapack_solve_lu_func_gen solveLU
\ingroup lapack_solve_mat
\brief Solve a system of equations
This function takes a co-efficient matrix **A** and an output matrix **B** as inputs to solve the following equation for **X**
\f$A * X = B\f$
This operation can be done in ArrayFire using the following code snippet.
\snippet test/solve_dense.cpp ex_solve_lu
This function along with \ref af::lu split up the task af::solve performs for square matrices.
\note This function is beneficial over \ref af::solve only in long running application where the coefficient matrix **A** stays the same, but the observed variables keep changing.
=======================================================================
\defgroup lapack_ops_func_inv inverse
\ingroup lapack_ops_mat
\brief Invert a matrix
This function inverts a square matrix **A**. The code snippet to demonstrate this can be seen below.
\snippet test/inverse_dense.cpp ex_inverse
The sample output can be seen below
\code
A [3 3 1 1]
0.0100 3.0100 6.0100
1.0100 4.0100 7.0000
2.0000 5.0000 8.0000
IA [3 3 1 1]
48.9076 -99.9927 50.7518
-99.1552 199.9852 -100.4968
49.7451 -99.9926 50.2475
I [3 3 1 1]
1.0000 0.0001 -0.0000
0.0000 1.0000 0.0000
0.0000 0.0000 1.0000
\endcode
==================================================================================
\defgroup lapack_ops_func_rank rank
\ingroup lapack_ops_mat
\brief Find the rank of the input matrix.
This function uses \ref af::qr to find the rank of the input matrix within the given tolerance.
=====================================================================================
\defgroup lapack_ops_func_det det
\ingroup lapack_ops_mat
\brief Find the determinant of the input matrix.
\note This function requires scratch space equal to the input array
===============================================================================
\defgroup lapack_ops_func_norm norm
\ingroup lapack_ops_mat
\brief Find the norm of the input matrix
This function can return the norm using various metrics based on the type paramter.
\note \ref AF_NORM_MATRIX_2 is currently not supported.
===============================================================================
\defgroup lapack_helper_func_available isLAPACKAvailable
\ingroup lapack_helper
\brief Returns true is ArrayFire is compiled with LAPACK support
===============================================================================
@}
*/
|