]> https://gitweb.dealii.org/ - dealii.git/commitdiff
move new blas extension template functions into a separate header 7502/head
authorDenis Davydov <davydden@gmail.com>
Thu, 6 Dec 2018 21:43:09 +0000 (22:43 +0100)
committerDenis Davydov <davydden@gmail.com>
Thu, 6 Dec 2018 22:43:19 +0000 (23:43 +0100)
include/deal.II/lac/blas_extension_templates.h [new file with mode: 0644]
include/deal.II/lac/lapack_templates.h
source/lac/lapack_full_matrix.cc

diff --git a/include/deal.II/lac/blas_extension_templates.h b/include/deal.II/lac/blas_extension_templates.h
new file mode 100644 (file)
index 0000000..a6b19cc
--- /dev/null
@@ -0,0 +1,174 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+#ifndef dealii_blas_extension_templates_h
+#define dealii_blas_extension_templates_h
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/lac/lapack_support.h>
+
+#ifdef DEAL_II_HAVE_FP_EXCEPTIONS
+#  include <cfenv>
+#endif
+
+// Intel-MKL specific functions
+#ifdef DEAL_II_LAPACK_WITH_MKL
+// see
+// https://software.intel.com/en-us/mkl-windows-developer-guide-using-complex-types-in-c-c
+#  define MKL_Complex8 std::complex<float>
+#  define MKL_Complex16 std::complex<double>
+#  include <mkl_trans.h>
+#endif
+
+
+DEAL_II_NAMESPACE_OPEN
+
+
+template <typename number1, typename number2, typename number3>
+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)
+{
+  Assert(false, ExcNotImplemented());
+}
+
+
+
+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)
+{
+#ifdef DEAL_II_LAPACK_WITH_MKL
+  mkl_somatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
+#else
+  (void)ordering;
+  (void)trans;
+  (void)rows;
+  (void)cols;
+  (void)alpha;
+  (void)A;
+  (void)lda;
+  (void)B;
+  (void)ldb;
+  Assert(false, LAPACKSupport::ExcMissing("mkl_somatcopy"));
+#endif // DEAL_II_LAPACK_WITH_MKL
+}
+
+
+
+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)
+{
+#ifdef DEAL_II_LAPACK_WITH_MKL
+  mkl_domatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
+#else
+  (void)ordering;
+  (void)trans;
+  (void)rows;
+  (void)cols;
+  (void)alpha;
+  (void)A;
+  (void)lda;
+  (void)B;
+  (void)ldb;
+  Assert(false, LAPACKSupport::ExcMissing("mkl_domatcopy"));
+#endif // DEAL_II_LAPACK_WITH_MKL
+}
+
+
+
+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)
+{
+#ifdef DEAL_II_LAPACK_WITH_MKL
+  mkl_comatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
+#else
+  (void)ordering;
+  (void)trans;
+  (void)rows;
+  (void)cols;
+  (void)alpha;
+  (void)A;
+  (void)lda;
+  (void)B;
+  (void)ldb;
+  Assert(false, LAPACKSupport::ExcMissing("mkl_comatcopy"));
+#endif // DEAL_II_LAPACK_WITH_MKL
+}
+
+
+
+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)
+{
+#ifdef DEAL_II_LAPACK_WITH_MKL
+  mkl_zomatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
+#else
+  (void)ordering;
+  (void)trans;
+  (void)rows;
+  (void)cols;
+  (void)alpha;
+  (void)A;
+  (void)lda;
+  (void)B;
+  (void)ldb;
+  Assert(false, LAPACKSupport::ExcMissing("mkl_zomatcopy"));
+#endif // DEAL_II_LAPACK_WITH_MKL
+}
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 62d5ebd86cf8d6257a27ec57729dc06ec1c38798..fde6c6ed31ede478bec6df2a7e889c588a76e89c 100644 (file)
 #  include <cfenv>
 #endif
 
