]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce type trait is_serial_vector<VectorType> 3784/head
authorUwe Köcher <koecher@hsu-hamburg.de>
Thu, 12 Jan 2017 15:06:33 +0000 (16:06 +0100)
committerUwe Köcher <koecher@hsu-hamburg.de>
Thu, 26 Jan 2017 19:20:59 +0000 (20:20 +0100)
renames type trait is_non_distributed_vector to is_serial_vector and put it into dealii::std_cxx11 namespace

intendation with astyle 2.04 / script

more obivious logical statements as requested

moves true_type and false_type to dealii/base/std_cxx11/type_traits.h

moves specialization of is_serial_vector for dealii::Vector classes to their implementation

correction of implementation

moves vector type traits from a single file to their declaration files

testsuite for vector_type_traits.h for current vectors and instanciations

change request for comments

declares the is_serial_vector template without defining it and minor work on documentation

squashes test cases for is_serial_vector and marks output files for trilinos, petsc, mpi

corrects output files for testsuite (results are correct)

corrects intendation

intendation of #include

23 files changed:
include/deal.II/base/std_cxx11/type_traits.h
include/deal.II/lac/block_vector.h
include/deal.II/lac/constraint_matrix.templates.h
include/deal.II/lac/la_parallel_block_vector.h
include/deal.II/lac/la_parallel_vector.h
include/deal.II/lac/la_vector.h
include/deal.II/lac/petsc_block_vector.h
include/deal.II/lac/petsc_parallel_block_vector.h
include/deal.II/lac/petsc_parallel_vector.h
include/deal.II/lac/petsc_vector.h
include/deal.II/lac/trilinos_block_vector.h
include/deal.II/lac/trilinos_vector.h
include/deal.II/lac/vector.h
include/deal.II/lac/vector_type_traits.h [new file with mode: 0644]
include/deal.II/numerics/vector_tools.templates.h
source/lac/vector.cc
source/lac/vector.inst.in
tests/lac/vector_type_traits_is_serial_01.cc [new file with mode: 0644]
tests/lac/vector_type_traits_is_serial_01.output [new file with mode: 0644]
tests/lac/vector_type_traits_is_serial_02.cc [new file with mode: 0644]
tests/lac/vector_type_traits_is_serial_02.with_trilinos=true.with_mpi=true.output [new file with mode: 0644]
tests/lac/vector_type_traits_is_serial_03.cc [new file with mode: 0644]
tests/lac/vector_type_traits_is_serial_03.with_petsc=true.with_mpi=true.output [new file with mode: 0644]

index d9692f7de6d77e430188bd55593e511399105a0e..a6f9a898788dd3234a1d50cb0baf4c5cce6a8fc8 100644 (file)
@@ -33,6 +33,8 @@ namespace std_cxx11
   using std::is_standard_layout;
   using std::is_trivial;
   using std::enable_if;
+  using std::true_type;
+  using std::false_type;
 }
 DEAL_II_NAMESPACE_CLOSE
 
@@ -77,6 +79,9 @@ namespace std_cxx11
                               boost::has_trivial_constructor<T>::value &&
                               boost::has_trivial_destructor<T>::value;
   };
+
+  using boost::true_type;
+  using boost::false_type;
 }
 DEAL_II_NAMESPACE_CLOSE
 
index 0186b56cef644c8aad8a7f0a6a1e176694bf48d2..e981223f3b9444d4b195bd325e7280ec60374e2d 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1999 - 2016 by the deal.II authors
+// Copyright (C) 1999 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -21,6 +21,7 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/lac/block_indices.h>
 #include <deal.II/lac/block_vector_base.h>
+#include <deal.II/lac/vector_type_traits.h>
 
 #include <cstdio>
 #include <vector>
@@ -517,6 +518,17 @@ namespace internal
   } /* namespace LinearOperator */
 } /* namespace internal */
 
