10 #undef NO_FTN_STRING_LEN_AT_END
16 #ifndef LAPACK_FTN_TRAILING_UNDERSCORE
19 #define _DGETRF_ dgetrf
20 #define _DGETRS_ dgetrs
21 #define _DGETRI_ dgetri
22 #define _DGELSS_ dgelss
23 #define _DGBCON_ dgbcon
25 #define _DGBTRF_ dgbtrf
26 #define _DGBTRS_ dgbtrs
27 #define _DGECON_ dgecon
28 #define _DLANGE_ dlange
30 #define _DGEQRF_ dgeqrf
31 #define _DORMQR_ dormqr
32 #define _DTRTRS_ dtrtrs
33 #define _DTRCON_ dtrcon
34 #define _DPOTRF_ dpotrf
35 #define _DPOTRS_ dpotrs
39 #define _DGEMV_ dgemv_
40 #define _DGETRF_ dgetrf_
41 #define _DGETRS_ dgetrs_
42 #define _DGETRI_ dgetri_
43 #define _DGELSS_ dgelss_
44 #define _DGBCON_ dgbcon_
45 #define _DGBSV_ dgbsv_
46 #define _DGBTRF_ dgbtrf_
47 #define _DGBTRS_ dgbtrs_
48 #define _DGECON_ dgecon_
49 #define _DLANGE_ dlange_
50 #define _DSCAL_ dscal_
51 #define _DGEQRF_ dgeqrf_
52 #define _DORMQR_ dormqr_
53 #define _DTRTRS_ dtrtrs_
54 #define _DTRCON_ dtrcon_
55 #define _DPOTRF_ dpotrf_
56 #define _DPOTRS_ dpotrs_
63 typedef enum {Transpose = 1, NoTranspose = 0} transpose_t;
64 typedef enum {ColMajor = 1, RowMajor = 0} storage_t;
65 typedef enum {UpperTriangular = 0, LowerTriangular = 1} upperlower_t;
66 typedef enum {Left = 0, Right = 1} side_t;
68 const char no_yes[2] = {
'N',
'T'};
69 const char upper_lower[2] = {
'U',
'L'};
70 const char left_right[2] = {
'L',
'R'};
75 #ifdef LAPACK_FTN_STRING_LEN_AT_END
76 int _DGEMV_(
const char* transpose,
77 const integer* m,
const integer* n,
const doublereal* alpha,
78 const doublereal* a,
const integer* lda,
const doublereal* x,
79 const integer* incX,
const doublereal* beta, doublereal* y,
80 const integer* incY, ftnlen trsize);
83 int _DGEMV_(
const char* transpose, ftnlen trsize,
84 const integer* m,
const integer* n,
const doublereal* alpha,
85 const doublereal* a,
const integer* lda,
const doublereal* x,
86 const integer* incX,
const doublereal* beta, doublereal* y,
90 int _DGETRF_(
const integer* m,
const integer* n,
91 doublereal* a, integer* lda, integer* ipiv,
94 #ifdef LAPACK_FTN_STRING_LEN_AT_END
96 int _DGETRS_(
const char* transpose,
const integer* n,
97 const integer* nrhs, doublereal* a,
const integer* lda,
98 integer* ipiv, doublereal* b,
const integer* ldb,
99 integer* info, ftnlen trsize);
103 int _DGETRS_(
const char* transpose, ftnlen trsize,
const integer* n,
104 const integer* nrhs,
const doublereal* a,
const integer* lda,
105 integer* ipiv, doublereal* b,
const integer* ldb, integer* info);
109 int _DGETRI_(
const integer* n, doublereal* a,
const integer* lda,
110 integer* ipiv, doublereal* work, integer* lwork, integer* info);
112 int _DGELSS_(
const integer* m,
const integer* n,
const integer* nrhs,
113 doublereal* a,
const integer* lda, doublereal* b,
114 const integer* ldb, doublereal* s,
const doublereal* rcond,
115 integer* rank, doublereal* work, integer* lwork, integer* info);
117 int _DGBSV_(integer* n, integer* kl, integer* ku, integer* nrhs,
118 doublereal* a, integer* lda, integer* ipiv, doublereal* b,
119 integer* ldb, integer* info);
121 int _DGBTRF_(integer* m, integer* n, integer* kl, integer* ku,
122 doublereal* a, integer* lda, integer* ipiv, integer* info);
124 #ifdef LAPACK_FTN_STRING_LEN_AT_END
125 int _DGBTRS_(
const char* trans, integer* n, integer* kl, integer* ku,
126 integer* nrhs, doublereal* a, integer* lda, integer* ipiv,
127 doublereal* b, integer* ldb, integer* info, ftnlen trsize);
129 int _DGBTRS_(
const char* trans, ftnlen trsize,
130 integer* n, integer* kl, integer* ku,
131 integer* nrhs, doublereal* a, integer* lda, integer* ipiv,
132 doublereal* b, integer* ldb, integer* info);
135 int _DSCAL_(integer* n, doublereal* da, doublereal* dx, integer* incx);
137 int _DGEQRF_(
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
138 doublereal* tau, doublereal* work,
const integer* lwork, integer* info);
140 #ifdef LAPACK_FTN_STRING_LEN_AT_END
141 int _DORMQR_(
const char* side,
const char* trans,
const integer* m,
const integer* n,
142 const integer* k, doublereal* a,
const integer* lda,
143 doublereal* tau, doublereal* c,
const integer* ldc,
144 doublereal* work,
const integer* lwork, integer* info, ftnlen sisize, ftnlen trsize);
146 int _DORMQR_(
const char* side, ftnlen sisize,
const char* trans, ftnlen trsize,
const integer* m,
147 const integer* n,
const integer* k, doublereal* a,
const integer* lda,
148 doublereal* tau,doublereal* c,
const integer* ldc,
149 doublereal* work,
const integer* lwork, integer* info);
152 #ifdef LAPACK_FTN_STRING_LEN_AT_END
153 int _DTRTRS_(
const char* uplo,
const char* trans,
const char* diag,
const integer* n,
154 const integer* nrhs, doublereal* a,
const integer* lda,
155 doublereal* b,
const integer* ldb, integer* info,
156 ftnlen upsize, ftnlen trsize, ftnlen disize);
158 int _DTRTRS_(
const char* uplo, ftnlen upsize,
const char* trans, ftnlen trsize,
const char* diag,
159 ftnlen disize,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
160 doublereal* b,
const integer* ldb, integer* info);
164 #ifdef LAPACK_FTN_STRING_LEN_AT_END
165 int _DTRCON_(
const char* norm,
const char* uplo,
const char* diag,
const integer* n,
166 doublereal* a,
const integer* lda,
const doublereal* rcond,
167 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize,
168 ftnlen upsize, ftnlen disize);
170 int _DTRCON_(
const char* norm, ftnlen nosize,
const char* uplo, ftnlen upsize,
const char* diag,
171 ftnlen disize,
const integer* n, doublereal* a,
const integer* lda,
const doublereal* rcond,
172 doublereal* work,
const integer* iwork, integer* info);
176 #ifdef LAPACK_FTN_STRING_LEN_AT_END
177 int _DPOTRF_(
const char* uplo,
const integer* n, doublereal* a,
const integer* lda, integer* info,
180 int _DPOTRF_(
const char* uplo, ftnlen upsize,
const integer* n, doublereal* a,
const integer* lda,
184 #ifdef LAPACK_FTN_STRING_LEN_AT_END
185 int _DPOTRS_(
const char* uplo,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
186 doublereal* b,
const integer* ldb, integer* info, ftnlen upsize);
188 int _DPOTRS_(
const char* uplo, ftnlen upsize,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
189 doublereal* b,
const integer* ldb, integer* info);
193 #ifdef LAPACK_FTN_STRING_LEN_AT_END
194 int _DGECON_(
const char* norm,
const integer* n, doublereal* a,
const integer* lda,
195 const doublereal* rnorm,
const doublereal* rcond,
196 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize);
198 int _DGECON_(
const char* norm, ftnlen nosize,
const integer* n, doublereal* a,
const integer* lda,
199 const doublereal* rnorm,
const doublereal* rcond,
200 doublereal* work,
const integer* iwork, integer* info);
204 #ifdef LAPACK_FTN_STRING_LEN_AT_END
205 int _DGBCON_(
const char* norm,
const integer* n, integer* kl, integer* ku, doublereal* ab,
const integer* ldab,
206 const integer* ipiv,
const doublereal* anorm,
const doublereal* rcond,
207 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize);
209 int _DGBCON_(
const char* norm, ftnlen nosize,
const integer* n, integer* kl, integer* ku, doublereal* ab,
const integer* ldab,
210 const integer* ipiv,
const doublereal* anorm,
const doublereal* rcond,
211 doublereal* work,
const integer* iwork, integer* info);
214 #ifdef LAPACK_FTN_STRING_LEN_AT_END
215 doublereal _DLANGE_(
const char* norm,
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
216 doublereal* work, ftnlen nosize);
218 doublereal _DLANGE_(
const char* norm, ftnlen nosize,
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
226 inline void ct_dgemv(ctlapack::storage_t storage,
227 ctlapack::transpose_t trans,
228 int m,
int n, doublereal alpha,
const doublereal* a,
int lda,
229 const doublereal* x,
int incX, doublereal beta,
230 doublereal* y,
int incY)
232 integer f_m = m, f_n = n, f_lda = lda, f_incX = incX, f_incY = incY;
233 doublereal f_alpha = alpha, f_beta = beta;
235 #ifdef NO_FTN_STRING_LEN_AT_END
236 _DGEMV_(&no_yes[trans], &f_m, &f_n, &f_alpha, a,
237 &f_lda, x, &f_incX, &f_beta, y, &f_incY);
239 #ifdef LAPACK_FTN_STRING_LEN_AT_END
240 _DGEMV_(&no_yes[trans], &f_m, &f_n, &f_alpha, a,
241 &f_lda, x, &f_incX, &f_beta, y, &f_incY, trsize);
243 _DGEMV_(&no_yes[trans], trsize, &f_m, &f_n, &f_alpha, a,
244 &f_lda, x, &f_incX, &f_beta, y, &f_incY);
249 inline void ct_dgbsv(
int n,
int kl,
int ku,
int nrhs,
250 doublereal* a,
int lda, integer* ipiv, doublereal* b,
int ldb,
253 integer f_n = n, f_kl = kl, f_ku = ku, f_nrhs = nrhs, f_lda = lda,
254 f_ldb = ldb, f_info = 0;
255 _DGBSV_(&f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
260 inline void ct_dgelss(
size_t m,
size_t n,
size_t nrhs, doublereal* a,
261 size_t lda, doublereal* b,
size_t ldb, doublereal* s,
262 doublereal rcond,
size_t& rank, doublereal* work,
263 int& lwork,
int& info)
265 integer f_m =
static_cast<integer
>(m);
266 integer f_n =
static_cast<integer
>(n);
267 integer f_nrhs =
static_cast<integer
>(nrhs);
268 integer f_lda =
static_cast<integer
>(lda);
269 integer f_ldb =
static_cast<integer
>(ldb);
270 integer f_lwork =
static_cast<integer
>(lwork);
271 integer f_rank, f_info;
273 _DGELSS_(&f_m, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, s, &rcond, &f_rank,
274 work, &f_lwork, &f_info);
276 info =
static_cast<int>(f_info);
277 rank =
static_cast<size_t>(f_rank);
278 lwork =
static_cast<int>(f_lwork);
281 inline void ct_dgbtrf(
size_t m,
size_t n,
size_t kl,
size_t ku,
282 doublereal* a,
size_t lda, integer* ipiv,
int& info)
284 integer f_m = (int) m;
285 integer f_n = (int) n;
286 integer f_kl = (int) kl;
287 integer f_ku = (int) ku;
288 integer f_lda = (int) lda;
290 _DGBTRF_(&f_m, &f_n, &f_kl, &f_ku, a, &f_lda, ipiv, &f_info);
294 inline void ct_dgbtrs(ctlapack::transpose_t trans,
size_t n,
295 size_t kl,
size_t ku,
size_t nrhs, doublereal* a,
size_t lda,
296 integer* ipiv, doublereal* b,
size_t ldb,
int& info)
298 integer f_n = (int) n;
299 integer f_kl = (int) kl;
300 integer f_ku = (int) ku;
301 integer f_nrhs = (int) nrhs;
302 integer f_lda = (int) lda;
303 integer f_ldb = (int) ldb;
305 char tr = no_yes[trans];
306 #ifdef NO_FTN_STRING_LEN_AT_END
307 _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
311 #ifdef LAPACK_FTN_STRING_LEN_AT_END
312 _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
313 b, &f_ldb, &f_info, trsize);
315 _DGBTRS_(&tr, trsize, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
322 inline void ct_dgetrf(
size_t m,
size_t n,
323 doublereal* a,
size_t lda, integer* ipiv,
int& info)
325 integer mm = (int) m;
326 integer nn = (int) n;
327 integer ldaa = (int) lda;
329 _DGETRF_(&mm, &nn, a, &ldaa, ipiv, &infoo);
333 inline void ct_dgetrs(ctlapack::transpose_t trans,
size_t n,
334 size_t nrhs, doublereal* a,
size_t lda,
335 integer* ipiv, doublereal* b,
size_t ldb,
int& info)
337 integer f_n = (int) n;
338 integer f_lda = (int) lda;
339 integer f_nrhs = (int) nrhs;
340 integer f_ldb = (int) ldb;
342 char tr = no_yes[trans];
344 #ifdef NO_FTN_STRING_LEN_AT_END
345 _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb,
349 #ifdef LAPACK_FTN_STRING_LEN_AT_END
350 _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info, trsize);
352 _DGETRS_(&tr, trsize, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info);
358 inline void ct_dgetri(
int n, doublereal* a,
int lda, integer* ipiv,
359 doublereal* work,
int lwork,
int& info)
361 integer f_n = n, f_lda = lda, f_lwork = lwork, f_info = 0;
362 _DGETRI_(&f_n, a, &f_lda, ipiv, work, &f_lwork, &f_info);
365 inline void ct_dscal(
int n, doublereal da, doublereal* dx,
int incx)
367 integer f_n = n, f_incx = incx;
368 _DSCAL_(&f_n, &da, dx, &f_incx);
371 inline void ct_dgeqrf(
size_t m,
size_t n, doublereal* a,
size_t lda, doublereal* tau,
372 doublereal* work,
size_t lwork,
int& info)
374 integer f_m =
static_cast<integer
>(m);
375 integer f_n =
static_cast<integer
>(n);
376 integer f_lda =
static_cast<integer
>(lda);
377 integer f_lwork =
static_cast<integer
>(lwork);
379 _DGEQRF_(&f_m, &f_n, a, &f_lda, tau, work, &f_lwork, &f_info);
383 inline void ct_dormqr(ctlapack::side_t rlside, ctlapack::transpose_t trans,
size_t m,
384 size_t n,
size_t k, doublereal* a,
size_t lda, doublereal* tau, doublereal* c,
size_t ldc,
385 doublereal* work,
size_t lwork,
int& info)
387 char side = left_right[rlside];
388 char tr = no_yes[trans];
389 integer f_m =
static_cast<integer
>(m);
390 integer f_n =
static_cast<integer
>(n);
391 integer f_k =
static_cast<integer
>(k);
392 integer f_lwork =
static_cast<integer
>(lwork);
393 integer f_lda =
static_cast<integer
>(lda);
394 integer f_ldc =
static_cast<integer
>(ldc);
396 #ifdef NO_FTN_STRING_LEN_AT_END
397 _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
400 #ifdef LAPACK_FTN_STRING_LEN_AT_END
401 _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info, trsize, trsize);
403 _DORMQR_(&side, trsize, &tr, trsize, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
409 inline void ct_dtrtrs(ctlapack::upperlower_t uplot, ctlapack::transpose_t trans,
const char* diag,
410 size_t n,
size_t nrhs, doublereal* a,
size_t lda, doublereal* b,
size_t ldb,
int& info)
412 char uplo = upper_lower[uplot];
413 char tr = no_yes[trans];
418 integer f_n =
static_cast<integer
>(n);
419 integer f_nrhs =
static_cast<integer
>(nrhs);
420 integer f_lda =
static_cast<integer
>(lda);
421 integer f_ldb =
static_cast<integer
>(ldb);
423 #ifdef NO_FTN_STRING_LEN_AT_END
424 _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
427 #ifdef LAPACK_FTN_STRING_LEN_AT_END
428 _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize, trsize, trsize);
430 _DTRTRS_(&uplo, trsize, &tr, trsize, &dd, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
441 inline doublereal
ct_dtrcon(
const char* norm, ctlapack::upperlower_t uplot,
const char* diag,
442 size_t n, doublereal* a,
size_t lda, doublereal* work,
int* iwork,
int& info)
444 char uplo = upper_lower[uplot];
453 integer f_n =
static_cast<integer
>(n);
454 integer f_lda =
static_cast<integer
>(lda);
457 #ifdef NO_FTN_STRING_LEN_AT_END
458 _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
461 #ifdef LAPACK_FTN_STRING_LEN_AT_END
462 _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info, trsize, trsize, trsize);
464 _DTRCON_(&nn, trsize, &uplo, trsize, &dd, trsize, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
471 inline void ct_dpotrf(ctlapack::upperlower_t uplot,
size_t n, doublereal* a,
size_t lda,
int& info)
473 char uplo = upper_lower[uplot];
474 integer f_n =
static_cast<integer
>(n);
475 integer f_lda =
static_cast<integer
>(lda);
478 #ifdef NO_FTN_STRING_LEN_AT_END
479 _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info);
482 #ifdef LAPACK_FTN_STRING_LEN_AT_END
483 _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info, trsize);
485 _DPOTRF_(&uplo, trsize, &f_n, a, &f_lda, &f_info);
492 inline void ct_dpotrs(ctlapack::upperlower_t uplot,
size_t n,
size_t nrhs, doublereal* a,
size_t lda,
493 doublereal* b,
size_t ldb,
int& info)
495 char uplo = upper_lower[uplot];
496 integer f_n =
static_cast<integer
>(n);
497 integer f_nrhs =
static_cast<integer
>(nrhs);
498 integer f_lda =
static_cast<integer
>(lda);
499 integer f_ldb =
static_cast<integer
>(ldb);
502 #ifdef NO_FTN_STRING_LEN_AT_END
503 _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
506 #ifdef LAPACK_FTN_STRING_LEN_AT_END
507 _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize);
509 _DPOTRS_(&uplo, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
516 inline doublereal ct_dgecon(
const char norm,
size_t n, doublereal* a,
size_t lda, doublereal anorm,
517 doublereal* work,
int* iwork,
int& info)
523 integer f_n =
static_cast<integer
>(n);
524 integer f_lda =
static_cast<integer
>(lda);
528 #ifdef NO_FTN_STRING_LEN_AT_END
529 _DGECON_(&cnorm, &f_n a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
532 #ifdef LAPACK_FTN_STRING_LEN_AT_END
533 _DGECON_(&cnorm, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info, trsize);
535 _DGECON_(&cnorm, trsize, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
542 inline doublereal ct_dgbcon(
const char norm,
size_t n,
size_t kl,
size_t ku,
543 doublereal* a,
size_t ldab,
int* ipiv, doublereal anorm,
544 doublereal* work,
int* iwork,
int& info)
550 integer f_n =
static_cast<integer
>(n);
551 integer f_kl =
static_cast<integer
>(kl);
552 integer f_ku =
static_cast<integer
>(ku);
553 integer f_ldab =
static_cast<integer
>(ldab);
557 #ifdef NO_FTN_STRING_LEN_AT_END
558 _DGBCON_(&cnorm, &f_n , &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
561 #ifdef LAPACK_FTN_STRING_LEN_AT_END
562 _DGBCON_(&cnorm, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info, trsize);
564 _DGBCON_(&cnorm, trsize, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
571 inline doublereal ct_dlange(
const char norm,
size_t m,
size_t n, doublereal* a,
size_t lda,
578 integer f_m =
static_cast<integer
>(m);
579 integer f_n =
static_cast<integer
>(n);
580 integer f_lda =
static_cast<integer
>(lda);
583 #ifdef NO_FTN_STRING_LEN_AT_END
584 anorm = _DLANGE_(&cnorm, &f_m, &f_n a, &f_lda, work);
587 #ifdef LAPACK_FTN_STRING_LEN_AT_END
588 anorm = _DLANGE_(&cnorm, &f_m, &f_n, a, &f_lda, work, trsize);
590 anorm = _DLANGE_(&cnorm, trsize, &f_m, &f_n, a, &f_lda, work);
This file contains definitions of terms that are used in internal routines and are unlikely to need m...
doublereal ct_dtrcon(const char *norm, ctlapack::upperlower_t uplot, const char *diag, size_t n, doublereal *a, size_t lda, doublereal *work, int *iwork, int &info)