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
38#include <boost/core/ignore_unused.hpp>
39
42
43// Translations for using Fortran version of blas
44namespace Blas
45{
46extern "C"
47{
48 // -- BLAS Level 1:
49 void F77NAME(dcopy)(const int &n, const double *x, const int &incx,
50 double *y, const int &incy);
51 void F77NAME(daxpy)(const int &n, const double &alpha, const double *x,
52 const int &incx, const double *y, const int &incy);
53 void F77NAME(dswap)(const int &n, double *x, const int &incx, double *y,
54 const int &incy);
55 void F77NAME(dscal)(const int &n, const double &alpha, double *x,
56 const int &incx);
57 void F77NAME(drot)(const int &n, double *x, const int &incx, double *y,
58 const int &incy, const double &c, const double &s);
59 double F77NAME(ddot)(const int &n, const double *x, const int &incx,
60 const double *y, const int &incy);
61 double F77NAME(dnrm2)(const int &n, const double *x, const int &incx);
62 double F77NAME(dasum)(const int &n, const double *x, const int &incx);
63 int F77NAME(idamax)(const int &n, const double *x, const int &incx);
64
65 // -- BLAS level 2
66 void F77NAME(dgemv)(const char &trans, const int &m, const int &n,
67 const double &alpha, const double *a, const int &lda,
68 const double *x, const int &incx, const double &beta,
69 double *y, const int &incy);
70 void F77NAME(sgemv)(const char &trans, const int &m, const int &n,
71 const float &alpha, const float *a, const int &lda,
72 const float *x, const int &incx, const float &beta,
73 float *y, const int &incy);
74
75 void F77NAME(dgbmv)(const char &trans, const int &m, const int &n,
76 const int &kl, const int &ku, const double &alpha,
77 const double *a, const int &lda, const double *x,
78 const int &incx, const double &beta, double *y,
79 const int &incy);
80 void F77NAME(sgbmv)(const char &trans, const int &m, const int &n,
81 const int &kl, const int &ku, const float &alpha,
82 const float *a, const int &lda, const float *x,
83 const int &incx, const float &beta, float *y,
84 const int &incy);
85
86 void F77NAME(dspmv)(const char &uplo, const int &n, const double &alpha,
87 const double *a, const double *x, const int &incx,
88 const double &beta, double *y, const int &incy);
89 void F77NAME(sspmv)(const char &uplo, const int &n, const float &alpha,
90 const float *a, const float *x, const int &incx,
91 const float &beta, float *y, const int &incy);
92
93 void F77NAME(dsbmv)(const char &uplo, const int &m, const int &k,
94 const double &alpha, const double *a, const int &lda,
95 const double *x, const int &incx, const double &beta,
96 double *y, const int &incy);
97 void F77NAME(ssbmv)(const char &uplo, const int &m, const int &k,
98 const float &alpha, const float *a, const int &lda,
99 const float *x, const int &incx, const float &beta,
100 float *y, const int &incy);
101
102 void F77NAME(dtpmv)(const char &uplo, const char &trans, const char &diag,
103 const int &n, const double *ap, double *x,
104 const int &incx);
105 void F77NAME(stpmv)(const char &uplo, const char &trans, const char &diag,
106 const int &n, const float *ap, float *x,
107 const int &incx);
108
109 void F77NAME(dger)(const int &m, const int &n, const double &alpha,
110 const double *x, const int &incx, const double *y,
111 const int &incy, double *a, const int &lda);
112 void F77NAME(sger)(const int &m, const int &n, const float &alpha,
113 const float *x, const int &incx, const float *y,
114 const int &incy, float *a, const int &lda);
115
116 // -- BLAS level 3:
117 void F77NAME(dgemm)(const char &trans, const char &transb, const int &m1,
118 const int &n, const int &k, const double &alpha,
119 const double *a, const int &lda, const double *b,
120 const int &ldb, const double &beta, double *c,
121 const int &ldc);
122 void F77NAME(sgemm)(const char &trans, const char &transb, const int &m1,
123 const int &n, const int &k, const float &alpha,
124 const float *a, const int &lda, const float *b,
125 const int &ldb, const float &beta, float *c,
126 const int &ldc);
127}
128
129/// \brief BLAS level 1: Copy \a x to \a y
130static inline void Dcopy(const int &n, const double *x, const int &incx,
131 double *y, const int &incy)
132{
133 F77NAME(dcopy)(n, x, incx, y, incy);
134}
135
136/// \brief BLAS level 1: y = alpha \a x plus \a y
137static inline void Daxpy(const int &n, const double &alpha, const double *x,
138 const int &incx, const double *y, const int &incy)
139{
140 F77NAME(daxpy)(n, alpha, x, incx, y, incy);
141}
142
143/// \brief BLAS level 1: Swap \a x with \a y
144static inline void Dswap(const int &n, double *x, const int &incx, double *y,
145 const int &incy)
146{
147 F77NAME(dswap)(n, x, incx, y, incy);
148}
149
150/// \brief BLAS level 1: x = alpha \a x
151static inline void Dscal(const int &n, const double &alpha, double *x,
152 const int &incx)
153{
154 F77NAME(dscal)(n, alpha, x, incx);
155}
156
157/// \brief BLAS level 1: Plane rotation by c = cos(theta), s = sin(theta)
158static inline void Drot(const int &n, double *x, const int &incx, double *y,
159 const int &incy, const double &c, const double &s)
160{
161 F77NAME(drot)(n, x, incx, y, incy, c, s);
162}
163
164/// \brief BLAS level 1: output = \f$ x^T y \f$
165static inline double Ddot(const int &n, const double *x, const int &incx,
166 const double *y, const int &incy)
167{
168 return F77NAME(ddot)(n, x, incx, y, incy);
169}
170
171// \brief BLAS level 1: output = \f$ ||x||_2 \f$
172static inline double Dnrm2(const int &n, const double *x, const int &incx)
173{
174 return F77NAME(dnrm2)(n, x, incx);
175}
176
177/// \brief BLAS level 1: output = \f$ ||x||_1 \f$
178static inline double Dasum(const int &n, const double *x, const int &incx)
179{
180 return F77NAME(dasum)(n, x, incx);
181}
182
183/// \brief BLAS level 1: output = 1st value where \f$ |x[i]| = max |x|_1 \f$
184/// Note it is modified to return a value between (0,n-1) as per
185/// the standard C convention
186static inline int Idamax(const int &n, const double *x, const int &incx)
187{
188 return F77NAME(idamax)(n, x, incx) - 1;
189}
190
191/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
192/// where A[m x n]
193static inline void Gemv(const char &trans, const int &m, const int &n,
194 const double &alpha, const double *a, const int &lda,
195 const double *x, const int &incx, const double &beta,
196 double *y, const int &incy)
197{
198 F77NAME(dgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
199}
200
201/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
202/// where A[m x n]
203static inline void Gemv(const char &trans, const int &m, const int &n,
204 const float &alpha, const float *a, const int &lda,
205 const float *x, const int &incx, const float &beta,
206 float *y, const int &incy)
207{
208 F77NAME(sgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
209}
210
211/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
212/// where A[m x n]
213static inline void Dgemv(const char &trans, const int &m, const int &n,
214 const double &alpha, const double *a, const int &lda,
215 const double *x, const int &incx, const double &beta,
216 double *y, const int &incy)
217{
218 F77NAME(dgemv)(trans, m, n, alpha, a, lda, x, incx, beta, y, incy);
219}
220
221/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
222/// where A[m x n] is banded
223static inline void Gbmv(const char &trans, const int &m, const int &n,
224 const int &kl, const int &ku, const double &alpha,
225 const double *a, const int &lda, const double *x,
226 const int &incx, const double &beta, double *y,
227 const int &incy)
228{
229 F77NAME(dgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
230}
231
232/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
233/// where A[m x n] is banded
234static inline void Gbmv(const char &trans, const int &m, const int &n,
235 const int &kl, const int &ku, const float &alpha,
236 const float *a, const int &lda, const float *x,
237 const int &incx, const float &beta, float *y,
238 const int &incy)
239{
240 F77NAME(sgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
241}
242
243/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
244/// where A[m x n] is banded
245static inline void Dgbmv(const char &trans, const int &m, const int &n,
246 const int &kl, const int &ku, const double &alpha,
247 const double *a, const int &lda, const double *x,
248 const int &incx, const double &beta, double *y,
249 const int &incy)
250{
251 F77NAME(dgbmv)(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy);
252}
253
254/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
255/// where A[n x n] is symmetric packed
256static inline void Spmv(const char &uplo, const int &n, const double &alpha,
257 const double *a, const double *x, const int &incx,
258 const double &beta, double *y, const int &incy)
259{
260 F77NAME(dspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
261}
262
263/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
264/// where A[n x n] is symmetric packed
265static inline void Spmv(const char &uplo, const int &n, const float &alpha,
266 const float *a, const float *x, const int &incx,
267 const float &beta, float *y, const int &incy)
268{
269 F77NAME(sspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
270}
271
272/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
273/// where A[n x n] is symmetric packed
274static inline void Dspmv(const char &uplo, const int &n, const double &alpha,
275 const double *a, const double *x, const int &incx,
276 const double &beta, double *y, const int &incy)
277{
278 F77NAME(dspmv)(uplo, n, alpha, a, x, incx, beta, y, incy);
279}
280
281/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
282/// where A[m x m] is symmetric banded
283static inline void Sbmv(const char &uplo, const int &m, const int &k,
284 const double &alpha, const double *a, const int &lda,
285 const double *x, const int &incx, const double &beta,
286 double *y, const int &incy)
287{
288 F77NAME(dsbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
289}
290
291/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
292/// where A[m x m] is symmetric banded
293static inline void Sbmv(const char &uplo, const int &m, const int &k,
294 const float &alpha, const float *a, const int &lda,
295 const float *x, const int &incx, const float &beta,
296 float *y, const int &incy)
297{
298 F77NAME(ssbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
299}
300
301/// \brief BLAS level 2: Matrix vector multiply y = alpha A \e x plus beta \a y
302/// where A[m x m] is symmetric banded
303static inline void Dsbmv(const char &uplo, const int &m, const int &k,
304 const double &alpha, const double *a, const int &lda,
305 const double *x, const int &incx, const double &beta,
306 double *y, const int &incy)
307{
308 F77NAME(dsbmv)(uplo, m, k, alpha, a, lda, x, incx, beta, y, incy);
309}
310
311/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
312static inline void Tpmv(const char &uplo, const char &trans, const char &diag,
313 const int &n, const double *ap, double *x,
314 const int &incx)
315{
316 F77NAME(dtpmv)(uplo, trans, diag, n, ap, x, incx);
317}
318
319/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
320static inline void Tpmv(const char &uplo, const char &trans, const char &diag,
321 const int &n, const float *ap, float *x,
322 const int &incx)
323{
324 F77NAME(stpmv)(uplo, trans, diag, n, ap, x, incx);
325}
326
327/// \brief BLAS level 2: Matrix vector multiply x = alpha A \e x where A[n x n]
328static inline void Dtpmv(const char &uplo, const char &trans, const char &diag,
329 const int &n, const double *ap, double *x,
330 const int &incx)
331{
332 F77NAME(dtpmv)(uplo, trans, diag, n, ap, x, incx);
333}
334
335/// \brief BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A
336/// where A[m x n]
337static inline void Ger(const int &m, const int &n, const double &alpha,
338 const double *x, const int &incx, const double *y,
339 const int &incy, double *a, const int &lda)
340{
341 F77NAME(dger)(m, n, alpha, x, incx, y, incy, a, lda);
342}
343
344/// \brief BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A
345/// where A[m x n]
346static inline void Ger(const int &m, const int &n, const float &alpha,
347 const float *x, const int &incx, const float *y,
348 const int &incy, float *a, const int &lda)
349{
350 F77NAME(sger)(m, n, alpha, x, incx, y, incy, a, lda);
351}
352
353/// \brief BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k],
354/// op(B)[k x n], C[m x n]
355/// DGEMM performs one of the matrix-matrix operations:
356/// C := alpha*op( A )*op( B ) + beta*C,
357static inline void Gemm(const char &transa, const char &transb, const int &m,
358 const int &n, const int &k, const double &alpha,
359 const double *a, const int &lda, const double *b,
360 const int &ldb, const double &beta, double *c,
361 const int &ldc)
362{
364 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
365}
366
367/// \brief BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k],
368/// op(B)[k x n], C[m x n]
369/// DGEMM performs one of the matrix-matrix operations:
370/// C := alpha*op( A )*op( B ) + beta*C,
371static inline void Gemm(const char &transa, const char &transb, const int &m,
372 const int &n, const int &k, const float &alpha,
373 const float *a, const int &lda, const float *b,
374 const int &ldb, const float &beta, float *c,
375 const int &ldc)
376{
378 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
379}
380
381/// \brief BLAS level 3: Matrix-matrix multiply C = A x B
382/// where op(A)[m x k], op(B)[k x n], C[m x n]
383/// DGEMM performs one of the matrix-matrix operations:
384/// C := alpha*op( A )*op( B ) + beta*C,
385static inline void Dgemm(const char &transa, const char &transb, const int &m,
386 const int &n, const int &k, const double &alpha,
387 const double *a, const int &lda, const double *b,
388 const int &ldb, const double &beta, double *c,
389 const int &ldc)
390{
392 (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc);
393}
394
395// \brief Wrapper to mutliply two (row major) matrices together C =
396// a*A*B + b*C
397static inline void Cdgemm(const int M, const int N, const int K, const double a,
398 const double *A, const int ldA, const double *B,
399 const int ldB, const double b, double *C,
400 const int ldC)
401{
402 boost::ignore_unused(ldA, ldB, ldC);
403 Dgemm('N', 'N', N, M, K, a, B, N, A, K, b, C, N);
404}
405} // namespace Blas
406#endif // NEKTAR_LIB_UTILITIES_LINEAR_ALGEBRA_BLAS_HPP
#define F77NAME(x)
Fortran routines need an underscore.
Definition: TransF77.hpp:46
Definition: Blas.hpp:45
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:213
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:245
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:144
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:158
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:151
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:283
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:357
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:193
static double Dasum(const int &n, const double *x, const int &incx)
BLAS level 1: output = .
Definition: Blas.hpp:178
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:223
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:303
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:186
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:130
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:172
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:165
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:397
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:385
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:312
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:274
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:137
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:337
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:256
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:328
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:61