+
+/**
+ * Declare dealii::BlockVector< Number > as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <typename Number>
+struct is_serial_vector< BlockVector< Number > > : std_cxx11::true_type
+{
+};
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index cadb380cc0a0cf114fa55352acf96021a41f980d..54949aa030bdfaf77945feb8b872d17947e038e5 100644 (file)
@@ -27,6 +27,7 @@
 #include <deal.II/lac/sparse_matrix.h>
 #include <deal.II/lac/block_sparsity_pattern.h>
 #include <deal.II/lac/block_sparse_matrix.h>
+
 #include <deal.II/lac/la_vector.h>
 #include <deal.II/lac/la_parallel_vector.h>
 #include <deal.II/lac/la_parallel_block_vector.h>
@@ -855,16 +856,7 @@ ConstraintMatrix::distribute (VectorType &vec) const
   // the last else is for the simple case (sequential vector)
   const IndexSet vec_owned_elements = vec.locally_owned_elements();
 
-  // FIXME: This has to be refactored into a typetrait
-  if ((typeid(vec) != typeid(Vector<double>)) &&
-      (typeid(vec) != typeid(Vector<float>)) &&
-      (typeid(vec) != typeid(Vector<std::complex<double> >)) &&
-      (typeid(vec) != typeid(BlockVector<double>)) &&
-      (typeid(vec) != typeid(BlockVector<float>)) &&
-      (typeid(vec) != typeid(BlockVector<std::complex<double> >)) &&
-      (typeid(vec) != typeid(LinearAlgebra::Vector<double>)) &&
-      (typeid(vec) != typeid(LinearAlgebra::Vector<float>)) &&
-      (typeid(vec) != typeid(LinearAlgebra::Vector<std::complex<double> >)))
+  if ( dealii::is_serial_vector< VectorType >::value == false )
     {
       // This processor owns only part of the vector. one may think that
       // every processor should be able to simply communicate those elements
index bbfeae8134e2337c579dd320c18f2b72f2719a52..04a3a786ceb80ecc15d576e5a453154bdcf20e6f 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1999 - 2016 by the deal.II authors
+// Copyright (C) 1999 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -22,6 +22,7 @@
 #include <deal.II/lac/block_indices.h>
 #include <deal.II/lac/block_vector_base.h>
 #include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/vector_type_traits.h>
 
 
 #include <cstdio>
@@ -578,7 +579,8 @@ namespace LinearAlgebra
 
   } // end of namespace distributed
 
-} // end of namespace parallel
+} // end of namespace LinearAlgebra
+
 
 /**
  * Global function which overloads the default implementation of the C++
@@ -596,6 +598,18 @@ void swap (LinearAlgebra::distributed::BlockVector<Number> &u,
   u.swap (v);
 }
 
+
+/**
+ * Declare dealii::LinearAlgebra::BlockVector< Number > as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <typename Number>
+struct is_serial_vector< LinearAlgebra::distributed::BlockVector< Number > > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #ifdef DEAL_II_MSVC
index b73fd4ea7ec61e588bebe01f381785ec96a14c67..95a02f9b804febc43b094a3d033755cac6673415 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2011 - 2016 by the deal.II authors
+// Copyright (C) 2011 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -23,6 +23,7 @@
 #include <deal.II/base/std_cxx11/shared_ptr.h>
 #include <deal.II/base/thread_management.h>
 #include <deal.II/lac/vector_space_vector.h>
+#include <deal.II/lac/vector_type_traits.h>
 
 #include <iomanip>
 
@@ -1449,7 +1450,6 @@ namespace LinearAlgebra
 }
 
 
-
 /**
  * Global function @p swap which overloads the default implementation of the
  * C++ standard library which uses a temporary object. The function simply
@@ -1466,6 +1466,18 @@ void swap (LinearAlgebra::distributed::Vector<Number> &u,
   u.swap (v);
 }
 
+
+/**
+ * Declare dealii::LinearAlgebra::Vector< Number > as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <typename Number>
+struct is_serial_vector< LinearAlgebra::distributed::Vector< Number > > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index 3894119dfd9a556afb3e0a4ed12f58c920768ff9..18dbcf5d0227952109f9322dc80ec8d3ac9ac30d 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C)  2015 by the deal.II authors
+// Copyright (C) 2015 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -23,6 +23,7 @@
 #include <deal.II/base/index_set.h>
 #include <deal.II/lac/read_write_vector.h>
 #include <deal.II/lac/vector_space_vector.h>
+#include <deal.II/lac/vector_type_traits.h>
 #include <boost/serialization/array.hpp>
 #include <boost/serialization/split_member.hpp>
 
@@ -376,6 +377,18 @@ namespace LinearAlgebra
   }
 } // end of namespace LinearAlgebra
 
+
+/**
+ * Declare dealii::LinearAlgebra::Vector< Number > as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <typename Number>
+struct is_serial_vector< LinearAlgebra::Vector< Number > > : std_cxx11::true_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #ifdef DEAL_II_MSVC
index c95b456cd1c4931a4209b27a1b1529cfa4e548ed..c71f4a94b219582273059a38a45e11e48b42ed75 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2004 - 2016 by the deal.II authors
+// Copyright (C) 2004 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -26,6 +26,7 @@
 #  include <deal.II/lac/block_indices.h>
 #  include <deal.II/lac/block_vector_base.h>
 #  include <deal.II/lac/exceptions.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -495,6 +496,18 @@ namespace internal
   } /* namespace LinearOperator */
 } /* namespace internal */
 
