mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
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 
29 namespace mx
30 {
31 namespace math
32 {
33 
34 template<>
35 void scal<float>( const int N,
36  const float & alpha,
37  float * X,
38  const int incX )
39 {
40  cblas_sscal(N, alpha, X, incX);
41 }
42 
43 template<>
44 void scal<double>( const int N,
45  const double & alpha,
46  double * X,
47  const int incX
48  )
49 {
50  cblas_dscal(N, alpha, X, incX);
51 }
52 
53 template<>
54 void scal<std::complex<float> >( const int N,
55  const std::complex<float> & alpha,
56  std::complex<float> * X,
57  const int incX
58  )
59 {
60  cblas_cscal(N, &alpha, X, incX);
61 }
62 
63 
64 template<>
65 void scal<std::complex<double> >( const int N,
66  const std::complex<double> & alpha,
67  std::complex<double> * X,
68  const int incX
69  )
70 {
71  cblas_zscal(N, &alpha, X, incX);
72 }
73 
74 template<>
75 void gemm<float>(const CBLAS_ORDER Order, const CBLAS_TRANSPOSE TransA,
76  const CBLAS_TRANSPOSE TransB, const int M, const int N,
77  const int K, const float & alpha, const float *A,
78  const int lda, const float *B, const int ldb,
79  const float & beta, float *C, const int ldc)
80 {
81  cblas_sgemm(Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc);
82 }
83 
84 template<>
85 void gemm<double>(const CBLAS_ORDER Order, const CBLAS_TRANSPOSE TransA,
86  const CBLAS_TRANSPOSE TransB, const int M, const int N,
87  const int K, const double & alpha, const double *A,
88  const int lda, const double *B, const int ldb,
89  const double & beta, double *C, const int ldc)
90 {
91  cblas_dgemm(Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc);
92 }
93 
94 
95 template<>
96 void gemm<std::complex<float> >(const CBLAS_ORDER Order, const CBLAS_TRANSPOSE TransA,
97  const CBLAS_TRANSPOSE TransB, const int M, const int N,
98  const int K, const std::complex<float> & alpha, const std::complex<float> *A,
99  const int lda, const std::complex<float> *B, const int ldb,
100  const std::complex<float> & beta, std::complex<float> *C, const int ldc)
101 {
102  cblas_cgemm(Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
103 }
104 
105 template<>
106 void gemm<std::complex<double> >(const CBLAS_ORDER Order, const CBLAS_TRANSPOSE TransA,
107  const CBLAS_TRANSPOSE TransB, const int M, const int N,
108  const int K, const std::complex<double> & alpha, const std::complex<double> *A,
109  const int lda, const std::complex<double> *B, const int ldb,
110  const std::complex<double> & beta, std::complex<double> *C, const int ldc)
111 {
112  cblas_zgemm(Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
113 }
114 
115 template<>
116 void syrk<float>(const CBLAS_ORDER Order, const CBLAS_UPLO Uplo,
117  const CBLAS_TRANSPOSE Trans, const int N, const int K,
118  const float & alpha, const float *A, const int lda,
119  const float & beta, float *C, const int ldc)
120 {
121  cblas_ssyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc);
122 }
123 
124 template<>
125 void syrk<double>(const CBLAS_ORDER Order, const CBLAS_UPLO Uplo,
126  const CBLAS_TRANSPOSE Trans, const int N, const int K,
127  const double & alpha, const double *A, const int lda,
128  const double & beta, double *C, const int ldc)
129 {
130  cblas_dsyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc);
131 }
132 
133 template<>
134 void syrk<std::complex<float> >(const CBLAS_ORDER Order, const CBLAS_UPLO Uplo,
135  const CBLAS_TRANSPOSE Trans, const int N, const int K,
136  const std::complex<float> & alpha, const std::complex<float> *A, const int lda,
137  const std::complex<float> & beta, std::complex<float> *C, const int ldc)
138 {
139  cblas_csyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc);
140 }
141 
142 template<>
143 void syrk<std::complex<double> >(const CBLAS_ORDER Order, const CBLAS_UPLO Uplo,
144  const CBLAS_TRANSPOSE Trans, const int N, const int K,
145  const std::complex<double> & alpha, const std::complex<double> *A, const int lda,
146  const std::complex<double> & beta, std::complex<double> *C, const int ldc)
147 {
148  cblas_zsyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc);
149 }
150 
151 } //namespace math
152 } //namespace mx
153 
154 
155 
The mxlib c++ namespace.
Definition: mxError.hpp:107
Declares and defines templatized wrappers for the BLAS.