This file is indexed.

/usr/include/af/lapack.h is in libarrayfire-dev 3.2.2+dfsg1-2.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  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
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
/*******************************************************
 * Copyright (c) 2014, ArrayFire
 * All rights reserved.
 *
 * This file is distributed under 3-clause BSD license.
 * The complete license agreement can be obtained at:
 * http://arrayfire.com/licenses/BSD-3-Clause
 ********************************************************/

#pragma once
#include <af/array.h>
#include <af/defines.h>

#ifdef __cplusplus
namespace af
{
#if AF_API_VERSION >= 31
    /**
       C++ Interface for SVD decomposition

       \param[out] u is the output array containing U
       \param[out] s is the output array containing the diagonal values of sigma, (singular values of the input matrix))
       \param[out] vt is the output array containing V^H
       \param[in] in is the input matrix

       \ingroup lapack_factor_func_svd
    */
    AFAPI void svd(array &u, array &s, array &vt, const array &in);
#endif

#if AF_API_VERSION >= 31
    /**
       C++ Interface for SVD decomposition

       \param[out] u is the output array containing U
       \param[out] s is the output array containing the diagonal values of sigma, (singular values of the input matrix))
       \param[out] vt is the output array containing V^H
       \param[inout] in is the input matrix and will contain random data after this operation

       \ingroup lapack_factor_func_svd
    */
    AFAPI void svdInPlace(array &u, array &s, array &vt, array &in);
#endif

    /**
       C++ Interface for LU decomposition in packed format

       \param[out] out is the output array containing the packed LU decomposition
       \param[out] pivot will contain the permutation indices to map the input to the decomposition
       \param[in] in is the input matrix
       \param[in] is_lapack_piv specifies if the pivot is returned in original LAPACK compliant format

       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_lu
    */
    AFAPI void lu(array &out, array &pivot, const array &in, const bool is_lapack_piv=true);

    /**
       C++ Interface for LU decomposition

       \param[out] lower will contain the lower triangular matrix of the LU decomposition
       \param[out] upper will contain the upper triangular matrix of the LU decomposition
       \param[out] pivot will contain the permutation indices to map the input to the decomposition
       \param[in] in is the input matrix

       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_lu
    */
    AFAPI void lu(array &lower, array &upper, array &pivot, const array &in);

    /**
      C++ Interface for in place LU decomposition

      \param[out] pivot will contain the permutation indices to map the input to the decomposition
      \param[inout] in contains the input on entry, the packed LU decomposition on exit
      \param[in] is_lapack_piv specifies if the pivot is returned in original LAPACK compliant format

      \note This function is not supported in GFOR

      \ingroup lapack_factor_func_lu
    */
    AFAPI void luInPlace(array &pivot, array &in, const bool is_lapack_piv=true);

    /**
       C++ Interface for QR decomposition in packed format

       \param[out] out is the output array containing the packed QR decomposition
       \param[out] tau will contain additional information needed for unpacking the data
       \param[in] in is the input matrix

       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_qr
    */
    AFAPI void qr(array &out, array &tau, const array &in);

    /**
       C++ Interface for QR decomposition

       \param[out] q is the orthogonal matrix from QR decomposition
       \param[out] r is the upper triangular matrix from QR decomposition
       \param[out] tau will contain additional information needed for solving a least squares problem using \p q and \p r
       \param[in] in is the input matrix

       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_qr
    */
    AFAPI void qr(array &q, array &r, array &tau, const array &in);

    /**
       C++ Interface for QR decomposition

       \param[out] tau will contain additional information needed for unpacking the data
       \param[inout] in is the input matrix on entry. It contains packed QR decomposition on exit

       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_qr
    */
    AFAPI void qrInPlace(array &tau, array &in);

    /**
       C++ Interface for cholesky decomposition

       \param[out] out contains the triangular matrix. Multiply \p out with its conjugate transpose reproduces the input \p in.
       \param[in] in is the input matrix
       \param[in] is_upper a boolean determining if \p out is upper or lower triangular

       \returns \p 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.

       \note The input matrix \b has to be a positive definite matrix, if it is not zero, the cholesky decomposition functions return a non-zero output.
       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_cholesky
    */
    AFAPI int cholesky(array &out, const array &in, const bool is_upper = true);

    /**
       C++ Interface for in place cholesky decomposition

       \param[inout] in is the input matrix on entry. It contains the triangular matrix on exit.
       \param[in] is_upper a boolean determining if \p in is upper or lower triangular

       \returns \p 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.

       \note The input matrix \b has to be a positive definite matrix, if it is not zero, the cholesky decomposition functions return a non-zero output.
       \note This function is not supported in GFOR

       \ingroup lapack_factor_func_cholesky
    */
    AFAPI int choleskyInPlace(array &in, const bool is_upper = true);

