cblas3_overloads.hpp

Go to the documentation of this file.
00001 /*
00002  * 
00003  * Copyright (c) Kresimir Fresl 2002 
00004  *
00005  * Distributed under the Boost Software License, Version 1.0.
00006  * (See accompanying file LICENSE_1_0.txt or copy at
00007  * http://www.boost.org/LICENSE_1_0.txt)
00008  *
00009  * Author acknowledges the support of the Faculty of Civil Engineering, 
00010  * University of Zagreb, Croatia.
00011  *
00012  */
00013 
00014 #ifndef BOOST_NUMERIC_BINDINGS_CBLAS3_OVERLOADS_HPP
00015 #define BOOST_NUMERIC_BINDINGS_CBLAS3_OVERLOADS_HPP
00016 
00017 #include <complex> 
00018 #include <boost/numeric/bindings/atlas/cblas_inc.hpp>
00019 #include <boost/numeric/bindings/traits/type.hpp>
00020 
00021 
00022 namespace boost { namespace numeric { namespace bindings { 
00023 
00024   namespace atlas { namespace detail {
00025 
00026     // C <- alpha * op (A) * op (B) + beta * C 
00027   
00028     inline 
00029     void gemm (CBLAS_ORDER const Order, 
00030                CBLAS_TRANSPOSE const TransA, CBLAS_TRANSPOSE const TransB, 
00031                int const M, int const N, int const K, 
00032                float const alpha, float const* A, int const lda,
00033                float const* B, int const ldb, 
00034                float const beta, float* C, int const ldc) 
00035     {
00036       cblas_sgemm (Order, TransA, TransB, M, N, K, 
00037                    alpha, A, lda, 
00038                    B, ldb,
00039                    beta, C, ldc); 
00040     }
00041 
00042     inline 
00043     void gemm (CBLAS_ORDER const Order, 
00044                CBLAS_TRANSPOSE const TransA, CBLAS_TRANSPOSE const TransB, 
00045                int const M, int const N, int const K, 
00046                double const alpha, double const* A, int const lda,
00047                double const* B, int const ldb, 
00048                double const beta, double* C, int const ldc) 
00049     {
00050       cblas_dgemm (Order, TransA, TransB, M, N, K, 
00051                    alpha, A, lda, 
00052                    B, ldb,
00053                    beta, C, ldc); 
00054     }
00055 
00056     inline 
00057     void gemm (CBLAS_ORDER const Order, 
00058                CBLAS_TRANSPOSE const TransA, CBLAS_TRANSPOSE const TransB, 
00059                int const M, int const N, int const K, 
00060                traits::complex_f const& alpha, 
00061                traits::complex_f const* A, int const lda,
00062                traits::complex_f const* B, int const ldb, 
00063                traits::complex_f const& beta, 
00064                traits::complex_f* C, int const ldc) 
00065     {
00066       cblas_cgemm (Order, TransA, TransB, M, N, K, 
00067                    static_cast<void const*> (&alpha), 
00068                    static_cast<void const*> (A), lda, 
00069                    static_cast<void const*> (B), ldb,
00070                    static_cast<void const*> (&beta), 
00071                    static_cast<void*> (C), ldc); 
00072     }
00073     
00074     inline 
00075     void gemm (CBLAS_ORDER const Order, 
00076                CBLAS_TRANSPOSE const TransA, CBLAS_TRANSPOSE const TransB, 
00077                int const M, int const N, int const K, 
00078                traits::complex_d const& alpha, 
00079                traits::complex_d const* A, int const lda,
00080                traits::complex_d const* B, int const ldb, 
00081                traits::complex_d const& beta, 
00082                traits::complex_d* C, int const ldc) 
00083     {
00084       cblas_zgemm (Order, TransA, TransB, M, N, K, 
00085                    static_cast<void const*> (&alpha), 
00086                    static_cast<void const*> (A), lda, 
00087                    static_cast<void const*> (B), ldb,
00088                    static_cast<void const*> (&beta), 
00089                    static_cast<void*> (C), ldc); 
00090     }
00091 
00092     
00093     // C <- alpha * A * B + beta * C 
00094     // C <- alpha * B * A + beta * C 
00095     // A == A^T
00096 
00097     inline 
00098     void symm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00099                CBLAS_UPLO const Uplo, int const M, int const N, 
00100                float const alpha, float const* A, int const lda,
00101                float const* B, int const ldb, 
00102                float const beta, float* C, int const ldc) 
00103     {
00104       cblas_ssymm (Order, Side, Uplo, M, N, 
00105                    alpha, A, lda, 
00106                    B, ldb,
00107                    beta, C, ldc); 
00108     }
00109   
00110     inline 
00111     void symm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00112                CBLAS_UPLO const Uplo, int const M, int const N, 
00113                double const alpha, double const* A, int const lda,
00114                double const* B, int const ldb, 
00115                double const beta, double* C, int const ldc) 
00116     {
00117       cblas_dsymm (Order, Side, Uplo, M, N, 
00118                    alpha, A, lda, 
00119                    B, ldb,
00120                    beta, C, ldc); 
00121     }
00122   
00123     inline 
00124     void symm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00125                CBLAS_UPLO const Uplo, int const M, int const N, 
00126                traits::complex_f const& alpha, 
00127                traits::complex_f const* A, int const lda,
00128                traits::complex_f const* B, int const ldb, 
00129                traits::complex_f const& beta, 
00130                traits::complex_f* C, int const ldc) 
00131     {
00132       cblas_csymm (Order, Side, Uplo, M, N, 
00133                    static_cast<void const*> (&alpha), 
00134                    static_cast<void const*> (A), lda, 
00135                    static_cast<void const*> (B), ldb,
00136                    static_cast<void const*> (&beta), 
00137                    static_cast<void*> (C), ldc); 
00138     }
00139   
00140     inline 
00141     void symm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00142                CBLAS_UPLO const Uplo, int const M, int const N, 
00143                traits::complex_d const& alpha, 
00144                traits::complex_d const* A, int const lda,
00145                traits::complex_d const* B, int const ldb, 
00146                traits::complex_d const& beta, 
00147                traits::complex_d* C, int const ldc) 
00148     {
00149       cblas_zsymm (Order, Side, Uplo, M, N, 
00150                    static_cast<void const*> (&alpha), 
00151                    static_cast<void const*> (A), lda, 
00152                    static_cast<void const*> (B), ldb,
00153                    static_cast<void const*> (&beta), 
00154                    static_cast<void*> (C), ldc); 
00155     }
00156   
00157 
00158     // C <- alpha * A * B + beta * C 
00159     // C <- alpha * B * A + beta * C 
00160     // A == A^H
00161   
00162     inline 
00163     void hemm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00164                CBLAS_UPLO const Uplo, int const M, int const N, 
00165                traits::complex_f const& alpha, 
00166                traits::complex_f const* A, int const lda,
00167                traits::complex_f const* B, int const ldb, 
00168                traits::complex_f const& beta, 
00169                traits::complex_f* C, int const ldc) 
00170     {
00171       cblas_chemm (Order, Side, Uplo, M, N, 
00172                    static_cast<void const*> (&alpha), 
00173                    static_cast<void const*> (A), lda, 
00174                    static_cast<void const*> (B), ldb,
00175                    static_cast<void const*> (&beta), 
00176                    static_cast<void*> (C), ldc); 
00177     }
00178   
00179     inline 
00180     void hemm (CBLAS_ORDER const Order, CBLAS_SIDE const Side,
00181                CBLAS_UPLO const Uplo, int const M, int const N, 
00182                traits::complex_d const& alpha, 
00183                traits::complex_d const* A, int const lda,
00184                traits::complex_d const* B, int const ldb, 
00185                traits::complex_d const& beta, 
00186                traits::complex_d* C, int const ldc) 
00187     {
00188       cblas_zhemm (Order, Side, Uplo, M, N, 
00189                    static_cast<void const*> (&alpha), 
00190                    static_cast<void const*> (A), lda, 
00191                    static_cast<void const*> (B), ldb,
00192                    static_cast<void const*> (&beta), 
00193                    static_cast<void*> (C), ldc); 
00194     }
00195 
00196 
00197     // C <- alpha * A * A^T + beta * C
00198     // C <- alpha * A^T * A + beta * C
00199     // C == C^T
00200 
00201     inline
00202     void syrk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00203                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00204                float const alpha, float const* A, int const lda,
00205                float const beta, float* C, int const ldc) 
00206     {
00207       cblas_ssyrk (Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc); 
00208     }
00209 
00210     inline
00211     void syrk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00212                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00213                double const alpha, double const* A, int const lda,
00214                double const beta, double* C, int const ldc) 
00215     {
00216       cblas_dsyrk (Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc); 
00217     }
00218 
00219     inline
00220     void syrk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00221                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00222                traits::complex_f const& alpha, 
00223                traits::complex_f const* A, int const lda,
00224                traits::complex_f const& beta, 
00225                traits::complex_f* C, int const ldc) 
00226     {
00227       cblas_csyrk (Order, Uplo, Trans, N, K, 
00228                    static_cast<void const*> (&alpha), 
00229                    static_cast<void const*> (A), lda, 
00230                    static_cast<void const*> (&beta), 
00231                    static_cast<void*> (C), ldc); 
00232     }
00233 
00234     inline
00235     void syrk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00236                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00237                traits::complex_d const& alpha, 
00238                traits::complex_d const* A, int const lda,
00239                traits::complex_d const& beta, 
00240                traits::complex_d* C, int const ldc) 
00241     {
00242       cblas_zsyrk (Order, Uplo, Trans, N, K, 
00243                    static_cast<void const*> (&alpha), 
00244                    static_cast<void const*> (A), lda, 
00245                    static_cast<void const*> (&beta), 
00246                    static_cast<void*> (C), ldc); 
00247     }
00248 
00249 
00250     // C <- alpha * A * B^T + conj(alpha) * B * A^T + beta * C
00251     // C <- alpha * A^T * B + conj(alpha) * B^T * A + beta * C
00252     // C == C^T
00253 
00254     inline
00255     void syr2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00256                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00257                 float const alpha, float const* A, int const lda,
00258                 float const* B, int const ldb,
00259                 float const beta, float* C, int const ldc) 
00260     {
00261       cblas_ssyr2k (Order, Uplo, Trans, N, K, 
00262                     alpha, A, lda, B, ldb, beta, C, ldc); 
00263     }
00264 
00265     inline
00266     void syr2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00267                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00268                 double const alpha, double const* A, int const lda,
00269                 double const* B, int const ldb,
00270                 double const beta, double* C, int const ldc) 
00271     {
00272       cblas_dsyr2k (Order, Uplo, Trans, N, K, 
00273                     alpha, A, lda, B, ldb, beta, C, ldc); 
00274     }
00275 
00276     inline
00277     void syr2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00278                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00279                 traits::complex_f const& alpha, 
00280                 traits::complex_f const* A, int const lda,
00281                 traits::complex_f const* B, int const ldb,
00282                 traits::complex_f const& beta, 
00283                 traits::complex_f* C, int const ldc) 
00284     {
00285       cblas_csyr2k (Order, Uplo, Trans, N, K, 
00286                     static_cast<void const*> (&alpha), 
00287                     static_cast<void const*> (A), lda, 
00288                     static_cast<void const*> (B), ldb, 
00289                     static_cast<void const*> (&beta), 
00290                     static_cast<void*> (C), ldc); 
00291     }
00292 
00293     inline
00294     void syr2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00295                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00296                 traits::complex_d const& alpha, 
00297                 traits::complex_d const* A, int const lda,
00298                 traits::complex_d const* B, int const ldb,
00299                 traits::complex_d const& beta, 
00300                 traits::complex_d* C, int const ldc) 
00301     {
00302       cblas_zsyr2k (Order, Uplo, Trans, N, K, 
00303                     static_cast<void const*> (&alpha), 
00304                     static_cast<void const*> (A), lda, 
00305                     static_cast<void const*> (B), ldb, 
00306                     static_cast<void const*> (&beta), 
00307                     static_cast<void*> (C), ldc); 
00308     }
00309 
00310 
00311     // C <- alpha * A * A^H + beta * C
00312     // C <- alpha * A^H * A + beta * C
00313     // C == C^H
00314 
00315     inline
00316     void herk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00317                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00318                float alpha, traits::complex_f const* A, int const lda,
00319                float beta, traits::complex_f* C, int const ldc) 
00320     {
00321       cblas_cherk (Order, Uplo, Trans, N, K, 
00322                    alpha, static_cast<void const*> (A), lda, 
00323                    beta, static_cast<void*> (C), ldc); 
00324     }
00325 
00326     inline
00327     void herk (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00328                CBLAS_TRANSPOSE const Trans, int const N, int const K,
00329                double alpha, traits::complex_d const* A, int const lda,
00330                double beta, traits::complex_d* C, int const ldc) 
00331     {
00332       cblas_zherk (Order, Uplo, Trans, N, K, 
00333                    alpha, static_cast<void const*> (A), lda, 
00334                    beta, static_cast<void*> (C), ldc); 
00335     }
00336 
00337 
00338     // C <- alpha * A * B^H + conj(alpha) * B * A^H + beta * C
00339     // C <- alpha * A^H * B + conj(alpha) * B^H * A + beta * C
00340     // C == C^H
00341 
00342     inline
00343     void her2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00344                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00345                 traits::complex_f const& alpha, 
00346                 traits::complex_f const* A, int const lda,
00347                 traits::complex_f const* B, int const ldb,
00348                 float beta, traits::complex_f* C, int const ldc) 
00349     {
00350       cblas_cher2k (Order, Uplo, Trans, N, K, 
00351                     static_cast<void const*> (&alpha), 
00352                     static_cast<void const*> (A), lda, 
00353                     static_cast<void const*> (B), ldb, 
00354                     beta, static_cast<void*> (C), ldc); 
00355     }
00356 
00357     inline
00358     void her2k (CBLAS_ORDER const Order, CBLAS_UPLO const Uplo,
00359                 CBLAS_TRANSPOSE const Trans, int const N, int const K,
00360                 traits::complex_d const& alpha, 
00361                 traits::complex_d const* A, int const lda,
00362                 traits::complex_d const* B, int const ldb,
00363                 double beta, traits::complex_d* C, int const ldc) 
00364     {
00365       cblas_zher2k (Order, Uplo, Trans, N, K, 
00366                     static_cast<void const*> (&alpha), 
00367                     static_cast<void const*> (A), lda, 
00368                     static_cast<void const*> (B), ldb, 
00369                     beta, static_cast<void*> (C), ldc); 
00370     }
00371 
00372   
00373   }} // namepaces detail & atlas
00374 
00375 }}} 
00376 
00377 
00378 #endif // BOOST_NUMERIC_BINDINGS_CBLAS3_OVERLOADS_HPP

Generated on Wed Nov 23 19:00:00 2011 for FreeCAD by  doxygen 1.6.1