Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ctlapack.h
Go to the documentation of this file.
1 /**
2  * @file ctlapack.h
3  */
4 // Copyright 2001 California Institute of Technology.
5 
6 #ifndef CT_CTLAPACK_H
7 #define CT_CTLAPACK_H
8 
9 #ifdef DARWIN
10 #undef NO_FTN_STRING_LEN_AT_END
11 #endif
12 
13 #include "cantera/base/ct_defs.h"
14 
15 // map BLAS names to names with or without a trailing underscore.
16 #ifndef LAPACK_FTN_TRAILING_UNDERSCORE
17 
18 #define _DGEMV_ dgemv
19 #define _DGETRF_ dgetrf
20 #define _DGETRS_ dgetrs
21 #define _DGETRI_ dgetri
22 #define _DGELSS_ dgelss
23 #define _DGBCON_ dgbcon
24 #define _DGBSV_ dgbsv
25 #define _DGBTRF_ dgbtrf
26 #define _DGBTRS_ dgbtrs
27 #define _DGECON_ dgecon
28 #define _DLANGE_ dlange
29 #define _DSCAL_ dscal
30 #define _DGEQRF_ dgeqrf
31 #define _DORMQR_ dormqr
32 #define _DTRTRS_ dtrtrs
33 #define _DTRCON_ dtrcon
34 #define _DPOTRF_ dpotrf
35 #define _DPOTRS_ dpotrs
36 
37 #else
38 
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_
57 
58 #endif
59 
60 
61 namespace ctlapack
62 {
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;
67 }
68 const char no_yes[2] = {'N', 'T'};
69 const char upper_lower[2] = {'U', 'L'};
70 const char left_right[2] = {'L', 'R'};
71 
72 // C interfaces for Fortran Lapack routines
73 extern "C" {
74 
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);
81 #else
82 
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,
87  const integer* incY);
88 #endif
89 
90  int _DGETRF_(const integer* m, const integer* n,
91  doublereal* a, integer* lda, integer* ipiv,
92  integer* info);
93 
94 #ifdef LAPACK_FTN_STRING_LEN_AT_END
95 
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);
100 
101 #else
102 
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);
106 
107 #endif
108 
109  int _DGETRI_(const integer* n, doublereal* a, const integer* lda,
110  integer* ipiv, doublereal* work, integer* lwork, integer* info);
111 
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);
116 
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);
120 
121  int _DGBTRF_(integer* m, integer* n, integer* kl, integer* ku,
122  doublereal* a, integer* lda, integer* ipiv, integer* info);
123 
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);
128 #else
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);
133 #endif
134 
135  int _DSCAL_(integer* n, doublereal* da, doublereal* dx, integer* incx);
136 
137  int _DGEQRF_(const integer* m, const integer* n, doublereal* a, const integer* lda,
138  doublereal* tau, doublereal* work, const integer* lwork, integer* info);
139 
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);
145 #else
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);
150 #endif
151 
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);
157 #else
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);
161 #endif
162 
163 
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);
169 #else
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);
173 #endif
174 
175 
176 #ifdef LAPACK_FTN_STRING_LEN_AT_END
177  int _DPOTRF_(const char* uplo, const integer* n, doublereal* a, const integer* lda, integer* info,
178  ftnlen upsize);
179 #else
180  int _DPOTRF_(const char* uplo, ftnlen upsize, const integer* n, doublereal* a, const integer* lda,
181  integer* info);
182 #endif
183 
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);
187 #else
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);
190 #endif
191 
192 
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);
197 #else
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);
201 #endif
202 
203 
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);
208 #else
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);
212 #endif
213 
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);
217 #else
218  doublereal _DLANGE_(const char* norm, ftnlen nosize, const integer* m, const integer* n, doublereal* a, const integer* lda,
219  doublereal* work);
220 #endif
221 }
222 
223 namespace Cantera
224 {
225 
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)
231 {
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;
234  ftnlen trsize = 1;
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);
238 #else
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);
242 #else
243  _DGEMV_(&no_yes[trans], trsize, &f_m, &f_n, &f_alpha, a,
244  &f_lda, x, &f_incX, &f_beta, y, &f_incY);
245 #endif
246 #endif
247 }
248 
249 inline void ct_dgbsv(int n, int kl, int ku, int nrhs,
250  doublereal* a, int lda, integer* ipiv, doublereal* b, int ldb,
251  int& info)
252 {
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,
256  b, &f_ldb, &f_info);
257  info = f_info;
258 }
259 
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)
264 {
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;
272 
273  _DGELSS_(&f_m, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, s, &rcond, &f_rank,
274  work, &f_lwork, &f_info);
275 
276  info = static_cast<int>(f_info);
277  rank = static_cast<size_t>(f_rank);
278  lwork = static_cast<int>(f_lwork);
279 }
280 
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)
283 {
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;
289  integer f_info = 0;
290  _DGBTRF_(&f_m, &f_n, &f_kl, &f_ku, a, &f_lda, ipiv, &f_info);
291  info = f_info;
292 }
293 
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)
297 {
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;
304  integer f_info = 0;
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,
308  b, &f_ldb, &f_info);
309 #else
310  ftnlen trsize = 1;
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);
314 #else
315  _DGBTRS_(&tr, trsize, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
316  b, &f_ldb, &f_info);
317 #endif
318 #endif
319  info = f_info;
320 }
321 
322 inline void ct_dgetrf(size_t m, size_t n,
323  doublereal* a, size_t lda, integer* ipiv, int& info)
324 {
325  integer mm = (int) m;
326  integer nn = (int) n;
327  integer ldaa = (int) lda;
328  integer infoo = 0;
329  _DGETRF_(&mm, &nn, a, &ldaa, ipiv, &infoo);
330  info = infoo;
331 }
332 
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)
336 {
337  integer f_n = (int) n;
338  integer f_lda = (int) lda;
339  integer f_nrhs = (int) nrhs;
340  integer f_ldb = (int) ldb;
341  integer f_info = 0;
342  char tr = no_yes[trans];
343 
344 #ifdef NO_FTN_STRING_LEN_AT_END
345  _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb,
346  &f_info);
347 #else
348  ftnlen trsize = 1;
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);
351 #else
352  _DGETRS_(&tr, trsize, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info);
353 #endif
354 #endif
355  info = f_info;
356 }
357 
358 inline void ct_dgetri(int n, doublereal* a, int lda, integer* ipiv,
359  doublereal* work, int lwork, int& info)
360 {
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);
363 }
364 
365 inline void ct_dscal(int n, doublereal da, doublereal* dx, int incx)
366 {
367  integer f_n = n, f_incx = incx;
368  _DSCAL_(&f_n, &da, dx, &f_incx);
369 }
370 
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)
373 {
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);
378  integer f_info = 0;
379  _DGEQRF_(&f_m, &f_n, a, &f_lda, tau, work, &f_lwork, &f_info);
380  info = f_info;
381 }
382 
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)
386 {
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);
395  integer f_info = 0;
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);
398 #else
399  ftnlen trsize = 1;
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);
402 #else
403  _DORMQR_(&side, trsize, &tr, trsize, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
404 #endif
405 #endif
406  info = f_info;
407 }
408 
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)
411 {
412  char uplo = upper_lower[uplot];
413  char tr = no_yes[trans];
414  char dd = 'N';
415  if (diag) {
416  dd = diag[0];
417  }
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);
422  integer f_info = 0;
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);
425 #else
426  ftnlen trsize = 1;
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);
429 #else
430  _DTRTRS_(&uplo, trsize, &tr, trsize, &dd, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
431 #endif
432 #endif
433  info = f_info;
434 }
435 
436 //!
437 /*!
438  * @param work Must be dimensioned equal to greater than 3N
439  * @param iwork Must be dimensioned equal to or greater than N
440  */
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)
443 {
444  char uplo = upper_lower[uplot];
445  char dd = 'N';
446  if (diag) {
447  dd = diag[0];
448  }
449  char nn = '1';
450  if (norm) {
451  nn = norm[0];
452  }
453  integer f_n = static_cast<integer>(n);
454  integer f_lda = static_cast<integer>(lda);
455  integer f_info = 0;
456  doublereal rcond;
457 #ifdef NO_FTN_STRING_LEN_AT_END
458  _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
459 #else
460  ftnlen trsize = 1;
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);
463 #else
464  _DTRCON_(&nn, trsize, &uplo, trsize, &dd, trsize, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
465 #endif
466 #endif
467  info = f_info;
468  return rcond;
469 }
470 
471 inline void ct_dpotrf(ctlapack::upperlower_t uplot, size_t n, doublereal* a, size_t lda, int& info)
472 {
473  char uplo = upper_lower[uplot];
474  integer f_n = static_cast<integer>(n);
475  integer f_lda = static_cast<integer>(lda);
476  integer f_info = 0;
477 
478 #ifdef NO_FTN_STRING_LEN_AT_END
479  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info);
480 #else
481  ftnlen trsize = 1;
482 #ifdef LAPACK_FTN_STRING_LEN_AT_END
483  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info, trsize);
484 #else
485  _DPOTRF_(&uplo, trsize, &f_n, a, &f_lda, &f_info);
486 #endif
487 #endif
488  info = f_info;
489  return;
490 }
491 
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)
494 {
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);
500  integer f_info = 0;
501 
502 #ifdef NO_FTN_STRING_LEN_AT_END
503  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
504 #else
505  ftnlen trsize = 1;
506 #ifdef LAPACK_FTN_STRING_LEN_AT_END
507  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize);
508 #else
509  _DPOTRS_(&uplo, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
510 #endif
511 #endif
512  info = f_info;
513  return;
514 }
515 
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)
518 {
519  char cnorm = '1';
520  if (norm) {
521  cnorm = norm;
522  }
523  integer f_n = static_cast<integer>(n);
524  integer f_lda = static_cast<integer>(lda);
525  integer f_info = 0;
526  doublereal rcond;
527 
528 #ifdef NO_FTN_STRING_LEN_AT_END
529  _DGECON_(&cnorm, &f_n a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
530 #else
531  ftnlen trsize = 1;
532 #ifdef LAPACK_FTN_STRING_LEN_AT_END
533  _DGECON_(&cnorm, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info, trsize);
534 #else
535  _DGECON_(&cnorm, trsize, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
536 #endif
537 #endif
538  info = f_info;
539  return rcond;
540 }
541 
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)
545 {
546  char cnorm = '1';
547  if (norm) {
548  cnorm = norm;
549  }
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);
554  integer f_info = 0;
555  doublereal rcond;
556 
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);
559 #else
560  ftnlen trsize = 1;
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);
563 #else
564  _DGBCON_(&cnorm, trsize, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
565 #endif
566 #endif
567  info = f_info;
568  return rcond;
569 }
570 
571 inline doublereal ct_dlange(const char norm, size_t m, size_t n, doublereal* a, size_t lda,
572  doublereal* work)
573 {
574  char cnorm = '1';
575  if (norm) {
576  cnorm = norm;
577  }
578  integer f_m = static_cast<integer>(m);
579  integer f_n = static_cast<integer>(n);
580  integer f_lda = static_cast<integer>(lda);
581  doublereal anorm;
582 
583 #ifdef NO_FTN_STRING_LEN_AT_END
584  anorm = _DLANGE_(&cnorm, &f_m, &f_n a, &f_lda, work);
585 #else
586  ftnlen trsize = 1;
587 #ifdef LAPACK_FTN_STRING_LEN_AT_END
588  anorm = _DLANGE_(&cnorm, &f_m, &f_n, a, &f_lda, work, trsize);
589 #else
590  anorm = _DLANGE_(&cnorm, trsize, &f_m, &f_n, a, &f_lda, work);
591 #endif
592 #endif
593  return anorm;
594 }
595 
596 }
597 
598 #endif
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)
Definition: ctlapack.h:441