+
+/**
+ * Declare dealii::PETScWrappers::BlockVector as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< PETScWrappers::BlockVector > : std_cxx11::true_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif  // DEAL_II_WITH_PETSC
index 03e5940c3322c7bbb0555d129f6015773993e9b9..66779bb586939d84853d163fc1874d34c740c664 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2004 - 2016 by the deal.II authors
+// Copyright (C) 2004 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -25,6 +25,7 @@
 #  include <deal.II/lac/block_indices.h>
 #  include <deal.II/lac/block_vector_base.h>
 #  include <deal.II/lac/exceptions.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -574,6 +575,18 @@ namespace internal
   } /* namespace LinearOperator */
 } /* namespace internal */
 
+
+/**
+ * Declare dealii::PETScWrappers::MPI::BlockVector as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< PETScWrappers::MPI::BlockVector > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif  // DEAL_II_WITH_PETSC
index 341511db40f3253335c922873c5b6d3bffbe404d..df36bb05d793d8e282241b32aa868e766ab155de 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2004 - 2016 by the deal.II authors
+// Copyright (C) 2004 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -26,6 +26,7 @@
 #  include <deal.II/lac/vector.h>
 #  include <deal.II/lac/petsc_vector_base.h>
 #  include <deal.II/base/index_set.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -619,6 +620,18 @@ namespace internal
 
 /**@}*/
 
+
+/**
+ * Declare dealii::PETScWrappers::MPI::Vector as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< PETScWrappers::MPI::Vector > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif // DEAL_II_WITH_PETSC
index 34287adf78d8f51bcce6cac690d963c8254d3178..c17688a56aee59653af735f59f648efa87f060d0 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2004 - 2016 by the deal.II authors
+// Copyright (C) 2004 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -26,6 +26,7 @@
 #  include <deal.II/lac/petsc_vector_base.h>
 #  include <deal.II/lac/petsc_parallel_vector.h>
 #  include <deal.II/lac/vector.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -429,6 +430,18 @@ namespace internal
   } /* namespace LinearOperator */
 } /* namespace internal */
 
+
+/**
+ * Declare dealii::PETScWrappers::Vector as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< PETScWrappers::Vector > : std_cxx11::true_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif // DEAL_II_WITH_PETSC
index 965bb147515fb777d37663dbe86c92d470533c6e..c3742e4a2dd3446733adc6dc905c9a7e273ba406 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2008 - 2016 by the deal.II authors
+// Copyright (C) 2008 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -26,6 +26,7 @@
 #  include <deal.II/lac/block_indices.h>
 #  include <deal.II/lac/block_vector_base.h>
 #  include <deal.II/lac/exceptions.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -494,6 +495,28 @@ namespace internal
 } /* namespace internal */
 
 
