]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add more reinit() to LA::distributed::BlockVector. 14486/head
authorMarc Fehling <mafehling.git@gmail.com>
Sat, 26 Nov 2022 07:13:25 +0000 (00:13 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Sun, 27 Nov 2022 20:38:35 +0000 (13:38 -0700)
include/deal.II/lac/la_parallel_block_vector.h
include/deal.II/lac/la_parallel_block_vector.templates.h
include/deal.II/lac/la_parallel_vector.h

index 2f1b00e8514b7088ed6dd2f96dfa778d327925ca..c250cf4004322fde7e104092526f72dc3131cce0 100644 (file)
@@ -299,6 +299,35 @@ namespace LinearAlgebra
       reinit(const BlockVector<Number2> &V,
              const bool                  omit_zeroing_entries = false);
 
+      /**
+       * Initialize the block vector. For each block, the local range is
+       * specified by the corresponding entry in @p local_ranges (note that this
+       * must be a contiguous interval, multiple intervals are not possible).
+       * The parameter @p ghost_indices specifies ghost indices for each block,
+       * i.e., indices which one might need to read data from or accumulate data
+       * from. It is allowed that the set of ghost indices also contains the
+       * local range, but it does not need to.
+       *
+       * This function involves global communication, so it should only be
+       * called once for a given layout. Use the @p reinit function with
+       * BlockVector<Number> argument to create additional vectors with the same
+       * parallel layout.
+       *
+       * @see
+       * @ref GlossGhostedVector "vectors with ghost elements"
+       */
+      void
+      reinit(const std::vector<IndexSet> &local_ranges,
+             const std::vector<IndexSet> &ghost_indices,
+             const MPI_Comm &             communicator);
+
+      /**
+       * Same as above, but without ghost entries.
+       */
+      void
+      reinit(const std::vector<IndexSet> &local_ranges,
+             const MPI_Comm &             communicator);
+
       /**
        * This function copies the data that has accumulated in the data buffer
        * for ghost indices to the owning processor. For the meaning of the
index 3a5f0c081600540804726e5fb7737dcc0cb43fd5..330a82adee3acc6d6120e6015261842c991059e6 100644 (file)
@@ -62,15 +62,7 @@ namespace LinearAlgebra
                                      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);
+      reinit(local_ranges, ghost_indices, communicator);
     }
 
 
@@ -78,15 +70,7 @@ namespace LinearAlgebra
     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);
+      reinit(local_ranges, communicator);
     }
 
 
@@ -95,11 +79,13 @@ namespace LinearAlgebra
     BlockVector<Number>::BlockVector(const BlockVector<Number> &v)
       : BlockVectorBase<Vector<Number>>()
     {
-      this->components.resize(v.n_blocks());
       this->block_indices = v.block_indices;
 
-      for (size_type i = 0; i < this->n_blocks(); ++i)
+      this->components.resize(this->n_blocks());
+      for (unsigned int i = 0; i < this->n_blocks(); ++i)
         this->components[i] = v.components[i];
+
+      this->collect_sizes();
     }
 
 
@@ -131,11 +117,12 @@ namespace LinearAlgebra
                                 const bool omit_zeroing_entries)
     {
       this->block_indices.reinit(block_sizes);
-      if (this->components.size() != this->n_blocks())
-        this->components.resize(this->n_blocks());
 
-      for (size_type i = 0; i < this->n_blocks(); ++i)
+      this->components.resize(this->n_blocks());
+      for (unsigned int i = 0; i < this->n_blocks(); ++i)
         this->components[i].reinit(block_sizes[i], omit_zeroing_entries);
+
+      this->collect_sizes();
     }
 
 
@@ -146,12 +133,57 @@ namespace LinearAlgebra
     BlockVector<Number>::reinit(const BlockVector<Number2> &v,
                                 const bool omit_zeroing_entries)
     {
-      this->block_indices = v.get_block_indices();
-      if (this->components.size() != this->n_blocks())
-        this->components.resize(this->n_blocks());
+      if (this->n_blocks() != v.n_blocks())
+        this->block_indices = v.get_block_indices();
+
+      this->components.resize(this->n_blocks());
+      for (unsigned int i = 0; i < this->n_blocks(); ++i)
+        this->components[i].reinit(v.block(i), omit_zeroing_entries);
+
+      this->collect_sizes();
+    }
+
+
+
+    template <typename Number>
+    void
+    BlockVector<Number>::reinit(const std::vector<IndexSet> &local_ranges,
+                                const std::vector<IndexSet> &ghost_indices,
+                                const MPI_Comm &             communicator)
+    {
+      AssertDimension(local_ranges.size(), ghost_indices.size());
+
+      // update the number of blocks
+      this->block_indices.reinit(local_ranges.size(), 0);
 
+      // initialize each block
+      this->components.resize(this->n_blocks());
       for (unsigned int i = 0; i < this->n_blocks(); ++i)
-        this->block(i).reinit(v.block(i), omit_zeroing_entries);
+        this->components[i].reinit(local_ranges[i],
+                                   ghost_indices[i],
+                                   communicator);
+
+      // update block_indices content
+      this->collect_sizes();
+    }
+
+
+
+    template <typename Number>
+    void
+    BlockVector<Number>::reinit(const std::vector<IndexSet> &local_ranges,
+                                const MPI_Comm &             communicator)
+    {
+      // update the number of blocks
+      this->block_indices.reinit(local_ranges.size(), 0);
+
+      // initialize each block
+      this->components.resize(this->n_blocks());
+      for (unsigned int i = 0; i < this->n_blocks(); ++i)
+        this->components[i].reinit(local_ranges[i], communicator);
+
+      // update block_indices content
+      this->collect_sizes();
     }
 
 
@@ -178,12 +210,14 @@ namespace LinearAlgebra
              ExcDimensionMismatch(this->n_blocks(), v.n_blocks()));
 
       if (this->n_blocks() != v.n_blocks())
-        reinit(v.n_blocks(), true);
+        this->block_indices = v.block_indices;
 
+      this->components.resize(this->n_blocks());
       for (size_type i = 0; i < this->n_blocks(); ++i)
-        this->components[i] = v.block(i);
+        this->components[i] = v.components[i];
 
       this->collect_sizes();
+
       return *this;
     }
 
index c8f09024f2ca9984e1a4813e33a31ec792ed920f..3ba62fd92788f2708befa3779d15805ff61bce43 100644 (file)
@@ -362,12 +362,12 @@ namespace LinearAlgebra
              const bool                          omit_zeroing_entries = false);
 
       /**
-       * Initialize the vector. The local range is specified by @p
-       * locally_owned_set (note that this must be a contiguous interval,
-       * multiple intervals are not possible). The IndexSet @p ghost_indices
-       * specifies ghost indices, i.e., indices which one might need to read
-       * data from or accumulate data from. It is allowed that the set of
-       * ghost indices also contains the local range, but it does not need to.
+       * Initialize the vector. The local range is specified by @p local_range
+       * (note that this must be a contiguous interval, multiple intervals are
+       * not possible). The IndexSet @p ghost_indices specifies ghost indices,
+       * i.e., indices which one might need to read data from or accumulate data
+       * from. It is allowed that the set of ghost indices also contains the
+       * local range, but it does not need to.
        *
        * This function involves global communication, so it should only be
        * called once for a given layout. Use the @p reinit function with

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.