// access to generic vectors that have operator ().
- template <typename VectorType>
+ // FIXME: this is wrong for Trilinos/Petsc MPI vectors
+ template <typename VectorType,
+ typename std::enable_if<!has_local_element<VectorType>::value,
+ VectorType>::type * = nullptr>
inline typename VectorType::value_type &
vector_access(VectorType &vec, const unsigned int entry)
{
// access to distributed MPI vectors that have a local_element(uint)
// method to access data in local index space, which is what we use in
// DoFInfo and hence in read_dof_values etc.
- template <typename Number>
- inline Number &
- vector_access(LinearAlgebra::distributed::Vector<Number> &vec,
- const unsigned int entry)
+ template <typename VectorType,
+ typename std::enable_if<has_local_element<VectorType>::value,
+ VectorType>::type * = nullptr>
+ inline typename VectorType::value_type &
+ vector_access(VectorType &vec, const unsigned int entry)
{
return vec.local_element(entry);
}
- // this is to make sure that the parallel partitioning in the
- // LinearAlgebra::distributed::Vector is really the same as stored in
- // MatrixFree
+ // this is to make sure that the parallel partitioning in VectorType
+ // is really the same as stored in MatrixFree
+ // FIXME: this is incorrect for PETSc/Trilinos MPI vectors
template <typename VectorType>
inline void
check_vector_compatibility(
AssertDimension(vec.size(), dof_info.vector_partitioner->size());
}
+
+ // this is to make sure that the parallel partitioning in the
+ // LinearAlgebra::distributed::Vector is really the same as stored in
+ // MatrixFree
template <typename Number>
inline void
check_vector_compatibility(
#include "../tests.h"
+// dummy class we use to check typetraits and internal function.
+// this one mimics LA::d::Vec
+template <typename Number>
+class Dummy
+{
+public:
+ using value_type = Number;
+
+ Number
+ local_element(const unsigned int) const
+ {
+ deallog << "Dummy::local_element() const" << std::endl;
+ return Number();
+ }
+
+ Number &
+ local_element(const unsigned int)
+ {
+ deallog << "Dummy::local_element()" << std::endl;
+ return dummy;
+ }
+
+ Number
+ operator()(const unsigned int) const
+ {
+ deallog << "Dummy::operator() const" << std::endl;
+ return Number();
+ }
+
+ Number &
+ operator()(const unsigned int)
+ {
+ deallog << "Dummy::operator()" << std::endl;
+ return dummy;
+ }
+
+private:
+ Number dummy;
+};
+
+
+template <typename Number>
+class Dummy2
+{
+public:
+ using value_type = Number;
+
+ Number
+ operator()(const unsigned int) const
+ {
+ deallog << "Dummy2::operator() const" << std::endl;
+ return Number();
+ }
+
+ Number &
+ operator()(const unsigned int)
+ {
+ deallog << "Dummy2::operator()" << std::endl;
+ return dummy;
+ }
+
+private:
+ Number dummy;
+};
+
int
main()
{
initlog();
+ Dummy<double> dummy;
+ Dummy2<double> dummy2;
+
deallog << "has_local_element:" << std::endl
<< "LinearAlgebra::distributed::Vector = "
<< internal::has_local_element<
<< std::endl
<< "TrilinosWrappers::MPI::Vector = "
<< internal::has_local_element<TrilinosWrappers::MPI::Vector>::value
+ << std::endl
+ << "Dummy = " << internal::has_local_element<Dummy<double>>::value
+ << std::endl
+ << "Dummy2 = " << internal::has_local_element<Dummy2<double>>::value
<< std::endl;
+ // now check internal::vector_access wrapper
+ // we expect local_element() to be called
+ deallog << "internal::vector_access:" << std::endl;
+ internal::vector_access(dummy, 0);
+ internal::vector_access(dummy2, 0);
+
+
deallog << "OK" << std::endl;
}