+/**
+ * Declare dealii::TrilinosWrappers::BlockVector as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< TrilinosWrappers::BlockVector > : std_cxx11::true_type
+{
+};
+
+
+/**
+ * Declare dealii::TrilinosWrappers::MPI::BlockVector as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< TrilinosWrappers::MPI::BlockVector > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif  // DEAL_II_WITH_TRILINOS
index 71f130ad18f9a147d67b5d4213b66c8861909c6f..dff8b6c9d635e0fc710255deb5a93c348664b1a4 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2008 - 2016 by the deal.II authors
+// Copyright (C) 2008 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -28,6 +28,7 @@
 #  include <deal.II/lac/exceptions.h>
 #  include <deal.II/lac/vector.h>
 #  include <deal.II/lac/trilinos_vector_base.h>
+#  include <deal.II/lac/vector_type_traits.h>
 
 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
 #  include "Epetra_Map.h"
@@ -928,7 +929,6 @@ namespace TrilinosWrappers
   };
 
 
-
 // ------------------- inline and template functions --------------
 
 
@@ -1082,6 +1082,28 @@ namespace internal
 } /* namespace internal */
 
 
+/**
+ * Declare dealii::TrilinosWrappers::Vector as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< TrilinosWrappers::Vector > : std_cxx11::true_type
+{
+};
+
+
+/**
+ * Declare dealii::TrilinosWrappers::MPI::Vector as distributed vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <>
+struct is_serial_vector< TrilinosWrappers::MPI::Vector > : std_cxx11::false_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif // DEAL_II_WITH_TRILINOS
index c6af741f76b3914fa342397d4a788e0c7da845aa..5b723715613069d7e9a086ce94a16916dc8c3060 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1999 - 2016 by the deal.II authors
+// Copyright (C) 1999 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -22,6 +22,7 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/subscriptor.h>
 #include <deal.II/base/index_set.h>
+#include <deal.II/lac/vector_type_traits.h>
 #include <boost/serialization/array.hpp>
 #include <boost/serialization/split_member.hpp>
 
@@ -1425,9 +1426,20 @@ operator << (LogStream &os, const Vector<number> &v)
   return os;
 }
 
-
 /*@}*/
 
