]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make variable names uniform. 3431/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 3 Nov 2016 16:39:47 +0000 (10:39 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 3 Nov 2016 16:39:47 +0000 (10:39 -0600)
Specifically, in DataPostprocessor-derived classes in the tutorial programs,
use the names of function arguments that we also use in DataPostprocessor.

examples/step-29/step-29.cc
examples/step-32/step-32.cc
examples/step-33/step-33.cc
examples/step-47/step-47.cc

index 25d9b4ab7e1236b01688ed4d504699480a115346..155a34c63e1b3faf8fb99b115c8b3f1a34a4d87e 100644 (file)
@@ -300,9 +300,9 @@ namespace Step29
 
     virtual
     void
-    compute_derived_quantities_vector (const std::vector<Vector<double> >               &uh,
-                                       const std::vector<std::vector<Tensor<1, dim> > > &duh,
-                                       const std::vector<std::vector<Tensor<2, dim> > > &dduh,
+    compute_derived_quantities_vector (const std::vector<Vector<double> >               &solution_values,
+                                       const std::vector<std::vector<Tensor<1, dim> > > &solution_gradients,
+                                       const std::vector<std::vector<Tensor<2, dim> > > &solution_hessians,
                                        const std::vector<Point<dim> >                   &normals,
                                        const std::vector<Point<dim> >                   &evaluation_points,
                                        std::vector<Vector<double> >                     &computed_quantities) const;
@@ -347,16 +347,16 @@ namespace Step29
   template <int dim>
   void
   ComputeIntensity<dim>::compute_derived_quantities_vector (
-    const std::vector<Vector<double> >                 &uh,
-    const std::vector<std::vector<Tensor<1, dim> > >   & /*duh*/,
-    const std::vector<std::vector<Tensor<2, dim> > >   & /*dduh*/,
+    const std::vector<Vector<double> >                 &solution_values,
+    const std::vector<std::vector<Tensor<1, dim> > >   & /*solution_gradients*/,
+    const std::vector<std::vector<Tensor<2, dim> > >   & /*solution_hessians*/,
     const std::vector<Point<dim> >                     & /*normals*/,
     const std::vector<Point<dim> >                     & /*evaluation_points*/,
     std::vector<Vector<double> >                       &computed_quantities
   ) const
   {
-    Assert(computed_quantities.size() == uh.size(),
-           ExcDimensionMismatch (computed_quantities.size(), uh.size()));
+    Assert(computed_quantities.size() == solution_values.size(),
+           ExcDimensionMismatch (computed_quantities.size(), solution_values.size()));
 
     // The computation itself is straightforward: We iterate over each entry
     // in the output vector and compute $|u|$ from the corresponding values of
@@ -365,9 +365,9 @@ namespace Step29
       {
         Assert(computed_quantities[i].size() == 1,
                ExcDimensionMismatch (computed_quantities[i].size(), 1));
-        Assert(uh[i].size() == 2, ExcDimensionMismatch (uh[i].size(), 2));
+        Assert(solution_values[i].size() == 2, ExcDimensionMismatch (solution_values[i].size(), 2));
 
-        computed_quantities[i](0) = std::sqrt(uh[i](0)*uh[i](0) + uh[i](1)*uh[i](1));
+        computed_quantities[i](0) = std::sqrt(solution_values[i](0)*solution_values[i](0) + solution_values[i](1)*solution_values[i](1));
       }
   }
 
index e1656ca381a86863fff53fd82bed279e49b06610..2bf00cf58e3e7e461c5003fa7968f464fda058f9 100644 (file)
@@ -3150,9 +3150,9 @@ namespace Step32
 
     virtual
     void
-    compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                       const std::vector<std::vector<Tensor<1,dim> > > &duh,
-                                       const std::vector<std::vector<Tensor<2,dim> > > &dduh,
+    compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                       const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
+                                       const std::vector<std::vector<Tensor<2,dim> > > &solution_hessians,
                                        const std::vector<Point<dim> >                  &normals,
                                        const std::vector<Point<dim> >                  &evaluation_points,
                                        std::vector<Vector<double> >                    &computed_quantities) const;
@@ -3242,33 +3242,33 @@ namespace Step32
   template <int dim>
   void
   BoussinesqFlowProblem<dim>::Postprocessor::
-  compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                     const std::vector<std::vector<Tensor<1,dim> > > &duh,
-                                     const std::vector<std::vector<Tensor<2,dim> > > &/*dduh*/,
+  compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                     const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
+                                     const std::vector<std::vector<Tensor<2,dim> > > &/*solution_hessians*/,
                                      const std::vector<Point<dim> >                  &/*normals*/,
                                      const std::vector<Point<dim> >                  &/*evaluation_points*/,
                                      std::vector<Vector<double> >                    &computed_quantities) const
   {
-    const unsigned int n_quadrature_points = uh.size();
-    Assert (duh.size() == n_quadrature_points,                  ExcInternalError());
+    const unsigned int n_quadrature_points = solution_values.size();
+    Assert (solution_gradients.size() == n_quadrature_points,                  ExcInternalError());
     Assert (computed_quantities.size() == n_quadrature_points,  ExcInternalError());
-    Assert (uh[0].size() == dim+2,                              ExcInternalError());
+    Assert (solution_values[0].size() == dim+2,                              ExcInternalError());
 
     for (unsigned int q=0; q<n_quadrature_points; ++q)
       {
         for (unsigned int d=0; d<dim; ++d)
           computed_quantities[q](d)
-            = (uh[q](d) *  EquationData::year_in_seconds * 100);
+            = (solution_values[q](d) *  EquationData::year_in_seconds * 100);
 
-        const double pressure = (uh[q](dim)-minimal_pressure);
+        const double pressure = (solution_values[q](dim)-minimal_pressure);
         computed_quantities[q](dim) = pressure;
 
-        const double temperature = uh[q](dim+1);
+        const double temperature = solution_values[q](dim+1);
         computed_quantities[q](dim+1) = temperature;
 
         Tensor<2,dim> grad_u;
         for (unsigned int d=0; d<dim; ++d)
-          grad_u[d] = duh[q][d];
+          grad_u[d] = solution_gradients[q][d];
         const SymmetricTensor<2,dim> strain_rate = symmetrize (grad_u);
         computed_quantities[q](dim+2) = 2 * EquationData::eta *
                                         strain_rate * strain_rate;
index b721d7ab867fdd0f81324c38cfd4d3181aa284f0..454934fa75357af121096a7da6cbaa5970ad6607 100644 (file)
@@ -543,9 +543,9 @@ namespace Step33
 
       virtual
       void
-      compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                         const std::vector<std::vector<Tensor<1,dim> > > &duh,
-                                         const std::vector<std::vector<Tensor<2,dim> > > &dduh,
+      compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                         const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
+                                         const std::vector<std::vector<Tensor<2,dim> > > &solution_hessians,
                                          const std::vector<Point<dim> >                  &normals,
                                          const std::vector<Point<dim> >                  &evaluation_points,
                                          std::vector<Vector<double> >                    &computed_quantities) const;
@@ -584,14 +584,14 @@ namespace Step33
   // quadrature point is itself vector-valued, namely the conserved
   // variables. What we're going to do here is to compute the quantities we're
   // interested in at each quadrature point. Note that for this we can ignore
-  // the Hessians ("dduh") and normal vectors; to avoid compiler warnings
+  // the Hessians ("solution_hessians") and normal vectors; to avoid compiler warnings
   // about unused variables, we comment out their names.
   template <int dim>
   void
   EulerEquations<dim>::Postprocessor::
-  compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                     const std::vector<std::vector<Tensor<1,dim> > > &duh,
-                                     const std::vector<std::vector<Tensor<2,dim> > > &/*dduh*/,
+  compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                     const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
+                                     const std::vector<std::vector<Tensor<2,dim> > > &/*solution_hessians*/,
                                      const std::vector<Point<dim> >                  &/*normals*/,
                                      const std::vector<Point<dim> >                  &/*evaluation_points*/,
                                      std::vector<Vector<double> >                    &computed_quantities) const
