Main Page
Modules
Namespaces
Classes
Files
File List
File Members
include
armadillo_bits
lapack_bones.hpp
Go to the documentation of this file.
1
// Copyright (C) 2008-2011 NICTA (www.nicta.com.au)
2
// Copyright (C) 2008-2011 Conrad Sanderson
3
// Copyright (C) 2009 Edmund Highcock
4
// Copyright (C) 2011 James Sanders
5
//
6
// This file is part of the Armadillo C++ library.
7
// It is provided without any warranty of fitness
8
// for any purpose. You can redistribute this file
9
// and/or modify it under the terms of the GNU
10
// Lesser General Public License (LGPL) as published
11
// by the Free Software Foundation, either version 3
12
// of the License or (at your option) any later version.
13
// (see http://www.opensource.org/licenses for more info)
14
15
16
17
#ifdef ARMA_USE_LAPACK
18
19
20
#if !defined(ARMA_BLAS_CAPITALS)
21
22
#define arma_sgetrf sgetrf
23
#define arma_dgetrf dgetrf
24
#define arma_cgetrf cgetrf
25
#define arma_zgetrf zgetrf
26
27
#define arma_sgetri sgetri
28
#define arma_dgetri dgetri
29
#define arma_cgetri cgetri
30
#define arma_zgetri zgetri
31
32
#define arma_strtri strtri
33
#define arma_dtrtri dtrtri
34
#define arma_ctrtri ctrtri
35
#define arma_ztrtri ztrtri
36
37
#define arma_ssyev ssyev
38
#define arma_dsyev dsyev
39
40
#define arma_cheev cheev
41
#define arma_zheev zheev
42
43
#define arma_sgeev sgeev
44
#define arma_dgeev dgeev
45
46
#define arma_cgeev cgeev
47
#define arma_zgeev zgeev
48
49
#define arma_spotrf spotrf
50
#define arma_dpotrf dpotrf
51
#define arma_cpotrf cpotrf
52
#define arma_zpotrf zpotrf
53
54
#define arma_spotri spotri
55
#define arma_dpotri dpotri
56
#define arma_cpotri cpotri
57
#define arma_zpotri zpotri
58
59
#define arma_sgeqrf sgeqrf
60
#define arma_dgeqrf dgeqrf
61
#define arma_cgeqrf cgeqrf
62
#define arma_zgeqrf zgeqrf
63
64
#define arma_sorgqr sorgqr
65
#define arma_dorgqr dorgqr
66
67
#define arma_cungqr cungqr
68
#define arma_zungqr zungqr
69
70
#define arma_sgesvd sgesvd
71
#define arma_dgesvd dgesvd
72
73
#define arma_cgesvd cgesvd
74
#define arma_zgesvd zgesvd
75
76
#define arma_sgesv sgesv
77
#define arma_dgesv dgesv
78
#define arma_cgesv cgesv
79
#define arma_zgesv zgesv
80
81
#define arma_sgels sgels
82
#define arma_dgels dgels
83
#define arma_cgels cgels
84
#define arma_zgels zgels
85
86
#define arma_strtrs strtrs
87
#define arma_dtrtrs dtrtrs
88
#define arma_ctrtrs ctrtrs
89
#define arma_ztrtrs ztrtrs
90
91
#define arma_sgees sgees
92
#define arma_dgees dgees
93
#define arma_cgees cgees
94
#define arma_zgees zgees
95
96
#define arma_strsyl strsyl
97
#define arma_dtrsyl dtrsyl
98
#define arma_ctrsyl ctrsyl
99
#define arma_ztrsyl ztrsyl
100
101
#define arma_ssytrf ssytrf
102
#define arma_dsytrf dsytrf
103
#define arma_csytrf csytrf
104
#define arma_zsytrf zsytrf
105
106
#define arma_ssytri ssytri
107
#define arma_dsytri dsytri
108
#define arma_csytri csytri
109
#define arma_zsytri zsytri
110
111
#else
112
113
#define arma_sgetrf SGETRF
114
#define arma_dgetrf DGETRF
115
#define arma_cgetrf CGETRF
116
#define arma_zgetrf ZGETRF
117
118
#define arma_sgetri SGETRI
119
#define arma_dgetri DGETRI
120
#define arma_cgetri CGETRI
121
#define arma_zgetri ZGETRI
122
123
#define arma_strtri STRTRI
124
#define arma_dtrtri DTRTRI
125
#define arma_ctrtri CTRTRI
126
#define arma_ztrtri ZTRTRI
127
128
#define arma_ssyev SSYEV
129
#define arma_dsyev DSYEV
130
131
#define arma_cheev CHEEV
132
#define arma_zheev ZHEEV
133
134
#define arma_sgeev SGEEV
135
#define arma_dgeev DGEEV
136
137
#define arma_cgeev CGEEV
138
#define arma_zgeev ZGEEV
139
140
#define arma_spotrf SPOTRF
141
#define arma_dpotrf DPOTRF
142
#define arma_cpotrf CPOTRF
143
#define arma_zpotrf ZPOTRF
144
145
#define arma_spotri SPOTRI
146
#define arma_dpotri DPOTRI
147
#define arma_cpotri CPOTRI
148
#define arma_zpotri ZPOTRI
149
150
#define arma_sgeqrf SGEQRF
151
#define arma_dgeqrf DGEQRF
152
#define arma_cgeqrf CGEQRF
153
#define arma_zgeqrf ZGEQRF
154
155
#define arma_sorgqr SORGQR
156
#define arma_dorgqr DORGQR
157
158
#define arma_cungqr CUNGQR
159
#define arma_zungqr ZUNGQR
160
161
#define arma_sgesvd SGESVD
162
#define arma_dgesvd DGESVD
163
164
#define arma_cgesvd CGESVD
165
#define arma_zgesvd ZGESVD
166
167
#define arma_sgesv SGESV
168
#define arma_dgesv DGESV
169
#define arma_cgesv CGESV
170
#define arma_zgesv ZGESV
171
172
#define arma_sgels SGELS
173
#define arma_dgels DGELS
174
#define arma_cgels CGELS
175
#define arma_zgels ZGELS
176
177
#define arma_strtrs STRTRS
178
#define arma_dtrtrs DTRTRS
179
#define arma_ctrtrs CTRTRS
180
#define arma_ztrtrs ZTRTRS
181
182
#define arma_sgees SGEES
183
#define arma_dgees DGEES
184
#define arma_cgees CGEES
185
#define arma_zgees ZGEES
186
187
#define arma_strsyl STRSYL
188
#define arma_dtrsyl DTRSYL
189
#define arma_ctrsyl CTRSYL
190
#define arma_ztrsyl ZTRSYL
191
192
#define arma_ssytrf SSYTRF
193
#define arma_dsytrf DSYTRF
194
#define arma_csytrf CSYTRF
195
#define arma_zsytrf ZSYTRF
196
197
#define arma_ssytri SSYTRI
198
#define arma_dsytri DSYTRI
199
#define arma_csytri CSYTRI
200
#define arma_zsytri ZSYTRI
201
202
#endif
203
204
205
206
extern
"C"
207
{
208
// LU factorisation
209
void
arma_fortran
(arma_sgetrf)(
blas_int
* m,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* ipiv,
blas_int
* info);
210
void
arma_fortran
(arma_dgetrf)(
blas_int
* m,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* ipiv,
blas_int
* info);
211
void
arma_fortran
(arma_cgetrf)(
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
blas_int
* info);
212
void
arma_fortran
(arma_zgetrf)(
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
blas_int
* info);
213
214
// matrix inversion (using LU factorisation result)
215
void
arma_fortran
(arma_sgetri)(
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* ipiv,
float
* work,
blas_int
* lwork,
blas_int
* info);
216
void
arma_fortran
(arma_dgetri)(
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* ipiv,
double
* work,
blas_int
* lwork,
blas_int
* info);
217
void
arma_fortran
(arma_cgetri)(
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* lwork,
blas_int
* info);
218
void
arma_fortran
(arma_zgetri)(
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* lwork,
blas_int
* info);
219
220
// matrix inversion (triangular matrices)
221
void
arma_fortran
(arma_strtri)(
char
* uplo,
char
* diag,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* info);
222
void
arma_fortran
(arma_dtrtri)(
char
* uplo,
char
* diag,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* info);
223
void
arma_fortran
(arma_ctrtri)(
char
* uplo,
char
* diag,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
224
void
arma_fortran
(arma_ztrtri)(
char
* uplo,
char
* diag,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
225
226
// eigenvector decomposition of symmetric real matrices
227
void
arma_fortran
(arma_ssyev)(
char
* jobz,
char
* uplo,
blas_int
* n,
float
* a,
blas_int
* lda,
float
* w,
float
* work,
blas_int
* lwork,
blas_int
* info);
228
void
arma_fortran
(arma_dsyev)(
char
* jobz,
char
* uplo,
blas_int
* n,
double
* a,
blas_int
* lda,
double
* w,
double
* work,
blas_int
* lwork,
blas_int
* info);
229
230
// eigenvector decomposition of hermitian matrices (complex)
231
void
arma_fortran
(arma_cheev)(
char
* jobz,
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
float
* w,
void
* work,
blas_int
* lwork,
float
* rwork,
blas_int
* info);
232
void
arma_fortran
(arma_zheev)(
char
* jobz,
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
double
* w,
void
* work,
blas_int
* lwork,
double
* rwork,
blas_int
* info);
233
234
// eigenvector decomposition of general real matrices
235
void
arma_fortran
(arma_sgeev)(
char
* jobvl,
char
* jobvr,
blas_int
* n,
float
* a,
blas_int
* lda,
float
* wr,
float
* wi,
float
* vl,
blas_int
* ldvl,
float
* vr,
blas_int
* ldvr,
float
* work,
blas_int
* lwork,
blas_int
* info);
236
void
arma_fortran
(arma_dgeev)(
char
* jobvl,
char
* jobvr,
blas_int
* n,
double
* a,
blas_int
* lda,
double
* wr,
double
* wi,
double
* vl,
blas_int
* ldvl,
double
* vr,
blas_int
* ldvr,
double
* work,
blas_int
* lwork,
blas_int
* info);
237
238
// eigenvector decomposition of general complex matrices
239
void
arma_fortran
(arma_cgeev)(
char
* jobvl,
char
* jobvr,
blas_int
* n,
void
* a,
blas_int
* lda,
void
* w,
void
* vl,
blas_int
* ldvl,
void
* vr,
blas_int
* ldvr,
void
* work,
blas_int
* lwork,
float
* rwork,
blas_int
* info);
240
void
arma_fortran
(arma_zgeev)(
char
* jobvl,
char
* jobvr,
blas_int
* n,
void
* a,
blas_int
* lda,
void
* w,
void
* vl,
blas_int
* ldvl,
void
* vr,
blas_int
* ldvr,
void
* work,
blas_int
* lwork,
double
* rwork,
blas_int
* info);
241
242
// Cholesky decomposition
243
void
arma_fortran
(arma_spotrf)(
char
* uplo,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* info);
244
void
arma_fortran
(arma_dpotrf)(
char
* uplo,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* info);
245
void
arma_fortran
(arma_cpotrf)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
246
void
arma_fortran
(arma_zpotrf)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
247
248
// matrix inversion (using Cholesky decomposition result)
249
void
arma_fortran
(arma_spotri)(
char
* uplo,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* info);
250
void
arma_fortran
(arma_dpotri)(
char
* uplo,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* info);
251
void
arma_fortran
(arma_cpotri)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
252
void
arma_fortran
(arma_zpotri)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* info);
253
254
// QR decomposition
255
void
arma_fortran
(arma_sgeqrf)(
blas_int
* m,
blas_int
* n,
float
* a,
blas_int
* lda,
float
* tau,
float
* work,
blas_int
* lwork,
blas_int
* info);
256
void
arma_fortran
(arma_dgeqrf)(
blas_int
* m,
blas_int
* n,
double
* a,
blas_int
* lda,
double
* tau,
double
* work,
blas_int
* lwork,
blas_int
* info);
257
void
arma_fortran
(arma_cgeqrf)(
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
void
* tau,
void
* work,
blas_int
* lwork,
blas_int
* info);
258
void
arma_fortran
(arma_zgeqrf)(
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
void
* tau,
void
* work,
blas_int
* lwork,
blas_int
* info);
259
260
// Q matrix calculation from QR decomposition (real matrices)
261
void
arma_fortran
(arma_sorgqr)(
blas_int
* m,
blas_int
* n,
blas_int
* k,
float
* a,
blas_int
* lda,
float
* tau,
float
* work,
blas_int
* lwork,
blas_int
* info);
262
void
arma_fortran
(arma_dorgqr)(
blas_int
* m,
blas_int
* n,
blas_int
* k,
double
* a,
blas_int
* lda,
double
* tau,
double
* work,
blas_int
* lwork,
blas_int
* info);
263
264
// Q matrix calculation from QR decomposition (complex matrices)
265
void
arma_fortran
(arma_cungqr)(
blas_int
* m,
blas_int
* n,
blas_int
* k,
void
* a,
blas_int
* lda,
void
* tau,
void
* work,
blas_int
* lwork,
blas_int
* info);
266
void
arma_fortran
(arma_zungqr)(
blas_int
* m,
blas_int
* n,
blas_int
* k,
void
* a,
blas_int
* lda,
void
* tau,
void
* work,
blas_int
* lwork,
blas_int
* info);
267
268
// SVD (real matrices)
269
void
arma_fortran
(arma_sgesvd)(
char
* jobu,
char
* jobvt,
blas_int
* m,
blas_int
* n,
float
* a,
blas_int
* lda,
float
* s,
float
* u,
blas_int
* ldu,
float
* vt,
blas_int
* ldvt,
float
* work,
blas_int
* lwork,
blas_int
* info);
270
void
arma_fortran
(arma_dgesvd)(
char
* jobu,
char
* jobvt,
blas_int
* m,
blas_int
* n,
double
* a,
blas_int
* lda,
double
* s,
double
* u,
blas_int
* ldu,
double
* vt,
blas_int
* ldvt,
double
* work,
blas_int
* lwork,
blas_int
* info);
271
272
// SVD (complex matrices)
273
void
arma_fortran
(arma_cgesvd)(
char
* jobu,
char
* jobvt,
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
float
* s,
void
* u,
blas_int
* ldu,
void
* vt,
blas_int
* ldvt,
void
* work,
blas_int
* lwork,
float
* rwork,
blas_int
* info);
274
void
arma_fortran
(arma_zgesvd)(
char
* jobu,
char
* jobvt,
blas_int
* m,
blas_int
* n,
void
* a,
blas_int
* lda,
double
* s,
void
* u,
blas_int
* ldu,
void
* vt,
blas_int
* ldvt,
void
* work,
blas_int
* lwork,
double
* rwork,
blas_int
* info);
275
276
// solve system of linear equations, using LU decomposition
277
void
arma_fortran
(arma_sgesv)(
blas_int
* n,
blas_int
* nrhs,
float
* a,
blas_int
* lda,
blas_int
* ipiv,
float
* b,
blas_int
* ldb,
blas_int
* info);
278
void
arma_fortran
(arma_dgesv)(
blas_int
* n,
blas_int
* nrhs,
double
* a,
blas_int
* lda,
blas_int
* ipiv,
double
* b,
blas_int
* ldb,
blas_int
* info);
279
void
arma_fortran
(arma_cgesv)(
blas_int
* n,
blas_int
* nrhs,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* b,
blas_int
* ldb,
blas_int
* info);
280
void
arma_fortran
(arma_zgesv)(
blas_int
* n,
blas_int
* nrhs,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* b,
blas_int
* ldb,
blas_int
* info);
281
282
// solve over/underdetermined system of linear equations
283
void
arma_fortran
(arma_sgels)(
char
*
trans
,
blas_int
* m,
blas_int
* n,
blas_int
* nrhs,
float
* a,
blas_int
* lda,
float
* b,
blas_int
* ldb,
float
* work,
blas_int
* lwork,
blas_int
* info);
284
void
arma_fortran
(arma_dgels)(
char
*
trans
,
blas_int
* m,
blas_int
* n,
blas_int
* nrhs,
double
* a,
blas_int
* lda,
double
* b,
blas_int
* ldb,
double
* work,
blas_int
* lwork,
blas_int
* info);
285
void
arma_fortran
(arma_cgels)(
char
*
trans
,
blas_int
* m,
blas_int
* n,
blas_int
* nrhs,
void
* a,
blas_int
* lda,
void
* b,
blas_int
* ldb,
void
* work,
blas_int
* lwork,
blas_int
* info);
286
void
arma_fortran
(arma_zgels)(
char
*
trans
,
blas_int
* m,
blas_int
* n,
blas_int
* nrhs,
void
* a,
blas_int
* lda,
void
* b,
blas_int
* ldb,
void
* work,
blas_int
* lwork,
blas_int
* info);
287
288
// solve a triangular system of linear equations
289
void
arma_fortran
(arma_strtrs)(
char
* uplo,
char
*
trans
,
char
* diag,
blas_int
* n,
blas_int
* nrhs,
const
float
* a,
blas_int
* lda,
float
* b,
blas_int
* ldb,
blas_int
* info);
290
void
arma_fortran
(arma_dtrtrs)(
char
* uplo,
char
*
trans
,
char
* diag,
blas_int
* n,
blas_int
* nrhs,
const
double
* a,
blas_int
* lda,
double
* b,
blas_int
* ldb,
blas_int
* info);
291
void
arma_fortran
(arma_ctrtrs)(
char
* uplo,
char
*
trans
,
char
* diag,
blas_int
* n,
blas_int
* nrhs,
const
void
* a,
blas_int
* lda,
void
* b,
blas_int
* ldb,
blas_int
* info);
292
void
arma_fortran
(arma_ztrtrs)(
char
* uplo,
char
*
trans
,
char
* diag,
blas_int
* n,
blas_int
* nrhs,
const
void
* a,
blas_int
* lda,
void
* b,
blas_int
* ldb,
blas_int
* info);
293
294
// Schur decomposition (real matrices)
295
void
arma_fortran
(arma_sgees)(
char
* jobvs,
char
*
sort
,
blas_int
* select,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* sdim,
float
* wr,
float
* wi,
float
* vs,
blas_int
* ldvs,
float
* work,
blas_int
* lwork,
blas_int
* bwork,
blas_int
* info);
296
void
arma_fortran
(arma_dgees)(
char
* jobvs,
char
*
sort
,
blas_int
* select,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* sdim,
double
* wr,
double
* wi,
double
* vs,
blas_int
* ldvs,
double
* work,
blas_int
* lwork,
blas_int
* bwork,
blas_int
* info);
297
298
// Schur decomposition (complex matrices)
299
void
arma_fortran
(arma_cgees)(
char
* jobvs,
char
*
sort
,
blas_int
* select,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* sdim,
void
* w,
void
* vs,
blas_int
* ldvs,
void
* work,
blas_int
* lwork,
float
* rwork,
blas_int
* bwork,
blas_int
* info);
300
void
arma_fortran
(arma_zgees)(
char
* jobvs,
char
*
sort
,
blas_int
* select,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* sdim,
void
* w,
void
* vs,
blas_int
* ldvs,
void
* work,
blas_int
* lwork,
double
* rwork,
blas_int
* bwork,
blas_int
* info);
301
302
// solve a Sylvester equation ax + xb = c, with a and b assumed to be in Schur form
303
void
arma_fortran
(arma_strsyl)(
char
* transa,
char
* transb,
blas_int
* isgn,
blas_int
* m,
blas_int
* n,
const
float
* a,
blas_int
* lda,
const
float
* b,
blas_int
* ldb,
float
* c,
blas_int
* ldc,
float
* scale,
blas_int
* info);
304
void
arma_fortran
(arma_dtrsyl)(
char
* transa,
char
* transb,
blas_int
* isgn,
blas_int
* m,
blas_int
* n,
const
double
* a,
blas_int
* lda,
const
double
* b,
blas_int
* ldb,
double
* c,
blas_int
* ldc,
double
* scale,
blas_int
* info);
305
void
arma_fortran
(arma_ctrsyl)(
char
* transa,
char
* transb,
blas_int
* isgn,
blas_int
* m,
blas_int
* n,
const
void
* a,
blas_int
* lda,
const
void
* b,
blas_int
* ldb,
void
* c,
blas_int
* ldc,
float
* scale,
blas_int
* info);
306
void
arma_fortran
(arma_ztrsyl)(
char
* transa,
char
* transb,
blas_int
* isgn,
blas_int
* m,
blas_int
* n,
const
void
* a,
blas_int
* lda,
const
void
* b,
blas_int
* ldb,
void
* c,
blas_int
* ldc,
double
* scale,
blas_int
* info);
307
308
void
arma_fortran
(arma_ssytrf)(
char
* uplo,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* ipiv,
float
* work,
blas_int
* lwork,
blas_int
* info);
309
void
arma_fortran
(arma_dsytrf)(
char
* uplo,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* ipiv,
double
* work,
blas_int
* lwork,
blas_int
* info);
310
void
arma_fortran
(arma_csytrf)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* lwork,
blas_int
* info);
311
void
arma_fortran
(arma_zsytrf)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* lwork,
blas_int
* info);
312
313
void
arma_fortran
(arma_ssytri)(
char
* uplo,
blas_int
* n,
float
* a,
blas_int
* lda,
blas_int
* ipiv,
float
* work,
blas_int
* info);
314
void
arma_fortran
(arma_dsytri)(
char
* uplo,
blas_int
* n,
double
* a,
blas_int
* lda,
blas_int
* ipiv,
double
* work,
blas_int
* info);
315
void
arma_fortran
(arma_csytri)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* info);
316
void
arma_fortran
(arma_zsytri)(
char
* uplo,
blas_int
* n,
void
* a,
blas_int
* lda,
blas_int
* ipiv,
void
* work,
blas_int
* info);
317
318
// void arma_fortran(arma_dgeqp3)(blas_int* m, blas_int* n, double* a, blas_int* lda, blas_int* jpvt, double* tau, double* work, blas_int* lwork, blas_int* info);
319
// void arma_fortran(arma_dormqr)(char* side, char* trans, blas_int* m, blas_int* n, blas_int* k, double* a, blas_int* lda, double* tau, double* c, blas_int* ldc, double* work, blas_int* lwork, blas_int* info);
320
// void arma_fortran(arma_dposv)(char* uplo, blas_int* n, blas_int* nrhs, double* a, blas_int* lda, double* b, blas_int* ldb, blas_int* info);
321
}
322
323
324
#endif
blas_int
int blas_int
Definition:
typedef_blas_int.hpp:23
trans
arma_inline const Op< T1, op_htrans > trans(const Base< typename T1::elem_type, T1 > &X)
Definition:
fn_trans.hpp:21
arma_fortran
#define arma_fortran(function)
Definition:
compiler_setup.hpp:38
sort
arma_inline const Op< T1, op_sort > sort(const Base< typename T1::elem_type, T1 > &X, const uword sort_type=0, const uword dim=0)
Definition:
fn_sort.hpp:21
armadillo_matrix
Author(s):
autogenerated on Fri Apr 16 2021 02:31:57