36 #ifndef MADNESS_LINALG_CBLAS_H__INCLUDED 
   37 #define MADNESS_LINALG_CBLAS_H__INCLUDED 
   47 #if defined(FORTRAN_LINKAGE_LC) 
   49 #   define F77_SGEMM sgemm 
   50 #   define F77_DGEMM dgemm 
   51 #   define F77_CGEMM cgemm 
   52 #   define F77_ZGEMM zgemm 
   53 #   define F77_SGEMV sgemv 
   54 #   define F77_DGEMV dgemv 
   55 #   define F77_CGEMV cgemv 
   56 #   define F77_ZGEMV zgemv 
   57 #   define F77_SGER sger 
   58 #   define F77_DGER dger 
   59 #   define F77_CGER cger 
   60 #   define F77_ZGER zger 
   61 #   define F77_SSCAL sscal 
   62 #   define F77_DSCAL dscal 
   63 #   define F77_CSCAL cscal 
   64 #   define F77_ZSCAL zscal 
   65 #   define F77_CSSCAL csscal 
   66 #   define F77_ZDSCAL zdscal 
   67 #   define F77_SDOT sdot 
   68 #   define F77_DDOT ddot 
   69 #   define F77_CDOTU cdotu 
   70 #   define F77_ZDOTU zdotu 
   72 #elif defined(FORTRAN_LINKAGE_LCU) 
   74 #   define F77_SGEMM sgemm_ 
   75 #   define F77_DGEMM dgemm_ 
   76 #   define F77_CGEMM cgemm_ 
   77 #   define F77_ZGEMM zgemm_ 
   78 #   define F77_SGEMV sgemv_ 
   79 #   define F77_DGEMV dgemv_ 
   80 #   define F77_CGEMV cgemv_ 
   81 #   define F77_ZGEMV zgemv_ 
   82 #   define F77_SGER sger_ 
   83 #   define F77_DGER dger_ 
   84 #   define F77_CGER cger_ 
   85 #   define F77_ZGER zger_ 
   86 #   define F77_SSCAL sscal_ 
   87 #   define F77_DSCAL dscal_ 
   88 #   define F77_CSCAL cscal_ 
   89 #   define F77_ZSCAL zscal_ 
   90 #   define F77_CSSCAL csscal_ 
   91 #   define F77_ZDSCAL zdscal_ 
   92 #   define F77_SDOT sdot_ 
   93 #   define F77_DDOT ddot_ 
   94 #   define F77_CDOTU cdotu_ 
   95 #   define F77_ZDOTU zdotu_ 
   97 #elif defined(FORTRAN_LINKAGE_LCUU) 
   99 #   define F77_SGEMM  sgemm__ 
  100 #   define F77_DGEMM  dgemm__ 
  101 #   define F77_CGEMM  cgemm__ 
  102 #   define F77_ZGEMM  zgemm__ 
  103 #   define F77_SGEMV  sgemv__ 
  104 #   define F77_DGEMV  dgemv__ 
  105 #   define F77_CGEMV  cgemv__ 
  106 #   define F77_ZGEMV  zgemv__ 
  107 #   define F77_SGER   sger__ 
  108 #   define F77_DGER   dger__ 
  109 #   define F77_CGER   cger__ 
  110 #   define F77_ZGER   zger__ 
  111 #   define F77_SSCAL  sscal__ 
  112 #   define F77_DSCAL  dscal__ 
  113 #   define F77_CSCAL  cscal__ 
  114 #   define F77_ZSCAL  zscal__ 
  115 #   define F77_CSSCAL csscal__ 
  116 #   define F77_ZDSCAL zdscal__ 
  117 #   define F77_SDOT   sdot__ 
  118 #   define F77_DDOT   ddot__ 
  119 #   define F77_CDOTU  cdotu__ 
  120 #   define F77_ZDOTU  zdotu__ 
  122 #elif defined(FORTRAN_LINKAGE_UC) 
  124 #   define F77_SGEMM  SGEMM 
  125 #   define F77_DGEMM  DGEMM 
  126 #   define F77_CGEMM  CGEMM 
  127 #   define F77_ZGEMM  ZGEMM 
  128 #   define F77_SGEMV  SGEMV 
  129 #   define F77_DGEMV  DGEMV 
  130 #   define F77_CGEMV  CGEMV 
  131 #   define F77_ZGEMV  ZGEMV 
  132 #   define F77_SGER   SGER 
  133 #   define F77_DGER   DGER 
  134 #   define F77_CGER   CGER 
  135 #   define F77_ZGER   ZGER 
  136 #   define F77_SSCAL  SSCAL 
  137 #   define F77_DSCAL  DSCAL 
  138 #   define F77_CSCAL  CSCAL 
  139 #   define F77_ZSCAL  ZSCAL 
  140 #   define F77_CSSCAL CSSCAL 
  141 #   define F77_ZDSCAL ZDSCAL 
  142 #   define F77_SDOT   SDOTU 
  143 #   define F77_DDOT   DDOTU 
  144 #   define F77_CDOTU  CDOTU 
  145 #   define F77_ZDOTU  ZDOTU 
  147 #elif defined(FORTRAN_LINKAGE_UCU) 
  149 #   define F77_SGEMM  SGEMM_ 
  150 #   define F77_DGEMM  DGEMM_ 
  151 #   define F77_CGEMM  CGEMM_ 
  152 #   define F77_ZGEMM  ZGEMM_ 
  153 #   define F77_SGEMV  SGEMV_ 
  154 #   define F77_DGEMV  DGEMV_ 
  155 #   define F77_CGEMV  CGEMV_ 
  156 #   define F77_ZGEMV  ZGEMV_ 
  157 #   define F77_SGER   SGER_ 
  158 #   define F77_DGER   DGER_ 
  159 #   define F77_CGER   CGER_ 
  160 #   define F77_ZGER   ZGER_ 
  161 #   define F77_SSCAL  SSCAL_ 
  162 #   define F77_DSCAL  DSCAL_ 
  163 #   define F77_CSCAL  CSCAL_ 
  164 #   define F77_ZSCAL  ZSCAL_ 
  165 #   define F77_CSSCAL CSSCAL_ 
  166 #   define F77_ZDSCAL ZDSCAL_ 
  167 #   define F77_SDOT   SDOT_ 
  168 #   define F77_DDOT   DDOTSUB_ 
  169 #   define F77_CDOTU  CDOTU_ 
  170 #   define F77_ZDOTU  ZDOTU_ 
  174 #   error "cblas.h does not support the current Fortran symbol convention -- please, edit and check in the changes." 
  182             const float*, 
