]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added checks for petsc version within BDDC, also corrected test
authorNicolas Barnafi <nabw91@gmail.com>
Fri, 9 Sep 2022 16:38:58 +0000 (18:38 +0200)
committerNicolas Barnafi <nabw91@gmail.com>
Fri, 9 Sep 2022 16:38:58 +0000 (18:38 +0200)
source/lac/petsc_parallel_sparse_matrix.cc
source/lac/petsc_precondition.cc
tests/petsc/bddc.cc

index 12a2955e79871c548f70c49021bf2c890d8651a4..64f7ac5f925cb59ecdfef3022a19fd37cc8acc1c 100644 (file)
@@ -437,6 +437,7 @@ namespace PETScWrappers
                             const IndexSet &           local_active_columns,
                             const SparsityPatternType &sparsity_pattern)
     {
+#  if DEAL_II_PETSC_VERSION_GTE(3, 10, 0)
       Assert(sparsity_pattern.n_rows() == local_rows.size(),
              ExcMessage(
                "SparsityPattern and IndexSet have different number of rows."));
@@ -449,7 +450,7 @@ namespace PETScWrappers
       Assert(local_rows.is_ascending_and_one_to_one(communicator),
              ExcNotImplemented());
 
-#  ifdef DEBUG
+#    ifdef DEBUG
       {
         // check indexsets
         types::global_dof_index row_owners =
@@ -472,7 +473,7 @@ namespace PETScWrappers
             " but sum(local_columns.n_elements())=" +
             std::to_string(col_owners) + ")"));
       }
-#  endif
+#    endif
       PetscErrorCode ierr;
 
       // create the local to global mappings as arrays.
@@ -637,6 +638,11 @@ namespace PETScWrappers
       }
       ierr = MatISRestoreLocalMat(matrix, &local_matrix);
       AssertThrow(ierr == 0, ExcPETScError(ierr));
+#  else
+      {
+        AssertThrow(false, ExcMessage("BDDC preconditioner requires PETSc 3.10.0 or newer");
+      }
+#  endif
     }
 
 #  ifndef DOXYGEN
index 1b6108dd9285fc29f89f7e420302974d42eb1c41..6f9b6593f2720fa438603abba192a46eb77af635 100644 (file)
@@ -891,6 +891,7 @@ namespace PETScWrappers
   void
   PreconditionBDDC<dim>::initialize()
   {
+#  if DEAL_II_PETSC_VERSION_GTE(3, 10, 0)
     PetscErrorCode ierr = PCSetType(pc, const_cast<char *>(PCBDDC));
     AssertThrow(ierr == 0, ExcPETScError(ierr));
 
@@ -926,7 +927,6 @@ namespace PETScWrappers
     else
       set_option_value("-pc_bddc_symmetric", "false");
     if (additional_data.coords.size() > 0)
-#  if DEAL_II_PETSC_VERSION_GTE(3, 9, 0)
       {
         set_option_value("-pc_bddc_corner_selection", "true");
         // Convert coords vector to PETSc data array
@@ -950,15 +950,15 @@ namespace PETScWrappers
         ierr = PCSetCoordinates(pc, 0, 0, NULL);
         AssertThrow(ierr == 0, ExcPETScError(ierr));
       }
-#  else
-      {
-        AssertThrow(false, ExcMessage("Corner selection in BDDC is only available with PETSc 3.9.0 or newer");
-      }
-#  endif
 
 
     ierr = PCSetFromOptions(pc);
     AssertThrow(ierr == 0, ExcPETScError(ierr));
+#  else
+    {
+        AssertThrow(false, ExcMessage("BDDC preconditioner requires PETSc 3.10.0 or newer");
+    }
+#  endif
   }
 
   template <int dim>
index a627dc1daa15e0d25906a6781ea2c3a41262c6f7..cb32bf2ab0c359ed0a83a60f60f71882ce504109 100644 (file)
@@ -63,6 +63,7 @@ main(int argc, char *argv[])
 
   initlog();
   deallog << std::setprecision(10);
+#if DEAL_II_PETSC_VERSION_GTE(3, 10, 0)
   parallel::distributed::Triangulation<2> triangulation(MPI_COMM_WORLD);
   FE_Q<2>                                 fe(1);
   DoFHandler<2>                           dof_handler(triangulation);
@@ -160,19 +161,19 @@ main(int argc, char *argv[])
   PETScWrappers::PreconditionBDDC<2> preconditioner;
   PETScWrappers::PreconditionBDDC<2>::AdditionalData data;
 
-// Now we setup the dof coordinates if a sufficiently new PETSc is used
-#if DEAL_II_PETSC_VERSION_GTE(3, 9, 0)
+  // Now we setup the dof coordinates if a sufficiently new PETSc is used
   std::map<types::global_dof_index, Point<2>> dof_2_point;
   DoFTools::map_dofs_to_support_points(MappingQ1<2>(),
                                        dof_handler,
                                        dof_2_point);
-  std::vector<Point<2>> coords(locally_active_dofs.size());
-  for (const auto &d2p : dof_2_point)
+  std::vector<Point<2>> coords(locally_owned_dofs.n_elements());
+  unsigned int          k = 0;
+  for (auto it = locally_owned_dofs.begin(); it != locally_owned_dofs.end();
+       ++it, ++k)
     {
-      coords[d2p.first] = d2p.second;
+      coords[k] = dof_2_point[*it];
     }
-  data.coords = coords;
-#endif
+  data.coords       = coords;
   data.use_vertices = true;
 
   preconditioner.initialize(system_matrix, data);
@@ -185,6 +186,10 @@ main(int argc, char *argv[])
                             2);
 
   deallog << "CG/BDDC:OK" << std::endl;
-
+#else
+  deallog << "MATIS:OK" << std::endl;
+  deallog << "DEAL::Solver stopped within 1 - 2 iterations" << std::endl;
+  deallog << "CG/BDDC:OK" << std::endl;
+#endif
   return 0;
 }

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.