10 #ifndef _TEUCHOS_BLAS_WRAPPERS_HPP_
11 #define _TEUCHOS_BLAS_WRAPPERS_HPP_
16 #pragma warning ( disable : 4190 )
26 #if defined(INTEL_CXML)
27 # define PREFIX __stdcall
28 # define Teuchos_fcd const char *, unsigned int
29 #elif defined(INTEL_MKL)
31 # define Teuchos_fcd const char *
34 # define Teuchos_fcd const char *
40 #if ( defined(_MSC_VER) )
41 #define Teuchos_Complex_double_type_name std::complex<double>
42 #define Teuchos_Complex_float_type_name std::complex<float>
44 #define Teuchos_Complex_double_type_name double _Complex
45 #define Teuchos_Complex_float_type_name float _Complex
52 #define DROTG_F77 F77_BLAS_MANGLE(drotg,DROTG)
53 #define DROT_F77 F77_BLAS_MANGLE(drot,DROT)
54 #define DASUM_F77 F77_BLAS_MANGLE(dasum,DASUM)
55 #define DAXPY_F77 F77_BLAS_MANGLE(daxpy,DAXPY)
56 #define DCOPY_F77 F77_BLAS_MANGLE(dcopy,DCOPY)
57 #define DDOT_F77 F77_BLAS_MANGLE(ddot,DDOT)
58 #define DNRM2_F77 F77_BLAS_MANGLE(dnrm2,DNRM2)
59 #define DSCAL_F77 F77_BLAS_MANGLE(dscal,DSCAL)
60 #define IDAMAX_F77 F77_BLAS_MANGLE(idamax,IDAMAX)
61 #define DGEMV_F77 F77_BLAS_MANGLE(dgemv,DGEMV)
62 #define DGER_F77 F77_BLAS_MANGLE(dger,DGER)
63 #define DTRMV_F77 F77_BLAS_MANGLE(dtrmv,DTRMV)
64 #define DGEMM_F77 F77_BLAS_MANGLE(dgemm,DGEMM)
65 #define DSWAP_F77 F77_BLAS_MANGLE(dswap,DSWAP)
66 #define DSYMM_F77 F77_BLAS_MANGLE(dsymm,DSYMM)
67 #define DSYRK_F77 F77_BLAS_MANGLE(dsyrk,DSYRK)
68 #define DTRMM_F77 F77_BLAS_MANGLE(dtrmm,DTRMM)
69 #define DTRSM_F77 F77_BLAS_MANGLE(dtrsm,DTRSM)
71 #ifdef HAVE_TEUCHOS_COMPLEX
73 #define ZROTG_F77 F77_BLAS_MANGLE(zrotg,ZROTG)
74 #define ZROT_F77 F77_BLAS_MANGLE(zrot,ZROT)
75 #define ZASUM_F77 F77_BLAS_MANGLE(dzasum,DZASUM)
76 #define ZAXPY_F77 F77_BLAS_MANGLE(zaxpy,ZAXPY)
77 #define ZCOPY_F77 F77_BLAS_MANGLE(zcopy,ZCOPY)
78 #define ZDOT_F77 F77_BLAS_MANGLE(zdotc,ZDOTC)
79 #define ZNRM2_F77 F77_BLAS_MANGLE(dznrm2,DZNRM2)
80 #define ZSCAL_F77 F77_BLAS_MANGLE(zscal,ZSCAL)
81 #define IZAMAX_F77 F77_BLAS_MANGLE(izamax,IZAMAX)
82 #define ZGEMV_F77 F77_BLAS_MANGLE(zgemv,ZGEMV)
83 #define ZGER_F77 F77_BLAS_MANGLE(zgeru,ZGERU)
84 #define ZTRMV_F77 F77_BLAS_MANGLE(ztrmv,ZTRMV)
85 #define ZGEMM_F77 F77_BLAS_MANGLE(zgemm,ZGEMM)
86 #define ZSWAP_F77 F77_BLAS_MANGLE(zswap,ZSWAP)
87 #define ZSYMM_F77 F77_BLAS_MANGLE(zsymm,ZSYMM)
88 #define ZSYRK_F77 F77_BLAS_MANGLE(zsyrk,ZSYRK)
89 #define ZHERK_F77 F77_BLAS_MANGLE(zherk,ZHERK)
90 #define ZTRMM_F77 F77_BLAS_MANGLE(ztrmm,ZTRMM)
91 #define ZTRSM_F77 F77_BLAS_MANGLE(ztrsm,ZTRSM)
95 #define SROTG_F77 F77_BLAS_MANGLE(srotg,SROTG)
96 #define SROT_F77 F77_BLAS_MANGLE(srot,SROT)
97 #define SSCAL_F77 F77_BLAS_MANGLE(sscal,SSCAL)
98 #define SCOPY_F77 F77_BLAS_MANGLE(scopy,SCOPY)
99 #define SAXPY_F77 F77_BLAS_MANGLE(saxpy,SAXPY)
100 #define SDOT_F77 F77_BLAS_MANGLE(sdot,SDOT)
101 #define SNRM2_F77 F77_BLAS_MANGLE(snrm2,SNRM2)
102 #define SASUM_F77 F77_BLAS_MANGLE(sasum,SASUM)
103 #define ISAMAX_F77 F77_BLAS_MANGLE(isamax,ISAMAX)
104 #define SGEMV_F77 F77_BLAS_MANGLE(sgemv,SGEMV)
105 #define SGER_F77 F77_BLAS_MANGLE(sger,SGER)
106 #define STRMV_F77 F77_BLAS_MANGLE(strmv,STRMV)
107 #define SGEMM_F77 F77_BLAS_MANGLE(sgemm,SGEMM)
108 #define SSWAP_F77 F77_BLAS_MANGLE(sswap,SSWAP)
109 #define SSYMM_F77 F77_BLAS_MANGLE(ssymm,SSYMM)
110 #define SSYRK_F77 F77_BLAS_MANGLE(ssyrk,SSYRK)
111 #define STRMM_F77 F77_BLAS_MANGLE(strmm,STRMM)
112 #define STRSM_F77 F77_BLAS_MANGLE(strsm,STRSM)
114 #ifdef HAVE_TEUCHOS_COMPLEX
116 #define CROTG_F77 F77_BLAS_MANGLE(crotg,CROTG)
117 #define CROT_F77 F77_BLAS_MANGLE(crot,CROT)
118 #define SCASUM_F77 F77_BLAS_MANGLE(scasum,SCASUM)
119 #define CAXPY_F77 F77_BLAS_MANGLE(caxpy,CAXPY)
120 #define CCOPY_F77 F77_BLAS_MANGLE(ccopy,CCOPY)
121 #define CDOT_F77 F77_BLAS_MANGLE(cdotc,CDOTC)
122 #define SCNRM2_F77 F77_BLAS_MANGLE(scnrm2,SCNRM2)
123 #define CSCAL_F77 F77_BLAS_MANGLE(cscal,CSCAL)
124 #define ICAMAX_F77 F77_BLAS_MANGLE(icamax,ICAMAX)
125 #define CGEMV_F77 F77_BLAS_MANGLE(cgemv,CGEMV)
126 #define CGER_F77 F77_BLAS_MANGLE(cgeru,CGERU)
127 #define CTRMV_F77 F77_BLAS_MANGLE(ctrmv,CTRMV)
128 #define CGEMM_F77 F77_BLAS_MANGLE(cgemm,CGEMM)
129 #define CSWAP_F77 F77_BLAS_MANGLE(cswap,CSWAP)
130 #define CSYMM_F77 F77_BLAS_MANGLE(csymm,CSYMM)
131 #define CSYRK_F77 F77_BLAS_MANGLE(csyrk,CSYRK)
132 #define CHERK_F77 F77_BLAS_MANGLE(cherk,CHERK)
133 #define CTRMM_F77 F77_BLAS_MANGLE(ctrmm,CTRMM)
134 #define CTRSM_F77 F77_BLAS_MANGLE(ctrsm,CTRSM)
135 #define TEUCHOS_BLAS_CONVERT_COMPLEX_FORTRAN_TO_CXX(TYPE, Z) \
136 reinterpret_cast<std::complex<TYPE>&>(Z);
150 void PREFIX DROTG_F77(
double* da,
double* db,
double* c,
double* s);
151 void PREFIX DROT_F77(
const int* n,
double* dx,
const int* incx,
double* dy,
const int* incy,
double* c,
double* s);
152 double PREFIX DASUM_F77(
const int* n,
const double x[],
const int* incx);
153 void PREFIX DAXPY_F77(
const int* n,
const double* alpha,
const double x[],
const int* incx,
double y[],
const int* incy);
154 void PREFIX DCOPY_F77(
const int* n,
const double *x,
const int* incx,
double *y,
const int* incy);
155 double PREFIX DDOT_F77(
const int* n,
const double x[],
const int* incx,
const double y[],
const int* incy);
156 double PREFIX DNRM2_F77(
const int* n,
const double x[],
const int* incx);
157 void PREFIX DSCAL_F77(
const int* n,
const double* alpha,
double *x,
const int* incx);
158 void PREFIX DSWAP_F77(
const int*
const n,
double*
const x,
const int*
const incx,
159 double*
const y,
const int*
const incy);
160 int PREFIX IDAMAX_F77(
const int* n,
const double *x,
const int* incx);
163 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
165 # if defined(HAVE_COMPLEX_BLAS_PROBLEM)
166 # if defined(HAVE_FIXABLE_COMPLEX_BLAS_PROBLEM)
167 void PREFIX ZDOT_F77(std::complex<double> *ret,
const int* n,
const std::complex<double> x[],
const int* incx,
const std::complex<double> y[],
const int* incy);
168 # elif defined(HAVE_VECLIB_COMPLEX_BLAS)
170 # include <vecLib/cblas.h>
176 # endif // HAVE_COMPLEX_BLAS_PROBLEM
178 Teuchos_Complex_double_type_name PREFIX ZDOT_F77(
const int* n,
const std::complex<double> x[],
const int* incx,
const std::complex<double> y[],
const int* incy);
179 # endif // defined(HAVE_COMPLEX_BLAS_PROBLEM)
181 double PREFIX ZNRM2_F77(
const int* n,
const std::complex<double> x[],
const int* incx);
182 double PREFIX ZASUM_F77(
const int* n,
const std::complex<double> x[],
const int* incx);
183 void PREFIX ZROTG_F77(std::complex<double>* da, std::complex<double>* db,
double* c, std::complex<double>* s);
184 void PREFIX ZROT_F77(
const int* n, std::complex<double>* dx,
const int* incx, std::complex<double>* dy,
const int* incy,
double* c, std::complex<double>* s);
185 void PREFIX ZAXPY_F77(
const int* n,
const std::complex<double>* alpha,
const std::complex<double> x[],
const int* incx, std::complex<double> y[],
const int* incy);
186 void PREFIX ZCOPY_F77(
const int* n,
const std::complex<double> *x,
const int* incx, std::complex<double> *y,
const int* incy);
187 void PREFIX ZSCAL_F77(
const int* n,
const std::complex<double>* alpha, std::complex<double> *x,
const int* incx);
188 void PREFIX ZSWAP_F77(
const int*
const n, std::complex<double>*
const x,
const int*
const incx,
189 std::complex<double>*
const y,
const int*
const incy);
190 int PREFIX IZAMAX_F77(
const int* n,
const std::complex<double> *x,
const int* incx);
195 #ifdef HAVE_TEUCHOS_BLASFLOAT
196 # ifdef HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX
197 # include <vecLib/cblas.h>
198 # elif defined(HAVE_TEUCHOS_BLASFLOAT_DOUBLE_RETURN)
199 double PREFIX SASUM_F77(
const int* n,
const float x[],
const int* incx);
200 double PREFIX SDOT_F77(
const int* n,
const float x[],
const int* incx,
const float y[],
const int* incy);
201 double PREFIX SNRM2_F77(
const int* n,
const float x[],
const int* incx);
203 float PREFIX SASUM_F77(
const int* n,
const float x[],
const int* incx);
204 float PREFIX SDOT_F77(
const int* n,
const float x[],
const int* incx,
const float y[],
const int* incy);
205 float PREFIX SNRM2_F77(
const int* n,
const float x[],
const int* incx);
206 # endif // which blasfloat
207 #endif // ifdef blasfloat
208 void PREFIX SROTG_F77(
float* da,
float* db,
float* c,
float* s);
209 void PREFIX SROT_F77(
const int* n,
float* dx,
const int* incx,
float* dy,
const int* incy,
float* c,
float* s);
210 void PREFIX SAXPY_F77(
const int* n,
const float* alpha,
const float x[],
const int* incx,
float y[],
const int* incy);
211 void PREFIX SCOPY_F77(
const int* n,
const float *x,
const int* incx,
float *y,
const int* incy);
212 void PREFIX SSCAL_F77(
const int* n,
const float* alpha,
float *x,
const int* incx);
213 void PREFIX SSWAP_F77(
const int*
const n,
float*
const x,
const int*
const incx,
214 float*
const y,
const int*
const incy);
215 int PREFIX ISAMAX_F77(
const int* n,
const float *x,
const int* incx);
218 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
219 # if defined(HAVE_TEUCHOS_BLASFLOAT)
220 # if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
222 # include <vecLib/cblas.h>
223 # elif defined(HAVE_TEUCHOS_BLASFLOAT_DOUBLE_RETURN)
224 double PREFIX SCASUM_F77(
const int* n,
const std::complex<float> x[],
const int* incx);
225 double PREFIX SCNRM2_F77(
const int* n,
const std::complex<float> x[],
const int* incx);
227 float PREFIX SCASUM_F77(
const int* n,
const std::complex<float> x[],
const int* incx);
228 float PREFIX SCNRM2_F77(
const int* n,
const std::complex<float> x[],
const int* incx);
229 # endif // Whether or not we have the veclib bugfix
230 #endif // defined(HAVE_TEUCHOS_BLASFLOAT)
232 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
234 #include <vecLib/cblas.h>
235 #elif defined(HAVE_COMPLEX_BLAS_PROBLEM) && defined(HAVE_FIXABLE_COMPLEX_BLAS_PROBLEM)
236 void PREFIX CDOT_F77(std::complex<float> *ret,
const int* n,
const std::complex<float> x[],
const int* incx,
const std::complex<float> y[],
const int* incy);
237 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
238 Teuchos_Complex_float_type_name PREFIX CDOT_F77(
const int* n,
const std::complex<float> x[],
const int* incx,
const std::complex<float> y[],
const int* incy);
243 void PREFIX CROTG_F77(std::complex<float>* da, std::complex<float>* db,
float* c, std::complex<float>* s);
244 void PREFIX CROT_F77(
const int* n, std::complex<float>* dx,
const int* incx, std::complex<float>* dy,
const int* incy,
float* c, std::complex<float>* s);
245 void PREFIX CAXPY_F77(
const int* n,
const std::complex<float>* alpha,
const std::complex<float> x[],
const int* incx, std::complex<float> y[],
const int* incy);
246 void PREFIX CCOPY_F77(
const int* n,
const std::complex<float> *x,
const int* incx, std::complex<float> *y,
const int* incy);
247 void PREFIX CSCAL_F77(
const int* n,
const std::complex<float>* alpha, std::complex<float> *x,
const int* incx);
248 void PREFIX CSWAP_F77(
const int*
const n, std::complex<float>*
const x,
const int*
const incx,
249 std::complex<float>*
const y,
const int*
const incy);
250 int PREFIX ICAMAX_F77(
const int* n,
const std::complex<float> *x,
const int* incx);
255 void PREFIX DGEMV_F77(Teuchos_fcd,
const int* m,
const int* n,
const double* alpha,
const double A[],
const int* lda,
256 const double x[],
const int* incx,
const double* beta,
double y[],
const int* incy);
257 void PREFIX DTRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
const int *n,
258 const double *a,
const int *lda,
double *x,
const int *incx);
259 void PREFIX DGER_F77(
const int *m,
const int *n,
const double *alpha,
const double *x,
const int *incx,
const double *y,
260 const int *incy,
double *a,
const int *lda);
263 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
265 void PREFIX ZGEMV_F77(Teuchos_fcd,
const int* m,
const int* n,
const std::complex<double>* alpha,
const std::complex<double> A[],
const int* lda,
266 const std::complex<double> x[],
const int* incx,
const std::complex<double>* beta, std::complex<double> y[],
const int* incy);
267 void PREFIX ZTRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
const int *n,
268 const std::complex<double> *a,
const int *lda, std::complex<double> *x,
const int *incx);
269 void PREFIX ZGER_F77(
const int *m,
const int *n,
const std::complex<double> *alpha,
const std::complex<double> *x,
const int *incx,
const std::complex<double> *y,
270 const int *incy, std::complex<double> *a,
const int *lda);
275 void PREFIX SGEMV_F77(Teuchos_fcd,
const int* m,
const int* n,
const float* alpha,
const float A[],
const int* lda,
276 const float x[],
const int* incx,
const float* beta,
float y[],
const int* incy);
277 void PREFIX STRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
const int *n,
278 const float *a,
const int *lda,
float *x,
const int *incx);
279 void PREFIX SGER_F77(
const int *m,
const int *n,
const float *alpha,
const float *x,
const int *incx,
const float *y,
280 const int *incy,
float *a,
const int *lda);
283 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
285 void PREFIX CGEMV_F77(Teuchos_fcd,
const int* m,
const int* n,
const std::complex<float>* alpha,
const std::complex<float> A[],
const int* lda,
286 const std::complex<float> x[],
const int* incx,
const std::complex<float>* beta, std::complex<float> y[],
const int* incy);
287 void PREFIX CTRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
const int *n,
288 const std::complex<float> *a,
const int *lda, std::complex<float> *x,
const int *incx);
289 void PREFIX CGER_F77(
const int *m,
const int *n,
const std::complex<float> *alpha,
const std::complex<float> *x,
const int *incx,
const std::complex<float> *y,
290 const int *incy, std::complex<float> *a,
const int *lda);
295 void PREFIX DGEMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int *
296 n,
const int *k,
const double *alpha,
const double *a,
const int *lda,
297 const double *b,
const int *ldb,
const double *beta,
double *c,
const int *ldc);
298 void PREFIX DSYMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int * n,
299 const double *alpha,
const double *a,
const int *lda,
300 const double *b,
const int *ldb,
const double *beta,
double *c,
const int *ldc);
301 void PREFIX DSYRK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
302 const double *alpha,
const double *a,
const int *lda,
303 const double *beta,
double *c,
const int *ldc);
304 void PREFIX DTRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
305 const int *m,
const int *n,
const double *alpha,
const double *a,
const int * lda,
double *b,
const int *ldb);
306 void PREFIX DTRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
307 const int *m,
const int *n,
const double *alpha,
const double *a,
const int *
308 lda,
double *b,
const int *ldb);
311 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
313 void PREFIX ZGEMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int *
314 n,
const int *k,
const std::complex<double> *alpha,
const std::complex<double> *a,
const int *lda,
315 const std::complex<double> *b,
const int *ldb,
const std::complex<double> *beta, std::complex<double> *c,
const int *ldc);
316 void PREFIX ZSYMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int * n,
317 const std::complex<double> *alpha,
const std::complex<double> *a,
const int *lda,
318 const std::complex<double> *b,
const int *ldb,
const std::complex<double> *beta, std::complex<double> *c,
const int *ldc);
319 void PREFIX ZSYRK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
320 const std::complex<double> *alpha,
const std::complex<double> *a,
const int *lda,
321 const std::complex<double> *beta, std::complex<double> *c,
const int *ldc);
322 void PREFIX ZHERK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
323 const std::complex<double> *alpha,
const std::complex<double> *a,
const int *lda,
324 const std::complex<double> *beta, std::complex<double> *c,
const int *ldc);
325 void PREFIX ZTRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
326 const int *m,
const int *n,
const std::complex<double> *alpha,
const std::complex<double> *a,
const int * lda, std::complex<double> *b,
const int *ldb);
327 void PREFIX ZTRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
328 const int *m,
const int *n,
const std::complex<double> *alpha,
const std::complex<double> *a,
const int *
329 lda, std::complex<double> *b,
const int *ldb);
334 void PREFIX SGEMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int *
335 n,
const int *k,
const float *alpha,
const float *a,
const int *lda,
336 const float *b,
const int *ldb,
const float *beta,
float *c,
const int *ldc);
337 void PREFIX SSYMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int * n,
338 const float *alpha,
const float *a,
const int *lda,
339 const float *b,
const int *ldb,
const float *beta,
float *c,
const int *ldc);
340 void PREFIX SSYRK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
341 const float *alpha,
const float *a,
const int *lda,
342 const float *beta,
float *c,
const int *ldc);
343 void PREFIX STRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
344 const int *m,
const int *n,
const float *alpha,
const float *a,
const int * lda,
float *b,
const int *ldb);
345 void PREFIX STRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
346 const int *m,
const int *n,
const float *alpha,
const float *a,
const int *
347 lda,
float *b,
const int *ldb);
351 #if defined(HAVE_TEUCHOS_COMPLEX) && defined(__cplusplus)
353 void PREFIX CGEMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int *
354 n,
const int *k,
const std::complex<float> *alpha,
const std::complex<float> *a,
const int *lda,
355 const std::complex<float> *b,
const int *ldb,
const std::complex<float> *beta, std::complex<float> *c,
const int *ldc);
356 void PREFIX CSYMM_F77(Teuchos_fcd, Teuchos_fcd,
const int *m,
const int * n,
357 const std::complex<float> *alpha,
const std::complex<float> *a,
const int *lda,
358 const std::complex<float> *b,
const int *ldb,
const std::complex<float> *beta, std::complex<float> *c,
const int *ldc);
359 void PREFIX CTRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
360 const int *m,
const int *n,
const std::complex<float> *alpha,
const std::complex<float> *a,
const int * lda, std::complex<float> *b,
const int *ldb);
361 void PREFIX CSYRK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
362 const std::complex<float> *alpha,
const std::complex<float> *a,
const int *lda,
363 const std::complex<float> *beta, std::complex<float> *c,
const int *ldc);
364 void PREFIX CHERK_F77(Teuchos_fcd, Teuchos_fcd,
const int *n,
const int * k,
365 const std::complex<float> *alpha,
const std::complex<float> *a,
const int *lda,
366 const std::complex<float> *beta, std::complex<float> *c,
const int *ldc);
367 void PREFIX CTRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd,
368 const int *m,
const int *n,
const std::complex<float> *alpha,
const std::complex<float> *a,
const int *
369 lda, std::complex<float> *b,
const int *ldb);
Teuchos header file which uses auto-configuration information to include necessary C++ headers...