14#ifndef LAPACK_FTN_TRAILING_UNDERSCORE
17#define _DGETRF_ dgetrf
18#define _DGETRS_ dgetrs
19#define _DGETRI_ dgetri
20#define _DGELSS_ dgelss
21#define _DGBCON_ dgbcon
23#define _DGBTRF_ dgbtrf
24#define _DGBTRS_ dgbtrs
25#define _DGECON_ dgecon
26#define _DLANGE_ dlange
28#define _DGEQRF_ dgeqrf
29#define _DORMQR_ dormqr
30#define _DTRTRS_ dtrtrs
31#define _DTRCON_ dtrcon
32#define _DPOTRF_ dpotrf
33#define _DPOTRS_ dpotrs
38#define _DGETRF_ dgetrf_
39#define _DGETRS_ dgetrs_
40#define _DGETRI_ dgetri_
41#define _DGELSS_ dgelss_
42#define _DGBCON_ dgbcon_
44#define _DGBTRF_ dgbtrf_
45#define _DGBTRS_ dgbtrs_
46#define _DGECON_ dgecon_
47#define _DLANGE_ dlange_
49#define _DGEQRF_ dgeqrf_
50#define _DORMQR_ dormqr_
51#define _DTRTRS_ dtrtrs_
52#define _DTRCON_ dtrcon_
53#define _DPOTRF_ dpotrf_
54#define _DPOTRS_ dpotrs_
61typedef enum {Transpose = 1, NoTranspose = 0} transpose_t;
62typedef enum {ColMajor = 1, RowMajor = 0} storage_t;
63typedef enum {UpperTriangular = 0, LowerTriangular = 1} upperlower_t;
64typedef enum {Left = 0, Right = 1} side_t;
66const char no_yes[2] = {
'N',
'T'};
67const char upper_lower[2] = {
'U',
'L'};
68const char left_right[2] = {
'L',
'R'};
73#ifdef LAPACK_FTN_STRING_LEN_AT_END
74 int _DGEMV_(
const char* transpose,
75 const integer* m,
const integer* n,
const doublereal* alpha,
76 const doublereal* a,
const integer* lda,
const doublereal* x,
77 const integer* incX,
const doublereal* beta, doublereal* y,
78 const integer* incY, ftnlen trsize);
80 int _DGEMV_(
const char* transpose, ftnlen trsize,
81 const integer* m,
const integer* n,
const doublereal* alpha,
82 const doublereal* a,
const integer* lda,
const doublereal* x,
83 const integer* incX,
const doublereal* beta, doublereal* y,
87 int _DGETRF_(
const integer* m,
const integer* n,
88 doublereal* a, integer* lda, integer* ipiv,
91#ifdef LAPACK_FTN_STRING_LEN_AT_END
92 int _DGETRS_(
const char* transpose,
const integer* n,
93 const integer* nrhs, doublereal* a,
const integer* lda,
94 integer* ipiv, doublereal* b,
const integer* ldb,
95 integer* info, ftnlen trsize);
97 int _DGETRS_(
const char* transpose, ftnlen trsize,
const integer* n,
98 const integer* nrhs,
const doublereal* a,
const integer* lda,
99 integer* ipiv, doublereal* b,
const integer* ldb, integer* info);
102 int _DGETRI_(
const integer* n, doublereal* a,
const integer* lda,
103 integer* ipiv, doublereal* work, integer* lwork, integer* info);
105 int _DGELSS_(
const integer* m,
const integer* n,
const integer* nrhs,
106 doublereal* a,
const integer* lda, doublereal* b,
107 const integer* ldb, doublereal* s,
const doublereal* rcond,
108 integer* rank, doublereal* work, integer* lwork, integer* info);
110 int _DGBSV_(integer* n, integer* kl, integer* ku, integer* nrhs,
111 doublereal* a, integer* lda, integer* ipiv, doublereal* b,
112 integer* ldb, integer* info);
114 int _DGBTRF_(integer* m, integer* n, integer* kl, integer* ku,
115 doublereal* a, integer* lda, integer* ipiv, integer* info);
117#ifdef LAPACK_FTN_STRING_LEN_AT_END
118 int _DGBTRS_(
const char* trans, integer* n, integer* kl, integer* ku,
119 integer* nrhs, doublereal* a, integer* lda, integer* ipiv,
120 doublereal* b, integer* ldb, integer* info, ftnlen trsize);
122 int _DGBTRS_(
const char* trans, ftnlen trsize,
123 integer* n, integer* kl, integer* ku,
124 integer* nrhs, doublereal* a, integer* lda, integer* ipiv,
125 doublereal* b, integer* ldb, integer* info);
128 int _DSCAL_(integer* n, doublereal* da, doublereal* dx, integer* incx);
130 int _DGEQRF_(
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
131 doublereal* tau, doublereal* work,
const integer* lwork, integer* info);
133#ifdef LAPACK_FTN_STRING_LEN_AT_END
134 int _DORMQR_(
const char* side,
const char* trans,
const integer* m,
const integer* n,
135 const integer* k, doublereal* a,
const integer* lda,
136 doublereal* tau, doublereal* c,
const integer* ldc,
137 doublereal* work,
const integer* lwork, integer* info, ftnlen sisize, ftnlen trsize);
139 int _DORMQR_(
const char* side, ftnlen sisize,
const char* trans, ftnlen trsize,
const integer* m,
140 const integer* n,
const integer* k, doublereal* a,
const integer* lda,
141 doublereal* tau,doublereal* c,
const integer* ldc,
142 doublereal* work,
const integer* lwork, integer* info);
145#ifdef LAPACK_FTN_STRING_LEN_AT_END
146 int _DTRTRS_(
const char* uplo,
const char* trans,
const char* diag,
const integer* n,
147 const integer* nrhs, doublereal* a,
const integer* lda,
148 doublereal* b,
const integer* ldb, integer* info,
149 ftnlen upsize, ftnlen trsize, ftnlen disize);
151 int _DTRTRS_(
const char* uplo, ftnlen upsize,
const char* trans, ftnlen trsize,
const char* diag,
152 ftnlen disize,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
153 doublereal* b,
const integer* ldb, integer* info);
157#ifdef LAPACK_FTN_STRING_LEN_AT_END
158 int _DTRCON_(
const char* norm,
const char* uplo,
const char* diag,
const integer* n,
159 doublereal* a,
const integer* lda,
const doublereal* rcond,
160 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize,
161 ftnlen upsize, ftnlen disize);
163 int _DTRCON_(
const char* norm, ftnlen nosize,
const char* uplo, ftnlen upsize,
const char* diag,
164 ftnlen disize,
const integer* n, doublereal* a,
const integer* lda,
const doublereal* rcond,
165 doublereal* work,
const integer* iwork, integer* info);
169#ifdef LAPACK_FTN_STRING_LEN_AT_END
170 int _DPOTRF_(
const char* uplo,
const integer* n, doublereal* a,
const integer* lda, integer* info,
173 int _DPOTRF_(
const char* uplo, ftnlen upsize,
const integer* n, doublereal* a,
const integer* lda,
177#ifdef LAPACK_FTN_STRING_LEN_AT_END
178 int _DPOTRS_(
const char* uplo,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
179 doublereal* b,
const integer* ldb, integer* info, ftnlen upsize);
181 int _DPOTRS_(
const char* uplo, ftnlen upsize,
const integer* n,
const integer* nrhs, doublereal* a,
const integer* lda,
182 doublereal* b,
const integer* ldb, integer* info);
185#ifdef LAPACK_FTN_STRING_LEN_AT_END
186 int _DGECON_(
const char* norm,
const integer* n, doublereal* a,
const integer* lda,
187 const doublereal* rnorm,
const doublereal* rcond,
188 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize);
190 int _DGECON_(
const char* norm, ftnlen nosize,
const integer* n, doublereal* a,
const integer* lda,
191 const doublereal* rnorm,
const doublereal* rcond,
192 doublereal* work,
const integer* iwork, integer* info);
195#ifdef LAPACK_FTN_STRING_LEN_AT_END
196 int _DGBCON_(
const char* norm,
const integer* n, integer* kl, integer* ku, doublereal* ab,
const integer* ldab,
197 const integer* ipiv,
const doublereal* anorm,
const doublereal* rcond,
198 doublereal* work,
const integer* iwork, integer* info, ftnlen nosize);
200 int _DGBCON_(
const char* norm, ftnlen nosize,
const integer* n, integer* kl, integer* ku, doublereal* ab,
const integer* ldab,
201 const integer* ipiv,
const doublereal* anorm,
const doublereal* rcond,
202 doublereal* work,
const integer* iwork, integer* info);
205#ifdef LAPACK_FTN_STRING_LEN_AT_END
206 doublereal _DLANGE_(
const char* norm,
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
207 doublereal* work, ftnlen nosize);
209 doublereal _DLANGE_(
const char* norm, ftnlen nosize,
const integer* m,
const integer* n, doublereal* a,
const integer* lda,
217inline void ct_dgemv(ctlapack::storage_t storage,
218 ctlapack::transpose_t trans,
219 int m,
int n, doublereal alpha,
const doublereal* a,
int lda,
220 const doublereal* x,
int incX, doublereal beta,
221 doublereal* y,
int incY)
223 integer f_m = m, f_n = n, f_lda = lda, f_incX = incX, f_incY = incY;
224 doublereal f_alpha = alpha, f_beta = beta;
226#ifdef LAPACK_FTN_STRING_LEN_AT_END
227 _DGEMV_(&no_yes[trans], &f_m, &f_n, &f_alpha, a,
228 &f_lda, x, &f_incX, &f_beta, y, &f_incY, trsize);
230 _DGEMV_(&no_yes[trans], trsize, &f_m, &f_n, &f_alpha, a,
231 &f_lda, x, &f_incX, &f_beta, y, &f_incY);
235inline void ct_dgbsv(
int n,
int kl,
int ku,
int nrhs,
236 doublereal* a,
int lda, integer* ipiv, doublereal* b,
int ldb,
239 integer f_n = n, f_kl = kl, f_ku = ku, f_nrhs = nrhs, f_lda = lda,
240 f_ldb = ldb, f_info = 0;
241 _DGBSV_(&f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
246inline void ct_dgelss(
size_t m,
size_t n,
size_t nrhs, doublereal* a,
247 size_t lda, doublereal* b,
size_t ldb, doublereal* s,
248 doublereal rcond,
size_t& rank, doublereal* work,
249 int& lwork,
int& info)
251 integer f_m =
static_cast<integer
>(m);
252 integer f_n =
static_cast<integer
>(n);
253 integer f_nrhs =
static_cast<integer
>(nrhs);
254 integer f_lda =
static_cast<integer
>(lda);
255 integer f_ldb =
static_cast<integer
>(ldb);
256 integer f_lwork =
static_cast<integer
>(lwork);
257 integer f_rank, f_info;
259 _DGELSS_(&f_m, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, s, &rcond, &f_rank,
260 work, &f_lwork, &f_info);
262 info =
static_cast<int>(f_info);
263 rank =
static_cast<size_t>(f_rank);
264 lwork =
static_cast<int>(f_lwork);
267inline void ct_dgbtrf(
size_t m,
size_t n,
size_t kl,
size_t ku,
268 doublereal* a,
size_t lda, integer* ipiv,
int& info)
270 integer f_m = (int) m;
271 integer f_n = (int) n;
272 integer f_kl = (int) kl;
273 integer f_ku = (int) ku;
274 integer f_lda = (int) lda;
276 _DGBTRF_(&f_m, &f_n, &f_kl, &f_ku, a, &f_lda, ipiv, &f_info);
280inline void ct_dgbtrs(ctlapack::transpose_t trans,
size_t n,
281 size_t kl,
size_t ku,
size_t nrhs, doublereal* a,
size_t lda,
282 integer* ipiv, doublereal* b,
size_t ldb,
int& info)
284 integer f_n = (int) n;
285 integer f_kl = (int) kl;
286 integer f_ku = (int) ku;
287 integer f_nrhs = (int) nrhs;
288 integer f_lda = (int) lda;
289 integer f_ldb = (int) ldb;
291 char tr = no_yes[trans];
293#ifdef LAPACK_FTN_STRING_LEN_AT_END
294 _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
295 b, &f_ldb, &f_info, trsize);
297 _DGBTRS_(&tr, trsize, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
303inline void ct_dgetrf(
size_t m,
size_t n,
304 doublereal* a,
size_t lda, integer* ipiv,
int& info)
306 integer mm = (int) m;
307 integer nn = (int) n;
308 integer ldaa = (int) lda;
310 _DGETRF_(&mm, &nn, a, &ldaa, ipiv, &infoo);
314inline void ct_dgetrs(ctlapack::transpose_t trans,
size_t n,
315 size_t nrhs, doublereal* a,
size_t lda,
316 integer* ipiv, doublereal* b,
size_t ldb,
int& info)
318 integer f_n = (int) n;
319 integer f_lda = (int) lda;
320 integer f_nrhs = (int) nrhs;
321 integer f_ldb = (int) ldb;
323 char tr = no_yes[trans];
325#ifdef LAPACK_FTN_STRING_LEN_AT_END
326 _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info, trsize);
328 _DGETRS_(&tr, trsize, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info);
333inline void ct_dgetri(
int n, doublereal* a,
int lda, integer* ipiv,
334 doublereal* work,
int lwork,
int& info)
336 integer f_n = n, f_lda = lda, f_lwork = lwork, f_info = 0;
337 _DGETRI_(&f_n, a, &f_lda, ipiv, work, &f_lwork, &f_info);
340inline void ct_dscal(
int n, doublereal da, doublereal* dx,
int incx)
342 integer f_n = n, f_incx = incx;
343 _DSCAL_(&f_n, &da, dx, &f_incx);
346inline void ct_dgeqrf(
size_t m,
size_t n, doublereal* a,
size_t lda, doublereal* tau,
347 doublereal* work,
size_t lwork,
int& info)
349 integer f_m =
static_cast<integer
>(m);
350 integer f_n =
static_cast<integer
>(n);
351 integer f_lda =
static_cast<integer
>(lda);
352 integer f_lwork =
static_cast<integer
>(lwork);
354 _DGEQRF_(&f_m, &f_n, a, &f_lda, tau, work, &f_lwork, &f_info);
358inline void ct_dormqr(ctlapack::side_t rlside, ctlapack::transpose_t trans,
size_t m,
359 size_t n,
size_t k, doublereal* a,
size_t lda, doublereal* tau, doublereal* c,
size_t ldc,
360 doublereal* work,
size_t lwork,
int& info)
362 char side = left_right[rlside];
363 char tr = no_yes[trans];
364 integer f_m =
static_cast<integer
>(m);
365 integer f_n =
static_cast<integer
>(n);
366 integer f_k =
static_cast<integer
>(k);
367 integer f_lwork =
static_cast<integer
>(lwork);
368 integer f_lda =
static_cast<integer
>(lda);
369 integer f_ldc =
static_cast<integer
>(ldc);
372#ifdef LAPACK_FTN_STRING_LEN_AT_END
373 _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info, trsize, trsize);
375 _DORMQR_(&side, trsize, &tr, trsize, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
380inline void ct_dtrtrs(ctlapack::upperlower_t uplot, ctlapack::transpose_t trans,
const char* diag,
381 size_t n,
size_t nrhs, doublereal* a,
size_t lda, doublereal* b,
size_t ldb,
int& info)
383 char uplo = upper_lower[uplot];
384 char tr = no_yes[trans];
389 integer f_n =
static_cast<integer
>(n);
390 integer f_nrhs =
static_cast<integer
>(nrhs);
391 integer f_lda =
static_cast<integer
>(lda);
392 integer f_ldb =
static_cast<integer
>(ldb);
395#ifdef LAPACK_FTN_STRING_LEN_AT_END
396 _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize, trsize, trsize);
398 _DTRTRS_(&uplo, trsize, &tr, trsize, &dd, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
407inline doublereal
ct_dtrcon(
const char* norm, ctlapack::upperlower_t uplot,
const char* diag,
408 size_t n, doublereal* a,
size_t lda, doublereal* work,
int* iwork,
int& info)
410 char uplo = upper_lower[uplot];
419 integer f_n =
static_cast<integer
>(n);
420 integer f_lda =
static_cast<integer
>(lda);
422 doublereal rcond = 0.0;
424#ifdef LAPACK_FTN_STRING_LEN_AT_END
425 _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info, trsize, trsize, trsize);
427 _DTRCON_(&nn, trsize, &uplo, trsize, &dd, trsize, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
433inline void ct_dpotrf(ctlapack::upperlower_t uplot,
size_t n, doublereal* a,
size_t lda,
int& info)
435 char uplo = upper_lower[uplot];
436 integer f_n =
static_cast<integer
>(n);
437 integer f_lda =
static_cast<integer
>(lda);
440#ifdef LAPACK_FTN_STRING_LEN_AT_END
441 _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info, trsize);
443 _DPOTRF_(&uplo, trsize, &f_n, a, &f_lda, &f_info);
449inline void ct_dpotrs(ctlapack::upperlower_t uplot,
size_t n,
size_t nrhs, doublereal* a,
size_t lda,
450 doublereal* b,
size_t ldb,
int& info)
452 char uplo = upper_lower[uplot];
453 integer f_n =
static_cast<integer
>(n);
454 integer f_nrhs =
static_cast<integer
>(nrhs);
455 integer f_lda =
static_cast<integer
>(lda);
456 integer f_ldb =
static_cast<integer
>(ldb);
459#ifdef LAPACK_FTN_STRING_LEN_AT_END
460 _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize);
462 _DPOTRS_(&uplo, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
468inline doublereal ct_dgecon(
const char norm,
size_t n, doublereal* a,
size_t lda, doublereal anorm,
469 doublereal* work,
int* iwork,
int& info)
475 integer f_n =
static_cast<integer
>(n);
476 integer f_lda =
static_cast<integer
>(lda);
478 doublereal rcond = 0.0;
480#ifdef LAPACK_FTN_STRING_LEN_AT_END
481 _DGECON_(&cnorm, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info, trsize);
483 _DGECON_(&cnorm, trsize, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
489inline doublereal ct_dgbcon(
const char norm,
size_t n,
size_t kl,
size_t ku,
490 doublereal* a,
size_t ldab,
int* ipiv, doublereal anorm,
491 doublereal* work,
int* iwork,
int& info)
497 integer f_n =
static_cast<integer
>(n);
498 integer f_kl =
static_cast<integer
>(kl);
499 integer f_ku =
static_cast<integer
>(ku);
500 integer f_ldab =
static_cast<integer
>(ldab);
502 doublereal rcond = 0.0;
504#ifdef LAPACK_FTN_STRING_LEN_AT_END
505 _DGBCON_(&cnorm, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info, trsize);
507 _DGBCON_(&cnorm, trsize, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
513inline doublereal ct_dlange(
const char norm,
size_t m,
size_t n, doublereal* a,
size_t lda,
520 integer f_m =
static_cast<integer
>(m);
521 integer f_n =
static_cast<integer
>(n);
522 integer f_lda =
static_cast<integer
>(lda);
525#ifdef LAPACK_FTN_STRING_LEN_AT_END
526 anorm = _DLANGE_(&cnorm, &f_m, &f_n, a, &f_lda, work, trsize);
528 anorm = _DLANGE_(&cnorm, trsize, &f_m, &f_n, a, &f_lda, work);
This file contains definitions of constants, types and terms that are used in internal routines and a...
Namespace for the Cantera kernel.
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)