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
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
boost::has_trivial_constructor<T>::value &&
boost::has_trivial_destructor<T>::value;
};
+
+ using boost::true_type;
+ using boost::false_type;
}
DEAL_II_NAMESPACE_CLOSE
// ---------------------------------------------------------------------
//
-// 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.
//
#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>
} /* 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
#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>
// 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
// ---------------------------------------------------------------------
//
-// 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.
//
#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>
} // end of namespace distributed
-} // end of namespace parallel
+} // end of namespace LinearAlgebra
+
/**
* Global function which overloads the default implementation of the C++
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
// ---------------------------------------------------------------------
//
-// 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.
//
#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>
}
-
/**
* Global function @p swap which overloads the default implementation of the
* C++ standard library which uses a temporary object. The function simply
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
// ---------------------------------------------------------------------
//
-// 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.
//
#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>
}
} // 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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
} /* 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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
} /* 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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
/**@}*/
+
+/**
+ * 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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
} /* 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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
} /* 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
// ---------------------------------------------------------------------
//
-// 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.
//
# 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"
};
-
// ------------------- inline and template functions --------------
} /* 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
// ---------------------------------------------------------------------
//
-// 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.
//
#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>
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
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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
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();
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
template
Vector<S2>& Vector<S2>::operator=<S1> (const Vector<S1> &);
}
-
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+ };
+}
--- /dev/null
+
+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
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+ };
+}
--- /dev/null
+
+DEAL::1
+DEAL::0
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::0
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// 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;
+ };
+}
--- /dev/null
+
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::1
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::OK
+DEAL::
+DEAL::0
+DEAL::OK