A high-performance general-purpose compute library
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
29 AFAPI void svd(array &u, array &s, array &vt, const array &in);
30#endif
31
32#if AF_API_VERSION >= 31
47 AFAPI void svdInPlace(array &u, array &s, array &vt, array &in);
48#endif
49
64 AFAPI void lu(array &out, array &pivot, const array &in, const bool is_lapack_piv=true);
65
79 AFAPI void lu(array &lower, array &upper, array &pivot, const array &in);
80
95 AFAPI void luInPlace(array &pivot, array &in, const bool is_lapack_piv=true);
96
108 AFAPI void qr(array &out, array &tau, const array &in);
109
123 AFAPI void qr(array &q, array &r, array &tau, const array &in);
124
135 AFAPI void qrInPlace(array &tau, array &in);
136
156 AFAPI int cholesky(array &out, const array &in, const bool is_upper = true);
157
173 AFAPI int choleskyInPlace(array &in, const bool is_upper = true);
174
190 AFAPI array solve(const array &a, const array &b, const matProp options = AF_MAT_NONE);
191
208 AFAPI array solveLU(const array &a, const array &piv,
209 const array &b, const matProp options = AF_MAT_NONE);
210
224 AFAPI array inverse(const array &in, const matProp options = AF_MAT_NONE);
225
226#if AF_API_VERSION >= 37
247 AFAPI array pinverse(const array &in, const double tol=1E-6,
248 const matProp options = AF_MAT_NONE);
249#endif
250
260 AFAPI unsigned rank(const array &in, const double tol=1E-5);
261
270 template<typename T> T det(const array &in);
271
285 AFAPI double norm(const array &in, const normType type=AF_NORM_EUCLID,
286 const double p=1, const double q=1);
287
288#if AF_API_VERSION >= 33
297#endif
298
299}
300#endif
301
302#ifdef __cplusplus
303extern "C" {
304#endif
305
306#if AF_API_VERSION >= 31
321#endif
322
323#if AF_API_VERSION >= 31
341#endif
342
356 AFAPI af_err af_lu(af_array *lower, af_array *upper, af_array *pivot, const af_array in);
357
374 AFAPI af_err af_lu_inplace(af_array *pivot, af_array in, const bool is_lapack_piv);
375
392
406
426 AFAPI af_err af_cholesky(af_array *out, int *info, const af_array in, const bool is_upper);
427
443 AFAPI af_err af_cholesky_inplace(int *info, af_array in, const bool is_upper);
444
463 const af_mat_prop options);
464
482 const af_array b, const af_mat_prop options);
483
497 AFAPI af_err af_inverse(af_array *out, const af_array in, const af_mat_prop options);
498
499#if AF_API_VERSION >= 37
523 AFAPI af_err af_pinverse(af_array *out, const af_array in, const double tol,
524 const af_mat_prop options);
525#endif
526
538 AFAPI af_err af_rank(unsigned *rank, const af_array in, const double tol);
539
551 AFAPI af_err af_det(double *det_real, double *det_imag, const af_array in);
552
568 AFAPI af_err af_norm(double *out, const af_array in, const af_norm_type type, const double p, const double q);
569
570#if AF_API_VERSION >= 33
582#endif
583
584
585#ifdef __cplusplus
586}
587#endif
A multi dimensional data container.
Definition: array.h:37
af_norm_type
Definition: defines.h:363
@ AF_NORM_EUCLID
The default. Same as AF_NORM_VECTOR_2.
Definition: defines.h:373
af_mat_prop
Definition: defines.h:348
@ AF_MAT_NONE
Default.
Definition: defines.h:349
af_err
Definition: defines.h:71
void * af_array
Definition: defines.h:240
#define AFAPI
Definition: defines.h:38
AFAPI array lower(const array &in, bool is_unit_diag=false)
C++ Interface to return the lower triangle array.
AFAPI array upper(const array &in, bool is_unit_diag=false)
C++ Interface to return the upper triangle array.
AFAPI int choleskyInPlace(array &in, const bool is_upper=true)
C++ Interface to perform in-place Cholesky decomposition.
AFAPI af_err af_cholesky(af_array *out, int *info, const af_array in, const bool is_upper)
C Interface to perform Cholesky decomposition.
AFAPI int cholesky(array &out, const array &in, const bool is_upper=true)
C++ Interface to perform Cholesky decomposition.
AFAPI af_err af_cholesky_inplace(int *info, af_array in, const bool is_upper)
C Interface to perform in-place Cholesky decomposition.
AFAPI af_err af_lu_inplace(af_array *pivot, af_array in, const bool is_lapack_piv)
C Interface to perform in-place LU decomposition.
AFAPI void luInPlace(array &pivot, array &in, const bool is_lapack_piv=true)
C++ Interface to perform in-place LU decomposition.
AFAPI void lu(array &out, array &pivot, const array &in, const bool is_lapack_piv=true)
C++ Interface to perform 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 to perform LU decomposition.
AFAPI void qrInPlace(array &tau, array &in)
C++ Interface to perform QR decomposition.
AFAPI af_err af_qr(af_array *q, af_array *r, af_array *tau, const af_array in)
C Interface to perform QR decomposition.
AFAPI af_err af_qr_inplace(af_array *tau, af_array in)
C Interface to perform QR decomposition.
AFAPI void qr(array &out, array &tau, const array &in)
C++ Interface to perform 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 to perform in-place singular value decomposition.
AFAPI void svdInPlace(array &u, array &s, array &vt, array &in)
C++ Interface to perform in-place singular value decomposition.
AFAPI af_err af_svd(af_array *u, af_array *s, af_array *vt, const af_array in)
C Interface to perform singular value decomposition.
AFAPI void svd(array &u, array &s, array &vt, const array &in)
C++ Interface to perform singular value decomposition.
AFAPI bool isLAPACKAvailable()
Returns true if ArrayFire is compiled with LAPACK support.
AFAPI af_err af_is_lapack_available(bool *out)
Returns true if ArrayFire is compiled with LAPACK support.
T det(const array &in)
C++ Interface to find the determinant of a matrix.
AFAPI af_err af_det(double *det_real, double *det_imag, const af_array in)
C Interface to find the determinant of a matrix.
AFAPI array inverse(const array &in, const matProp options=AF_MAT_NONE)
C++ Interface to invert a matrix.
AFAPI af_err af_inverse(af_array *out, const af_array in, const af_mat_prop options)
C Interface to invert a matrix.
AFAPI double norm(const array &in, const normType type=AF_NORM_EUCLID, const double p=1, const double q=1)
C++ Interface to find the 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 to find the norm of a matrix.
AFAPI af_err af_pinverse(af_array *out, const af_array in, const double tol, const af_mat_prop options)
C Interface to pseudo-invert (Moore-Penrose) a matrix.
AFAPI array pinverse(const array &in, const double tol=1E-6, const matProp options=AF_MAT_NONE)
C++ Interface to pseudo-invert (Moore-Penrose) a matrix.
AFAPI unsigned rank(const array &in, const double tol=1E-5)
C++ Interface to find the rank of a matrix.
AFAPI af_err af_rank(unsigned *rank, const af_array in, const double tol)
C Interface to find the rank of a matrix.
AFAPI array solve(const array &a, const array &b, const matProp options=AF_MAT_NONE)
C++ Interface to solve 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 to solve a system of equations.
AFAPI array solveLU(const array &a, const array &piv, const array &b, const matProp options=AF_MAT_NONE)
C++ Interface to solve 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 to solve a system of equations.
Definition: algorithm.h:15