1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii_blas_extension_templates_h
18 #define dealii_blas_extension_templates_h
19 
20 #include <deal.II/base/config.h>
21 
22 #include <deal.II/lac/lapack_support.h>
23 
24 #ifdef DEAL_II_HAVE_FP_EXCEPTIONS
25 #  include <cfenv>
26 #endif
27 
28 // Intel-MKL specific functions
29 #ifdef DEAL_II_LAPACK_WITH_MKL
30 // see
31 // https://software.intel.com/en-us/mkl-windows-developer-guide-using-complex-types-in-c-c
32 #  define MKL_Complex8 std::complex<float>
33 #  define MKL_Complex16 std::complex<double>
34 #  include <mkl_trans.h>
35 #endif
36 
37 
38 DEAL_II_NAMESPACE_OPEN
39 
40 
41 template <typename number1, typename number2, typename number3>
42 inline void
omatcopy(char,char,dealii::types::blas_int,dealii::types::blas_int,const number1,const number2 *,dealii::types::blas_int,number3 *,dealii::types::blas_int)43 omatcopy(char,
44          char,
45          dealii::types::blas_int,
46          dealii::types::blas_int,
47          const number1,
48          const number2 *,
49          dealii::types::blas_int,
50          number3 *,
51          dealii::types::blas_int)
52 {
53   Assert(false, ExcNotImplemented());
54 }
55 
56 
57 
58 inline void
omatcopy(char ordering,char trans,dealii::types::blas_int rows,dealii::types::blas_int cols,const float alpha,const float * A,dealii::types::blas_int lda,float * B,dealii::types::blas_int ldb)59 omatcopy(char                    ordering,
60          char                    trans,
61          dealii::types::blas_int rows,
62          dealii::types::blas_int cols,
63          const float             alpha,
64          const float *           A,
65          dealii::types::blas_int lda,
66          float *                 B,
67          dealii::types::blas_int ldb)
68 {
69 #ifdef DEAL_II_LAPACK_WITH_MKL
70   mkl_somatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
71 #else
72   (void)ordering;
73   (void)trans;
74   (void)rows;
75   (void)cols;
76   (void)alpha;
77   (void)A;
78   (void)lda;
79   (void)B;
80   (void)ldb;
81   Assert(false, LAPACKSupport::ExcMissing("mkl_somatcopy"));
82 #endif // DEAL_II_LAPACK_WITH_MKL
83 }
84 
85 
86 
87 inline void
omatcopy(char ordering,char trans,dealii::types::blas_int rows,dealii::types::blas_int cols,const double alpha,const double * A,dealii::types::blas_int lda,double * B,dealii::types::blas_int ldb)88 omatcopy(char                    ordering,
89          char                    trans,
90          dealii::types::blas_int rows,
91          dealii::types::blas_int cols,
92          const double            alpha,
93          const double *          A,
94          dealii::types::blas_int lda,
95          double *                B,
96          dealii::types::blas_int ldb)
97 {
98 #ifdef DEAL_II_LAPACK_WITH_MKL
99   mkl_domatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
100 #else
101   (void)ordering;
102   (void)trans;
103   (void)rows;
104   (void)cols;
105   (void)alpha;
106   (void)A;
107   (void)lda;
108   (void)B;
109   (void)ldb;
110   Assert(false, LAPACKSupport::ExcMissing("mkl_domatcopy"));
111 #endif // DEAL_II_LAPACK_WITH_MKL
112 }
113 
114 
115 
116 inline void
omatcopy(char ordering,char trans,dealii::types::blas_int rows,dealii::types::blas_int cols,const std::complex<float> alpha,const std::complex<float> * A,dealii::types::blas_int lda,std::complex<float> * B,dealii::types::blas_int ldb)117 omatcopy(char                       ordering,
118          char                       trans,
119          dealii::types::blas_int    rows,
120          dealii::types::blas_int    cols,
121          const std::complex<float>  alpha,
122          const std::complex<float> *A,
123          dealii::types::blas_int    lda,
124          std::complex<float> *      B,
125          dealii::types::blas_int    ldb)
126 {
127 #ifdef DEAL_II_LAPACK_WITH_MKL
128   mkl_comatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
129 #else
130   (void)ordering;
131   (void)trans;
132   (void)rows;
133   (void)cols;
134   (void)alpha;
135   (void)A;
136   (void)lda;
137   (void)B;
138   (void)ldb;
139   Assert(false, LAPACKSupport::ExcMissing("mkl_comatcopy"));
140 #endif // DEAL_II_LAPACK_WITH_MKL
141 }
142 
143 
144 
145 inline void
omatcopy(char ordering,char trans,dealii::types::blas_int rows,dealii::types::blas_int cols,const std::complex<double> alpha,const std::complex<double> * A,dealii::types::blas_int lda,std::complex<double> * B,dealii::types::blas_int ldb)146 omatcopy(char                        ordering,
147          char                        trans,
148          dealii::types::blas_int     rows,
149          dealii::types::blas_int     cols,
150          const std::complex<double>  alpha,
151          const std::complex<double> *A,
152          dealii::types::blas_int     lda,
153          std::complex<double> *      B,
154          dealii::types::blas_int     ldb)
155 {
156 #ifdef DEAL_II_LAPACK_WITH_MKL
157   mkl_zomatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
158 #else
159   (void)ordering;
160   (void)trans;
161   (void)rows;
162   (void)cols;
163   (void)alpha;
164   (void)A;
165   (void)lda;
166   (void)B;
167   (void)ldb;
168   Assert(false, LAPACKSupport::ExcMissing("mkl_zomatcopy"));
169 #endif // DEAL_II_LAPACK_WITH_MKL
170 }
171 
172 DEAL_II_NAMESPACE_CLOSE
173 
174 #endif
175