@@ -599,21 +599,21 @@ namespace Step33
     // At the beginning of the function, let us make sure that all variables
     // have the correct sizes, so that we can access individual vector
     // elements without having to wonder whether we might read or write
-    // invalid elements; we also check that the <code>duh</code> vector only
+    // invalid elements; we also check that the <code>solution_gradients</code> vector only
     // contains data if we really need it (the system knows about this because
     // we say so in the <code>get_needed_update_flags()</code> function
     // below). For the inner vectors, we check that at least the first element
     // of the outer vector has the correct inner size:
-    const unsigned int n_quadrature_points = uh.size();
+    const unsigned int n_quadrature_points = solution_values.size();
 
     if (do_schlieren_plot == true)
-      Assert (duh.size() == n_quadrature_points,
+      Assert (solution_gradients.size() == n_quadrature_points,
               ExcInternalError());
 
     Assert (computed_quantities.size() == n_quadrature_points,
             ExcInternalError());
 
-    Assert (uh[0].size() == n_components,
+    Assert (solution_values[0].size() == n_components,
             ExcInternalError());
 
     if (do_schlieren_plot == true)
@@ -630,17 +630,17 @@ namespace Step33
     // <code>density_component</code> information:
     for (unsigned int q=0; q<n_quadrature_points; ++q)
       {
-        const double density = uh[q](density_component);
+        const double density = solution_values[q](density_component);
 
         for (unsigned int d=0; d<dim; ++d)
           computed_quantities[q](d)
-            = uh[q](first_momentum_component+d) / density;
+            = solution_values[q](first_momentum_component+d) / density;
 
-        computed_quantities[q](dim) = compute_pressure (uh[q]);
+        computed_quantities[q](dim) = compute_pressure (solution_values[q]);
 
         if (do_schlieren_plot == true)
-          computed_quantities[q](dim+1) = duh[q][density_component] *
-                                          duh[q][density_component];
+          computed_quantities[q](dim+1) = solution_gradients[q][density_component] *
+                                          solution_gradients[q][density_component];
       }
   }
 
index 6959040a4c5c395c7fb4cff42a8d354886233793..eceafa87f715eefbfbd69d653388eba59afb85f0 100644 (file)
@@ -905,9 +905,9 @@ namespace Step47
   public:
     virtual
     void
-    compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                       const std::vector<std::vector<Tensor<1,dim> > > &duh,
-                                       const std::vector<std::vector<Tensor<2,dim> > > &dduh,
+    compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                       const std::vector<std::vector<Tensor<1,dim> > > &solution_gradients,
+                                       const std::vector<std::vector<Tensor<2,dim> > > &solution_hessians,
                                        const std::vector<Point<dim> >                  &normals,
                                        const std::vector<Point<dim> >                  &evaluation_points,
                                        std::vector<Vector<double> >                    &computed_quantities) const;
