]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix according to review
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sat, 11 Jun 2016 20:24:52 +0000 (22:24 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Tue, 14 Jun 2016 15:10:12 +0000 (17:10 +0200)
cmake/config/template-arguments.in
include/deal.II/lac/la_parallel_vector.h
include/deal.II/lac/la_parallel_vector.templates.h
include/deal.II/lac/parallel_vector.h
include/deal.II/lac/precondition.h
include/deal.II/lac/read_write_vector.h
include/deal.II/lac/read_write_vector.templates.h
include/deal.II/matrix_free/fe_evaluation.h
source/lac/la_parallel_vector.cc
source/lac/la_parallel_vector.inst.in

index adf41b2ab8217e0860e1cb2c8c84f8f859f27d79..ec553ea34de8085f1b8194b1070319be4cada08a 100644 (file)
@@ -23,8 +23,8 @@ SERIAL_VECTORS := { Vector<double>;
                     BlockVector<double>;
                     BlockVector<float>;
 
-                    parallel::distributed::Vector<double>;
-                    parallel::distributed::Vector<float> ;
+                    LinearAlgebra::distributed::Vector<double>;
+                    LinearAlgebra::distributed::Vector<float> ;
 
                     parallel::distributed::BlockVector<double>;
                     parallel::distributed::BlockVector<float> ;
index e9bf996890be4caadafc93d3da099210310841ae..2cc4e0999dab049497080f805d54dcfd577dde02 100644 (file)
@@ -39,9 +39,29 @@ namespace LinearAlgebra
   template <typename> class ReadWriteVector;
 }
 
+#ifdef DEAL_II_WITH_PETSC
+namespace PETScWrappers
+{
+  namespace MPI
+  {
+    class Vector;
+  }
+}
+#endif
+
+#ifdef DEAL_II_WITH_TRILINOS
+namespace TrilinosWrappers
+{
+  namespace MPI
+  {
+    class Vector;
+  }
+}
+#endif
+
 namespace LinearAlgebra
 {
-  namespace parallel
+  namespace distributed
   {
     /*! @addtogroup Vectors
      *@{
@@ -335,6 +355,37 @@ namespace LinearAlgebra
       template <typename Number2>
       Vector<Number> &
       operator = (const Vector<Number2> &in_vector);
+
+#ifdef DEAL_II_WITH_PETSC
+      /**
+       * Copy the content of a PETSc vector into the calling vector. This
+       * function assumes that the vectors layouts have already been
+       * initialized to match.
+       *
+       * This operator is only available if deal.II was configured with PETSc.
+       *
+       * This function is deprecated. Use the interface through
+       * ReadWriteVector instead.
+       */
+      Vector<Number> &
+      operator = (const PETScWrappers::MPI::Vector &petsc_vec) DEAL_II_DEPRECATED;
+#endif
+
+#ifdef DEAL_II_WITH_TRILINOS
+      /**
+       * Copy the content of a Trilinos vector into the calling vector. This
+       * function assumes that the vectors layouts have already been
+       * initialized to match.
+       *
+       * This operator is only available if deal.II was configured with
+       * Trilinos.
+       *
+       * This function is deprecated. Use the interface through
+       * ReadWriteVector instead.
+       */
+      Vector<Number> &
+      operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec) DEAL_II_DEPRECATED;
+#endif
       //@}
 
       /**
@@ -906,7 +957,7 @@ namespace LinearAlgebra
        * Exception
        */
       DeclException3 (ExcNonMatchingElements,
-                      double, double, unsigned int,
+                      Number, Number, unsigned int,
                       << "Called compress(VectorOperation::insert), but"
                       << " the element received from a remote processor, value "
                       << std::setprecision(16) << arg1
@@ -1372,8 +1423,8 @@ namespace LinearAlgebra
  */
 template <typename Number>
 inline
-void swap (LinearAlgebra::parallel::Vector<Number> &u,
-           LinearAlgebra::parallel::Vector<Number> &v)
+void swap (LinearAlgebra::distributed::Vector<Number> &u,
+           LinearAlgebra::distributed::Vector<Number> &v)
 {
   u.swap (v);
 }
index 2e6b7cf7ea531bcbb401d7a6cdd77f7490be8e42..7c18e5b419b49d2403d3ef19a1d3ceb6cef2b4d0 100644 (file)
@@ -21,6 +21,8 @@
 #include <deal.II/lac/la_parallel_vector.h>
 #include <deal.II/lac/vector_operations_internal.h>
 #include <deal.II/lac/read_write_vector.h>
+#include <deal.II/lac/petsc_parallel_vector.h>
+#include <deal.II/lac/trilinos_vector.h>
 
 
 DEAL_II_NAMESPACE_OPEN
@@ -28,7 +30,7 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace LinearAlgebra
 {
-  namespace parallel
+  namespace distributed
   {
 
     template <typename Number>
@@ -384,6 +386,50 @@ namespace LinearAlgebra
 
 
 
+#ifdef DEAL_II_WITH_PETSC
+
+    template <typename Number>
+    Vector<Number> &
+    Vector<Number>::operator = (const PETScWrappers::MPI::Vector &petsc_vec)
+    {
+      IndexSet combined_set = partitioner->locally_owned_range();
+      combined_set.add_indices(partitioner->ghost_indices());
+      ReadWriteVector<Number> rw_vector(combined_set);
+      rw_vector.import(petsc_vec, VectorOperation::insert);
+      import(rw_vector, VectorOperation::insert);
+
+      if (vector_is_ghosted || petsc_vec.has_ghost_elements())
+        update_ghost_values();
+
+      return *this;
+    }
+
+#endif
+
+
+
+#ifdef DEAL_II_WITH_TRILINOS
+
+    template <typename Number>
+    Vector<Number> &
+    Vector<Number>::operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec)
+    {
+      IndexSet combined_set = partitioner->locally_owned_range();
+      combined_set.add_indices(partitioner->ghost_indices());
+      ReadWriteVector<Number> rw_vector(combined_set);
+      rw_vector.import(trilinos_vec, VectorOperation::insert);
+      import(rw_vector, VectorOperation::insert);
+
+      if (vector_is_ghosted || trilinos_vec.has_ghost_elements())
+        update_ghost_values();
+
+      return *this;
+    }
+
+#endif
+
+
+
     template <typename Number>
     void
     Vector<Number>::compress (::dealii::VectorOperation::values operation)
@@ -576,10 +622,10 @@ namespace LinearAlgebra
             for ( ; my_imports!=part.import_indices().end(); ++my_imports)
               for (unsigned int j=my_imports->first; j<my_imports->second;
                    j++, read_position++)
-                Assert(*read_position == 0. ||
+                Assert(*read_position == Number() ||
                        std::abs(local_element(j) - *read_position) <=
                        std::abs(local_element(j)) * 1000. *
-                       std::numeric_limits<Number>::epsilon(),
+                       std::numeric_limits<real_type>::epsilon(),
                        ExcNonMatchingElements(*read_position, local_element(j),
                                               part.this_mpi_process()));
           AssertDimension(read_position-import_data,part.n_import_indices());
@@ -738,14 +784,15 @@ namespace LinearAlgebra
           comm_pattern =
             std_cxx11::dynamic_pointer_cast<const Utilities::MPI::Partitioner> (communication_pattern);
           AssertThrow(comm_pattern != NULL,
-                      ExcMessage(std::string("The communication pattern is not of type ") +
+                      ExcMessage("The communication pattern is not of type "
                                  "Utilities::MPI::Partitioner."));
         }
-      LinearAlgebra::parallel::Vector<Number> tmp_vector(comm_pattern);
+      Vector<Number> tmp_vector(comm_pattern);
 
       // fill entries from ReadWriteVector into the distributed vector,
       // including ghost entries. this is not really efficient right now
-      // because indices are translated twice, once for
+      // because indices are translated twice, once by nth_index_in_set(i) and
+      // once for operator() of tmp_vector
       const IndexSet &v_stored = V.get_stored_elements();
       for (size_type i=0; i<v_stored.n_elements(); ++i)
         tmp_vector(v_stored.nth_index_in_set(i)) = V.local_element(i);
@@ -1033,7 +1080,7 @@ namespace LinearAlgebra
     Vector<Number> &
     Vector<Number>::operator /= (const Number factor)
     {
-      operator *= (1./factor);
+      operator *= (static_cast<Number>(1.)/factor);
       return *this;
     }
 
index 8cb7113f4258182b635c8e594ef0f55a6bffbf46..146e35085d84f02a47edc8c640d37ed7b13fcbe9 100644 (file)
@@ -139,7 +139,7 @@ namespace parallel
      *
      * @author Katharina Kormann, Martin Kronbichler, 2010, 2011
      */
-    using LinearAlgebra::parallel::Vector;
+    using LinearAlgebra::distributed::Vector;
   }
 }
 
index b633050b42a9a5b24417c541492dbd4e050693a6..1439a66ae4d9f3084afbf63dfe03056de2545028 100644 (file)
@@ -1021,7 +1021,7 @@ template<class VectorType>
 inline void
 PreconditionIdentity::vmult_add (VectorType &dst, const VectorType &src) const
 {
-  dst.add(src);
+  dst += src;
 }
 
 
@@ -1030,7 +1030,7 @@ template<class VectorType>
 inline void
 PreconditionIdentity::Tvmult_add (VectorType &dst, const VectorType &src) const
 {
-  dst.add(src);
+  dst += src;
 }
 
 inline PreconditionIdentity::size_type
index d42be66c2f18759110c9c79372676c2c9648eff6..51746c395a4a66a70e5b41beff7aa2c09a6bf006 100644 (file)
@@ -39,6 +39,10 @@ DEAL_II_NAMESPACE_OPEN
 namespace LinearAlgebra
 {
   class CommunicationPatternBase;
+  namespace distributed
+  {
+    template <typename> class Vector;
+  }
 }
 
 #ifdef DEAL_II_WITH_PETSC
@@ -246,6 +250,19 @@ namespace LinearAlgebra
      */
     ReadWriteVector<Number> &operator = (const Number s);
 
+    /**
+     * Imports all the elements present in the vector's IndexSet from the
+     * input vector @p vec. VectorOperation::values @p operation
+     * is used to decide if the elements in @p V should be added to the
+     * current vector or replace the current elements. The last parameter can
+     * be used if the same communication pattern is used multiple times. This
+     * can be used to improve performance.
+     */
+    void import(const distributed::Vector<Number> &vec,
+                VectorOperation::values operation,
+                std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern =
+                  std_cxx11::shared_ptr<const CommunicationPatternBase> ());
+
 #ifdef DEAL_II_WITH_PETSC
     /**
      * Imports all the elements present in the vector's IndexSet from the input
index 11d5dfbe5540ddc3ec467926a587ef54975d9ee6..3a745de51326644a3d395ba0e6303882f801ddad 100644 (file)
@@ -20,6 +20,8 @@
 #include <deal.II/base/config.h>
 #include <deal.II/lac/read_write_vector.h>
 #include <deal.II/lac/vector_operations_internal.h>
+#include <deal.II/base/partitioner.h>
+#include <deal.II/lac/la_parallel_vector.h>
 
 #ifdef DEAL_II_WITH_PETSC
 #  include <deal.II/lac/petsc_parallel_vector.h>
@@ -188,6 +190,45 @@ namespace LinearAlgebra
 
 
 
+  template <typename Number>
+  void
+  ReadWriteVector<Number>::import(const distributed::Vector<Number> &vec,
+                                  VectorOperation::values operation,
+                                  std_cxx11::shared_ptr<const CommunicationPatternBase> communication_pattern)
+  {
+    // If no communication pattern is given, create one. Otherwise, use the
+    // given one.
+    std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> comm_pattern;
+    if (communication_pattern.get() == NULL)
+      {
+        comm_pattern.reset(new Utilities::MPI::Partitioner(vec.locally_owned_elements(),
+                                                           get_stored_elements(),
+                                                           vec.get_mpi_communicator()));
+      }
+    else
+      {
+        comm_pattern =
+          std_cxx11::dynamic_pointer_cast<const Utilities::MPI::Partitioner> (communication_pattern);
+        AssertThrow(comm_pattern != NULL,
+                    ExcMessage("The communication pattern is not of type "
+                               "Utilities::MPI::Partitioner."));
+      }
+    distributed::Vector<Number> tmp_vector(comm_pattern);
+
+    std::copy(vec.begin(), vec.end(), tmp_vector.begin());
+    tmp_vector.update_ghost_values();
+
+    const IndexSet &stored = get_stored_elements();
+    if (operation == VectorOperation::add)
+      for (size_type i=0; i<stored.n_elements(); ++i)
+        local_element(i) += tmp_vector(stored.nth_index_in_set(i));
+    else
+      for (size_type i=0; i<stored.n_elements(); ++i)
+        local_element(i) = tmp_vector(stored.nth_index_in_set(i));
+  }
+
+
+
 #ifdef DEAL_II_WITH_PETSC
   namespace internal
   {
index 0903bc2f377a324502f8d3df95182e7531bb644c..d6490e2f57773945fce7ddcf4797a81cbebbb7de 100644 (file)
@@ -36,7 +36,7 @@ DEAL_II_NAMESPACE_OPEN
 // forward declarations
 namespace LinearAlgebra
 {
-  namespace parallel
+  namespace distributed
   {
     template <typename> class Vector;
   }
@@ -2117,8 +2117,8 @@ namespace internal
   template <typename Number>
   inline
   Number &
-  vector_access (LinearAlgebra::parallel::Vector<Number> &vec,
-                 const unsigned int                     entry)
+  vector_access (LinearAlgebra::distributed::Vector<Number> &vec,
+                 const unsigned int                          entry)
   {
     return vec.local_element(entry);
   }
@@ -2131,8 +2131,8 @@ namespace internal
   template <typename Number>
   inline
   Number
-  vector_access (const LinearAlgebra::parallel::Vector<Number> &vec,
-                 const unsigned int                           entry)
+  vector_access (const LinearAlgebra::distributed::Vector<Number> &vec,
+                 const unsigned int                                entry)
   {
     return vec.local_element(entry);
   }
@@ -2140,7 +2140,8 @@ namespace internal
 
 
   // this is to make sure that the parallel partitioning in the
-  // parallel::distributed::Vector is really the same as stored in MatrixFree
+  // LinearAlgebra::distributed::Vector is really the same as stored in
+  // MatrixFree
   template <typename VectorType>
   inline
   void check_vector_compatibility (const VectorType                             &vec,
@@ -2155,8 +2156,8 @@ namespace internal
 
   template <typename Number>
   inline
-  void check_vector_compatibility (const LinearAlgebra::parallel::Vector<Number>  &vec,
-                                   const internal::MatrixFreeFunctions::DoFInfo &dof_info)
+  void check_vector_compatibility (const LinearAlgebra::distributed::Vector<Number> &vec,
+                                   const internal::MatrixFreeFunctions::DoFInfo     &dof_info)
   {
     Assert (vec.partitioners_are_compatible(*dof_info.vector_partitioner),
             ExcMessage("The parallel layout of the given vector is not "
index f62ac75049af23a0a63bbeaccf617eea5deb888c..e34776129644f325f616dff93ddff3f3cade48d5 100644 (file)
@@ -27,7 +27,7 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace LinearAlgebra
 {
-  namespace parallel
+  namespace distributed
   {
 #define TEMPL_COPY_CONSTRUCTOR(S1,S2)                                   \
   template Vector<S1>& Vector<S1>::operator=<S2> (const Vector<S2> &)
index 99b7760bd4fb2545bdce676db619450dce5354d9..e06292f4e9baea3c1f001ffff9d80fd0027fd807 100644 (file)
@@ -19,7 +19,7 @@ for (SCALAR : REAL_SCALARS)
 {
   namespace LinearAlgebra
   \{
-    namespace parallel
+    namespace distributed
     \{
       template class Vector<SCALAR>;
     \}
@@ -30,7 +30,30 @@ for (S1, S2 : REAL_SCALARS)
 {
   namespace LinearAlgebra
   \{
-    namespace parallel
+    namespace distributed
+    \{
+      template void Vector<S1>::reinit<S2> (const Vector<S2>&,
+                                            const bool);
+    \}
+  \}
+}
+
+for (SCALAR : COMPLEX_SCALARS)
+{
+  namespace LinearAlgebra
+  \{
+    namespace distributed
+    \{
+      template class Vector<SCALAR>;
+    \}
+  \}
+}
+
+for (S1, S2 : COMPLEX_SCALARS)
+{
+  namespace LinearAlgebra
+  \{
+    namespace distributed
     \{
       template void Vector<S1>::reinit<S2> (const Vector<S2>&,
                                             const bool);

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.