-// Intel-MKL specific functions
-#ifdef DEAL_II_LAPACK_WITH_MKL
-// see
-// https://software.intel.com/en-us/mkl-windows-developer-guide-using-complex-types-in-c-c
-#  define MKL_Complex8 std::complex<float>
-#  define MKL_Complex16 std::complex<double>
-#  include "mkl_trans.h"
-#endif
-
 extern "C"
 {
   void
@@ -1418,138 +1409,6 @@ extern "C"
 DEAL_II_NAMESPACE_OPEN
 
 
-template <typename number1, typename number2, typename number3>
-inline void
-mkl_omatcopy(char,
-             char,
-             dealii::types::blas_int,
-             dealii::types::blas_int,
-             const number1,
-             const number2 *,
-             dealii::types::blas_int,
-             number3 *,
-             dealii::types::blas_int)
-{
-  Assert(false, ExcNotImplemented());
-}
-
-
-
-inline void
-mkl_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)
-{
-#ifdef DEAL_II_LAPACK_WITH_MKL
-  mkl_somatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
-#else
-  (void)ordering;
-  (void)trans;
-  (void)rows;
-  (void)cols;
-  (void)alpha;
-  (void)A;
-  (void)lda;
-  (void)B;
-  (void)ldb;
-  Assert(false, LAPACKSupport::ExcMissing("mkl_somatcopy"));
-#endif // DEAL_II_LAPACK_WITH_MKL
-}
-
-
-
-inline void
-mkl_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)
-{
-#ifdef DEAL_II_LAPACK_WITH_MKL
-  mkl_domatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
-#else
-  (void)ordering;
-  (void)trans;
-  (void)rows;
-  (void)cols;
-  (void)alpha;
-  (void)A;
-  (void)lda;
-  (void)B;
-  (void)ldb;
-  Assert(false, LAPACKSupport::ExcMissing("mkl_domatcopy"));
-#endif // DEAL_II_LAPACK_WITH_MKL
-}
-
-
-
-inline void
-mkl_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)
-{
-#ifdef DEAL_II_LAPACK_WITH_MKL
-  mkl_comatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
-#else
-  (void)ordering;
-  (void)trans;
-  (void)rows;
-  (void)cols;
-  (void)alpha;
-  (void)A;
-  (void)lda;
-  (void)B;
-  (void)ldb;
-  Assert(false, LAPACKSupport::ExcMissing("mkl_comatcopy"));
-#endif // DEAL_II_LAPACK_WITH_MKL
-}
-
-
-
-inline void
-mkl_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)
-{
-#ifdef DEAL_II_LAPACK_WITH_MKL
-  mkl_zomatcopy(ordering, trans, rows, cols, alpha, A, lda, B, ldb);
-#else
-  (void)ordering;
-  (void)trans;
-  (void)rows;
-  (void)cols;
-  (void)alpha;
-  (void)A;
-  (void)lda;
-  (void)B;
-  (void)ldb;
-  Assert(false, LAPACKSupport::ExcMissing("mkl_zomatcopy"));
-#endif // DEAL_II_LAPACK_WITH_MKL
-}
-
-
 
 template <typename number1, typename number2, typename number3>
 inline void
index 837191d34d99f19ab84d04157478d575e5884a84..90c303a12b3d91ec154063dbd93b2c8dea9851f3 100644 (file)
@@ -15,6 +15,7 @@
 
 #include <deal.II/base/numbers.h>
 
+#include <deal.II/lac/blas_extension_templates.h>
 #include <deal.II/lac/block_vector.h>
 #include <deal.II/lac/full_matrix.h>
 #include <deal.II/lac/lapack_full_matrix.h>
@@ -1042,7 +1043,7 @@ LAPACKFullMatrix<number>::transpose(LAPACKFullMatrix<number> &B) const
   const types::blas_int n = B.n();
 #ifdef DEAL_II_LAPACK_WITH_MKL
   const number one = 1.;
-  mkl_omatcopy('C', 'C', n, m, one, &A.values[0], n, &B.values[0], m);
+  omatcopy('C', 'C', n, m, one, &A.values[0], n, &B.values[0], m);
 #else
   for (types::blas_int i = 0; i < m; ++i)
     for (types::blas_int j = 0; j < n; ++j)

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.