mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
templateBLAS.cpp
Go to the documentation of this file.
1/** \file templateBLAS.cpp
2 * \brief Implementation of templatized wrappers for the BLAS
3 * \ingroup gen_math_files
4 * \author Jared R. Males (jaredmales@gmail.com)
5 *
6 */
7
8//***********************************************************************//
9// Copyright 2015-2020 Jared R. Males (jaredmales@gmail.com)
10//
11// This file is part of mxlib.
12//
13// mxlib is free software: you can redistribute it and/or modify
14// it under the terms of the GNU General Public License as published by
15// the Free Software Foundation, either version 3 of the License, or
16// (at your option) any later version.
17//
18// mxlib is distributed in the hope that it will be useful,
19// but WITHOUT ANY WARRANTY; without even the implied warranty of
20// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21// GNU General Public License for more details.
22//
23// You should have received a copy of the GNU General Public License
24// along with mxlib. If not, see <http://www.gnu.org/licenses/>.
25//***********************************************************************//
26
27#include "math/templateBLAS.hpp"
28
29namespace mx
30{
31namespace math
32{
33
34template <>
35void scal<float>( const int N, const float &alpha, float *X, const int incX )
36{
37 cblas_sscal( N, alpha, X, incX );
38}
39
40template <>
41void scal<double>( const int N, const double &alpha, double *X, const int incX )
42{
43 cblas_dscal( N, alpha, X, incX );
44}
45
46template <>
47void scal<std::complex<float>>( const int N, const std::complex<float> &alpha, std::complex<float> *X, const int incX )
48{
49 cblas_cscal( N, &alpha, X, incX );
50}
51
52template <>
53void scal<std::complex<double>>( const int N,
54 const std::complex<double> &alpha,
55 std::complex<double> *X,
56 const int incX )
57{
58 cblas_zscal( N, &alpha, X, incX );
59}
60
61template <>
62void gemm<float>( const CBLAS_ORDER Order,
65 const int M,
66 const int N,
67 const int K,
68 const float &alpha,
69 const float *A,
70 const int lda,
71 const float *B,
72 const int ldb,
73 const float &beta,
74 float *C,
75 const int ldc )
76{
77 cblas_sgemm( Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc );
78}
79
80template <>
81void gemm<double>( const CBLAS_ORDER Order,
84 const int M,
85 const int N,
86 const int K,
87 const double &alpha,
88 const double *A,
89 const int lda,
90 const double *B,
91 const int ldb,
92 const double &beta,
93 double *C,
94 const int ldc )
95{
96 cblas_dgemm( Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc );
97}
98
99template <>
103 const int M,
104 const int N,
105 const int K,
106 const std::complex<float> &alpha,
107 const std::complex<float> *A,
108 const int lda,
109 const std::complex<float> *B,
110 const int ldb,
111 const std::complex<float> &beta,
112 std::complex<float> *C,
113 const int ldc )
114{
115 cblas_cgemm( Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc );
116}
117
118template <>
122 const int M,
123 const int N,
124 const int K,
125 const std::complex<double> &alpha,
126 const std::complex<double> *A,
127 const int lda,
128 const std::complex<double> *B,
129 const int ldb,
130 const std::complex<double> &beta,
131 std::complex<double> *C,
132 const int ldc )
133{
134 cblas_zgemm( Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc );
135}
136
137template <>
138void syrk<float>( const CBLAS_ORDER Order,
139 const CBLAS_UPLO Uplo,
141 const int N,
142 const int K,
143 const float &alpha,
144 const float *A,
145 const int lda,
146 const float &beta,
147 float *C,
148 const int ldc )
149{
150 cblas_ssyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc );
151}
152
153template <>
154void syrk<double>( const CBLAS_ORDER Order,
155 const CBLAS_UPLO Uplo,
157 const int N,
158 const int K,
159 const double &alpha,
160 const double *A,
161 const int lda,
162 const double &beta,
163 double *C,
164 const int ldc )
165{
166 cblas_dsyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc );
167}
168
169template <>
171 const CBLAS_UPLO Uplo,
173 const int N,
174 const int K,
175 const std::complex<float> &alpha,
176 const std::complex<float> *A,
177 const int lda,
178 const std::complex<float> &beta,
179 std::complex<float> *C,
180 const int ldc )
181{
182 cblas_csyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc );
183}
184
185template <>
187 const CBLAS_UPLO Uplo,
189 const int N,
190 const int K,
191 const std::complex<double> &alpha,
192 const std::complex<double> *A,
193 const int lda,
194 const std::complex<double> &beta,
195 std::complex<double> *C,
196 const int ldc )
197{
198 cblas_zsyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc );
199}
200
201} // namespace math
202} // namespace mx
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
The mxlib c++ namespace.
Definition mxError.hpp:106
Declares and defines templatized wrappers for the BLAS.