    /**
       C++ Interface for solving a system of equations

       \param[in] a is the coefficient matrix
       \param[in] b is the measured values
       \param[in] options determining various properties of matrix \p a
       \returns \p x, the matrix of unknown variables

       \note \p options needs to be one of \ref AF_MAT_NONE, \ref AF_MAT_LOWER or \ref AF_MAT_UPPER
       \note This function is not supported in GFOR

       \ingroup lapack_solve_func_gen
    */
    AFAPI array solve(const array &a, const array &b, const matProp options = AF_MAT_NONE);


    /**
       C++ Interface for solving a system of equations

       \param[in] a is the output matrix from packed LU decomposition of the coefficient matrix
       \param[in] piv is the pivot array from packed LU decomposition of the coefficient matrix
       \param[in] b is the matrix of measured values
       \param[in] options determining various properties of matrix \p a
       \returns \p x, the matrix of unknown variables

       \ingroup lapack_solve_lu_func_gen

       \note \p options currently needs to be \ref AF_MAT_NONE
       \note This function is not supported in GFOR
    */
    AFAPI array solveLU(const array &a, const array &piv,
                        const array &b, const matProp options = AF_MAT_NONE);

    /**
       C++ Interface for inverting a matrix

       \param[in] in is input matrix
       \param[in] options determining various properties of matrix \p in
       \returns \p x, the inverse of the input matrix

       \note \p options currently needs to be \ref AF_MAT_NONE
       \note This function is not supported in GFOR

       \ingroup lapack_ops_func_inv
    */
    AFAPI array inverse(const array &in, const matProp options = AF_MAT_NONE);

    /**
       C++ Interface for finding the rank of a matrix

       \param[in] in is input matrix
       \param[in] tol is the tolerance value

       \returns the rank of the matrix

       \ingroup lapack_ops_func_rank
    */
    AFAPI unsigned rank(const array &in, const double tol=1E-5);

    /**
       C++ Interface for finding the determinant of a matrix

       \param[in] in is input matrix

       \returns the determinant of the matrix

       \ingroup lapack_ops_func_det
    */
    template<typename T> T det(const array &in);

    /**
       C++ Interface for norm of a matrix

       \param[in] in is the input matrix
       \param[in] type specifies the \ref af::normType. Default: \ref AF_NORM_VECTOR_1
       \param[in] p specifies the value of P when \p type is one of \ref AF_NORM_VECTOR_P, AF_NORM_MATRIX_L_PQ is used. It is ignored for other values of \p type
       \param[in] q specifies the value of Q when \p type is AF_NORM_MATRIX_L_PQ. This parameter is ignored if \p type is anything else

       \returns the norm of \p inbased on \p type

       \ingroup lapack_ops_func_norm
    */
    AFAPI double norm(const array &in, const normType type=AF_NORM_EUCLID,
                      const double p=1, const double q=1);
}
#endif

