]> https://gitweb.dealii.org/ - dealii.git/commitdiff
step-31: C++ modernization
authorDavid Wells <drwells@email.unc.edu>
Thu, 9 May 2019 15:58:36 +0000 (11:58 -0400)
committerDavid Wells <drwells@email.unc.edu>
Thu, 9 May 2019 17:43:03 +0000 (13:43 -0400)
examples/step-31/step-31.cc

index 857dee792a2dc2af53649a1b70946d2fb97c73ce..1b3ee08a84dcd86b88db307201a150d665583fd4 100644 (file)
@@ -119,10 +119,10 @@ namespace Step31
   // thermal expansion coefficient $\beta$):
   namespace EquationData
   {
-    const double eta     = 1;
-    const double kappa   = 1e-6;
-    const double beta    = 10;
-    const double density = 1;
+    constexpr double eta     = 1;
+    constexpr double kappa   = 1e-6;
+    constexpr double beta    = 10;
+    constexpr double density = 1;
 
 
     template <int dim>
@@ -639,11 +639,7 @@ namespace Step31
 
     const FEValuesExtractors::Vector velocities(0);
 
-    typename DoFHandler<dim>::active_cell_iterator cell = stokes_dof_handler
-                                                            .begin_active(),
-                                                   endc =
-                                                     stokes_dof_handler.end();
-    for (; cell != endc; ++cell)
+    for (const auto &cell : stokes_dof_handler.active_cell_iterators())
       {
         fe_values.reinit(cell);
         fe_values[velocities].get_function_values(stokes_solution,
@@ -699,10 +695,7 @@ namespace Step31
         double min_temperature = std::numeric_limits<double>::max(),
                max_temperature = -std::numeric_limits<double>::max();
 
-        typename DoFHandler<dim>::active_cell_iterator
-          cell = temperature_dof_handler.begin_active(),
-          endc = temperature_dof_handler.end();
-        for (; cell != endc; ++cell)
+        for (const auto &cell : temperature_dof_handler.active_cell_iterators())
           {
             fe_values.reinit(cell);
             fe_values.get_function_values(old_temperature_solution,
@@ -728,10 +721,7 @@ namespace Step31
         double min_temperature = std::numeric_limits<double>::max(),
                max_temperature = -std::numeric_limits<double>::max();
 
-        typename DoFHandler<dim>::active_cell_iterator
-          cell = temperature_dof_handler.begin_active(),
-          endc = temperature_dof_handler.end();
-        for (; cell != endc; ++cell)
+        for (const auto &cell : temperature_dof_handler.active_cell_iterators())
           {
             fe_values.reinit(cell);
             fe_values.get_function_values(old_temperature_solution,
@@ -790,8 +780,8 @@ namespace Step31
     const double                       global_T_variation,
     const double                       cell_diameter) const
   {
-    const double beta  = 0.017 * dim;
-    const double alpha = 1;
+    constexpr double beta  = 0.017 * dim;
+    constexpr double alpha = 1.0;
 
     if (global_u_infty == 0)
       return 5e-3 * cell_diameter;
@@ -1091,11 +1081,7 @@ namespace Step31
     const FEValuesExtractors::Vector velocities(0);
     const FEValuesExtractors::Scalar pressure(dim);
 
-    typename DoFHandler<dim>::active_cell_iterator cell = stokes_dof_handler
-                                                            .begin_active(),
-                                                   endc =
-                                                     stokes_dof_handler.end();
-    for (; cell != endc; ++cell)
+    for (const auto &cell : stokes_dof_handler.active_cell_iterators())
       {
         stokes_fe_values.reinit(cell);
         local_matrix = 0;
@@ -1330,12 +1316,9 @@ namespace Step31
     // the update flags, zeroing out the local arrays and getting the values
     // of the old solution at the quadrature points. Then we are ready to loop
     // over the quadrature points on the cell.
-    typename DoFHandler<dim>::active_cell_iterator cell = stokes_dof_handler
-                                                            .begin_active(),
-                                                   endc =
-                                                     stokes_dof_handler.end();
-    typename DoFHandler<dim>::active_cell_iterator temperature_cell =
-      temperature_dof_handler.begin_active();
+    auto cell             = stokes_dof_handler.begin_active(),
+         endc             = stokes_dof_handler.end();
+    auto temperature_cell = temperature_dof_handler.begin_active();
 
     for (; cell != endc; ++cell, ++temperature_cell)
       {
@@ -1474,10 +1457,7 @@ namespace Step31
     // <code>EquationData::kappa</code>. Finally, we let the constraints
     // object insert these values into the global matrix, and directly
     // condense the constraints into the matrix.
-    typename DoFHandler<dim>::active_cell_iterator
-      cell = temperature_dof_handler.begin_active(),
-      endc = temperature_dof_handler.end();
-    for (; cell != endc; ++cell)
+    for (const auto &cell : temperature_dof_handler.active_cell_iterators())
       {
         local_mass_matrix      = 0;
         local_stiffness_matrix = 0;
@@ -1615,11 +1595,9 @@ namespace Step31
     // Stokes part we restrict ourselves to extracting the velocity part (and
     // ignoring the pressure part) by using
     // <code>stokes_fe_values[velocities].get_function_values</code>.
-    typename DoFHandler<dim>::active_cell_iterator
-      cell = temperature_dof_handler.begin_active(),
-      endc = temperature_dof_handler.end();
-    typename DoFHandler<dim>::active_cell_iterator stokes_cell =
-      stokes_dof_handler.begin_active();
+    auto cell        = temperature_dof_handler.begin_active(),
+         endc        = temperature_dof_handler.end();
+    auto stokes_cell = stokes_dof_handler.begin_active();
 
     for (; cell != endc; ++cell, ++stokes_cell)
       {
@@ -1972,22 +1950,19 @@ namespace Step31
   {
     Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
 
-    KellyErrorEstimator<dim>::estimate(
-      temperature_dof_handler,
-      QGauss<dim - 1>(temperature_degree + 1),
-      std::map<types::boundary_id, const Function<dim> *>(),
-      temperature_solution,
-      estimated_error_per_cell);
+    KellyErrorEstimator<dim>::estimate(temperature_dof_handler,
+                                       QGauss<dim - 1>(temperature_degree + 1),
+                                       {},
+                                       temperature_solution,
+                                       estimated_error_per_cell);
 
     GridRefinement::refine_and_coarsen_fixed_fraction(triangulation,
                                                       estimated_error_per_cell,
                                                       0.8,
                                                       0.1);
     if (triangulation.n_levels() > max_grid_level)
-      for (typename Triangulation<dim>::active_cell_iterator cell =
-             triangulation.begin_active(max_grid_level);
-           cell != triangulation.end();
-           ++cell)
+      for (auto &cell :
+           triangulation.active_cell_iterators_on_level(max_grid_level))
         cell->clear_refine_flag();
 
     // As part of mesh refinement we need to transfer the solution vectors

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.