51 #if defined (INTEL_CXML)
59 #if defined (INTEL_CXML)
60 #define CHAR_MACRO(char_var) &char_var, one
62 #define CHAR_MACRO(char_var) &char_var
85 template<
typename Scalar>
86 Scalar generic_dot(
const int&
n,
const Scalar* x,
const int& incx,
87 const Scalar* y,
const int& incy)
91 if (incx==1 && incy==1) {
92 for (
int i = 0; i <
n; ++i)
93 dot += (*x++)*ST::conjugate(*y++);
100 for (
int i = 0; i <
n; ++i, x+=incx, y+=incy)
101 dot += (*x)*ST::conjugate(*y);
115 # ifdef HAVE_TEUCHOS_COMPLEX
116 template BLAS<long int, std::complex<float> >;
117 template BLAS<long int, std::complex<double> >;
119 template BLAS<long int, float>;
120 template BLAS<long int, double>;
128 void BLAS<int, float>::ROT(
const int& n,
float* dx,
const int& incx,
float* dy,
const int& incy,
float* c,
float* s)
const
129 {
SROT_F77(&n, dx, &incx, dy, &incy, c, s); }
134 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
135 return cblas_sasum(n, x, incx);
136 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
143 for (
int i = 0; i <
n; ++i)
144 sum += ST::magnitude(*x++);
147 for (
int i = 0; i <
n; ++i, x+=incx)
148 sum += ST::magnitude(*x);
154 void BLAS<int, float>::AXPY(
const int& n,
const float& alpha,
const float* x,
const int& incx,
float* y,
const int& incy)
const
155 {
SAXPY_F77(&n, &alpha, x, &incx, y, &incy); }
162 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
163 return cblas_sdot(n, x, incx, y, incy);
164 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
165 return SDOT_F77(&n, x, &incx, y, &incy);
167 return generic_dot(n, x, incx, y, incy);
176 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
177 return cblas_snrm2(n, x, incx);
178 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
188 void BLAS<int, float>::GEMV(
ETransp trans,
const int& m,
const int& n,
const float& alpha,
const float*
A,
const int& lda,
const float* x,
const int& incx,
const float& beta,
float* y,
const int& incy)
const
189 {
SGEMV_F77(
CHAR_MACRO(
ETranspChar[trans]), &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy); }
191 void BLAS<int, float>::GER(
const int& m,
const int& n,
const float& alpha,
const float* x,
const int& incx,
const float* y,
const int& incy,
float*
A,
const int& lda)
const
192 {
SGER_F77(&m, &n, &alpha, x, &incx, y, &incy, A, &lda); }
197 void BLAS<int, float>::GEMM(
ETransp transa,
ETransp transb,
const int& m,
const int& n,
const int& k,
const float& alpha,
const float*
A,
const int& lda,
const float*
B,
const int& ldb,
const float& beta,
float*
C,
const int& ldc)
const
198 {
SGEMM_F77(
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
ETranspChar[transb]), &m, &n, &k, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
205 void BLAS<int, float>::SYMM(
ESide side,
EUplo uplo,
const int& m,
const int& n,
const float& alpha,
const float*
A,
const int& lda,
const float*
B,
const int& ldb,
const float& beta,
float*
C,
const int& ldc)
const
206 {
SSYMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]), &m, &n, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
208 void BLAS<int, float>::SYRK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const float& alpha,
const float*
A,
const int& lda,
const float& beta,
float*
C,
const int& ldc)
const
211 void BLAS<int, float>::HERK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const float& alpha,
const float*
A,
const int& lda,
const float& beta,
float*
C,
const int& ldc)
const
214 void BLAS<int, float>::TRMM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
const int& m,
const int& n,
const float& alpha,
const float*
A,
const int& lda,
float*
B,
const int& ldb)
const
215 {
STRMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
217 void BLAS<int, float>::TRSM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
const int& m,
const int& n,
const float& alpha,
const float*
A,
const int& lda,
float*
B,
const int& ldb)
const
218 {
STRSM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
225 void BLAS<int, double>::ROT(
const int& n,
double* dx,
const int& incx,
double* dy,
const int& incy,
double* c,
double* s)
const
226 {
DROT_F77(&n, dx, &incx, dy, &incy, c, s); }
231 void BLAS<int, double>::AXPY(
const int& n,
const double& alpha,
const double* x,
const int& incx,
double* y,
const int& incy)
const
232 {
DAXPY_F77(&n, &alpha, x, &incx, y, &incy); }
239 return DDOT_F77(&n, x, &incx, y, &incy);
251 void BLAS<int, double>::GEMV(
ETransp trans,
const int& m,
const int& n,
const double& alpha,
const double*
A,
const int& lda,
const double* x,
const int& incx,
const double& beta,
double* y,
const int& incy)
const
252 {
DGEMV_F77(
CHAR_MACRO(
ETranspChar[trans]), &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy); }
254 void BLAS<int, double>::GER(
const int& m,
const int& n,
const double& alpha,
const double* x,
const int& incx,
const double* y,
const int& incy,
double*
A,
const int& lda)
const
255 {
DGER_F77(&m, &n, &alpha, x, &incx, y, &incy, A, &lda); }
260 void BLAS<int, double>::GEMM(
ETransp transa,
ETransp transb,
const int& m,
const int& n,
const int& k,
const double& alpha,
const double*
A,
const int& lda,
const double*
B,
const int& ldb,
const double& beta,
double*
C,
const int& ldc)
const
261 {
DGEMM_F77(
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
ETranspChar[transb]), &m, &n, &k, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
268 void BLAS<int, double>::SYMM(
ESide side,
EUplo uplo,
const int& m,
const int& n,
const double& alpha,
const double*
A,
const int& lda,
const double*
B,
const int& ldb,
const double& beta,
double*
C,
const int& ldc)
const
269 {
DSYMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]), &m, &n, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
271 void BLAS<int, double>::SYRK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const double& alpha,
const double*
A,
const int& lda,
const double& beta,
double*
C,
const int& ldc)
const
274 void BLAS<int, double>::HERK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const double& alpha,
const double*
A,
const int& lda,
const double& beta,
double*
C,
const int& ldc)
const
277 void BLAS<int, double>::TRMM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
const int& m,
const int& n,
const double& alpha,
const double*
A,
const int& lda,
double*
B,
const int& ldb)
const
278 {
DTRMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
280 void BLAS<int, double>::TRSM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
const int& m,
const int& n,
const double& alpha,
const double*
A,
const int& lda,
double*
B,
const int& ldb)
const
281 {
DTRSM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
283 #ifdef HAVE_TEUCHOS_COMPLEX
288 { CROTG_F77(da, db, c, s ); }
290 void BLAS<int, std::complex<float> >::ROT(
const int& n, std::complex<float>* dx,
const int& incx, std::complex<float>* dy,
const int& incy,
float* c, std::complex<float>* s)
const
291 { CROT_F77(&n, dx, &incx, dy, &incy, c, s); }
293 float BLAS<int, std::complex<float> >::ASUM(
const int& n,
const std::complex<float>* x,
const int& incx)
const
295 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
296 return cblas_scasum(n, x, incx);
297 #elif defined(HAVE_TEUCHOS_BLASFLOAT_DOUBLE_RETURN)
298 return (
float) SCASUM_F77(&n, x, &incx);
299 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
300 return SCASUM_F77(&n, x, &incx);
301 #else // Wow, you just plain don't have this routine.
306 for (
int i = 0; i <
n; ++i) {
307 result += std::abs (std::real (x[i])) + std::abs (std::imag (x[i]));
310 const int nincx = n * incx;
311 for (
int i = 0; i < nincx; i += incx) {
312 result += std::abs (std::real (x[i])) + std::abs (std::imag (x[i]));
315 return static_cast<float> (result);
319 void BLAS<int, std::complex<float> >::AXPY(
const int& n,
const std::complex<float> alpha,
const std::complex<float>* x,
const int& incx, std::complex<float>* y,
const int& incy)
const
320 { CAXPY_F77(&n, &alpha, x, &incx, y, &incy); }
322 void BLAS<int, std::complex<float> >::COPY(
const int& n,
const std::complex<float>* x,
const int& incx, std::complex<float>* y,
const int& incy)
const
323 { CCOPY_F77(&n, x, &incx, y, &incy); }
325 std::complex<float> BLAS<int, std::complex<float> >::DOT(
const int& n,
const std::complex<float>* x,
const int& incx,
const std::complex<float>* y,
const int& incy)
const
327 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
328 std::complex<float> z;
329 cblas_cdotc_sub(n,x,incx,y,incy,&z);
331 #elif defined(HAVE_COMPLEX_BLAS_PROBLEM) && defined(HAVE_FIXABLE_COMPLEX_BLAS_PROBLEM)
332 std::complex<float> z;
333 CDOT_F77(&z, &n, x, &incx, y, &incy);
335 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
336 return CDOT_F77(&n, x, &incx, y, &incy);
337 #else // Wow, you just plain don't have this routine.
340 std::complex<double> result (0, 0);
342 if (incx == 1 && incy == 1) {
343 for (
int i = 0; i <
n; ++i) {
344 result += std::conj (x[i]) * y[i];
355 for (
int i = 0; i <
n; ++i) {
356 result += std::conj (x[ix]) * y[iy];
362 return static_cast<std::complex<float>
> (result);
366 int BLAS<int, std::complex<float> >::IAMAX(
const int& n,
const std::complex<float>* x,
const int& incx)
const
367 {
return ICAMAX_F77(&n, x, &incx); }
369 float BLAS<int, std::complex<float> >::NRM2(
const int& n,
const std::complex<float>* x,
const int& incx)
const
371 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
372 return cblas_scnrm2(n, x, incx);
373 #elif defined(HAVE_TEUCHOS_BLASFLOAT_DOUBLE_RETURN)
374 return (
float) SCNRM2_F77(&n, x, &incx);
375 #elif defined(HAVE_TEUCHOS_BLASFLOAT)
376 return SCNRM2_F77(&n, x, &incx);
377 #else // Wow, you just plain don't have this routine.
380 if (n < 1 || incx < 1) {
386 const int upper = 1 + (n-1)*incx;
387 for (
int ix = 0; ix < upper; ix += incx) {
392 if (std::real (x[ix]) != 0) {
393 const double temp = std::abs (std::real (x[ix]));
395 const double scale_over_temp = scale / temp;
396 ssq = 1 + ssq * scale_over_temp*scale_over_temp;
400 const double temp_over_scale = temp / scale;
401 ssq = ssq + temp_over_scale*temp_over_scale;
404 if (std::imag (x[ix]) != 0) {
405 const double temp = std::abs (std::imag (x[ix]));
407 const double scale_over_temp = scale / temp;
408 ssq = 1 + ssq * scale_over_temp*scale_over_temp;
412 const double temp_over_scale = temp / scale;
413 ssq = ssq + temp_over_scale*temp_over_scale;
417 return static_cast<float> (scale * std::sqrt (ssq));
422 void BLAS<int, std::complex<float> >::SCAL(
const int& n,
const std::complex<float> alpha, std::complex<float>* x,
const int& incx)
const
423 { CSCAL_F77(&n, &alpha, x, &incx); }
425 void BLAS<int, std::complex<float> >::GEMV(
ETransp trans,
const int& m,
const int& n,
const std::complex<float> alpha,
const std::complex<float>*
A,
const int& lda,
const std::complex<float>* x,
const int& incx,
const std::complex<float> beta, std::complex<float>* y,
const int& incy)
const
426 { CGEMV_F77(
CHAR_MACRO(
ETranspChar[trans]), &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy); }
428 void BLAS<int, std::complex<float> >::GER(
const int& m,
const int& n,
const std::complex<float> alpha,
const std::complex<float>* x,
const int& incx,
const std::complex<float>* y,
const int& incy, std::complex<float>* A,
const int& lda)
const
429 { CGER_F77(&m, &n, &alpha, x, &incx, y, &incy, A, &lda); }
431 void BLAS<int, std::complex<float> >::TRMV(
EUplo uplo,
ETransp trans,
EDiag diag,
const int& n,
const std::complex<float>* A,
const int& lda, std::complex<float>* x,
const int& incx)
const
434 void BLAS<int, std::complex<float> >::GEMM(
ETransp transa,
ETransp transb,
const int& m,
const int& n,
const int& k,
const std::complex<float> alpha,
const std::complex<float>* A,
const int& lda,
const std::complex<float>*
B,
const int& ldb,
const std::complex<float> beta, std::complex<float>*
C,
const int& ldc)
const
435 { CGEMM_F77(
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
ETranspChar[transb]), &m, &n, &k, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
437 void BLAS<int, std::complex<float> >::SWAP(
const int& n, std::complex<float>*
const x,
const int& incx, std::complex<float>*
const y,
const int& incy)
const
439 CSWAP_F77 (&n, x, &incx, y, &incy);
442 void BLAS<int, std::complex<float> >::SYMM(
ESide side,
EUplo uplo,
const int& m,
const int& n,
const std::complex<float> alpha,
const std::complex<float>* A,
const int& lda,
const std::complex<float>* B,
const int& ldb,
const std::complex<float> beta, std::complex<float>* C,
const int& ldc)
const
443 { CSYMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]), &m, &n, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
445 void BLAS<int, std::complex<float> >::SYRK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const std::complex<float> alpha,
const std::complex<float>* A,
const int& lda,
const std::complex<float> beta, std::complex<float>* C,
const int& ldc)
const
448 void BLAS<int, std::complex<float> >::HERK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const std::complex<float> alpha,
const std::complex<float>* A,
const int& lda,
const std::complex<float> beta, std::complex<float>* C,
const int& ldc)
const
451 void BLAS<int, std::complex<float> >::TRMM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
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)
const
452 { CTRMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
454 void BLAS<int, std::complex<float> >::TRSM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
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)
const
455 { CTRSM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
459 void BLAS<int, std::complex<double> >::ROTG(std::complex<double>* da, std::complex<double>* db,
double* c, std::complex<double>* s)
const
460 { ZROTG_F77(da, db, c, s); }
462 void BLAS<int, std::complex<double> >::ROT(
const int& n, std::complex<double>* dx,
const int& incx, std::complex<double>* dy,
const int& incy,
double* c, std::complex<double>* s)
const
463 { ZROT_F77(&n, dx, &incx, dy, &incy, c, s); }
465 double BLAS<int, std::complex<double> >::ASUM(
const int& n,
const std::complex<double>* x,
const int& incx)
const
466 {
return ZASUM_F77(&n, x, &incx); }
468 void BLAS<int, std::complex<double> >::AXPY(
const int& n,
const std::complex<double> alpha,
const std::complex<double>* x,
const int& incx, std::complex<double>* y,
const int& incy)
const
469 { ZAXPY_F77(&n, &alpha, x, &incx, y, &incy); }
471 void BLAS<int, std::complex<double> >::COPY(
const int& n,
const std::complex<double>* x,
const int& incx, std::complex<double>* y,
const int& incy)
const
472 { ZCOPY_F77(&n, x, &incx, y, &incy); }
474 std::complex<double> BLAS<int, std::complex<double> >::DOT(
const int& n,
const std::complex<double>* x,
const int& incx,
const std::complex<double>* y,
const int& incy)
const
476 #if defined(HAVE_TEUCHOS_BLASFLOAT_APPLE_VECLIB_BUGFIX)
477 std::complex<double> z;
478 cblas_zdotc_sub(n,x,incx,y,incy,&z);
480 #elif defined(HAVE_COMPLEX_BLAS_PROBLEM)
481 # if defined(HAVE_FIXABLE_COMPLEX_BLAS_PROBLEM)
482 std::complex<double> z;
483 ZDOT_F77(&z, &n, x, &incx, y, &incy);
489 std::complex<double> ztemp (0, 0);
491 if (incx == 1 && incy == 1) {
492 for (
int i = 0; i <
n; ++i) {
493 ztemp += std::conj (x[i]) * y[i];
504 for (
int i = 0; i <
n; ++i) {
505 ztemp += std::conj (x[ix]) * y[iy];
513 # endif // defined(HAVE_FIXABLE_COMPLEX_BLAS_PROBLEM)
515 return ZDOT_F77(&n, x, &incx, y, &incy);
519 int BLAS<int, std::complex<double> >::IAMAX(
const int& n,
const std::complex<double>* x,
const int& incx)
const
520 {
return IZAMAX_F77(&n, x, &incx); }
522 double BLAS<int, std::complex<double> >::NRM2(
const int& n,
const std::complex<double>* x,
const int& incx)
const
523 {
return ZNRM2_F77(&n, x, &incx); }
525 void BLAS<int, std::complex<double> >::SCAL(
const int& n,
const std::complex<double> alpha, std::complex<double>* x,
const int& incx)
const
526 { ZSCAL_F77(&n, &alpha, x, &incx); }
528 void BLAS<int, std::complex<double> >::GEMV(
ETransp trans,
const int& m,
const int& n,
const std::complex<double> alpha,
const std::complex<double>* A,
const int& lda,
const std::complex<double>* x,
const int& incx,
const std::complex<double> beta, std::complex<double>* y,
const int& incy)
const
529 { ZGEMV_F77(
CHAR_MACRO(
ETranspChar[trans]), &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy); }
531 void BLAS<int, std::complex<double> >::GER(
const int& m,
const int& n,
const std::complex<double> alpha,
const std::complex<double>* x,
const int& incx,
const std::complex<double>* y,
const int& incy, std::complex<double>* A,
const int& lda)
const
532 { ZGER_F77(&m, &n, &alpha, x, &incx, y, &incy, A, &lda); }
534 void BLAS<int, std::complex<double> >::TRMV(
EUplo uplo,
ETransp trans,
EDiag diag,
const int& n,
const std::complex<double>* A,
const int& lda, std::complex<double>* x,
const int& incx)
const
537 void BLAS<int, std::complex<double> >::GEMM(
ETransp transa,
ETransp transb,
const int& m,
const int& n,
const int& k,
const std::complex<double> alpha,
const std::complex<double>* A,
const int& lda,
const std::complex<double>* B,
const int& ldb,
const std::complex<double> beta, std::complex<double>* C,
const int& ldc)
const
538 { ZGEMM_F77(
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
ETranspChar[transb]), &m, &n, &k, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
540 void BLAS<int, std::complex<double> >::SWAP(
const int& n, std::complex<double>*
const x,
const int& incx, std::complex<double>*
const y,
const int& incy)
const
542 ZSWAP_F77 (&n, x, &incx, y, &incy);
545 void BLAS<int, std::complex<double> >::SYMM(
ESide side,
EUplo uplo,
const int& m,
const int& n,
const std::complex<double> alpha,
const std::complex<double>* A,
const int& lda,
const std::complex<double> *B,
const int& ldb,
const std::complex<double> beta, std::complex<double> *C,
const int& ldc)
const
546 { ZSYMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]), &m, &n, &alpha, A, &lda, B, &ldb, &beta, C, &ldc); }
548 void BLAS<int, std::complex<double> >::SYRK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const std::complex<double> alpha,
const std::complex<double>* A,
const int& lda,
const std::complex<double> beta, std::complex<double>* C,
const int& ldc)
const
551 void BLAS<int, std::complex<double> >::HERK(
EUplo uplo,
ETransp trans,
const int& n,
const int& k,
const std::complex<double> alpha,
const std::complex<double>* A,
const int& lda,
const std::complex<double> beta, std::complex<double>* C,
const int& ldc)
const
554 void BLAS<int, std::complex<double> >::TRMM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
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)
const
555 { ZTRMM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
557 void BLAS<int, std::complex<double> >::TRSM(
ESide side,
EUplo uplo,
ETransp transa,
EDiag diag,
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)
const
558 { ZTRSM_F77(
CHAR_MACRO(
ESideChar[side]),
CHAR_MACRO(
EUploChar[uplo]),
CHAR_MACRO(
ETranspChar[transa]),
CHAR_MACRO(
EDiagChar[diag]), &m, &n, &alpha, A, &lda, B, &ldb); }
560 #endif // HAVE_TEUCHOS_COMPLEX
void TRSM(ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
Solves the matrix equations: op(A)*X=alpha*B or X*op(A)=alpha*B where X and B are m by n matrices...
void GER(const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy, ScalarType *A, const OrdinalType &lda) const
Performs the rank 1 operation: A <- alpha*x*y'+A.
void AXPY(const OrdinalType &n, const alpha_type alpha, const x_type *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
Perform the operation: y <- y+alpha*x.
void PREFIX DTRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *n, const double *a, const int *lda, double *x, const int *incx)
void PREFIX DROTG_F77(double *da, double *db, double *c, double *s)
void TRMV(EUplo uplo, ETransp trans, EDiag diag, const OrdinalType &n, const A_type *A, const OrdinalType &lda, ScalarType *x, const OrdinalType &incx) const
Performs the matrix-vector operation: x <- A*x or x <- A'*x where A is a unit/non-unit n by n upper/low...
void PREFIX DGEMM_F77(Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const int *k, const double *alpha, const double *a, const int *lda, const double *b, const int *ldb, const double *beta, double *c, const int *ldc)
void PREFIX DGEMV_F77(Teuchos_fcd, const int *m, const int *n, const double *alpha, const double A[], const int *lda, const double x[], const int *incx, const double *beta, double y[], const int *incy)
int PREFIX ISAMAX_F77(const int *n, const float *x, const int *incx)
Templated interface class to BLAS routines.
void GEMV(ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const
Performs the matrix-vector operation: y <- alpha*A*x+beta*y or y <- alpha*A'*x+beta*y where A is a gene...
void ROT(const OrdinalType &n, ScalarType *dx, const OrdinalType &incx, ScalarType *dy, const OrdinalType &incy, MagnitudeType *c, ScalarType *s) const
Applies a Givens plane rotation.
void PREFIX DROT_F77(const int *n, double *dx, const int *incx, double *dy, const int *incy, double *c, double *s)
void PREFIX SGEMV_F77(Teuchos_fcd, const int *m, const int *n, const float *alpha, const float A[], const int *lda, const float x[], const int *incx, const float *beta, float y[], const int *incy)
void PREFIX SCOPY_F77(const int *n, const float *x, const int *incx, float *y, const int *incy)
static T squareroot(T x)
Returns a number of magnitudeType that is the square root of this scalar type x.
void PREFIX SROTG_F77(float *da, float *db, float *c, float *s)
void PREFIX DSCAL_F77(const int *n, const double *alpha, double *x, const int *incx)
void PREFIX DSWAP_F77(const int *const n, double *const x, const int *const incx, double *const y, const int *const incy)
void PREFIX DSYRK_F77(Teuchos_fcd, Teuchos_fcd, const int *n, const int *k, const double *alpha, const double *a, const int *lda, const double *beta, double *c, const int *ldc)
void PREFIX SSYMM_F77(Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const float *alpha, const float *a, const int *lda, const float *b, const int *ldb, const float *beta, float *c, const int *ldc)
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Compute the 2-norm of the vector x.
OrdinalType IAMAX(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Return the index of the element of x with the maximum magnitude.
void GEMM(ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
General matrix-matrix multiply.
double PREFIX DNRM2_F77(const int *n, const double x[], const int *incx)
The Templated BLAS wrappers.
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char EDiagChar[]
void PREFIX DAXPY_F77(const int *n, const double *alpha, const double x[], const int *incx, double y[], const int *incy)
This structure defines some basic traits for a scalar field type.
ScalarTraits< ScalarType >::magnitudeType ASUM(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Sum the absolute values of the entries of x.
void PREFIX SSWAP_F77(const int *const n, float *const x, const int *const incx, float *const y, const int *const incy)
void COPY(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
Copy the vector x to the vector y.
void PREFIX SGER_F77(const int *m, const int *n, const float *alpha, const float *x, const int *incx, const float *y, const int *incy, float *a, const int *lda)
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char EUploChar[]
void PREFIX DTRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const double *alpha, const double *a, const int *lda, double *b, const int *ldb)
double PREFIX DASUM_F77(const int *n, const double x[], const int *incx)
void PREFIX SSYRK_F77(Teuchos_fcd, Teuchos_fcd, const int *n, const int *k, const float *alpha, const float *a, const int *lda, const float *beta, float *c, const int *ldc)
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ETypeChar[]
ScalarType DOT(const OrdinalType &n, const x_type *x, const OrdinalType &incx, const y_type *y, const OrdinalType &incy) const
Form the dot product of the vectors x and y.
void PREFIX STRMV_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *n, const float *a, const int *lda, float *x, const int *incx)
void PREFIX SSCAL_F77(const int *n, const float *alpha, float *x, const int *incx)
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ESideChar[]
int PREFIX IDAMAX_F77(const int *n, const double *x, const int *incx)
void PREFIX SGEMM_F77(Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const int *k, const float *alpha, const float *a, const int *lda, const float *b, const int *ldb, const float *beta, float *c, const int *ldc)
#define CHAR_MACRO(char_var)
void SYMM(ESide side, EUplo uplo, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
Performs the matrix-matrix operation: C <- alpha*A*B+beta*C or C <- alpha*B*A+beta*C where A is an m ...
void TRMM(ESide side, EUplo uplo, ETransp transa, EDiag diag, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, ScalarType *B, const OrdinalType &ldb) const
Performs the matrix-matrix operation: B <- alpha*op(A)*B or B <- alpha*B*op(A) where op(A) is an unit...
void SCAL(const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const
Scale the vector x by the constant alpha.
void ROTG(ScalarType *da, ScalarType *db, rotg_c_type *c, ScalarType *s) const
Computes a Givens plane rotation.
void PREFIX SAXPY_F77(const int *n, const float *alpha, const float x[], const int *incx, float y[], const int *incy)
TEUCHOSNUMERICS_LIB_DLL_EXPORT const char ETranspChar[]
void SYRK(EUplo uplo, ETransp trans, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
Performs the symmetric rank k operation: C <- alpha*A*A'+beta*C or C <- alpha*A'*A+beta*C, where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case or k by n matrix in the second case.
void PREFIX STRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const float *alpha, const float *a, const int *lda, float *b, const int *ldb)
void PREFIX DTRMM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const double *alpha, const double *a, const int *lda, double *b, const int *ldb)
void PREFIX SROT_F77(const int *n, float *dx, const int *incx, float *dy, const int *incy, float *c, float *s)
void PREFIX STRSM_F77(Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const float *alpha, const float *a, const int *lda, float *b, const int *ldb)
void PREFIX DGER_F77(const int *m, const int *n, const double *alpha, const double *x, const int *incx, const double *y, const int *incy, double *a, const int *lda)
void PREFIX DCOPY_F77(const int *n, const double *x, const int *incx, double *y, const int *incy)
double PREFIX DDOT_F77(const int *n, const double x[], const int *incx, const double y[], const int *incy)
void SWAP(const OrdinalType &n, ScalarType *const x, const OrdinalType &incx, ScalarType *const y, const OrdinalType &incy) const
Swap the entries of x and y.
void PREFIX DSYMM_F77(Teuchos_fcd, Teuchos_fcd, const int *m, const int *n, const double *alpha, const double *a, const int *lda, const double *b, const int *ldb, const double *beta, double *c, const int *ldc)