Nektar++
Blas.hpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Blas.hpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: wrapper of functions around standard BLAS routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_BLAS_HPP
36#define NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_BLAS_HPP
37
40
41// Translations for using Fortran version of blas
42namespace Blas
43{
44extern "C"
45{
46 // -- BLAS Level 1:
47 void F77NAME(dcopy)(const int &n, const double *x, const int &incx,
48 double *y, const int &incy);
49 void F77NAME(daxpy)(const int &n, const double &alpha, const double *x,
50 const int &incx, const double *y, const int &incy);
51 void F77NAME(dswap)(const int &n, double *x, const int &incx, double *y,
52 const int &incy);
53 void F77NAME(dscal)(const int &n, const double &alpha, double *x,
54 const int &incx);
55 void F77NAME(drot)(const int &n, double *x, const int &incx, double *y,
56 const int &incy, const double &c, const double &s);
57 double F77NAME(ddot)(const int &n, const double *x, const int &incx,
58 const double *y, const int &incy);
59 double F77NAME(dnrm2)(const int &n, const double *x, const int &incx);
60 double F77NAME(dasum)(const int &n, const double *x, const int &incx);
61 int F77NAME(idamax)(const int &n, const double *x, const int &incx);
62
63 // -- BLAS level 2
64 void F77NAME(dgemv)(const char &trans, const int &m, const int &n,
65 const double &alpha, const double *a, const int &lda,
66 const double *x, const int &incx, const double &beta,
67 double *y, const int &incy);
68 void F77NAME(sgemv)(const char &trans, const int &m, const int &n,
69 const float &alpha, const float *a, const int &lda,
70 const float *x, const int &incx, const float &beta,
71 float *y, const int &incy);
72
73 void F77NAME(dgbmv)(const char &trans, const int &m, const int &n,
74 const int &kl, const int &ku, const double &alpha,
75 const double *a, const int &lda, const double *x,
76 const int &incx, const double &beta, double *y,
77 const int &incy);
78 void F77NAME(sgbmv)(const char &trans, const int &m, const int &n,
79 const int &kl, const int &ku, const float &alpha,
80 const float *a, const int &lda, const float *x,
81 const int &incx, const float &beta, float *y,
82 const int &incy);
83
84 void F77NAME(dspmv)(const char &uplo, const int &n, const double &alpha,
85 const double *a, const double *x, const int &incx,
86 const double &beta, double *y, const int &incy);
87 void F77NAME(sspmv)(const char &uplo, const int &n, const float &alpha,
88 const float *a, const float *x, const int &incx,
89 const float &beta, float *y, const int &incy);
90
91 void F77NAME(dsbmv)(const char &uplo, const int &m, const int &k,
92 const double &alpha, const double *a, const int &lda,
93 const double *x, const int &incx, const double &beta,
94 double *y, const int &incy);
95 void F77NAME(ssbmv)(const char &uplo, const int &m, const int &k,
96 const float &alpha, const float *a, const int &lda,
97 const float *x, const int &incx, const float &beta,
98 float *y, const int &incy);
99
100 void F77NAME(dtpmv)(const char &uplo, const char &trans, const char &diag,
101 const int &n, const double *ap, double *x,
102 const int &incx);
103 void F77NAME(stpmv)(const char &uplo, const char &trans, const char &diag,
104 const int &n, const float *ap, float *x,
105 const int &incx);
106
107 void F77NAME(dger)(const int &m, const int &n, const double &alpha,
108 const double *x, const int &incx, const double *y,
109 const int &incy, double *a, const int &lda);
110 void F77NAME(sger)(const int &m, const int &n, const float &alpha,
111 const float *x, const int &incx, const float *y,
112 const int &incy, float *a, const int &lda);
113
114 // -- BLAS level 3:
115 void F77NAME(dgemm)(const char &trans, const char &transb, const int &m1,
116 const int &n, const int &k, const double &alpha,
117 const double *a, const int &lda, const double *b,
118 const int &ldb, const double &beta, double *c,
119 const int &ldc);
120 void F77NAME(sgemm)(const char &trans, const char &transb, const int &m1,
121 const int &n, const int &k, const float &alpha,
122 const float *a, const int &lda, const float *b,
123 const int &ldb, const float &beta, float *c,
124 const int &ldc);
125}
126
127/// \brief BLAS level 1: Copy \a x to \a y
128static inline void Dcopy(const int &n, const double *x, const int &incx,
129 double *y, const int &incy)
130{
131 F77NAME(dcopy)(n, x, incx, y, incy);
132}
133
134/// \brief BLAS level 1: y = alpha \a x plus \a y
135static inline void Daxpy(const int &n, const double &alpha, const double *x,
136 const int &incx, const double *y, const int &incy)
137{
138 F77NAME(daxpy)(n, alpha, x, incx, y, incy);
139}
140
141/// \brief BLAS level 1: Swap \a x with \a y
142static inline void Dswap(const int &n, double *x, const int &incx, double *y,
143 const int &incy)
144{
145 F77NAME(dswap)(n, x, incx, y, incy);
146}
147
148/// \brief BLAS level 1: x = alpha \a x
149static inline void Dscal(const int &n, const double &alpha, double *x,
150 const int &incx)
151{
152 F77NAME(dscal)(n, alpha, x, incx);
153}
154
155/// \brief BLAS level 1: Plane rotation by c = cos(theta), s = sin(theta)
156static inline void Drot(const int &n, double *x, const int &incx, double *y,
157 const int &incy, const double &c, const double &s)
158{
159 F77NAME(drot)(n, x, incx, y, incy, c, s);
160}
161
162/// \brief BLAS level 1: output = \f$ x^T y \f$
163static inline double Ddot(const int &n, const double *x, const int &incx,
164 const double *y, const int &incy)
165{
166 return F77NAME(ddot)(n, x, incx, y, incy);
167}
168
169// \brief BLAS level 1: output = \f$ ||x||_2 \f$
170static inline double Dnrm2(const int &n, const double *x, const int &incx)
171{
172 return F77NAME(dnrm2)(n, x, incx);
173}
174
175/// \brief BLAS level 1: output = \f$ ||x||_1 \f$
176static inline double Dasum(const int &n, const double *x, const int &incx)
177{
178 return F77NAME(dasum)(n, x, incx);
179}
180
181/// \brief BLAS level 1: output = 1st value where \f$ |x[i]| = max |x|_1 \f$
182/// Note it is modified to return a value between (0,n-1) as per
183/// the standard C convention
184static inline int Idamax(const int &n, const double *x, const int &incx)
185{
186 return F77NAME(idamax)(n, x, incx) - 1;
187}
188
189/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
190/// where A[m x n]
191static inline void Gemv(const char &trans, const int &m, const int &n,
192 const double &alpha, const double *a, const int &lda,
193 const double *x, const int &incx, const double &beta,
194 double *y, const int &incy)
195{
196 F77NAME(dgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
197}
198
199/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
200/// where A[m x n]
201static inline void Gemv(const char &trans, const int &m, const int &n,
202 const float &alpha, const float *a, const int &lda,
203 const float *x, const int &incx, const float &beta,
204 float *y, const int &incy)
205{
206 F77NAME(sgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
207}
208
209/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
210/// where A[m x n]
211static inline void Dgemv(const char &trans, const int &m, const int &n,
212 const double &alpha, const double *a, const int &lda,
213 const double *x, const int &incx, const double &beta,
214 double *y, const int &incy)
215{
216 F77NAME(dgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
217}
218
219/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
220/// where A[m x n] is banded
221static inline void Gbmv(const char &trans, const int &m, const int &n,
222 const int &kl, const int &ku, const double &alpha,
223 const double *a, const int &lda, const double *x,
224 const int &incx, const double &beta, double *y,
225 const int &incy)
226{
227 F77NAME(dgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
228}
229
230/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
231/// where A[m x n] is banded
232static inline void Gbmv(const char &trans, const int &m, const int &n,
233 const int &kl, const int &ku, const float &alpha,
234 const float *a, const int &lda, const float *x,
235 const int &incx, const float &beta, float *y,
236 const int &incy)
237{
238 F77NAME(sgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
239}
240
241/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
242/// where A[m x n] is banded
243static inline void Dgbmv(const char &trans, const int &m, const int &n,
244 const int &kl, const int &ku, const double &alpha,
245 const double *a, const int &lda, const double *x,
246 const int &incx, const double &beta, double *y,
247 const int &incy)
248{
249 F77NAME(dgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
250}
251
252/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
253/// where A[n x n] is symmetric packed
254static inline void Spmv(const char &uplo, const int &n, const double &alpha,
255 const double *a, const double *x, const int &incx,
256 const double &beta, double *y, const int &incy)
257{
258 F77NAME(dspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
259}
260
261/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
262/// where A[n x n] is symmetric packed
263static inline void Spmv(const char &uplo, const int &n, const float &alpha,
264 const float *a, const float *x, const int &incx,
265 const float &beta, float *y, const int &incy)
266{
267 F77NAME(sspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
268}
269
270/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
271/// where A[n x n] is symmetric packed
272static inline void Dspmv(const char &uplo, const int &n, const double &alpha,
273 const double *a, const double *x, const int &incx,
274 const double &beta, double *y, const int &incy)
275{
276 F77NAME(dspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
277}
278
279/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
280/// where A[m x m] is symmetric banded
281static inline void Sbmv(const char &uplo, const int &m, const int &k,
282 const double &alpha, const double *a, const int &lda,
283 const double *x, const int &incx, const double &beta,
284 double *y, const int &incy)
285{
286 F77NAME(dsbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
287}
288
289/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
290/// where A[m x m] is symmetric banded
291static inline void Sbmv(const char &uplo, const int &m, const int &k,
292 const float &alpha, const float *a, const int &lda,
293 const float *x, const int &incx, const float &beta,
294 float *y, const int &incy)
295{
296 F77NAME(ssbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
297}
298
299/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
300/// where A[m x m] is symmetric banded
301static inline void Dsbmv(const char &uplo, const int &m, const int &k,
302 const double &alpha, const double *a, const int &lda,
303 const double *x, const int &incx, const double &beta,
304 double *y, const int &incy)
305{
306 F77NAME(dsbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
307}
308
309/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
310static inline void Tpmv(const char &uplo, const char &trans, const char &diag,
311 const int &n, const double *ap, double *x,
312 const int &incx)
313{
314 F77NAME(dtpmv)(uplo, trans, diag, n, ap, x, incx);
315}
316
317/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
318static inline void Tpmv(const char &uplo, const char &trans, const char &diag,
319 const int &n, const float *ap, float *x,
320 const int &incx)
321{
322 F77NAME(stpmv)(uplo, trans, diag, n, ap, x, incx);
323}
324
325/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
326static inline void Dtpmv(const char &uplo, const char &trans, const char &diag,
327 const int &n, const double *ap, double *x,
328 const int &incx)
329{
330 F77NAME(dtpmv)(uplo, trans, diag, n, ap, x, incx);
331}
332
333/// \brief BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A
334/// where A[m x n]
335static inline void Ger(const int &m, const int &n, const double &alpha,
336 const double *x, const int &incx, const double *y,
337 const int &incy, double *a, const int &lda)
338{
339 F77NAME(dger)(m, n, alpha, x, incx, y, incy, a, lda);
340}
341
342/// \brief BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A
343/// where A[m x n]
344static inline void Ger(const int &m, const int &n, const float &alpha,
345 const float *x, const int &incx, const float *y,
346 const int &incy, float *a, const int &lda)
347{
348 F77NAME(sger)(m, n, alpha, x, incx, y, incy, a, lda);
349}
350
351/// \brief BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k],
352/// op(B)[k x n], C[m x n]
353/// DGEMM performs one of the matrix-matrix operations:
354/// C := alpha*op( A )*op( B ) + beta*C,
355static inline void Gemm(const char &transa, const char &transb, const int &m,
356 const int &n, const int &k, const double &alpha,
357 const double *a, const int &lda, const double *b,
358 const int &ldb, const double &beta, double *c,
359 const int &ldc)
360{
362 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
363}
364
365/// \brief BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k],
366/// op(B)[k x n], C[m x n]
367/// DGEMM performs one of the matrix-matrix operations:
368/// C := alpha*op( A )*op( B ) + beta*C,
369static inline void Gemm(const char &transa, const char &transb, const int &m,
370 const int &n, const int &k, const float &alpha,
371 const float *a, const int &lda, const float *b,
372 const int &ldb, const float &beta, float *c,
373 const int &ldc)
374{
376 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
377}
378
379/// \brief BLAS level 3: Matrix-matrix multiply C = A x B
380/// where op(A)[m x k], op(B)[k x n], C[m x n]
381/// DGEMM performs one of the matrix-matrix operations:
382/// C := alpha*op( A )*op( B ) + beta*C,
383static inline void Dgemm(const char &transa, const char &transb, const int &m,
384 const int &n, const int &k, const double &alpha,
385 const double *a, const int &lda, const double *b,
386 const int &ldb, const double &beta, double *c,
387 const int &ldc)
388{
390 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
391}
392
393// \brief Wrapper to mutliply two (row major) matrices together C =
394// a*A*B + b*C
395static inline void Cdgemm(const int M, const int N, const int K, const double a,
396 const double *A, [[maybe_unused]] const int ldA,
397 const double *B, [[maybe_unused]] const int ldB,
398 const double b, double *C,
399 [[maybe_unused]] const int ldC)
400{
401 Dgemm('N', 'N', N, M, K, a, B, N, A, K, b, C, N);
402}
403} // namespace Blas
404#endif // NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_BLAS_HPP
#define F77NAME(x)
Fortran routines need an underscore.
Definition: TransF77.hpp:46
Definition: Blas.hpp:43
int F77NAME() idamax(const int &n, const double *x, const int &incx)
static void Dgemv(const char &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)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:211
static void Dgbmv(const char &trans, const int &m, const int &n, const int &kl, const int &ku, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n] is banded.
Definition: Blas.hpp:243
static void Dswap(const int &n, double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Swap x with y.
Definition: Blas.hpp:142
static void Drot(const int &n, double *x, const int &incx, double *y, const int &incy, const double &c, const double &s)
BLAS level 1: Plane rotation by c = cos(theta), s = sin(theta)
Definition: Blas.hpp:156
void F77NAME() daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:149
static void Sbmv(const char &uplo, const int &m, const int &k, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x m] is symmetric banded.
Definition: Blas.hpp:281
void F77NAME() dswap(const int &n, double *x, const int &incx, double *y, const int &incy)
void F77NAME() drot(const int &n, double *x, const int &incx, double *y, const int &incy, const double &c, const double &s)
static void Gemm(const char &transa, const char &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)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:355
double F77NAME() dasum(const int &n, const double *x, const int &incx)
void F77NAME() dtpmv(const char &uplo, const char &trans, const char &diag, const int &n, const double *ap, double *x, const int &incx)
void F77NAME() stpmv(const char &uplo, const char &trans, const char &diag, const int &n, const float *ap, float *x, const int &incx)
static void Gemv(const char &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)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:191
static double Dasum(const int &n, const double *x, const int &incx)
BLAS level 1: output = .
Definition: Blas.hpp:176
void F77NAME() dspmv(const char &uplo, const int &n, const double &alpha, const double *a, const double *x, const int &incx, const double &beta, double *y, const int &incy)
void F77NAME() dscal(const int &n, const double &alpha, double *x, const int &incx)
double F77NAME() ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
static void Gbmv(const char &trans, const int &m, const int &n, const int &kl, const int &ku, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x n] is banded.
Definition: Blas.hpp:221
static void Dsbmv(const char &uplo, const int &m, const int &k, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[m x m] is symmetric banded.
Definition: Blas.hpp:301
static int Idamax(const int &n, const double *x, const int &incx)
BLAS level 1: output = 1st value where Note it is modified to return a value between (0,...
Definition: Blas.hpp:184
void F77NAME() ssbmv(const char &uplo, const int &m, const int &k, const float &alpha, const float *a, const int &lda, const float *x, const int &incx, const float &beta, float *y, const int &incy)
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
Definition: Blas.hpp:128
void F77NAME() dsbmv(const char &uplo, const int &m, const int &k, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
void F77NAME() sspmv(const char &uplo, const int &n, const float &alpha, const float *a, const float *x, const int &incx, const float &beta, float *y, const int &incy)
void F77NAME() sgbmv(const char &trans, const int &m, const int &n, const int &kl, const int &ku, const float &alpha, const float *a, const int &lda, const float *x, const int &incx, const float &beta, float *y, const int &incy)
static double Dnrm2(const int &n, const double *x, const int &incx)
Definition: Blas.hpp:170
void F77NAME() dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
void F77NAME() sger(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)
void F77NAME() sgemv(const char &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)
void F77NAME() dgbmv(const char &trans, const int &m, const int &n, const int &kl, const int &ku, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
void F77NAME() dgemm(const char &trans, const char &transb, const int &m1, 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 F77NAME() dgemv(const char &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)
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:163
static void Cdgemm(const int M, const int N, const int K, const double a, const double *A, const int ldA, const double *B, const int ldB, const double b, double *C, const int ldC)
Definition: Blas.hpp:395
static void Dgemm(const char &transa, const char &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)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:383
void F77NAME() sgemm(const char &trans, const char &transb, const int &m1, 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)
double F77NAME() dnrm2(const int &n, const double *x, const int &incx)
static void Tpmv(const char &uplo, const char &trans, const char &diag, const int &n, const double *ap, double *x, const int &incx)
BLAS level 2: Matrix vector multiply x = alpha A x where A[n x n].
Definition: Blas.hpp:310
static void Dspmv(const char &uplo, const int &n, const double &alpha, const double *a, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[n x n] is symmetric packed.
Definition: Blas.hpp:272
void F77NAME() dger(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)
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:135
static void 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)
BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A where A[m x n].
Definition: Blas.hpp:335
static void Spmv(const char &uplo, const int &n, const double &alpha, const double *a, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = alpha A x plus beta y where A[n x n] is symmetric packed.
Definition: Blas.hpp:254
static void Dtpmv(const char &uplo, const char &trans, const char &diag, const int &n, const double *ap, double *x, const int &incx)
BLAS level 2: Matrix vector multiply x = alpha A x where A[n x n].
Definition: Blas.hpp:326
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59