@@ -955,24 +955,24 @@ namespace Step47
   template <int dim>
   void
   Postprocessor<dim>::
-  compute_derived_quantities_vector (const std::vector<Vector<double> >              &uh,
-                                     const std::vector<std::vector<Tensor<1,dim> > > &/*duh*/,
-                                     const std::vector<std::vector<Tensor<2,dim> > > &/*dduh*/,
+  compute_derived_quantities_vector (const std::vector<Vector<double> >              &solution_values,
+                                     const std::vector<std::vector<Tensor<1,dim> > > &/*solution_gradients*/,
+                                     const std::vector<std::vector<Tensor<2,dim> > > &/*solution_hessians*/,
                                      const std::vector<Point<dim> >                  &/*normals*/,
                                      const std::vector<Point<dim> >                  &evaluation_points,
                                      std::vector<Vector<double> >                    &computed_quantities) const
   {
-    const unsigned int n_quadrature_points = uh.size();
+    const unsigned int n_quadrature_points = solution_values.size();
     Assert (computed_quantities.size() == n_quadrature_points,  ExcInternalError());
-    Assert (uh[0].size() == 2,                                  ExcInternalError());
+    Assert (solution_values[0].size() == 2,                                  ExcInternalError());
 
     for (unsigned int q=0; q<n_quadrature_points; ++q)
       {
         computed_quantities[q](0)
-          = (uh[q](0)
+          = (solution_values[q](0)
              +
 //TODO: shift in weight function is missing!
-             uh[q](1) * std::fabs(level_set(evaluation_points[q])));
+             solution_values[q](1) * std::fabs(level_set(evaluation_points[q])));
         computed_quantities[q](1)
           = (computed_quantities[q](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.