]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix a failing PETSc test. 12320/head
authorDavid Wells <drwells@email.unc.edu>
Wed, 26 May 2021 21:00:29 +0000 (17:00 -0400)
committerDavid Wells <drwells@email.unc.edu>
Wed, 26 May 2021 21:20:37 +0000 (17:20 -0400)
This test depended on a now-removed constructor.

tests/lac/linear_operator_09.cc

index ebeac49686dd1fe3392637e6de6cdfd15cb1fa4c..4840b0528d10d57836d58c4f59f5963905826abe 100644 (file)
@@ -17,6 +17,9 @@
 // different kinds of PETSc matrices and vectors
 // TODO: A bit more tests...
 
+#include <deal.II/base/index_set.h>
+
+#include <deal.II/lac/dynamic_sparsity_pattern.h>
 #include <deal.II/lac/linear_operator.h>
 
 #include "../tests.h"
@@ -37,18 +40,29 @@ main(int argc, char *argv[])
   using size_type = PETScWrappers::MPI::SparseMatrix::size_type;
 
   Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+  const auto rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
 
   initlog();
   deallog << std::setprecision(10);
 
   {
-    unsigned int np = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+    const unsigned int np     = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
+    const unsigned int n_dofs = 4;
     if (4 % np == 0 && np <= 4)
       {
-        PETScWrappers::MPI::SparseMatrix a(
-          MPI_COMM_WORLD, 4, 4, 4 / np, 4 / np, 1);
-        for (unsigned int i = 0; i < 4; ++i)
-          for (unsigned int j = 0; j < 4; ++j)
+        const auto dofs_per_processor = n_dofs / np;
+        IndexSet   locally_owned_dofs(n_dofs);
+        locally_owned_dofs.add_range(rank * dofs_per_processor,
+                                     (rank + 1) * dofs_per_processor);
+        locally_owned_dofs.compress();
+        DynamicSparsityPattern dsp(n_dofs, n_dofs);
+        for (const auto &index : locally_owned_dofs)
+          dsp.add(index, index);
+
+        PETScWrappers::MPI::SparseMatrix a;
+        a.reinit(locally_owned_dofs, locally_owned_dofs, dsp, MPI_COMM_WORLD);
+        for (const auto &i : locally_owned_dofs)
+          for (const auto &j : locally_owned_dofs)
             a.add(i, i, 1);
         a.compress(VectorOperation::add);
         auto op_a = linear_operator<PETScWrappers::MPI::Vector>(a);

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.