]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add similar interface for copying from a Trilinos MPI vector as for Petsc.
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Fri, 31 Jan 2014 19:03:37 +0000 (19:03 +0000)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Fri, 31 Jan 2014 19:03:37 +0000 (19:03 +0000)
git-svn-id: https://svn.dealii.org/trunk@32358 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/lac/parallel_block_vector.h
deal.II/include/deal.II/lac/parallel_vector.h
deal.II/include/deal.II/lac/parallel_vector.templates.h

index 7e14d6133b7a40fc2cdf2104e1ba2f3b2ba4bd0c..f935735aa875e4c925c92229ccbd7ffd06e2a00c 100644 (file)
@@ -103,9 +103,10 @@ inconvenience this causes.
 
 <ol>
   <li>New: There is now a method to copy the content from a
-  PETScWrappers::MPI::Vector to deal.II's parallel distributed vector.
+  PETScWrappers::MPI::Vector and TrilinosWrappers::MPI::Vector to
+  deal.II's parallel distributed vector.
   <br>
-  (Ben Thompson, 2014/01/31)
+  (Ben Thompson, Martin Kronbichler, 2014/01/31)
 
   <li>Fixed: The SolutionTransfer class had all sorts of problems when
   used with hp::DoFHandler that made its results at least questionable.
index 7ffd04094bdb794776f2dd4df6aa4ca7ef02be1d..d7a68d4c3b1962d3346062af4288970326cc8c4b 100644 (file)
 #include <deal.II/lac/block_vector_base.h>
 #include <deal.II/lac/parallel_vector.h>
 
+#include <deal.II/lac/petsc_parallel_block_vector.h>
+#include <deal.II/lac/trilinos_parallel_block_vector.h>
+
+
 #include <cstdio>
 #include <vector>
 
 DEAL_II_NAMESPACE_OPEN
 
 
-
 namespace parallel
 {
   namespace distributed
@@ -171,6 +174,31 @@ namespace parallel
       BlockVector &
       operator= (const Vector<Number> &V);
 
+#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.
+       */
+      BlockVector<Number> &
+      operator = (const PETScWrappers::MPI::BlockVector &petsc_vec);
+#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.
+       */
+      BlockVector<Number> &
+      operator = (const TrilinosWrappers::MPI::BlockVector &trilinos_vec);
+#endif
+
       /**
        * Reinitialize the BlockVector to contain <tt>num_blocks</tt> blocks of
        * size <tt>block_size</tt> each.
@@ -594,6 +622,40 @@ namespace parallel
 
 
 
+#ifdef DEAL_II_WITH_PETSC
+
+    template <typename Number>
+    BlockVector<Number> &
+    BlockVector<Number>::operator = (const PETScWrappers::MPI::BlockVector &petsc_vec)
+    {
+      AssertDimension(this->n_blocks(), petsc_vec.n_blocks());
+      for (unsigned int i=0; i<this->n_blocks(); ++i)
+        this->block(i) = petsc_vec.block(i);
+
+      return *this;
+    }
+
+#endif
+
+
+
+#ifdef DEAL_II_WITH_TRILINOS
+
+    template <typename Number>
+    BlockVector<Number> &
+    BlockVector<Number>::operator = (const TrilinosWrappers::MPI::BlockVector &trilinos_vec)
+    {
+      AssertDimension(this->n_blocks(), trilinos_vec.n_blocks());
+      for (unsigned int i=0; i<this->n_blocks(); ++i)
+        this->block(i) = trilinos_vec.block(i);
+
+      return *this;
+    }
+
+#endif
+
+
+
     template <typename Number>
     inline
     void
index 9029db91c5fc1c636a4bbee47a2cbad9c0efa04c..126737eba779078c1d80dd0be70a4c916f6f9401 100644 (file)
@@ -43,6 +43,16 @@ namespace PETScWrappers
 }
 #endif
 
+#ifdef DEAL_II_WITH_TRILINOS
+namespace TrilinosWrappers
+{
+  namespace MPI
+  {
+    class Vector;
+  }
+}
+#endif
+
 
 namespace parallel
 {
@@ -288,9 +298,9 @@ namespace parallel
 
 #ifdef DEAL_II_WITH_PETSC
       /**
-       * Copy the content of a PETSc vector into a vector. This function
-       * assumes that the vectors layouts have already been initialized to
-       * match.
+       * 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.
        */
@@ -298,6 +308,19 @@ namespace parallel
       operator = (const PETScWrappers::MPI::Vector &petsc_vec);
 #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.
+       */
+      Vector<Number> &
+      operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec);
+#endif
+
       /**
        * This method copies the local range from another vector with the same
        * local range, but possibly different layout of ghost indices.
index fbce1fb85dfd37fb07b8d7177fc138be4c6839de..1b60485822fdf1907207d02de203b27dccc4a6d2 100644 (file)
 
 #include <deal.II/base/config.h>
 #include <deal.II/lac/parallel_vector.h>
+#include <deal.II/lac/vector_view.h>
 
 #include <deal.II/lac/petsc_parallel_vector.h>
+#include <deal.II/lac/trilinos_vector.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -209,7 +211,6 @@ namespace parallel
 #ifdef DEAL_II_WITH_PETSC
 
     template <typename Number>
-    inline
     Vector<Number> &
     Vector<Number>::operator = (const PETScWrappers::MPI::Vector &petsc_vec)
     {
@@ -229,9 +230,34 @@ namespace parallel
       AssertThrow (ierr == 0, ExcPETScError(ierr));
 
       // spread ghost values between processes?
-      bool must_update_ghost_values = vector_is_ghosted || 
-        petsc_vec.has_ghost_elements();
-      if (must_update_ghost_values)
+      if (vector_is_ghosted || petsc_vec.has_ghost_elements())
+        update_ghost_values();
+
+      // return a pointer to this object per normal c++ operator overloading
+      // semantics
+      return *this;
+    }
+
+#endif
+
+
+
+#ifdef DEAL_II_WITH_TRILINOS
+
+    template <typename Number>
+    Vector<Number> &
+    Vector<Number>::operator = (const TrilinosWrappers::MPI::Vector &trilinos_vec)
+    {
+      Assert(trilinos_vec.locally_owned_elements() == locally_owned_elements(),
+             StandardExceptions::ExcInvalidState());
+
+      // create on trilinos data
+      const VectorView<double> in_view (local_size(), trilinos_vec.begin());
+      static_cast<::dealii::Vector<Number>&>(vector_view) =
+        static_cast<const ::dealii::Vector<double>&>(in_view);
+
+      // spread ghost values between processes?
+      if (vector_is_ghosted || trilinos_vec.has_ghost_elements())
         update_ghost_values();
 
       // return a pointer to this object per normal c++ operator overloading

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.