Cantera  2.0
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 _DGBCON_ dgbcon
23 #define _DGBSV_ dgbsv
24 #define _DGBTRF_ dgbtrf
25 #define _DGBTRS_ dgbtrs
26 #define _DGECON_ dgecon
27 #define _DLANGE_ dlange
28 #define _DSCAL_ dscal
29 #define _DGEQRF_ dgeqrf
30 #define _DORMQR_ dormqr
31 #define _DTRTRS_ dtrtrs
32 #define _DTRCON_ dtrcon
33 #define _DPOTRF_ dpotrf
34 #define _DPOTRS_ dpotrs
35 
36 #else
37 
38 #define _DGEMV_ dgemv_
39 #define _DGETRF_ dgetrf_
40 #define _DGETRS_ dgetrs_
41 #define _DGETRI_ dgetri_
42 #define _DGBCON_ dgbcon_
43 #define _DGBSV_ dgbsv_
44 #define _DGBTRF_ dgbtrf_
45 #define _DGBTRS_ dgbtrs_
46 #define _DGECON_ dgecon_
47 #define _DLANGE_ dlange_
48 #define _DSCAL_ dscal_
49 #define _DGEQRF_ dgeqrf_
50 #define _DORMQR_ dormqr_
51 #define _DTRTRS_ dtrtrs_
52 #define _DTRCON_ dtrcon_
53 #define _DPOTRF_ dpotrf_
54 #define _DPOTRS_ dpotrs_
55 
56 #endif
57 
58 
59 namespace ctlapack
60 {
61 typedef enum {Transpose = 1, NoTranspose = 0} transpose_t;
62 typedef enum {ColMajor = 1, RowMajor = 0} storage_t;
63 typedef enum {UpperTriangular = 0, LowerTriangular = 1} upperlower_t;
64 typedef enum {Left = 0, Right = 1} side_t;
65 }
66 const char no_yes[2] = {'N', 'T'};
67 const char upper_lower[2] = {'U', 'L'};
68 const char left_right[2] = {'L', 'R'};
69 
70 // C interfaces for Fortran Lapack routines
71 extern "C" {
72 
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);
79 #else
80 
81  int _DGEMV_(const char* transpose, ftnlen trsize,
82  const integer* m, const integer* n, const doublereal* alpha,
83  const doublereal* a, const integer* lda, const doublereal* x,
84  const integer* incX, const doublereal* beta, doublereal* y,
85  const integer* incY);
86 #endif
87 
88  int _DGETRF_(const integer* m, const integer* n,
89  doublereal* a, integer* lda, integer* ipiv,
90  integer* info);
91 
92 #ifdef LAPACK_FTN_STRING_LEN_AT_END
93 
94  int _DGETRS_(const char* transpose, const integer* n,
95  const integer* nrhs, doublereal* a, const integer* lda,
96  integer* ipiv, doublereal* b, const integer* ldb,
97  integer* info, ftnlen trsize);
98 
99 #else
100 
101  int _DGETRS_(const char* transpose, ftnlen trsize, const integer* n,
102  const integer* nrhs, const doublereal* a, const integer* lda,
103  integer* ipiv, doublereal* b, const integer* ldb, integer* info);
104 
105 #endif
106 
107  int _DGETRI_(const integer* n, doublereal* a, const integer* lda,
108  integer* ipiv, doublereal* work, integer* lwork, integer* info);
109 
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);
113 
114  int _DGBTRF_(integer* m, integer* n, integer* kl, integer* ku,
115  doublereal* a, integer* lda, integer* ipiv, integer* info);
116 
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);
121 #else
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);
126 #endif
127 
128  int _DSCAL_(integer* n, doublereal* da, doublereal* dx, integer* incx);
129 
130  int _DGEQRF_(const integer* m, const integer* n, doublereal* a, const integer* lda,
131  doublereal* tau, doublereal* work, const integer* lwork, integer* info);
132 
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);
138 #else
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);
143 #endif
144 
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);
150 #else
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);
154 #endif
155 
156 
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);
162 #else
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);
166 #endif
167 
168 
169 #ifdef LAPACK_FTN_STRING_LEN_AT_END
170  int _DPOTRF_(const char* uplo, const integer* n, doublereal* a, const integer* lda, integer* info,
171  ftnlen upsize);
172 #else
173  int _DPOTRF_(const char* uplo, ftnlen upsize, const integer* n, doublereal* a, const integer* lda,
174  integer* info);
175 #endif
176 
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);
180 #else
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);
183 #endif
184 
185 
186 #ifdef LAPACK_FTN_STRING_LEN_AT_END
187  int _DGECON_(const char* norm, const integer* n, doublereal* a, const integer* lda,
188  const doublereal* rnorm, const doublereal* rcond,
189  doublereal* work, const integer* iwork, integer* info, ftnlen nosize);
190 #else
191  int _DGECON_(const char* norm, ftnlen nosize, const integer* n, doublereal* a, const integer* lda,
192  const doublereal* rnorm, const doublereal* rcond,
193  doublereal* work, const integer* iwork, integer* info);
194 #endif
195 
196 
197 #ifdef LAPACK_FTN_STRING_LEN_AT_END
198  int _DGBCON_(const char* norm, const integer* n, integer* kl, integer* ku, doublereal* ab, const integer* ldab,
199  const integer* ipiv, const doublereal* anorm, const doublereal* rcond,
200  doublereal* work, const integer* iwork, integer* info, ftnlen nosize);
201 #else
202  int _DGBCON_(const char* norm, ftnlen nosize, const integer* n, integer* kl, integer* ku, doublereal* ab, const integer* ldab,
203  const integer* ipiv, const doublereal* anorm, const doublereal* rcond,
204  doublereal* work, const integer* iwork, integer* info);
205 #endif
206 
207 #ifdef LAPACK_FTN_STRING_LEN_AT_END
208  doublereal _DLANGE_(const char* norm, const integer* m, const integer* n, doublereal* a, const integer* lda,
209  doublereal* work, ftnlen nosize);
210 #else
211  doublereal _DLANGE_(const char* norm, ftnlen nosize, const integer* m, const integer* n, doublereal* a, const integer* lda,
212  doublereal* work);
213 #endif
214 }
215 
216 namespace Cantera
217 {
218 //====================================================================================================================
219 inline void ct_dgemv(ctlapack::storage_t storage,
220  ctlapack::transpose_t trans,
221  int m, int n, doublereal alpha, const doublereal* a, int lda,
222  const doublereal* x, int incX, doublereal beta,
223  doublereal* y, int incY)
224 {
225  integer f_m = m, f_n = n, f_lda = lda, f_incX = incX, f_incY = incY;
226  doublereal f_alpha = alpha, f_beta = beta;
227  ftnlen trsize = 1;
228 #ifdef NO_FTN_STRING_LEN_AT_END
229  _DGEMV_(&no_yes[trans], &f_m, &f_n, &f_alpha, a,
230  &f_lda, x, &f_incX, &f_beta, y, &f_incY);
231 #else
232 #ifdef LAPACK_FTN_STRING_LEN_AT_END
233  _DGEMV_(&no_yes[trans], &f_m, &f_n, &f_alpha, a,
234  &f_lda, x, &f_incX, &f_beta, y, &f_incY, trsize);
235 #else
236  _DGEMV_(&no_yes[trans], trsize, &f_m, &f_n, &f_alpha, a,
237  &f_lda, x, &f_incX, &f_beta, y, &f_incY);
238 #endif
239 #endif
240 }
241 
242 //====================================================================================================================
243 inline void ct_dgbsv(int n, int kl, int ku, int nrhs,
244  doublereal* a, int lda, integer* ipiv, doublereal* b, int ldb,
245  int& info)
246 {
247  integer f_n = n, f_kl = kl, f_ku = ku, f_nrhs = nrhs, f_lda = lda,
248  f_ldb = ldb, f_info = info;
249  _DGBSV_(&f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
250  b, &f_ldb, &f_info);
251  info = f_info;
252 }
253 
254 inline void ct_dgbtrf(size_t m, size_t n, size_t kl, size_t ku,
255  doublereal* a, size_t lda, integer* ipiv, int& info)
256 {
257  integer f_m = (int) m;
258  integer f_n = (int) n;
259  integer f_kl = (int) kl;
260  integer f_ku = (int) ku;
261  integer f_lda = (int) lda;
262  integer f_info = info;
263  _DGBTRF_(&f_m, &f_n, &f_kl, &f_ku, a, &f_lda, ipiv, &f_info);
264  info = f_info;
265 }
266 
267 inline void ct_dgbtrs(ctlapack::transpose_t trans, size_t n,
268  size_t kl, size_t ku, size_t nrhs, doublereal* a, size_t lda,
269  integer* ipiv, doublereal* b, size_t ldb, int& info)
270 {
271  integer f_n = (int) n;
272  integer f_kl = (int) kl;
273  integer f_ku = (int) ku;
274  integer f_nrhs = (int) nrhs;
275  integer f_lda = (int) lda;
276  integer f_ldb = (int) ldb;
277  integer f_info = info;
278  char tr = no_yes[trans];
279 #ifdef NO_FTN_STRING_LEN_AT_END
280  _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
281  b, &f_ldb, &f_info);
282 #else
283  ftnlen trsize = 1;
284 #ifdef LAPACK_FTN_STRING_LEN_AT_END
285  _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
286  b, &f_ldb, &f_info, trsize);
287 #else
288  _DGBTRS_(&tr, trsize, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
289  b, &f_ldb, &f_info);
290 #endif
291 #endif
292  info = f_info;
293 }
294 
295 inline void ct_dgetrf(size_t m, size_t n,
296  doublereal* a, size_t lda, integer* ipiv, int& info)
297 {
298  integer mm = (int) m;
299  integer nn = (int) n;
300  integer ldaa = (int) lda;
301  integer infoo = info;
302  _DGETRF_(&mm, &nn, a, &ldaa, ipiv, &infoo);
303  info = infoo;
304 }
305 
306 inline void ct_dgetrs(ctlapack::transpose_t trans, size_t n,
307  size_t nrhs, doublereal* a, size_t lda,
308  integer* ipiv, doublereal* b, size_t ldb, int& info)
309 {
310  integer f_n = (int) n;
311  integer f_lda = (int) lda;
312  integer f_nrhs = (int) nrhs;
313  integer f_ldb = (int) ldb;
314  integer f_info = info;
315  char tr = no_yes[trans];
316 
317 #ifdef NO_FTN_STRING_LEN_AT_END
318  _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb,
319  &f_info);
320 #else
321  ftnlen trsize = 1;
322 #ifdef LAPACK_FTN_STRING_LEN_AT_END
323  _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info, trsize);
324 #else
325  _DGETRS_(&tr, trsize, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info);
326 #endif
327 #endif
328  info = f_info;
329 }
330 //====================================================================================================================
331 inline void ct_dgetri(int n, doublereal* a, int lda, integer* ipiv,
332  doublereal* work, int lwork, int& info)
333 {
334  integer f_n = n, f_lda = lda, f_lwork = lwork, f_info = info;
335  _DGETRI_(&f_n, a, &f_lda, ipiv, work, &f_lwork, &f_info);
336 }
337 
338 inline void ct_dscal(int n, doublereal da, doublereal* dx, int incx)
339 {
340  integer f_n = n, f_incx = incx;
341  _DSCAL_(&f_n, &da, dx, &f_incx);
342 }
343 //====================================================================================================================
344 inline void ct_dgeqrf(size_t m, size_t n, doublereal* a, size_t lda, doublereal* tau,
345  doublereal* work, size_t lwork, int& info)
346 {
347  integer f_m = static_cast<integer>(m);
348  integer f_n = static_cast<integer>(n);
349  integer f_lda = static_cast<integer>(lda);
350  integer f_lwork = static_cast<integer>(lwork);
351  integer f_info = info;
352  _DGEQRF_(&f_m, &f_n, a, &f_lda, tau, work, &f_lwork, &f_info);
353  info = f_info;
354 }
355 //====================================================================================================================
356 inline void ct_dormqr(ctlapack::side_t rlside, ctlapack::transpose_t trans, size_t m,
357  size_t n, size_t k, doublereal* a, size_t lda, doublereal* tau, doublereal* c, size_t ldc,
358  doublereal* work, size_t lwork, int& info)
359 {
360  char side = left_right[rlside];
361  char tr = no_yes[trans];
362  integer f_m = static_cast<integer>(m);
363  integer f_n = static_cast<integer>(n);
364  integer f_k = static_cast<integer>(k);
365  integer f_lwork = static_cast<integer>(lwork);
366  integer f_lda = static_cast<integer>(lda);
367  integer f_ldc = static_cast<integer>(ldc);
368  integer f_info = info;
369 #ifdef NO_FTN_STRING_LEN_AT_END
370  _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
371 #else
372  ftnlen trsize = 1;
373 #ifdef LAPACK_FTN_STRING_LEN_AT_END
374  _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info, trsize, trsize);
375 #else
376  _DORMQR_(&side, trsize, &tr, trsize, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
377 #endif
378 #endif
379  info = f_info;
380 }
381 //====================================================================================================================
382 inline void ct_dtrtrs(ctlapack::upperlower_t uplot, ctlapack::transpose_t trans, const char* diag,
383  size_t n, size_t nrhs, doublereal* a, size_t lda, doublereal* b, size_t ldb, int& info)
384 {
385  char uplo = upper_lower[uplot];
386  char tr = no_yes[trans];
387  char dd = 'N';
388  if (diag) {
389  dd = diag[0];
390  }
391  integer f_n = static_cast<integer>(n);
392  integer f_nrhs = static_cast<integer>(nrhs);
393  integer f_lda = static_cast<integer>(lda);
394  integer f_ldb = static_cast<integer>(ldb);
395  integer f_info = info;
396 #ifdef NO_FTN_STRING_LEN_AT_END
397  _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
398 #else
399  ftnlen trsize = 1;
400 #ifdef LAPACK_FTN_STRING_LEN_AT_END
401  _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize, trsize, trsize);
402 #else
403  _DTRTRS_(&uplo, trsize, &tr, trsize, &dd, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
404 #endif
405 #endif
406  info = f_info;
407 }
408 //====================================================================================================================
409 //!
410 /*!
411  * @param work Must be dimensioned equal to greater than 3N
412  * @param iwork Must be dimensioned equal to or greater than N
413  */
414 inline doublereal ct_dtrcon(const char* norm, ctlapack::upperlower_t uplot, const char* diag,
415  size_t n, doublereal* a, size_t lda, doublereal* work, int* iwork, int& info)
416 {
417  char uplo = upper_lower[uplot];
418  char dd = 'N';
419  if (diag) {
420  dd = diag[0];
421  }
422  char nn = '1';
423  if (norm) {
424  nn = norm[0];
425  }
426  integer f_n = static_cast<integer>(n);
427  integer f_lda = static_cast<integer>(lda);
428  integer f_info = info;
429  doublereal rcond;
430 #ifdef NO_FTN_STRING_LEN_AT_END
431  _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
432 #else
433  ftnlen trsize = 1;
434 #ifdef LAPACK_FTN_STRING_LEN_AT_END
435  _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info, trsize, trsize, trsize);
436 #else
437  _DTRCON_(&nn, trsize, &uplo, trsize, &dd, trsize, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
438 #endif
439 #endif
440  info = f_info;
441  return rcond;
442 }
443 //====================================================================================================================
444 
445 inline void ct_dpotrf(ctlapack::upperlower_t uplot, size_t n, doublereal* a, size_t lda, int& info)
446 {
447  char uplo = upper_lower[uplot];
448  integer f_n = static_cast<integer>(n);
449  integer f_lda = static_cast<integer>(lda);
450  integer f_info = info;
451 
452 #ifdef NO_FTN_STRING_LEN_AT_END
453  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info);
454 #else
455  ftnlen trsize = 1;
456 #ifdef LAPACK_FTN_STRING_LEN_AT_END
457  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info, trsize);
458 #else
459  _DPOTRF_(&uplo, trsize, &f_n, a, &f_lda, &f_info);
460 #endif
461 #endif
462  info = f_info;
463  return;
464 }
465 //====================================================================================================================
466 //!
467 /*!
468  */
469 inline void ct_dpotrs(ctlapack::upperlower_t uplot, size_t n, size_t nrhs, doublereal* a, size_t lda,
470  doublereal* b, size_t ldb, int& info)
471 {
472  char uplo = upper_lower[uplot];
473  integer f_n = static_cast<integer>(n);
474  integer f_nrhs = static_cast<integer>(nrhs);
475  integer f_lda = static_cast<integer>(lda);
476  integer f_ldb = static_cast<integer>(ldb);
477  integer f_info = info;
478 
479 #ifdef NO_FTN_STRING_LEN_AT_END
480  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
481 #else
482  ftnlen trsize = 1;
483 #ifdef LAPACK_FTN_STRING_LEN_AT_END
484  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize);
485 #else
486  _DPOTRS_(&uplo, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
487 #endif
488 #endif
489  info = f_info;
490  return;
491 }
492 
493 //====================================================================================================================
494 //!
495 /*!
496  */
497 inline doublereal ct_dgecon(const char norm, size_t n, doublereal* a, size_t lda, doublereal anorm,
498  doublereal* work, int* iwork, int& info)
499 {
500  char cnorm = '1';
501  if (norm) {
502  cnorm = norm;
503  }
504  integer f_n = static_cast<integer>(n);
505  integer f_lda = static_cast<integer>(lda);
506  integer f_info = info;
507  doublereal rcond;
508 
509 #ifdef NO_FTN_STRING_LEN_AT_END
510  _DGECON_(&cnorm, &f_n a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
511 #else
512  ftnlen trsize = 1;
513 #ifdef LAPACK_FTN_STRING_LEN_AT_END
514  _DGECON_(&cnorm, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info, trsize);
515 #else
516  _DGECON_(&cnorm, trsize, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
517 #endif
518 #endif
519  info = f_info;
520  return rcond;
521 }
522 
523 //====================================================================================================================
524 //!
525 /*!
526  */
527 inline doublereal ct_dgbcon(const char norm, size_t n, size_t kl, size_t ku,
528  doublereal* a, size_t ldab, int* ipiv, doublereal anorm,
529  doublereal* work, int* iwork, int& info)
530 {
531  char cnorm = '1';
532  if (norm) {
533  cnorm = norm;
534  }
535  integer f_n = static_cast<integer>(n);
536  integer f_kl = static_cast<integer>(kl);
537  integer f_ku = static_cast<integer>(ku);
538  integer f_ldab = static_cast<integer>(ldab);
539  integer f_info = info;
540  doublereal rcond;
541 
542 #ifdef NO_FTN_STRING_LEN_AT_END
543  _DGBCON_(&cnorm, &f_n , &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
544 #else
545  ftnlen trsize = 1;
546 #ifdef LAPACK_FTN_STRING_LEN_AT_END
547  _DGBCON_(&cnorm, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info, trsize);
548 #else
549  _DGBCON_(&cnorm, trsize, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
550 #endif
551 #endif
552  info = f_info;
553  return rcond;
554 }
555 
556 //====================================================================================================================
557 //!
558 /*!
559  */
560 inline doublereal ct_dlange(const char norm, size_t m, size_t n, doublereal* a, size_t lda,
561  doublereal* work)
562 {
563  char cnorm = '1';
564  if (norm) {
565  cnorm = norm;
566  }
567  integer f_m = static_cast<integer>(m);
568  integer f_n = static_cast<integer>(n);
569  integer f_lda = static_cast<integer>(lda);
570  doublereal anorm;
571 
572 #ifdef NO_FTN_STRING_LEN_AT_END
573  anorm = _DLANGE_(&cnorm, &f_m, &f_n a, &f_lda, work);
574 #else
575  ftnlen trsize = 1;
576 #ifdef LAPACK_FTN_STRING_LEN_AT_END
577  anorm = _DLANGE_(&cnorm, &f_m, &f_n, a, &f_lda, work, trsize);
578 #else
579  anorm = _DLANGE_(&cnorm, trsize, &f_m, &f_n, a, &f_lda, work);
580 #endif
581 #endif
582  return anorm;
583 }
584 //====================================================================================================================
585 }
586 
587 #endif