const integer*, 
const float*, 
float*, 
const integer*);
 
  185             const double*, 
const integer*, 
const double*, 
double*, 
const integer*);
 
  198             const float*, 
float*, 
const integer*);
 
  201             const double*, 
double*, 
const integer*);
 
  233             const double *, 
const integer*);
 
  277         static const char *
op[] = { 
"n",
"t" };
 
  278         F77_SGEMM(op[OpA], op[OpB], &m, &n, &k, &alpha, a, &lda, b, &ldb, &beta, c, &ldc);
 
  287         static const char *
op[] = { 
"n",
"t" };
 
  288         F77_DGEMM(op[OpA], op[OpB], &m, &n, &k, &alpha, a, &lda, b, &ldb, &beta, c, &ldc);
 
  297       static const char *
op[] = { 
"n",
"t",
"c" };
 
  298       F77_CGEMM(op[OpA], op[OpB], &m, &n, &k, &alpha, a, &lda, b, &ldb, &beta, c, &ldc);
 
  306       static const char *
op[] = { 
"n",
"t",
"c" };
 
  307       F77_ZGEMM(op[OpA], op[OpB], &m, &n, &k, &alpha, a, &lda, b, &ldb, &beta, c, &ldc);
 
  329        const float alpha, 
const float *A, 
const integer lda, 
const float *x,
 
  333         static const char *
op[] = { 
"n",
"t" };
 
  334         F77_SGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
 
  338        const double alpha, 
const double *A, 
const integer lda, 
const double *x,
 
  342         static const char *
op[] = { 
"n",
"t" };
 
  343         F77_DGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
 
  351         static const char *
op[] = { 
"n",
"t",
"c" };
 
  352         F77_CGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
 
  360         static const char *
op[] = { 
"n",
"t",
"c" };
 
  361         F77_ZGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
 
  381         const float *x, 
const integer incx, 
const float *y, 
const integer incy,
 
  384         F77_SGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
 
  388         const double *x, 
const integer incx, 
const double *y, 
const integer incy,
 
  391         F77_DGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
 
  398         F77_CGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
 
  405         F77_ZGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
 
  422         const float* y, 
const integer incy)
 
  424         return F77_SDOT(&n, x, &incx, y, &incy);
 
  428         const double* y, 
const integer incy)
 
  430         return F77_DDOT(&n, x, &incx, y, &incy);
 
  437         F77_CDOTU(&n, x, &incx, y, &incy, &result);
 
  445         F77_ZDOTU(&n, x, &incx, y, &incy, &result);
 
  488 #endif // MADNESS_LINALG_CBLAS_H__INCLUDED 
