]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fixed bddc test, added points as well
authorNicolas Barnafi <nabw91@gmail.com>
Wed, 7 Sep 2022 07:51:48 +0000 (09:51 +0200)
committerNicolas Barnafi <nabw91@gmail.com>
Wed, 7 Sep 2022 07:51:48 +0000 (09:51 +0200)
include/deal.II/lac/petsc_precondition.h
source/lac/petsc_precondition.cc
tests/petsc/bddc.cc

index 70416e60f08bf640f716caa40b054218fbe8af99..6ac0b8a1f9ada22b1ec274dcee50d7aeccb23338 100644 (file)
@@ -989,11 +989,11 @@ namespace PETScWrappers
        * Constructor. Note that BDDC offers a lot more options to set
        * than what is exposed here.
        */
-      AdditionalData(const bool              use_vertices = true,
-                     const bool              use_edges    = false,
-                     const bool              use_faces    = false,
-                     const bool              symmetric    = false,
-                     std::vector<Point<dim>> coords       = {});
+      AdditionalData(const bool                    use_vertices = true,
+                     const bool                    use_edges    = false,
+                     const bool                    use_faces    = false,
+                     const bool                    symmetric    = false,
+                     const std::vector<Point<dim>> coords       = {});
 
       /**
        * This flag sets the use of degrees of freedom in the vertices of the
@@ -1024,7 +1024,7 @@ namespace PETScWrappers
        * Set the location of each DoF. This helps in improving the definition of
        * the vertices for unstructured meshes.
        */
-      const std::vector<Point<dim>> coords;
+      std::vector<Point<dim>> coords;
     };
 
     /**
@@ -1079,6 +1079,8 @@ namespace PETScWrappers
   using PreconditionerBase DEAL_II_DEPRECATED = PreconditionBase;
 } // namespace PETScWrappers
 
+template class PETScWrappers::PreconditionBDDC<2>;
+template class PETScWrappers::PreconditionBDDC<3>;
 
 DEAL_II_NAMESPACE_CLOSE
 
index 7950579cc64687b021fdb79b1a7a268120772ec6..467fa4ed43753d9d412d09167217a682474648a3 100644 (file)
@@ -912,8 +912,7 @@ namespace PETScWrappers
       set_option_value("-pc_bddc_symmetric", "true");
     else
       set_option_value("-pc_bddc_symmetric", "false");
-    if (additional_data.coords_cdim && additional_data.coords_n &&
-        additional_data.coords_data)
+    if (additional_data.coords.size() > 0)
       {
         set_option_value("-pc_bddc_corner_selection", "true");
         // Convert coords vector to PETSc data array
index b5778ad23e5d8ac0503674dcf9b3cac008ab3f46..1dfc9fb76e34971f0334400508858661cd40d8f9 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2016 - 2021 by the deal.II authors
+// Copyright (C) 2016 - 2022 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
 //
 // ---------------------------------------------------------------------
 
-// Test that it is possible to instantiate a LinearOperator object for all
-// different kinds of PETSc matrices and vectors
-// TODO: A bit more tests...
+// This test provides a minimal Laplace equation in 2D for testing
+// the BDDC preconditioner. It cannot be simpler as it fails to setup
+// in simpler matrices (i.e. diagonal), and the preset matrix factories
+// create non-IS matrices, which are not supported by BDDC.
 
 #include <deal.II/base/index_set.h>
 
@@ -44,6 +45,7 @@
 
 #include <deal.II/fe/fe_q.h>
 #include <deal.II/fe/fe_values.h>
+#include <deal.II/fe/mapping_q1.h>
 
 #include <deal.II/grid/grid_generator.h>
 
@@ -87,12 +89,12 @@ main(int argc, char *argv[])
                                              locally_relevant_dofs);
 
   PETScWrappers::MPI::SparseMatrix system_matrix;
-  system_matrix.reinit_IS(locally_owned_dofs,
-                          locally_active_dofs,
-                          locally_owned_dofs,
-                          locally_active_dofs,
-                          dsp,
-                          MPI_COMM_WORLD);
+  system_matrix.reinit(locally_owned_dofs,
+                       locally_active_dofs,
+                       locally_owned_dofs,
+                       locally_active_dofs,
+                       dsp,
+                       MPI_COMM_WORLD);
   deallog << "MATIS:OK" << std::endl;
   PETScWrappers::MPI::Vector locally_relevant_solution(locally_owned_dofs,
                                                        locally_relevant_dofs,
@@ -154,9 +156,22 @@ main(int argc, char *argv[])
 
   SolverControl solver_control(dof_handler.n_dofs(), 1e-12);
 
-  PETScWrappers::SolverCG         solver(solver_control, MPI_COMM_WORLD);
-  PETScWrappers::PreconditionBDDC preconditioner(system_matrix);
-
+  PETScWrappers::SolverCG            solver(solver_control, MPI_COMM_WORLD);
+  PETScWrappers::PreconditionBDDC<2> preconditioner;
+  PETScWrappers::PreconditionBDDC<2>::AdditionalData data;
+
+  // Now we setup the dof coordinates
+  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)
+    {
+      coords[d2p.first] = d2p.second;
+    }
+
+  preconditioner.initialize(system_matrix, data);
   check_solver_within_range(solver.solve(system_matrix,
                                          completely_distributed_solution,
                                          system_rhs,

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.