]> https://gitweb.dealii.org/ - dealii.git/commitdiff
ghost import for distributed::parallel:BlockVector::operator= and add new constructors
authorTimo Heister <timo.heister@gmail.com>
Wed, 4 Sep 2013 17:02:38 +0000 (17:02 +0000)
committerTimo Heister <timo.heister@gmail.com>
Wed, 4 Sep 2013 17:02:38 +0000 (17:02 +0000)
git-svn-id: https://svn.dealii.org/trunk@30597 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 9baf5e8c27cc3814dccc38dff0ffbbea88525d8d..5a3549443d6a5c29f226183b1e0a55920cb40f82 100644 (file)
@@ -68,6 +68,15 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li>
+  Changed: distributed::parallel:BlockVector::operator= now allows importing
+  of ghost values like all other vector types. Also added some new constructors
+  for BlockVector and Vector using IndexSets to mirror the other linear algebra
+  classes.
+  <br>
+  (Timo Heister, 2013/09/04)
+  </li>
+
   <li>
   Fixed: VectorTools::compute_no_normal_flux_constraints had a bug that
   only manifested on complex meshes. This is now fixed.
index d0f538785167c006a337b38decebe52cd9f87344..747d12a33eb1f00561054fc44c507392e4cfd81b 100644 (file)
@@ -125,6 +125,20 @@ namespace parallel
        */
       BlockVector (const std::vector<size_type> &block_sizes);
 
+      /**
+       * Construct a block vector with an IndexSet for the local range
+       * and ghost entries for each block.
+       */
+      BlockVector (const std::vector<IndexSet> &local_ranges,
+                   const std::vector<IndexSet> &ghost_indices,
+                   const MPI_Comm  communicator);
+
+      /**
+       * Same as above but the ghost indicies are assumed to be empty.
+       */
+      BlockVector (const std::vector<IndexSet> &local_ranges,
+                   const MPI_Comm  communicator);
+
       /**
        * Destructor. Clears memory.
        */
@@ -356,6 +370,41 @@ namespace parallel
     }
 
 
+    template <typename Number>
+    inline
+    BlockVector<Number>::BlockVector (const std::vector<IndexSet> &local_ranges,
+        const std::vector<IndexSet> &ghost_indices,
+        const MPI_Comm  communicator)
+    {
+      std::vector<size_type> sizes(local_ranges.size());
+      for (unsigned int i=0; i<local_ranges.size(); ++i)
+        sizes[i] = local_ranges[i].size();
+
+      this->block_indices.reinit(sizes);
+      this->components.resize(this->n_blocks());
+
+      for (unsigned int i=0; i<this->n_blocks(); ++i)
+        this->block(i).reinit(local_ranges[i], ghost_indices[i], communicator);
+    }
+
+
+    template <typename Number>
+    inline
+    BlockVector<Number>::BlockVector (const std::vector<IndexSet> &local_ranges,
+        const MPI_Comm  communicator)
+    {
+      std::vector<size_type> sizes(local_ranges.size());
+      for (unsigned int i=0; i<local_ranges.size(); ++i)
+        sizes[i] = local_ranges[i].size();
+
+      this->block_indices.reinit(sizes);
+      this->components.resize(this->n_blocks());
+
+      for (unsigned int i=0; i<this->n_blocks(); ++i)
+        this->block(i).reinit(local_ranges[i], communicator);
+    }
+
+
 
     template <typename Number>
     inline
@@ -454,8 +503,18 @@ namespace parallel
     BlockVector<Number> &
     BlockVector<Number>::operator = (const BlockVector &v)
     {
-      reinit (v, true);
-      BaseClass::operator = (v);
+      // we only allow assignment to vectors with the same number of blocks
+      // or to an empty BlockVector
+      Assert (this->n_blocks() == 0 || this->n_blocks() == v.n_blocks(),
+                    ExcDimensionMismatch(this->n_blocks(), v.n_blocks()));
+
+      if (this->n_blocks() != v.n_blocks())
+        reinit(v.n_blocks(), true);
+
+      for (size_type i=0; i<this->n_blocks(); ++i)
+        this->components[i] = v.block(i);
+
+      this->collect_sizes();
       return *this;
     }
 
index 6d966154a07b9740f801f1f8d6b7e4efcbd3388d..ffea6ae8435b6d716efb3999e7820962465e06bf 100644 (file)
@@ -154,6 +154,12 @@ namespace parallel
               const IndexSet &ghost_indices,
               const MPI_Comm  communicator);
 
+      /**
+       * Same constructor as above but without any ghost indices.
+       */
+      Vector (const IndexSet &local_range,
+              const MPI_Comm  communicator);
+
       /**
        * Create the vector based on the parallel partitioning described in @p
        * partitioner. The input argument is a shared pointer, which store the
@@ -204,6 +210,12 @@ namespace parallel
                    const IndexSet &ghost_indices,
                    const MPI_Comm  communicator);
 
+      /**
+       * Same as above, but without ghost entries.
+       */
+      void reinit (const IndexSet &local_range,
+                   const MPI_Comm  communicator);
+
       /**
        * Initialize the vector given to the parallel partitioning described in
        * @p partitioner. The input argument is a shared pointer, which store
@@ -1002,6 +1014,22 @@ namespace parallel
 
 
 
+    template <typename Number>
+    inline
+    Vector<Number>::Vector (const IndexSet &local_range,
+                             const MPI_Comm  communicator)
+      :
+      allocated_size (0),
+      val (0),
+      import_data (0),
+      vector_view (0, static_cast<Number *>(0))
+    {
+      IndexSet ghost_indices(local_range.size());
+      reinit (local_range, ghost_indices, communicator);
+    }
+
+
+
     template <typename Number>
     inline
     Vector<Number>::Vector (const size_type size)
index 23f0734d0329f17ea3829a14f380e336803d21cd..334966c5f099f57c4c842d5bc689c229ee984886 100644 (file)
@@ -144,8 +144,22 @@ namespace parallel
                             const IndexSet &ghost_indices,
                             const MPI_Comm  communicator)
     {
-      // set up parallel partitioner with index sets
-      // and communicator
+      // set up parallel partitioner with index sets and communicator
+      std_cxx1x::shared_ptr<const Utilities::MPI::Partitioner> new_partitioner
+      (new Utilities::MPI::Partitioner (locally_owned_indices,
+                                        ghost_indices, communicator));
+      reinit (new_partitioner);
+    }
+
+
+
+    template <typename Number>
+    void
+    Vector<Number>::reinit (const IndexSet &locally_owned_indices,
+                            const MPI_Comm  communicator)
+    {
+      // set up parallel partitioner with index sets and communicator
+      IndexSet ghost_indices(locally_owned_indices.size());
       std_cxx1x::shared_ptr<const Utilities::MPI::Partitioner> new_partitioner
       (new Utilities::MPI::Partitioner (locally_owned_indices,
                                         ghost_indices, communicator));

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.