]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Generalize template specialization
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 18 Feb 2019 13:59:28 +0000 (14:59 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 20 Feb 2019 12:46:19 +0000 (13:46 +0100)
include/deal.II/lac/vector_element_access.h

index 236203e1984061315987fccce95f1cc37cfe285a..9e00b8de58078416f22cb093c0d22387bbe74c94 100644 (file)
@@ -126,147 +126,96 @@ namespace internal
 
 
 #  ifdef DEAL_II_TRILINOS_WITH_TPETRA
-  template <>
-  inline void
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<double>>::add(
-    const double                                   value,
-    const types::global_dof_index                  i,
-    LinearAlgebra::TpetraWrappers::Vector<double> &V)
+  template <typename NumberType>
+  struct ElementAccess<LinearAlgebra::TpetraWrappers::Vector<NumberType>>
   {
-    // Extract local indices in the vector.
-    Tpetra::Vector<double, int, types::global_dof_index> vector =
-      V.trilinos_vector();
-    TrilinosWrappers::types::int_type trilinos_i =
-      vector.getMap()->getLocalElement(
-        static_cast<TrilinosWrappers::types::int_type>(i));
-
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
-    auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
-    // We're going to modify the data on host.
-    vector.modify<Kokkos::HostSpace>();
-    vector_1d(trilinos_i) += value;
-    vector.sync<Tpetra::Vector<double, int, types::global_dof_index>::
-                  device_type::memory_space>();
-  }
-
-
+  public:
+    using VectorType = LinearAlgebra::TpetraWrappers::Vector<NumberType>;
+    static void
+    add(const typename VectorType::value_type value,
+        const types::global_dof_index         i,
+        VectorType &                          V);
 
-  template <>
-  inline void
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<float>>::add(
-    const float                                   value,
-    const types::global_dof_index                 i,
-    LinearAlgebra::TpetraWrappers::Vector<float> &V)
-  {
-    // Extract local indices in the vector.
-    Tpetra::Vector<float, int, types::global_dof_index> vector =
-      V.trilinos_vector();
-    TrilinosWrappers::types::int_type trilinos_i =
-      vector.getMap()->getLocalElement(
-        static_cast<TrilinosWrappers::types::int_type>(i));
+    static void
+    set(typename VectorType::value_type value,
+        const types::global_dof_index   i,
+        VectorType &                    V);
 
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
-    auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
-    // We're going to modify the data on host.
-    vector.modify<Kokkos::HostSpace>();
-    vector_1d(trilinos_i) += value;
-    vector.sync<Tpetra::Vector<float, int, types::global_dof_index>::
-                  device_type::memory_space>();
-  }
+    static typename VectorType::value_type
+    get(const VectorType &V, const types::global_dof_index i);
+  };
 
 
 
-  template <>
+  template <typename NumberType>
   inline void
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<double>>::set(
-    const double                                   value,
-    const types::global_dof_index                  i,
-    LinearAlgebra::TpetraWrappers::Vector<double> &V)
+  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<NumberType>>::add(
+    const typename VectorType::value_type              value,
+    const types::global_dof_index                      i,
+    LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<double, int, types::global_dof_index> vector =
+    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
       V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
 
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
+    vector.template sync<Kokkos::HostSpace>();
+    auto vector_2d = vector.template getLocalView<Kokkos::HostSpace>();
     auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
     // We're going to modify the data on host.
-    vector.modify<Kokkos::HostSpace>();
-    vector_1d(trilinos_i) = value;
-    vector.sync<Tpetra::Vector<double, int, types::global_dof_index>::
-                  device_type::memory_space>();
+    vector.template modify<Kokkos::HostSpace>();
+    vector_1d(trilinos_i) += value;
+    vector.template sync<
+      typename Tpetra::Vector<NumberType, int, types::global_dof_index>::
+        device_type::memory_space>();
   }
 
 
 
-  template <>
+  template <typename NumberType>
   inline void
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<float>>::set(
-    const float                                   value,
-    const types::global_dof_index                 i,
-    LinearAlgebra::TpetraWrappers::Vector<float> &V)
+  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<NumberType>>::set(
+    const typename VectorType::value_type              value,
+    const types::global_dof_index                      i,
+    LinearAlgebra::TpetraWrappers::Vector<NumberType> &V)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<float, int, types::global_dof_index> vector =
+    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
       V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
 
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
+    vector.template sync<Kokkos::HostSpace>();
+    auto vector_2d = vector.template getLocalView<Kokkos::HostSpace>();
     auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
     // We're going to modify the data on host.
-    vector.modify<Kokkos::HostSpace>();
+    vector.template modify<Kokkos::HostSpace>();
     vector_1d(trilinos_i) = value;
-    vector.sync<Tpetra::Vector<float, int, types::global_dof_index>::
-                  device_type::memory_space>();
+    vector.template sync<
+      typename Tpetra::Vector<NumberType, int, types::global_dof_index>::
+        device_type::memory_space>();
   }
 
 
 
-  template <>
-  inline double
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<double>>::get(
-    const LinearAlgebra::TpetraWrappers::Vector<double> &V,
-    const types::global_dof_index                        i)
-  {
-    // Extract local indices in the vector.
-    Tpetra::Vector<double, int, types::global_dof_index> vector =
-      V.trilinos_vector();
-    TrilinosWrappers::types::int_type trilinos_i =
-      vector.getMap()->getLocalElement(
-        static_cast<TrilinosWrappers::types::int_type>(i));
-
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
-    auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
-    // We're going to modify the data on host.
-    return vector_1d(trilinos_i);
-  }
-
-
-
-  template <>
-  inline float
-  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<float>>::get(
-    const LinearAlgebra::TpetraWrappers::Vector<float> &V,
-    const types::global_dof_index                       i)
+  template <typename NumberType>
+  inline typename LinearAlgebra::TpetraWrappers::Vector<NumberType>::value_type
+  ElementAccess<LinearAlgebra::TpetraWrappers::Vector<NumberType>>::get(
+    const LinearAlgebra::TpetraWrappers::Vector<NumberType> &V,
+    const types::global_dof_index                            i)
   {
     // Extract local indices in the vector.
-    Tpetra::Vector<float, int, types::global_dof_index> vector =
+    Tpetra::Vector<NumberType, int, types::global_dof_index> vector =
       V.trilinos_vector();
     TrilinosWrappers::types::int_type trilinos_i =
       vector.getMap()->getLocalElement(
         static_cast<TrilinosWrappers::types::int_type>(i));
 
-    vector.sync<Kokkos::HostSpace>();
-    auto vector_2d = vector.getLocalView<Kokkos::HostSpace>();
+    vector.template sync<Kokkos::HostSpace>();
+    auto vector_2d = vector.template getLocalView<Kokkos::HostSpace>();
     auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0);
     // We're going to modify the data on host.
     return vector_1d(trilinos_i);

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.