]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a test for the conversion BlockVector->Vector
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Tue, 16 Aug 2016 16:48:42 +0000 (18:48 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 16 Sep 2016 12:43:08 +0000 (14:43 +0200)
include/deal.II/lac/trilinos_vector.h
source/lac/trilinos_vector.cc
tests/mpi/blockvec_02.cc [new file with mode: 0644]
tests/mpi/blockvec_02.with_trilinos=true.mpirun=2.output [new file with mode: 0644]

index 8e57c44931fb9c9f857b3bc52f26d095d5979b2d..e2225accc161c549ea8abcef1596a1fba6363f25 100644 (file)
@@ -583,6 +583,10 @@ namespace TrilinosWrappers
        * or may not have ghost elements. See the general documentation of this
        * class for more information.
        *
+       * In case @p parallel_partitioning is overlapping, it is not clear which
+       * process should own which elements. Hence, locally_owned_elements()
+       * returns an empty IndexSet in this case.
+       *
        * @see
        * @ref GlossGhostedVector "vectors with ghost elements"
        */
index 9ea961c8d1d118451e4036d511ce25d7a12f7b4b..cfae13069ac6b731996336fc7d96ac8357ea1dc5 100644 (file)
@@ -262,10 +262,10 @@ namespace TrilinosWrappers
       // which process owns what. So we decide that no process
       // owns anything in that case.
       if (has_ghosts)
-      {
-        owned_elements.clear();
-        owned_elements.set_size(parallel_partitioner.size());
-      }
+        {
+          owned_elements.clear();
+          owned_elements.set_size(parallel_partitioner.size());
+        }
       else
         owned_elements = parallel_partitioner;
 
diff --git a/tests/mpi/blockvec_02.cc b/tests/mpi/blockvec_02.cc
new file mode 100644 (file)
index 0000000..4544a49
--- /dev/null
@@ -0,0 +1,73 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2004 - 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// Test constructor/reinit of BlockVector with IndexSets and the conversion
+// to a Vector
+
+#include "../tests.h"
+#include "coarse_grid_common.h"
+#include <deal.II/lac/petsc_parallel_block_vector.h>
+#include <deal.II/lac/trilinos_parallel_block_vector.h>
+#include <deal.II/base/index_set.h>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#include <vector>
+
+void test ()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+  unsigned int numproc = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
+
+  std::vector<IndexSet> local_active(2);
+
+  // block 0:
+  local_active[0].set_size(numproc);
+  local_active[0].add_range(myid,myid+1);
+
+  local_active[1].set_size(2*numproc);
+  local_active[1].add_range(myid*2,myid*2+2);
+
+  TrilinosWrappers::MPI::BlockVector v_block(local_active, MPI_COMM_WORLD);
+
+  v_block(myid) = 100.0 + myid;
+
+  v_block.block(1)(myid*2)=myid*2.0;
+  v_block.block(1)(myid*2+1)=myid*2.0+1.0;
+
+  v_block.compress(VectorOperation::insert);
+
+  deallog << "block size: " << v_block.size() << std::endl;
+  deallog << "block size[0]: " << v_block.block(0).size() << std::endl;
+  deallog << "block size[1]: " << v_block.block(1).size() << std::endl;
+  v_block.locally_owned_elements().print(deallog.get_file_stream());
+
+  TrilinosWrappers::MPI::Vector v;
+  v.reinit (v_block);
+
+  deallog << "size: " << v.size() << std::endl;
+  v.locally_owned_elements().print(deallog.get_file_stream());
+}
+
+
+
+int main (int argc, char **argv)
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1);
+  MPILogInitAll log;
+  test();
+}
diff --git a/tests/mpi/blockvec_02.with_trilinos=true.mpirun=2.output b/tests/mpi/blockvec_02.with_trilinos=true.mpirun=2.output
new file mode 100644 (file)
index 0000000..b2bc951
--- /dev/null
@@ -0,0 +1,15 @@
+
+DEAL:0::block size: 6
+DEAL:0::block size[0]: 2
+DEAL:0::block size[1]: 4
+{0, [2,3]}
+DEAL:0::size: 6
+{0, [2,3]}
+
+DEAL:1::block size: 6
+DEAL:1::block size[0]: 2
+DEAL:1::block size[1]: 4
+{1, [4,5]}
+DEAL:1::size: 6
+{1, [4,5]}
+

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.