+
+/**
+ * Declare dealii::Vector< Number > as serial vector.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template <typename Number>
+struct is_serial_vector< Vector<Number> > : std_cxx11::true_type
+{
+};
+
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
diff --git a/include/deal.II/lac/vector_type_traits.h b/include/deal.II/lac/vector_type_traits.h
new file mode 100644 (file)
index 0000000..db8a45f
--- /dev/null
@@ -0,0 +1,50 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii__vector_type_traits_h
+#define dealii__vector_type_traits_h
+
+#include <deal.II/base/std_cxx11/type_traits.h>
+
+
+DEAL_II_NAMESPACE_OPEN
+
+
+/**
+ * Type trait for a serial vector, i.e. a vector class for which storage is not
+ * supported to be distributed over processes.
+ *
+ * The specialization
+ * @code
+ *   template <>
+ *   struct is_serial_vector< VectorType > : std_cxx11::true_type {};
+ * @endcode
+ * for a serial vector type, respectively,
+ * @code
+ *   template <>
+ *   struct is_serial_vector< VectorType > : std_cxx11::false_type {};
+ * @endcode
+ * for a vector type with support of distributed storage,
+ * must be done in a header file of a vector declaration.
+ *
+ * @author Uwe Koecher, 2017
+ */
+template<typename T>
+struct is_serial_vector;
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 0b64ae5dc32bf94aa75662895e9c7f56871f6f62..9c311ddf0ee3ae1bec8aefe4ca132508f78bdab3 100644 (file)
@@ -7782,12 +7782,7 @@ namespace VectorTools
     else
       {
         // This function is not implemented for distributed vectors.
-        Assert(typeid(v) == typeid(Vector<double>) ||
-               typeid(v) == typeid(Vector<float>) ||
-               typeid(v) == typeid(BlockVector<double>) ||
-               typeid(v) == typeid(BlockVector<float>) ||
-               typeid(v) == typeid(LinearAlgebra::Vector<double>) ||
-               typeid(v) == typeid(LinearAlgebra::Vector<float>),
+        Assert((dealii::is_serial_vector< VectorType >::value == true),
                ExcNotImplemented());
 
         const unsigned int n = v.size();
index e7f0a283ed28d4c2b7c83b8e1ac09182f2275665..4f24cdbbc30c820769e4deb0e9fd8128303ac1ea 100644 (file)
@@ -30,7 +30,6 @@ void Vector<int>::reinit<double>(const Vector<double> &, const bool);
 template class Vector<long double>;
 template long double Vector<long double>::operator *<long double>(const Vector<long double> &) const;
 
-
 // do a few functions that currently don't fit the scheme because they have
 // two template arguments that need to be different (the case of same
 // arguments is covered by the default copy constructor and copy operator that
index e9c29f58562ffc69705c3106734a3db3caf21fa0..7965594f036a001f1e3ff3c76e74b59de09baa4c 100644 (file)
@@ -56,4 +56,3 @@ for (S1: REAL_SCALARS; S2: COMPLEX_SCALARS)
     template
     Vector<S2>& Vector<S2>::operator=<S1> (const Vector<S1> &);
 }
-
diff --git a/tests/lac/vector_type_traits_is_serial_01.cc b/tests/lac/vector_type_traits_is_serial_01.cc
new file mode 100644 (file)
index 0000000..a868652
--- /dev/null
@@ -0,0 +1,128 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check is_serial_vector type trait for internal deal.II Vector types
+
+#include "../tests.h"
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/block_vector.h>
+#include <deal.II/lac/la_vector.h>
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/la_parallel_block_vector.h>
+#include <fstream>
+#include <iomanip>
+
+void test ()
+{
+  // make sure that is_serial_vector< dealii::Vector<Number> > is working
+  Assert (is_serial_vector< dealii::Vector<double> >::value == true,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::Vector<float> >::value << std::endl;
+  deallog << is_serial_vector< dealii::Vector<double> >::value << std::endl;
+  deallog << is_serial_vector< dealii::Vector<long double> >::value << std::endl;
+  deallog << is_serial_vector< dealii::Vector<int> >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::Vector< std::complex<float> > >::value << std::endl;
+  deallog << is_serial_vector< dealii::Vector< std::complex<double> > >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that is_serial_vector< dealii::BlockVector<Number> > is working
+  Assert (is_serial_vector< dealii::BlockVector<double> >::value == true,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::BlockVector<float> >::value << std::endl;
+  deallog << is_serial_vector< dealii::BlockVector<double> >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::BlockVector< std::complex<float> > >::value << std::endl;
+  deallog << is_serial_vector< dealii::BlockVector< std::complex<double> > >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+
+  // make sure that is_serial_vector< dealii::LinearAlgebra::Vector<Number> > is working
+  Assert (is_serial_vector< dealii::LinearAlgebra::Vector<double> >::value == true,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::Vector<float> >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::Vector<double> >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::Vector< std::complex<float> > >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::Vector< std::complex<double> > >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that dealii::LinearAlgebra::distributed::Vector<Number> > is working
+  Assert (is_serial_vector< dealii::LinearAlgebra::distributed::Vector<double> >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::Vector<float> >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::Vector<double> >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::Vector< std::complex<float> > >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::Vector< std::complex<double> > >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that dealii::LinearAlgebra::distributed::BlockVector<Number> > is working
+  Assert (is_serial_vector< dealii::LinearAlgebra::distributed::BlockVector<double> >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::BlockVector<float> >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::BlockVector<double> >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::BlockVector< std::complex<float> > >::value << std::endl;
+  deallog << is_serial_vector< dealii::LinearAlgebra::distributed::BlockVector< std::complex<double> > >::value << std::endl;
+
+  deallog << "OK" << std::endl;
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+
+  try
+    {
+      test ();
+    }
+  catch (std::exception &exc)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Exception on processing: " << std::endl
+              << exc.what() << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+
+      return 1;
+    }
+  catch (...)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Unknown exception!" << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      return 1;
+    };
+}
diff --git a/tests/lac/vector_type_traits_is_serial_01.output b/tests/lac/vector_type_traits_is_serial_01.output
new file mode 100644 (file)
index 0000000..9aaeb30
--- /dev/null
@@ -0,0 +1,32 @@
+
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::0
+DEAL::0
+DEAL::0
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::0
+DEAL::0
+DEAL::0
+DEAL::OK
diff --git a/tests/lac/vector_type_traits_is_serial_02.cc b/tests/lac/vector_type_traits_is_serial_02.cc
new file mode 100644 (file)
index 0000000..e263ef9
--- /dev/null
@@ -0,0 +1,89 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check is_serial_vector type trait
+
+#include "../tests.h"
+#include <deal.II/lac/trilinos_vector.h>
+#include <deal.II/lac/trilinos_block_vector.h>
+#include <fstream>
+#include <iomanip>
+
+void test ()
+{
+  // make sure that is_serial_vector< dealii::TrilinosWrappers::Vector > is working
+  Assert (is_serial_vector< dealii::TrilinosWrappers::Vector >::value == true,
+          ExcInternalError());
+
+  Assert (is_serial_vector< dealii::TrilinosWrappers::MPI::Vector >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::TrilinosWrappers::Vector >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::TrilinosWrappers::MPI::Vector >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+
+  // make sure that dealii::TrilinosWrappers::BlockVector > is working
+  Assert (is_serial_vector< dealii::TrilinosWrappers::BlockVector >::value == true,
+          ExcInternalError());
+
+  Assert (is_serial_vector< dealii::TrilinosWrappers::MPI::BlockVector >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::TrilinosWrappers::BlockVector >::value << std::endl;
+
+  deallog << is_serial_vector< dealii::TrilinosWrappers::MPI::BlockVector >::value << std::endl;
+
+  deallog << "OK" << std::endl;
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+
+  try
+    {
+      test ();
+    }
+  catch (std::exception &exc)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Exception on processing: " << std::endl
+              << exc.what() << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+
+      return 1;
+    }
+  catch (...)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Unknown exception!" << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      return 1;
+    };
+}
diff --git a/tests/lac/vector_type_traits_is_serial_02.with_trilinos=true.with_mpi=true.output b/tests/lac/vector_type_traits_is_serial_02.with_trilinos=true.with_mpi=true.output
new file mode 100644 (file)
index 0000000..9536f3b
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::1
+DEAL::0
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::0
+DEAL::OK
diff --git a/tests/lac/vector_type_traits_is_serial_03.cc b/tests/lac/vector_type_traits_is_serial_03.cc
new file mode 100644 (file)
index 0000000..5b6536a
--- /dev/null
@@ -0,0 +1,96 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check is_serial_vector type trait
+
+#include "../tests.h"
+#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/petsc_block_vector.h>
+#include <deal.II/lac/petsc_parallel_vector.h>
+#include <deal.II/lac/petsc_parallel_block_vector.h>
+#include <fstream>
+#include <iomanip>
+
+void test ()
+{
+  // make sure that is_serial_vector< dealii::PETScWrappers::Vector > is working
+  Assert (is_serial_vector< dealii::PETScWrappers::Vector >::value == true,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::PETScWrappers::Vector >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that is_serial_vector< dealii::PETScWrappers::BlockVector > is working
+  Assert (is_serial_vector< dealii::PETScWrappers::BlockVector >::value == true,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::PETScWrappers::BlockVector >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that is_serial_vector< dealii::PETScWrappers::MPI::Vector > is working
+  Assert (is_serial_vector< dealii::PETScWrappers::MPI::Vector >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::PETScWrappers::MPI::Vector >::value << std::endl;
+
+  deallog << "OK" << std::endl << std::endl;
+
+  // make sure that is_serial_vector< dealii::PETScWrappers::MPI::BlockVector > is working
+  Assert (is_serial_vector< dealii::PETScWrappers::MPI::BlockVector >::value == false,
+          ExcInternalError());
+
+  deallog << is_serial_vector< dealii::PETScWrappers::MPI::BlockVector >::value << std::endl;
+
+  deallog << "OK" << std::endl;
+}
+
+int main ()
+{
+  std::ofstream logfile("output");
+  deallog.attach(logfile);
+
+  try
+    {
+      test ();
+    }
+  catch (std::exception &exc)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Exception on processing: " << std::endl
+              << exc.what() << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+
+      return 1;
+    }
+  catch (...)
+    {
+      deallog << std::endl << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      deallog << "Unknown exception!" << std::endl
+              << "Aborting!" << std::endl
+              << "----------------------------------------------------"
+              << std::endl;
+      return 1;
+    };
+}
diff --git a/tests/lac/vector_type_traits_is_serial_03.with_petsc=true.with_mpi=true.output b/tests/lac/vector_type_traits_is_serial_03.with_petsc=true.with_mpi=true.output
new file mode 100644 (file)
index 0000000..dee23be
--- /dev/null
@@ -0,0 +1,12 @@
+
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::OK

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.