]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Modified reinit function to check if parallel_partitioner is_ascending_and_one_to_one... 6725/head
authorjthano <josh.hanophy@gmail.com>
Mon, 4 Jun 2018 05:56:00 +0000 (22:56 -0700)
committerjthano <josh.hanophy@gmail.com>
Tue, 19 Jun 2018 01:30:05 +0000 (18:30 -0700)
added draft test

Added test. Test passes

Modified reinit function to check if parallel_partitioner is_ascending_and_one_to_one and then pass this result to make_trilinos_map. The reason for this check is that if the check returns true, then make_trilinos_map may be able to create a linear epetra map. Also added test and test output.

Fixed issues related to test and added change log

Renamed test files so that mpi test will happen. Also made minor edits to test file

doc/news/changes/minor/20180616JohsuaHanophy [new file with mode: 0644]
source/lac/trilinos_vector.cc
tests/trilinos/vector_reinit_to_linear_map.cc [new file with mode: 0644]
tests/trilinos/vector_reinit_to_linear_map.mpirun=2.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20180616JohsuaHanophy b/doc/news/changes/minor/20180616JohsuaHanophy
new file mode 100644 (file)
index 0000000..3732977
--- /dev/null
@@ -0,0 +1,5 @@
+Fixed: TrilinosWrappers::MPI::Vector::reinit check  now checks if
+parallel_partitioner is_ascending_and_one_to_one before calling
+make_trilinos_map because make_trilinos_map may be able to make
+a linear map if is_ascending_and_one_to_one is true.
+(Joshua Hanophy, 2018/06/16)
index 3ad0b5472fb871b131f03c39b70ad99fcc1d519c..f2acfcec4bfce13e0ba33761d74c89fe2dfbc07a 100644 (file)
@@ -157,8 +157,11 @@ namespace TrilinosWrappers
     {
       nonlocal_vector.reset();
 
+      const bool overlapping =
+        !parallel_partitioner.is_ascending_and_one_to_one(communicator);
+
       Epetra_Map map =
-        parallel_partitioner.make_trilinos_map(communicator, true);
+        parallel_partitioner.make_trilinos_map(communicator, overlapping);
 
       vector = std_cxx14::make_unique<Epetra_FEVector>(map);
 
diff --git a/tests/trilinos/vector_reinit_to_linear_map.cc b/tests/trilinos/vector_reinit_to_linear_map.cc
new file mode 100644 (file)
index 0000000..b1964d4
--- /dev/null
@@ -0,0 +1,80 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 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.
+//
+// ---------------------------------------------------------------------
+
+// run TrilinosWrappers::MPI::Vector::reinit() in MPI environment
+// and with an IndexSet that is_ascending_and_one_to_one to check
+// that make_trilinos_map will construct the object with a Map that
+// is linear
+
+#include <deal.II/base/index_set.h>
+
+#include <deal.II/distributed/grid_refinement.h>
+#include <deal.II/distributed/tria.h>
+
+#include <deal.II/dofs/dof_handler.h>
+#include <deal.II/dofs/dof_tools.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/lac/trilinos_vector.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+void
+test()
+{
+  MPI_Comm mpi_communicator(MPI_COMM_WORLD);
+
+  parallel::distributed::Triangulation<2> tria(mpi_communicator);
+
+  GridGenerator::hyper_cube(tria, -1, 0);
+  tria.refine_global(2);
+
+  const unsigned int poly_degree = 1;
+  FE_Q<2>            fe(poly_degree);
+
+  DoFHandler<2> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  IndexSet locally_owned_dofs = dof_handler.locally_owned_dofs();
+
+  TrilinosWrappers::MPI::Vector vector_linear;
+  vector_linear.reinit(locally_owned_dofs, mpi_communicator);
+
+  const int is_linear_index_set =
+    locally_owned_dofs.is_ascending_and_one_to_one(mpi_communicator);
+
+  const int is_linear_map = vector_linear.trilinos_vector().Map().LinearMap();
+
+  if (is_linear_index_set == 1 && is_linear_map == 1)
+    {
+      deallog << "OK" << std::endl;
+    }
+}
+
+
+int
+main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  mpi_initlog();
+
+  test();
+}
diff --git a/tests/trilinos/vector_reinit_to_linear_map.mpirun=2.output b/tests/trilinos/vector_reinit_to_linear_map.mpirun=2.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK

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.