Cantera  2.1.2
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 //====================================================================================================================
250 inline void ct_dgbsv(int n, int kl, int ku, int nrhs,
251  doublereal* a, int lda, integer* ipiv, doublereal* b, int ldb,
252  int& info)
253 {
254  integer f_n = n, f_kl = kl, f_ku = ku, f_nrhs = nrhs, f_lda = lda,
255  f_ldb = ldb, f_info = info;
256  _DGBSV_(&f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
257  b, &f_ldb, &f_info);
258  info = f_info;
259 }
260 
261 inline void ct_dgelss(size_t m, size_t n, size_t nrhs, doublereal* a,
262  size_t lda, doublereal* b, size_t ldb, doublereal* s,
263  doublereal rcond, size_t& rank, doublereal* work,
264  int& lwork, int& info)
265 {
266  integer f_m = static_cast<integer>(m);
267  integer f_n = static_cast<integer>(n);
268  integer f_nrhs = static_cast<integer>(nrhs);
269  integer f_lda = static_cast<integer>(lda);
270  integer f_ldb = static_cast<integer>(ldb);
271  integer f_lwork = static_cast<integer>(lwork);
272  integer f_rank, f_info;
273 
274  _DGELSS_(&f_m, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, s, &rcond, &f_rank,
275  work, &f_lwork, &f_info);
276 
277  info = static_cast<int>(f_info);
278  rank = static_cast<size_t>(f_rank);
279  lwork = static_cast<int>(f_lwork);
280 }
281 
282 inline void ct_dgbtrf(size_t m, size_t n, size_t kl, size_t ku,
283  doublereal* a, size_t lda, integer* ipiv, int& info)
284 {
285  integer f_m = (int) m;
286  integer f_n = (int) n;
287  integer f_kl = (int) kl;
288  integer f_ku = (int) ku;
289  integer f_lda = (int) lda;
290  integer f_info = info;
291  _DGBTRF_(&f_m, &f_n, &f_kl, &f_ku, a, &f_lda, ipiv, &f_info);
292  info = f_info;
293 }
294 
295 inline void ct_dgbtrs(ctlapack::transpose_t trans, size_t n,
296  size_t kl, size_t ku, size_t nrhs, doublereal* a, size_t lda,
297  integer* ipiv, doublereal* b, size_t ldb, int& info)
298 {
299  integer f_n = (int) n;
300  integer f_kl = (int) kl;
301  integer f_ku = (int) ku;
302  integer f_nrhs = (int) nrhs;
303  integer f_lda = (int) lda;
304  integer f_ldb = (int) ldb;
305  integer f_info = info;
306  char tr = no_yes[trans];
307 #ifdef NO_FTN_STRING_LEN_AT_END
308  _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
309  b, &f_ldb, &f_info);
310 #else
311  ftnlen trsize = 1;
312 #ifdef LAPACK_FTN_STRING_LEN_AT_END
313  _DGBTRS_(&tr, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
314  b, &f_ldb, &f_info, trsize);
315 #else
316  _DGBTRS_(&tr, trsize, &f_n, &f_kl, &f_ku, &f_nrhs, a, &f_lda, ipiv,
317  b, &f_ldb, &f_info);
318 #endif
319 #endif
320  info = f_info;
321 }
322 
323 inline void ct_dgetrf(size_t m, size_t n,
324  doublereal* a, size_t lda, integer* ipiv, int& info)
325 {
326  integer mm = (int) m;
327  integer nn = (int) n;
328  integer ldaa = (int) lda;
329  integer infoo = info;
330  _DGETRF_(&mm, &nn, a, &ldaa, ipiv, &infoo);
331  info = infoo;
332 }
333 
334 inline void ct_dgetrs(ctlapack::transpose_t trans, size_t n,
335  size_t nrhs, doublereal* a, size_t lda,
336  integer* ipiv, doublereal* b, size_t ldb, int& info)
337 {
338  integer f_n = (int) n;
339  integer f_lda = (int) lda;
340  integer f_nrhs = (int) nrhs;
341  integer f_ldb = (int) ldb;
342  integer f_info = info;
343  char tr = no_yes[trans];
344 
345 #ifdef NO_FTN_STRING_LEN_AT_END
346  _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb,
347  &f_info);
348 #else
349  ftnlen trsize = 1;
350 #ifdef LAPACK_FTN_STRING_LEN_AT_END
351  _DGETRS_(&tr, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info, trsize);
352 #else
353  _DGETRS_(&tr, trsize, &f_n, &f_nrhs, a, &f_lda, ipiv, b, &f_ldb, &f_info);
354 #endif
355 #endif
356  info = f_info;
357 }
358 //====================================================================================================================
359 inline void ct_dgetri(int n, doublereal* a, int lda, integer* ipiv,
360  doublereal* work, int lwork, int& info)
361 {
362  integer f_n = n, f_lda = lda, f_lwork = lwork, f_info = info;
363  _DGETRI_(&f_n, a, &f_lda, ipiv, work, &f_lwork, &f_info);
364 }
365 
366 inline void ct_dscal(int n, doublereal da, doublereal* dx, int incx)
367 {
368  integer f_n = n, f_incx = incx;
369  _DSCAL_(&f_n, &da, dx, &f_incx);
370 }
371 //====================================================================================================================
372 inline void ct_dgeqrf(size_t m, size_t n, doublereal* a, size_t lda, doublereal* tau,
373  doublereal* work, size_t lwork, int& info)
374 {
375  integer f_m = static_cast<integer>(m);
376  integer f_n = static_cast<integer>(n);
377  integer f_lda = static_cast<integer>(lda);
378  integer f_lwork = static_cast<integer>(lwork);
379  integer f_info = info;
380  _DGEQRF_(&f_m, &f_n, a, &f_lda, tau, work, &f_lwork, &f_info);
381  info = f_info;
382 }
383 //====================================================================================================================
384 inline void ct_dormqr(ctlapack::side_t rlside, ctlapack::transpose_t trans, size_t m,
385  size_t n, size_t k, doublereal* a, size_t lda, doublereal* tau, doublereal* c, size_t ldc,
386  doublereal* work, size_t lwork, int& info)
387 {
388  char side = left_right[rlside];
389  char tr = no_yes[trans];
390  integer f_m = static_cast<integer>(m);
391  integer f_n = static_cast<integer>(n);
392  integer f_k = static_cast<integer>(k);
393  integer f_lwork = static_cast<integer>(lwork);
394  integer f_lda = static_cast<integer>(lda);
395  integer f_ldc = static_cast<integer>(ldc);
396  integer f_info = info;
397 #ifdef NO_FTN_STRING_LEN_AT_END
398  _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
399 #else
400  ftnlen trsize = 1;
401 #ifdef LAPACK_FTN_STRING_LEN_AT_END
402  _DORMQR_(&side, &tr, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info, trsize, trsize);
403 #else
404  _DORMQR_(&side, trsize, &tr, trsize, &f_m, &f_n, &f_k, a, &f_lda, tau, c, &f_ldc, work, &f_lwork, &f_info);
405 #endif
406 #endif
407  info = f_info;
408 }
409 //====================================================================================================================
410 inline void ct_dtrtrs(ctlapack::upperlower_t uplot, ctlapack::transpose_t trans, const char* diag,
411  size_t n, size_t nrhs, doublereal* a, size_t lda, doublereal* b, size_t ldb, int& info)
412 {
413  char uplo = upper_lower[uplot];
414  char tr = no_yes[trans];
415  char dd = 'N';
416  if (diag) {
417  dd = diag[0];
418  }
419  integer f_n = static_cast<integer>(n);
420  integer f_nrhs = static_cast<integer>(nrhs);
421  integer f_lda = static_cast<integer>(lda);
422  integer f_ldb = static_cast<integer>(ldb);
423  integer f_info = info;
424 #ifdef NO_FTN_STRING_LEN_AT_END
425  _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
426 #else
427  ftnlen trsize = 1;
428 #ifdef LAPACK_FTN_STRING_LEN_AT_END
429  _DTRTRS_(&uplo, &tr, &dd, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize, trsize, trsize);
430 #else
431  _DTRTRS_(&uplo, trsize, &tr, trsize, &dd, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
432 #endif
433 #endif
434  info = f_info;
435 }
436 //====================================================================================================================
437 //!
438 /*!
439  * @param work Must be dimensioned equal to greater than 3N
440  * @param iwork Must be dimensioned equal to or greater than N
441  */
442 inline doublereal ct_dtrcon(const char* norm, ctlapack::upperlower_t uplot, const char* diag,
443  size_t n, doublereal* a, size_t lda, doublereal* work, int* iwork, int& info)
444 {
445  char uplo = upper_lower[uplot];
446  char dd = 'N';
447  if (diag) {
448  dd = diag[0];
449  }
450  char nn = '1';
451  if (norm) {
452  nn = norm[0];
453  }
454  integer f_n = static_cast<integer>(n);
455  integer f_lda = static_cast<integer>(lda);
456  integer f_info = info;
457  doublereal rcond;
458 #ifdef NO_FTN_STRING_LEN_AT_END
459  _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
460 #else
461  ftnlen trsize = 1;
462 #ifdef LAPACK_FTN_STRING_LEN_AT_END
463  _DTRCON_(&nn, &uplo, &dd, &f_n, a, &f_lda, &rcond, work, iwork, &f_info, trsize, trsize, trsize);
464 #else
465  _DTRCON_(&nn, trsize, &uplo, trsize, &dd, trsize, &f_n, a, &f_lda, &rcond, work, iwork, &f_info);
466 #endif
467 #endif
468  info = f_info;
469  return rcond;
470 }
471 //====================================================================================================================
472 
473 inline void ct_dpotrf(ctlapack::upperlower_t uplot, size_t n, doublereal* a, size_t lda, int& info)
474 {
475  char uplo = upper_lower[uplot];
476  integer f_n = static_cast<integer>(n);
477  integer f_lda = static_cast<integer>(lda);
478  integer f_info = info;
479 
480 #ifdef NO_FTN_STRING_LEN_AT_END
481  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info);
482 #else
483  ftnlen trsize = 1;
484 #ifdef LAPACK_FTN_STRING_LEN_AT_END
485  _DPOTRF_(&uplo, &f_n, a, &f_lda, &f_info, trsize);
486 #else
487  _DPOTRF_(&uplo, trsize, &f_n, a, &f_lda, &f_info);
488 #endif
489 #endif
490  info = f_info;
491  return;
492 }
493 //====================================================================================================================
494 //!
495 /*!
496  */
497 inline void ct_dpotrs(ctlapack::upperlower_t uplot, size_t n, size_t nrhs, doublereal* a, size_t lda,
498  doublereal* b, size_t ldb, int& info)
499 {
500  char uplo = upper_lower[uplot];
501  integer f_n = static_cast<integer>(n);
502  integer f_nrhs = static_cast<integer>(nrhs);
503  integer f_lda = static_cast<integer>(lda);
504  integer f_ldb = static_cast<integer>(ldb);
505  integer f_info = info;
506 
507 #ifdef NO_FTN_STRING_LEN_AT_END
508  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
509 #else
510  ftnlen trsize = 1;
511 #ifdef LAPACK_FTN_STRING_LEN_AT_END
512  _DPOTRS_(&uplo, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info, trsize);
513 #else
514  _DPOTRS_(&uplo, trsize, &f_n, &f_nrhs, a, &f_lda, b, &f_ldb, &f_info);
515 #endif
516 #endif
517  info = f_info;
518  return;
519 }
520 
521 //====================================================================================================================
522 //!
523 /*!
524  */
525 inline doublereal ct_dgecon(const char norm, size_t n, doublereal* a, size_t lda, doublereal anorm,
526  doublereal* work, int* iwork, int& info)
527 {
528  char cnorm = '1';
529  if (norm) {
530  cnorm = norm;
531  }
532  integer f_n = static_cast<integer>(n);
533  integer f_lda = static_cast<integer>(lda);
534  integer f_info = info;
535  doublereal rcond;
536 
537 #ifdef NO_FTN_STRING_LEN_AT_END
538  _DGECON_(&cnorm, &f_n a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
539 #else
540  ftnlen trsize = 1;
541 #ifdef LAPACK_FTN_STRING_LEN_AT_END
542  _DGECON_(&cnorm, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info, trsize);
543 #else
544  _DGECON_(&cnorm, trsize, &f_n, a, &f_lda, &anorm, &rcond, work, iwork, &f_info);
545 #endif
546 #endif
547  info = f_info;
548  return rcond;
549 }
550 
551 //====================================================================================================================
552 //!
553 /*!
554  */
555 inline doublereal ct_dgbcon(const char norm, size_t n, size_t kl, size_t ku,
556  doublereal* a, size_t ldab, int* ipiv, doublereal anorm,
557  doublereal* work, int* iwork, int& info)
558 {
559  char cnorm = '1';
560  if (norm) {
561  cnorm = norm;
562  }
563  integer f_n = static_cast<integer>(n);
564  integer f_kl = static_cast<integer>(kl);
565  integer f_ku = static_cast<integer>(ku);
566  integer f_ldab = static_cast<integer>(ldab);
567  integer f_info = info;
568  doublereal rcond;
569 
570 #ifdef NO_FTN_STRING_LEN_AT_END
571  _DGBCON_(&cnorm, &f_n , &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
572 #else
573  ftnlen trsize = 1;
574 #ifdef LAPACK_FTN_STRING_LEN_AT_END
575  _DGBCON_(&cnorm, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info, trsize);
576 #else
577  _DGBCON_(&cnorm, trsize, &f_n, &f_kl, &f_ku, a, &f_ldab, ipiv, &anorm, &rcond, work, iwork, &f_info);
578 #endif
579 #endif
580  info = f_info;
581  return rcond;
582 }
583 
584 //====================================================================================================================
585 //!
586 /*!
587  */
588 inline doublereal ct_dlange(const char norm, size_t m, size_t n, doublereal* a, size_t lda,
589  doublereal* work)
590 {
591  char cnorm = '1';
592  if (norm) {
593  cnorm = norm;
594  }
595  integer f_m = static_cast<integer>(m);
596  integer f_n = static_cast<integer>(n);
597  integer f_lda = static_cast<integer>(lda);
598  doublereal anorm;
599 
600 #ifdef NO_FTN_STRING_LEN_AT_END
601  anorm = _DLANGE_(&cnorm, &f_m, &f_n a, &f_lda, work);
602 #else
603  ftnlen trsize = 1;
604 #ifdef LAPACK_FTN_STRING_LEN_AT_END
605  anorm = _DLANGE_(&cnorm, &f_m, &f_n, a, &f_lda, work, trsize);
606 #else
607  anorm = _DLANGE_(&cnorm, trsize, &f_m, &f_n, a, &f_lda, work);
608 #endif
609 #endif
610  return anorm;
611 }
612 //====================================================================================================================
613 }
614 
615 #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:442