]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate MatrixFree::reinit() functions with default Mapping argument 11169/head
authorPeter Munch <peterrmuench@gmail.com>
Wed, 11 Nov 2020 22:11:17 +0000 (23:11 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 13 Nov 2020 13:24:44 +0000 (14:24 +0100)
doc/news/changes/incompatibilities/20201111Munch [new file with mode: 0644]
examples/step-37/step-37.cc
examples/step-48/step-48.cc
examples/step-50/step-50.cc
examples/step-59/step-59.cc
include/deal.II/matrix_free/matrix_free.h

diff --git a/doc/news/changes/incompatibilities/20201111Munch b/doc/news/changes/incompatibilities/20201111Munch
new file mode 100644 (file)
index 0000000..8acc0cf
--- /dev/null
@@ -0,0 +1,4 @@
+Deprecated: All `MatrixFree::reinit()` functions without `Mapping` as argument 
+have been deprecated. Use the functions that explicit provide the Mapping instead.
+<br>
+(Peter Munch, 2020/11/11)
index 7a996300583d5db923de7f60ae00b93ea33ccead..663bd0f16073d00475298da85cf8759dfa484e8c 100644 (file)
@@ -688,6 +688,8 @@ namespace Step37
     FE_Q<dim>       fe;
     DoFHandler<dim> dof_handler;
 
+    MappingQ1<dim> mapping;
+
     AffineConstraints<double> constraints;
     using SystemMatrixType =
       LaplaceOperator<dim, degree_finite_element, double>;
@@ -796,10 +798,8 @@ namespace Step37
     constraints.clear();
     constraints.reinit(locally_relevant_dofs);
     DoFTools::make_hanging_node_constraints(dof_handler, constraints);
-    VectorTools::interpolate_boundary_values(dof_handler,
-                                             0,
-                                             Functions::ZeroFunction<dim>(),
-                                             constraints);
+    VectorTools::interpolate_boundary_values(
+      mapping, dof_handler, 0, Functions::ZeroFunction<dim>(), constraints);
     constraints.close();
     setup_time += time.wall_time();
     time_details << "Distribute DoFs & B.C.     (CPU/wall) " << time.cpu_time()
@@ -814,7 +814,8 @@ namespace Step37
         (update_gradients | update_JxW_values | update_quadrature_points);
       std::shared_ptr<MatrixFree<dim, double>> system_mf_storage(
         new MatrixFree<dim, double>());
-      system_mf_storage->reinit(dof_handler,
+      system_mf_storage->reinit(mapping,
+                                dof_handler,
                                 constraints,
                                 QGauss<1>(fe.degree + 1),
                                 additional_data);
@@ -869,7 +870,8 @@ namespace Step37
         additional_data.mg_level = level;
         std::shared_ptr<MatrixFree<dim, float>> mg_mf_storage_level(
           new MatrixFree<dim, float>());
-        mg_mf_storage_level->reinit(dof_handler,
+        mg_mf_storage_level->reinit(mapping,
+                                    dof_handler,
                                     level_constraints,
                                     QGauss<1>(fe.degree + 1),
                                     additional_data);
@@ -1126,7 +1128,7 @@ namespace Step37
     solution.update_ghost_values();
     data_out.attach_dof_handler(dof_handler);
     data_out.add_data_vector(solution, "solution");
-    data_out.build_patches();
+    data_out.build_patches(mapping);
 
     DataOutBase::VtkFlags flags;
     flags.compression_level = DataOutBase::VtkFlags::best_speed;
index 5006e9433ab516b5f0ae5dc9d9c61394cbfcc814..7d455339a1ca3fec817150e012cf94a31df0fa40 100644 (file)
@@ -308,8 +308,11 @@ namespace Step48
 #else
     Triangulation<dim> triangulation;
 #endif
-    FE_Q<dim>                 fe;
-    DoFHandler<dim>           dof_handler;
+    FE_Q<dim>       fe;
+    DoFHandler<dim> dof_handler;
+
+    MappingQ1<dim> mapping;
+
     AffineConstraints<double> constraints;
     IndexSet                  locally_relevant_dofs;
 
@@ -434,7 +437,8 @@ namespace Step48
     additional_data.tasks_parallel_scheme =
       MatrixFree<dim>::AdditionalData::TasksParallelScheme::partition_partition;
 
-    matrix_free_data.reinit(dof_handler,
+    matrix_free_data.reinit(mapping,
+                            dof_handler,
                             constraints,
                             QGaussLobatto<1>(fe_degree + 1),
                             additional_data);
@@ -479,7 +483,8 @@ namespace Step48
 
     Vector<float> norm_per_cell(triangulation.n_active_cells());
     solution.update_ghost_values();
-    VectorTools::integrate_difference(dof_handler,
+    VectorTools::integrate_difference(mapping,
+                                      dof_handler,
                                       solution,
                                       Functions::ZeroFunction<dim>(),
                                       norm_per_cell,
@@ -498,7 +503,7 @@ namespace Step48
 
     data_out.attach_dof_handler(dof_handler);
     data_out.add_data_vector(solution, "solution");
-    data_out.build_patches();
+    data_out.build_patches(mapping);
 
     data_out.write_vtu_with_pvtu_record(
       "./", "solution", timestep_number, MPI_COMM_WORLD, 3);
@@ -563,10 +568,12 @@ namespace Step48
     // get later consumed by the SineGordonOperation::apply() function. Next,
     // an instance of the <code> SineGordonOperation class </code> based on
     // the finite element degree specified at the top of this file is set up.
-    VectorTools::interpolate(dof_handler,
+    VectorTools::interpolate(mapping,
+                             dof_handler,
                              InitialCondition<dim>(1, time),
                              solution);
-    VectorTools::interpolate(dof_handler,
+    VectorTools::interpolate(mapping,
+                             dof_handler,
                              InitialCondition<dim>(1, time - time_step),
                              old_solution);
     output_results(0);
index a949ce81c92e2068d8ef75969d9da9136ed7656d..349b8d93a9717a62c9c47a64baef2b867635abdb 100644 (file)
@@ -513,7 +513,8 @@ void LaplaceProblem<dim, degree>::setup_system()
             (update_gradients | update_JxW_values | update_quadrature_points);
           std::shared_ptr<MatrixFree<dim, double>> mf_storage =
             std::make_shared<MatrixFree<dim, double>>();
-          mf_storage->reinit(dof_handler,
+          mf_storage->reinit(mapping,
+                             dof_handler,
                              constraints,
                              QGauss<1>(degree + 1),
                              additional_data);
@@ -615,7 +616,8 @@ void LaplaceProblem<dim, degree>::setup_multigrid()
               additional_data.mg_level = level;
               std::shared_ptr<MatrixFree<dim, float>> mf_storage_level(
                 new MatrixFree<dim, float>());
-              mf_storage_level->reinit(dof_handler,
+              mf_storage_level->reinit(mapping,
+                                       dof_handler,
                                        level_constraints,
                                        QGauss<1>(degree + 1),
                                        additional_data);
index 4df20ebcd7f4309dd1652112deb53485053d7a0a..8d03e3c96c49b18ee79910ad1ff1fc653315cf64 100644 (file)
@@ -931,6 +931,8 @@ namespace Step59
     FE_DGQHermite<dim> fe;
     DoFHandler<dim>    dof_handler;
 
+    MappingQ1<dim> mapping;
+
     using SystemMatrixType = LaplaceOperator<dim, fe_degree, double>;
     SystemMatrixType system_matrix;
 
@@ -1022,10 +1024,8 @@ namespace Step59
          update_quadrature_points);
       const auto system_mf_storage =
         std::make_shared<MatrixFree<dim, double>>();
-      system_mf_storage->reinit(dof_handler,
-                                dummy,
-                                QGauss<1>(fe.degree + 1),
-                                additional_data);
+      system_mf_storage->reinit(
+        mapping, dof_handler, dummy, QGauss<1>(fe.degree + 1), additional_data);
       system_matrix.initialize(system_mf_storage);
     }
 
@@ -1054,7 +1054,8 @@ namespace Step59
         additional_data.mg_level = level;
         const auto mg_mf_storage_level =
           std::make_shared<MatrixFree<dim, float>>();
-        mg_mf_storage_level->reinit(dof_handler,
+        mg_mf_storage_level->reinit(mapping,
+                                    dof_handler,
                                     dummy,
                                     QGauss<1>(fe.degree + 1),
                                     additional_data);
@@ -1290,7 +1291,8 @@ namespace Step59
   void LaplaceProblem<dim, fe_degree>::analyze_results() const
   {
     Vector<float> error_per_cell(triangulation.n_active_cells());
-    VectorTools::integrate_difference(dof_handler,
+    VectorTools::integrate_difference(mapping,
+                                      dof_handler,
                                       solution,
                                       Solution<dim>(),
                                       error_per_cell,
index 7dadca5886381be04e40b4f6baf6056b940ad603..dce6e4b812a3df15e1b33e3cadbdc7a7eac66495 100644 (file)
@@ -591,9 +591,11 @@ public:
   /**
    * Initializes the data structures. Same as above, but using a $Q_1$
    * mapping.
+   *
+   * @deprecated Use the overload taking a Mapping object instead.
    */
   template <typename QuadratureType, typename number2>
-  void
+  DEAL_II_DEPRECATED void
   reinit(const DoFHandler<dim> &           dof_handler,
          const AffineConstraints<number2> &constraint,
          const QuadratureType &            quad,
@@ -644,9 +646,11 @@ public:
   /**
    * Initializes the data structures. Same as above, but  using a $Q_1$
    * mapping.
+   *
+   * @deprecated Use the overload taking a Mapping object instead.
    */
   template <typename QuadratureType, typename number2>
-  void
+  DEAL_II_DEPRECATED void
   reinit(const std::vector<const DoFHandler<dim> *> &           dof_handler,
          const std::vector<const AffineConstraints<number2> *> &constraint,
          const std::vector<QuadratureType> &                    quad,
@@ -658,7 +662,7 @@ public:
    * @deprecated Use the overload taking a DoFHandler object instead.
    */
   template <typename QuadratureType, typename number2>
-  void
+  DEAL_II_DEPRECATED void
   reinit(const std::vector<const hp::DoFHandler<dim> *> &       dof_handler_hp,
          const std::vector<const AffineConstraints<number2> *> &constraint,
          const std::vector<QuadratureType> &                    quad,
@@ -695,9 +699,11 @@ public:
   /**
    * Initializes the data structures. Same as above, but  using a $Q_1$
    * mapping.
+   *
+   * @deprecated Use the overload taking a Mapping object instead.
    */
   template <typename QuadratureType, typename number2>
-  void
+  DEAL_II_DEPRECATED void
   reinit(const std::vector<const DoFHandler<dim> *> &           dof_handler,
          const std::vector<const AffineConstraints<number2> *> &constraint,
          const QuadratureType &                                 quad,

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.