#ifdef __cplusplus
extern "C" {
#endif

#if AF_API_VERSION >= 31
    /**
       C Interface for SVD decomposition

       \param[out] u is the output array containing U
       \param[out] s is the output array containing the diagonal values of sigma, (singular values of the input matrix))
       \param[out] vt is the output array containing V^H
       \param[in] in is the input matrix

       \ingroup lapack_factor_func_svd
    */
    AFAPI af_err af_svd(af_array *u, af_array *s, af_array *vt, const af_array in);
#endif

#if AF_API_VERSION >= 31
    /**
       C Interface for SVD decomposition

       \param[out] u is the output array containing U
       \param[out] s is the output array containing the diagonal values of sigma, (singular values of the input matrix))
       \param[out] vt is the output array containing V^H
       \param[inout] in is the input matrix that will contain random data after this operation

       \ingroup lapack_factor_func_svd
    */
    AFAPI af_err af_svd_inplace(af_array *u, af_array *s, af_array *vt, af_array in);
#endif

    /**
       C Interface for LU decomposition

       \param[out] lower will contain the lower triangular matrix of the LU decomposition
       \param[out] upper will contain the upper triangular matrix of the LU decomposition
       \param[out] pivot will contain the permutation indices to map the input to the decomposition
       \param[in] in is the input matrix

       \ingroup lapack_factor_func_lu
    */
    AFAPI af_err af_lu(af_array *lower, af_array *upper, af_array *pivot, const af_array in);

    /**
       C Interface for in place LU decomposition

       \param[out] pivot will contain the permutation indices to map the input to the decomposition
       \param[inout] in contains the input on entry, the packed LU decomposition on exit
       \param[in] is_lapack_piv specifies if the pivot is returned in original LAPACK compliant format

       \ingroup lapack_factor_func_lu
    */
    AFAPI af_err af_lu_inplace(af_array *pivot, af_array in, const bool is_lapack_piv);

    /**
       C Interface for QR decomposition

       \param[out] q is the orthogonal matrix from QR decomposition
       \param[out] r is the upper triangular matrix from QR decomposition
       \param[out] tau will contain additional information needed for solving a least squares problem using \p q and \p r
       \param[in] in is the input matrix

       \ingroup lapack_factor_func_qr
    */
    AFAPI af_err af_qr(af_array *q, af_array *r, af_array *tau, const af_array in);

    /**
       C Interface for QR decomposition

       \param[out] tau will contain additional information needed for unpacking the data
       \param[inout] in is the input matrix on entry. It contains packed QR decomposition on exit

       \ingroup lapack_factor_func_qr
    */
    AFAPI af_err af_qr_inplace(af_array *tau, af_array in);

    /**
       C++ Interface for cholesky decomposition

       \param[out] out contains the triangular matrix. Multiply \p out with it conjugate transpose reproduces the input \p in.
       \param[out] info is \p 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.
       \param[in] in is the input matrix
       \param[in] is_upper a boolean determining if \p out is upper or lower triangular

       \note The input matrix \b has to be a positive definite matrix, if it is not zero, the cholesky decomposition functions return a non zero output.

       \ingroup lapack_factor_func_cholesky
    */
    AFAPI af_err af_cholesky(af_array *out, int *info, const af_array in, const bool is_upper);

    /**
       C Interface for in place cholesky decomposition

       \param[out] info is \p 0 if cholesky decomposition passes, if not it returns the rank at which the decomposition failed.
       \param[inout] in is the input matrix on entry. It contains the triangular matrix on exit.
       \param[in] is_upper a boolean determining if \p in is upper or lower triangular

       \note The input matrix \b has to be a positive definite matrix, if it is not zero, the cholesky decomposition functions return a non zero output.

       \ingroup lapack_factor_func_cholesky
    */
    AFAPI af_err af_cholesky_inplace(int *info, af_array in, const bool is_upper);

    /**
       C Interface for solving a system of equations

       \param[out] x is the matrix of unknown variables
       \param[in] a is the coefficient matrix
       \param[in] b is the measured values
       \param[in] options determining various properties of matrix \p a

       \ingroup lapack_solve_func_gen

       \note \p options needs to be one of \ref AF_MAT_NONE, \ref AF_MAT_LOWER or \ref AF_MAT_UPPER
    */
    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

       \param[out] x will contain the matrix of unknown variables
       \param[in] a is the output matrix from packed LU decomposition of the coefficient matrix
       \param[in] piv is the pivot array from packed LU decomposition of the coefficient matrix
       \param[in] b is the matrix of measured values
       \param[in] options determining various properties of matrix \p a

       \ingroup lapack_solve_lu_func_gen

       \note \p options currently needs to be \ref AF_MAT_NONE
       \note This function is not supported in GFOR
    */
    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 inverting a matrix

       \param[out] out will contain the inverse of matrix \p in
       \param[in] in is input matrix
       \param[in] options determining various properties of matrix \p in

       \ingroup lapack_ops_func_inv

       \note currently options needs to be \ref AF_MAT_NONE
    */
    AFAPI af_err af_inverse(af_array *out, const af_array in, const af_mat_prop options);

    /**
       C Interface for finding the rank of a matrix

       \param[out] rank will contain the rank of \p in
       \param[in] in is input matrix
       \param[in] tol is the tolerance value

       \ingroup lapack_ops_func_rank
    */
    AFAPI af_err af_rank(unsigned *rank, const af_array in, const double tol);

    /**
       C Interface for finding the determinant of a matrix

       \param[out] det_real will contain the real part of the determinant of \p in
       \param[out] det_imag will contain the imaginary part of the determinant of \p in
       \param[in] in is input matrix

       \ingroup lapack_ops_func_det
    */
    AFAPI af_err af_det(double *det_real, double *det_imag, const af_array in);

    /**
       C Interface for norm of a matrix

       \param[out] out will contain the norm of \p in
       \param[in] in is the input matrix
       \param[in] type specifies the \ref af::normType. Default: \ref AF_NORM_VECTOR_1
       \param[in] p specifies the value of P when \p type is one of \ref AF_NORM_VECTOR_P,  AF_NORM_MATRIX_L_PQ is used. It is ignored for other values of \p type
       \param[in] q specifies the value of Q when \p type is AF_NORM_MATRIX_L_PQ. This parameter is ignored if \p type is anything else


       \ingroup lapack_ops_func_norm
    */
    AFAPI af_err af_norm(double *out, const af_array in, const af_norm_type type, const double p, const double q);


#ifdef __cplusplus
}
#endif