My Project
lapack.h
Go to the documentation of this file.
1/*******************************************************
2 * Copyright (c) 2014, ArrayFire
3 * All rights reserved.
4 *
5 * This file is distributed under 3-clause BSD license.
6 * The complete license agreement can be obtained at:
7 * http://arrayfire.com/licenses/BSD-3-Clause
8 ********************************************************/
9
10#pragma once
11#include <af/array.h>
12#include <af/defines.h>
13
14#ifdef __cplusplus
15namespace af
16{
17#if AF_API_VERSION >= 31
28 AFAPI void svd(array &u, array &s, array &vt, const array &in);
29#endif
30
31#if AF_API_VERSION >= 31
42 AFAPI void svdInPlace(array &u, array &s, array &vt, array &in);
43#endif
44
57 AFAPI void lu(array &out, array &pivot, const array &in, const bool is_lapack_piv=true);
58
71 AFAPI void lu(array &lower, array &upper, array &pivot, const array &in);
72
84 AFAPI void luInPlace(array &pivot, array &in, const bool is_lapack_piv=true);
85
97 AFAPI void qr(array &out, array &tau, const array &in);
98
111 AFAPI void qr(array &q, array &r, array &tau, const array &in);
112
123 AFAPI void qrInPlace(array &tau, array &in);
124
139 AFAPI int cholesky(array &out, const array &in, const bool is_upper = true);
140
154 AFAPI int choleskyInPlace(array &in, const bool is_upper = true);
155
169 AFAPI array solve(const array &a, const array &b, const matProp options = AF_MAT_NONE);
170
171
186 AFAPI array solveLU(const array &a, const array &piv,
187 const array &b, const matProp options = AF_MAT_NONE);
188
201 AFAPI array inverse(const array &in, const matProp options = AF_MAT_NONE);
202
213 AFAPI unsigned rank(const array &in, const double tol=1E-5);
214
224 template<typename T> T det(const array &in);
225
238 AFAPI double norm(const array &in, const normType type=AF_NORM_EUCLID,
239 const double p=1, const double q=1);
240
241#if AF_API_VERSION >= 33
250#endif
251
252}
253#endif
254
255#ifdef __cplusplus
256extern "C" {
257#endif
258
259#if AF_API_VERSION >= 31
271#endif
272
273#if AF_API_VERSION >= 31
285#endif
286
298
308 AFAPI af_err af_lu_inplace(af_array *pivot, af_array in, const bool is_lapack_piv);
309
321
331
344 AFAPI af_err af_cholesky(af_array *out, int *info, const af_array in, const bool is_upper);
345
357 AFAPI af_err af_cholesky_inplace(int *info, af_array in, const bool is_upper);
358
372 const af_mat_prop options);
373
389 const af_array b, const af_mat_prop options);
390
402 AFAPI af_err af_inverse(af_array *out, const af_array in, const af_mat_prop options);
403
413 AFAPI af_err af_rank(unsigned *rank, const af_array in, const double tol);
414
424 AFAPI af_err af_det(double *det_real, double *det_imag, const af_array in);
425
438 AFAPI af_err af_norm(double *out, const af_array in, const af_norm_type type, const double p, const double q);
439
440#if AF_API_VERSION >= 33
451#endif
452
453
454#ifdef __cplusplus
455}
456#endif
A multi dimensional data container.
Definition: array.h:27
af_norm_type
Definition: defines.h:319
@ AF_NORM_EUCLID
The default. Same as AF_NORM_VECTOR_2.
Definition: defines.h:329
af_mat_prop
Definition: defines.h:304
@ AF_MAT_NONE
Default.
Definition: defines.h:305
af_err
Definition: defines.h:63
void * af_array
Definition: defines.h:222
#define AFAPI
Definition: defines.h:31
AFAPI array lower(const array &in, bool is_unit_diag=false)
AFAPI array upper(const array &in, bool is_unit_diag=false)
AFAPI void info()
AFAPI int choleskyInPlace(array &in, const bool is_upper=true)
C++ Interface for in place cholesky decomposition.
AFAPI af_err af_cholesky(af_array *out, int *info, const af_array in, const bool is_upper)
C++ Interface for cholesky decomposition.
AFAPI int cholesky(array &out, const array &in, const bool is_upper=true)
C++ Interface for cholesky decomposition.
AFAPI af_err af_cholesky_inplace(int *info, af_array in, const bool is_upper)
C Interface for in place cholesky decomposition.
AFAPI af_err af_lu_inplace(af_array *pivot, af_array in, const bool is_lapack_piv)
C Interface for in place LU decomposition.
AFAPI void luInPlace(array &pivot, array &in, const bool is_lapack_piv=true)
C++ Interface for in place LU decomposition.
AFAPI void lu(array &out, array &pivot, const array &in, const bool is_lapack_piv=true)
C++ Interface for LU decomposition in packed format.
AFAPI af_err af_lu(af_array *lower, af_array *upper, af_array *pivot, const af_array in)
C Interface for LU decomposition.
AFAPI void qrInPlace(array &tau, array &in)
C++ Interface for QR decomposition.
AFAPI af_err af_qr(af_array *q, af_array *r, af_array *tau, const af_array in)
C Interface for QR decomposition.
AFAPI af_err af_qr_inplace(af_array *tau, af_array in)
C Interface for QR decomposition.
AFAPI void qr(array &out, array &tau, const array &in)
C++ Interface for QR decomposition in packed format.
AFAPI af_err af_svd_inplace(af_array *u, af_array *s, af_array *vt, af_array in)
C Interface for SVD decomposition.
AFAPI void svdInPlace(array &u, array &s, array &vt, array &in)
C++ Interface for SVD decomposition.
AFAPI af_err af_svd(af_array *u, af_array *s, af_array *vt, const af_array in)
C Interface for SVD decomposition.
AFAPI void svd(array &u, array &s, array &vt, const array &in)
C++ Interface for SVD decomposition.
T det(const array &in)
C++ Interface for finding the determinant of a matrix.
AFAPI af_err af_det(double *det_real, double *det_imag, const af_array in)
C Interface for finding the determinant of a matrix.
AFAPI array inverse(const array &in, const matProp options=AF_MAT_NONE)
C++ Interface for inverting a matrix.
AFAPI af_err af_inverse(af_array *out, const af_array in, const af_mat_prop options)
C Interface for inverting a matrix.
AFAPI double norm(const array &in, const normType type=AF_NORM_EUCLID, const double p=1, const double q=1)
C++ Interface for norm of a matrix.
AFAPI af_err af_norm(double *out, const af_array in, const af_norm_type type, const double p, const double q)
C Interface for norm of a matrix.
AFAPI bool isLAPACKAvailable()
Returns true is ArrayFire is compiled with LAPACK support.
AFAPI af_err af_is_lapack_available(bool *out)
Returns true is ArrayFire is compiled with LAPACK support.
AFAPI unsigned rank(const array &in, const double tol=1E-5)
C++ Interface for finding the rank of a matrix.
AFAPI af_err af_rank(unsigned *rank, const af_array in, const double tol)
C Interface for finding the rank of a matrix.
AFAPI array solve(const array &a, const array &b, const matProp options=AF_MAT_NONE)
C++ Interface for solving a system of equations.
AFAPI af_err af_solve(af_array *x, const af_array a, const af_array b, const af_mat_prop options)
C Interface for solving a system of equations.
AFAPI array solveLU(const array &a, const array &piv, const array &b, const matProp options=AF_MAT_NONE)
C++ Interface for solving a system of equations.
AFAPI af_err af_solve_lu(af_array *x, const af_array a, const af_array piv, const af_array b, const af_mat_prop options)
C Interface for solving a system of equations.
Definition: algorithm.h:15