MADNESS  version 0.9
cblas.h
Go to the documentation of this file.
1 /*
2  This file is part of MADNESS.
3 
4  Copyright (C) 2007,2010 Oak Ridge National Laboratory
5 
6  This program is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2 of the License, or
9  (at your option) any later version.
10 
11  This program is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with this program; if not, write to the Free Software
18  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20  For more information please contact:
21 
22  Robert J. Harrison
23  Oak Ridge National Laboratory
24  One Bethel Valley Road
25  P.O. Box 2008, MS-6367
26 
27  email: harrisonrj@ornl.gov
28  tel: 865-241-3937
29  fax: 865-572-0680
30 
31 
32  $Id$
33 */
34 
35 
36 #ifndef MADNESS_LINALG_CBLAS_H__INCLUDED
37 #define MADNESS_LINALG_CBLAS_H__INCLUDED
38 
41 
42 
43 #include <madness/fortran_ctypes.h>
44 #include <madness/madness_config.h>
45 #include <madness/world/worldexc.h>
46 
47 #if defined(FORTRAN_LINKAGE_LC)
48 
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
71 
72 #elif defined(FORTRAN_LINKAGE_LCU)
73 
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_
96 
97 #elif defined(FORTRAN_LINKAGE_LCUU)
98 
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__
121 
122 #elif defined(FORTRAN_LINKAGE_UC)
123 
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
146 
147 #elif defined(FORTRAN_LINKAGE_UCU)
148 
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_
171 
172 #else
173 // If detected another convention complain loudly.
174 # error "cblas.h does not support the current Fortran symbol convention -- please, edit and check in the changes."
175 #endif
176 
177 extern "C" {
178 
179  // BLAS _GEMM declarations
180  void F77_SGEMM(const char*, const char*, const integer*, const integer*,
181  const integer*, const float*, const float*, const integer*,
182  const float*, const integer*, const float*, float*, const integer*);
183  void F77_DGEMM(const char*, const char*, const integer*, const integer*,
184  const integer*, const double*, const double*, const integer*,
185  const double*, const integer*, const double*, double*, const integer*);
186  void F77_CGEMM(const char*, const char*, const integer*, const integer*,
187  const integer*, const complex_real4*, const complex_real4*,
188  const integer*, const complex_real4*, const integer*,
189  const complex_real4*, complex_real4*, const integer*);
190  void F77_ZGEMM(const char*, const char*, const integer*, const integer*,
191  const integer*, const complex_real8*, const complex_real8*,
192  const integer*, const complex_real8*, const integer*,
193  const complex_real8*, complex_real8*, const integer*);
194 
195  // BLAS _GEMV declarations
196  void F77_SGEMV(const char*, const integer*, const integer*, const float*,
197  const float*, const integer*, const float*, const integer*,
198  const float*, float*, const integer*);
199  void F77_DGEMV(const char*, const integer*, const integer*, const double*,
200  const double*, const integer*, const double*, const integer*,
201  const double*, double*, const integer*);
202  void F77_CGEMV(const char*, const integer*, const integer*, const complex_real4*,
203  const complex_real4*, const integer*, const complex_real4*,
204  const integer*, const complex_real4*, complex_real4*, const integer*);
205  void F77_ZGEMV(const char*, const integer*, const integer*, const complex_real8*,
206  const complex_real8*, const integer*, const complex_real8*,
207  const integer*, const complex_real8*, complex_real8*, const integer*);
208 
209  // BLAS _GER declarations
210  void F77_SGER(const integer*, const integer*, const float*, const float*,
211  const integer*, const float*, const integer*, float*, const integer*);
212  void F77_DGER(const integer*, const integer*, const double*, const double*,
213  const integer*, const double*, const integer*, double*, const integer*);
214  void F77_CGER(const integer*, const integer*, const complex_real4*,
215  const complex_real4*, const integer*, const complex_real4*,
216  const integer*, complex_real4*, const integer*);
217  void F77_ZGER(const integer*, const integer*, const complex_real8*,
218  const complex_real8*, const integer*, const complex_real8*,
219  const integer*, complex_real8*, const integer*);
220 
221  // BLAS _SCAL declarations
222  void F77_SSCAL(const integer*, const float*, float*, const integer*);
223  void F77_DSCAL(const integer*, const double*, double*, const integer*);
224  void F77_CSCAL(const integer*, const complex_real4*, complex_real4*, const integer*);
225  void F77_CSSCAL(const integer*, const float*, complex_real4*, const integer*);
226  void F77_ZSCAL(const integer*, const complex_real8*, complex_real8*, const integer*);
227  void F77_ZDSCAL(const integer*, const double*, complex_real8*, const integer*);
228 
229  // BLAS _DOT declarations
230  float F77_SDOT(const integer*, const float*, const integer*, const float*,
231  const integer*);
232  double F77_DDOT(const integer*, const double *, const integer*,
233  const double *, const integer*);
234  void F77_CDOTU(const integer*, const complex_real4*, const integer*,
235  const complex_real4*, const integer*, complex_real4*);
236  void F77_ZDOTU(const integer*, const complex_real8*, const integer*,
237  const complex_real8*, const integer*, complex_real8*);
238 }
239 
240 
241 namespace madness {
242 namespace cblas {
243 
245  typedef enum {
247  Trans=1,
249  } CBLAS_TRANSPOSE;
250 
252 
270  inline void gemm(const CBLAS_TRANSPOSE OpA, const CBLAS_TRANSPOSE OpB,
271  const integer m, const integer n, const integer k, const float alpha,
272  const float* a, const integer lda, const float* b, const integer ldb,
273  const float beta, float* c, const integer ldc)
274  {
275  MADNESS_ASSERT(OpA != ConjTrans);
276  MADNESS_ASSERT(OpB != ConjTrans);
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);
279  }
280 
281  inline void gemm(const CBLAS_TRANSPOSE OpA, const CBLAS_TRANSPOSE OpB,
282  const integer m, const integer n, const integer k, const double alpha,
283  const double* a, const integer lda, const double* b, const integer ldb,
284  const double beta, double* c, const integer ldc) {
285  MADNESS_ASSERT(OpA != ConjTrans);
286  MADNESS_ASSERT(OpB != ConjTrans);
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);
289  }
290 
291  inline void gemm(const CBLAS_TRANSPOSE OpA, const CBLAS_TRANSPOSE OpB,
292  const integer m, const integer n, const integer k,
293  const complex_real4 alpha, const complex_real4* a, const integer lda,
294  const complex_real4* b, const integer ldb, const complex_real4 beta,
295  complex_real4* c, const integer ldc)
296  {
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);
299  }
300 
301  inline void gemm(const CBLAS_TRANSPOSE OpA, const CBLAS_TRANSPOSE OpB,
302  const integer m, const integer n, const integer k,
303  const complex_real8 alpha, const complex_real8* a, const integer lda,
304  const complex_real8* b, const integer ldb, const complex_real8 beta,
305  complex_real8* c, const integer 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);
308  }
310 
312 
328  inline void gemv(const CBLAS_TRANSPOSE OpA, const integer m, const integer n,
329  const float alpha, const float *A, const integer lda, const float *x,
330  const integer incx, const float beta, float *y, const integer incy)
331  {
332  MADNESS_ASSERT(OpA != ConjTrans);
333  static const char *op[] = { "n","t" };
334  F77_SGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
335  }
336 
337  inline void gemv(const CBLAS_TRANSPOSE OpA, const integer m, const integer n,
338  const double alpha, const double *A, const integer lda, const double *x,
339  const integer incx, const double beta, double *y, const integer incy)
340  {
341  MADNESS_ASSERT(OpA != ConjTrans);
342  static const char *op[] = { "n","t" };
343  F77_DGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
344  }
345 
346  inline void gemv(const CBLAS_TRANSPOSE OpA, const integer m, const integer n,
347  const complex_real4 alpha, const complex_real4 *A, const integer lda,
348  const complex_real4 *x, const integer incx, const complex_real4 beta,
349  complex_real4 *y, const integer incy)
350  {
351  static const char *op[] = { "n","t","c" };
352  F77_CGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
353  }
354 
355  inline void gemv(const CBLAS_TRANSPOSE OpA, const integer m, const integer n,
356  const complex_real8 alpha, const complex_real8 *A, const integer lda,
357  const complex_real8 *x, const integer incx, const complex_real8 beta,
358  complex_real8 *y, const integer incy)
359  {
360  static const char *op[] = { "n","t","c" };
361  F77_ZGEMV(op[OpA], &m, &n, &alpha, A, &lda, x, &incx, &beta, y, &incy);
362  }
364 
366 
380  inline void ger(const integer m, const integer n, const float alpha,
381  const float *x, const integer incx, const float *y, const integer incy,
382  float *A, const integer lda)
383  {
384  F77_SGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
385  }
386 
387  inline void ger(const integer m, const integer n, const double alpha,
388  const double *x, const integer incx, const double *y, const integer incy,
389  double *A, const integer lda)
390  {
391  F77_DGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
392  }
393 
394  inline void ger(const integer m, const integer n, const complex_real4 alpha,
395  const complex_real4 *x, const integer incx, const complex_real4 *y,
396  const integer incy, complex_real4 *A, const integer lda)
397  {
398  F77_CGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
399  }
400 
401  inline void ger(const integer m, const integer n, const complex_real8 alpha,
402  const complex_real8 *x, const integer incx, const complex_real8 *y,
403  const integer incy, complex_real8 *A, const integer lda)
404  {
405  F77_ZGER(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
406  }
408 
410 
421  inline float dot(const integer n, const float* x, const integer incx,
422  const float* y, const integer incy)
423  {
424  return F77_SDOT(&n, x, &incx, y, &incy);
425  }
426 
427  inline double dot(const integer n, const double* x, const integer incx,
428  const double* y, const integer incy)
429  {
430  return F77_DDOT(&n, x, &incx, y, &incy);
431  }
432 
433  inline complex_real4 dot(const integer n, const complex_real4* x,
434  const integer incx, const complex_real4* y, const integer incy)
435  {
436  complex_real4 result(0.0, 0.0);
437  F77_CDOTU(&n, x, &incx, y, &incy, &result);
438  return result;
439  }
440 
441  inline complex_real8 dot(const integer n, const complex_real8* x,
442  const integer incx, const complex_real8* y, const integer incy)
443  {
444  complex_real8 result(0.0, 0.0);
445  F77_ZDOTU(&n, x, &incx, y, &incy, &result);
446  return result;
447  }
449 
451 
460  inline void scal(const integer n, const float alpha, float* x, const integer incx) {
461  F77_SSCAL(&n, &alpha, x, &incx);
462  }
463 
464  inline void scal(const integer n, const double alpha, double* x, const integer incx) {
465  F77_DSCAL(&n, &alpha, x, &incx);
466  }
467 
468  inline void scal(const integer n, const complex_real4 alpha, complex_real4* x, const integer incx) {
469  F77_CSCAL(&n, &alpha, x, &incx);
470  }
471 
472  inline void scal(const integer n, const complex_real8 alpha, complex_real8* x, const integer incx) {
473  F77_ZSCAL(&n, &alpha, x, &incx);
474  }
475 
476  inline void scal(const integer n, const float alpha, complex_real4* x, const integer incx) {
477  F77_CSSCAL(&n, &alpha, x, &incx);
478  }
479 
480  inline void scal(const integer n, const double alpha, complex_real8* x, const integer incx) {
481  F77_ZDSCAL(&n, &alpha, x, &incx);
482  }
484 
485 } // namespace cblas
486 } // namespace madness
487 
488 #endif // MADNESS_LINALG_CBLAS_H__INCLUDED
489 
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.
Definition: cblas.h:248
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
Definition: cblas.h:246
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 *)
Definition: cblas.h:247
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