void ger(const integer m, const integer n, const float alpha, const float *x, const integer incx, const float *y, const integer incy, float *A, const integer lda)
Multiplies vector  by the transform of vector . 
Definition: cblas.h:380
void gemv(const CBLAS_TRANSPOSE OpA, const integer m, const integer n, const float alpha, const float *A, const integer lda, const float *x, const integer incx, const float beta, float *y, const integer incy)
Multiplies a matrix by a vector. 
Definition: cblas.h:328
Corresponding C and Fortran types. 
void F77_ZGER(const integer *, const integer *, const complex_real8 *, const complex_real8 *, const integer *, const complex_real8 *, const integer *, complex_real8 *, const integer *)
void F77_DSCAL(const integer *, const double *, double *, const integer *)
void F77_CGEMM(const char *, const char *, const integer *, const integer *, const integer *, const complex_real4 *, const complex_real4 *, const integer *, const complex_real4 *, const integer *, const complex_real4 *, complex_real4 *, const integer *)
void F77_ZGEMM(const char *, const char *, const integer *, const integer *, const integer *, const complex_real8 *, const complex_real8 *, const integer *, const complex_real8 *, const integer *, const complex_real8 *, complex_real8 *, const integer *)
void F77_CGEMV(const char *, const integer *, const integer *, const complex_real4 *, const complex_real4 *, const integer *, const complex_real4 *, const integer *, const complex_real4 *, complex_real4 *, const integer *)
void F77_CSSCAL(const integer *, const float *, complex_real4 *, const integer *)
void scal(const integer n, const float alpha, float *x, const integer incx)
Scale a vector. 
Definition: cblas.h:460
const double beta
Definition: gygi_soltion.cc:63
std::complex< float > complex_real4
Fortran single complex. 
Definition: fortran_ctypes.h:86
void F77_ZGEMV(const char *, const integer *, const integer *, const complex_real8 *, const complex_real8 *, const integer *, const complex_real8 *, const integer *, const complex_real8 *, complex_real8 *, const integer *)
void F77_CDOTU(const integer *, const complex_real4 *, const integer *, const complex_real4 *, const integer *, complex_real4 *)
void F77_DGEMM(const char *, const char *, const integer *, const integer *, const integer *, const double *, const double *, const integer *, const double *, const integer *, const double *, double *, const integer *)
void F77_SGEMM(const char *, const char *, const integer *, const integer *, const integer *, const float *, const float *, const integer *, const float *, const integer *, const float *, float *, const integer *)
void F77_SGER(const integer *, const integer *, const float *, const float *, const integer *, const float *, const integer *, float *, const integer *)
void F77_CGER(const integer *, const integer *, const complex_real4 *, const complex_real4 *, const integer *, const complex_real4 *, const integer *, complex_real4 *, const integer *)
void F77_ZSCAL(const integer *, const complex_real8 *, complex_real8 *, const integer *)
CBLAS_TRANSPOSE
Matrix operations for BLAS function calls. 
Definition: cblas.h:245
float dot(const integer n, const float *x, const integer incx, const float *y, const integer incy)
Compute the dot product of vectors  and . 
Definition: cblas.h:421
FLOAT a(int j, FLOAT z)
Definition: y1.cc:86
float F77_SDOT(const integer *, const float *, const integer *, const float *, const integer *)
std::complex< double > complex_real8
Fortran double complex. 
Definition: fortran_ctypes.h:82
void F77_SSCAL(const integer *, const float *, float *, const integer *)
int integer
Definition: DFcode/fci/crayio.c:25
void F77_ZDOTU(const integer *, const complex_real8 *, const integer *, const complex_real8 *, const integer *, complex_real8 *)
const double m
Definition: gfit.cc:199
double F77_DDOT(const integer *, const double *, const integer *, const double *, const integer *)
void F77_CSCAL(const integer *, const complex_real4 *, complex_real4 *, const integer *)
void F77_ZDSCAL(const integer *, const double *, complex_real8 *, const integer *)
void F77_DGEMV(const char *, const integer *, const integer *, const double *, const double *, const integer *, const double *, const integer *, const double *, double *, const integer *)
Tensor< double > op(const Tensor< double > &x)
Definition: kain.cc:508
Implements MadnessException. 
Holds machinery to set up Functions/FuncImpls using various Factories and Interfaces. 
Definition: chem/atomutil.cc:45
const double c
Definition: gfit.cc:200
void F77_SGEMV(const char *, const integer *, const integer *, const float *, const float *, const integer *, const float *, const integer *, const float *, float *, const integer *)
void F77_DGER(const integer *, const integer *, const double *, const double *, const integer *, const double *, const integer *, double *, const integer *)
FLOAT b(int j, FLOAT z)
Definition: y1.cc:79
void gemm(const CBLAS_TRANSPOSE OpA, const CBLAS_TRANSPOSE OpB, const integer m, const integer n, const integer k, const float alpha, const float *a, const integer lda, const float *b, const integer ldb, const float beta, float *c, const integer ldc)
Multiplies a matrix by a vector. 
Definition: cblas.h:270