]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Better document FunctionMap::type.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 11 Feb 2014 02:31:17 +0000 (02:31 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 11 Feb 2014 02:31:17 +0000 (02:31 +0000)
git-svn-id: https://svn.dealii.org/trunk@32452 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/examples/step-16/step-16.cc
deal.II/examples/step-6/step-6.cc
deal.II/examples/step-7/step-7.cc
deal.II/examples/step-8/step-8.cc
deal.II/include/deal.II/dofs/function_map.h

index 67b1a5b4138395c8c318c0c639ae2cc96bbcd69e..1eb7c93fac7bfb2b19acb717c7a3bb31fffbf31f 100644 (file)
@@ -273,7 +273,7 @@ namespace Step16
     // the global matrix from local contributions. This works, but the same
     // can be done in a slightly simpler way if we already take care of these
     // constraints at the time of copying local contributions into the global
-    // matrix. To this end, we here do not just compute the constraints do to
+    // matrix. To this end, we here do not just compute the constraints due to
     // hanging nodes, but also due to zero boundary conditions. We will use
     // this set of constraints later on to help us copy local contributions
     // correctly into the global linear system right away, without the need
@@ -283,11 +283,11 @@ namespace Step16
     DoFTools::make_hanging_node_constraints (mg_dof_handler, hanging_node_constraints);
     DoFTools::make_hanging_node_constraints (mg_dof_handler, constraints);
 
-    typename FunctionMap<dim>::type      dirichlet_boundary;
+    typename FunctionMap<dim>::type      dirichlet_boundary_functions;
     ZeroFunction<dim>                    homogeneous_dirichlet_bc (1);
-    dirichlet_boundary[0] = &homogeneous_dirichlet_bc;
+    dirichlet_boundary_functions[0] = &homogeneous_dirichlet_bc;
     VectorTools::interpolate_boundary_values (static_cast<const DoFHandler<dim>&>(mg_dof_handler),
-                                              dirichlet_boundary,
+                                              dirichlet_boundary_functions,
                                               constraints);
     constraints.close ();
     hanging_node_constraints.close ();
@@ -299,7 +299,7 @@ namespace Step16
     // about the boundary values as well, so we pass the
     // <code>dirichlet_boundary</code> here as well.
     mg_constrained_dofs.clear();
-    mg_constrained_dofs.initialize(mg_dof_handler, dirichlet_boundary);
+    mg_constrained_dofs.initialize(mg_dof_handler, dirichlet_boundary_functions);
 
 
     // Now for the things that concern the multigrid data structures. First,
index 52dda1eb453ab2250fae4de3d96de1ca73e4c68c..a703376eaa92dd7b6f9c5d3cda8c81840ece8839 100644 (file)
@@ -530,20 +530,23 @@ void Step6<dim>::solve ()
 // quartic function, for which a 3 point Gauss formula is sufficient since it
 // integrates polynomials up to order 5 exactly.)
 //
-// Secondly, the function wants a list of boundaries where we have imposed
-// Neumann value, and the corresponding Neumann values. This information is
+// Secondly, the function wants a list of boundary indicators for those
+// boundaries where we have imposed Neumann values of the kind
+// $\partial_n u(\mathbf x) = h(\mathbf x)$, along with a function $h(\mathbf x)$
+// for each such boundary. This information is
 // represented by an object of type <code>FunctionMap::type</code> that is
-// essentially a map from boundary indicators to function objects describing
-// Neumann boundary values (in the present example program, we do not use
+// a typedef to a map from boundary indicators to function objects describing
+// the Neumann boundary values. In the present example program, we do not use
 // Neumann boundary values, so this map is empty, and in fact constructed
 // using the default constructor of the map in the place where the function
-// call expects the respective function argument).
+// call expects the respective function argument.
 //
 // The output, as mentioned is a vector of values for all cells. While it may
-// make sense to compute the *value* of a degree of freedom very accurately,
-// it is usually not helpful to compute the *error indicator* corresponding to
-// a cell particularly accurately. We therefore typically use a vector of
-// floats instead of a vector of doubles to represent error indicators.
+// make sense to compute the <b>value</b> of a solution degree of freedom
+// very accurately, it is usually not necessary to compute the <b>error indicator</b>
+// corresponding to the solution on a cell particularly accurately. We therefore
+// typically use a vector of floats instead of a vector of doubles to represent
+// error indicators.
 template <int dim>
 void Step6<dim>::refine_grid ()
 {
index 058e7cf1cdc05da7cf7e0b254d910de74a075622..dee41ef7e93980f2bc7982c14b2135287282972c 100644 (file)
@@ -775,7 +775,8 @@ namespace Step7
   // since we have Neumann boundary conditions on part of the boundaries, but
   // since we don't have a function here that describes the Neumann values (we
   // only construct these values from the exact solution when assembling the
-  // matrix), we omit this detail even though they would not be hard to add.
+  // matrix), we omit this detail even though doing this in a strictly correct
+  // way would not be hard to add.
   //
   // At the end of the switch, we have a default case that looks slightly
   // strange: an <code>Assert</code> statement with a <code>false</code>
@@ -809,10 +810,9 @@ namespace Step7
       {
         Vector<float> estimated_error_per_cell (triangulation.n_active_cells());
 
-        typename FunctionMap<dim>::type neumann_boundary;
         KellyErrorEstimator<dim>::estimate (dof_handler,
                                             QGauss<dim-1>(3),
-                                            neumann_boundary,
+                                            typename FunctionMap<dim>::type(),
                                             solution,
                                             estimated_error_per_cell);
 
index fbc9051c67387cb14690c9c09fbdd5e1d5d0dcf6..6fcd742e809e8d41fa182a2993dadb5bfab5d2b4 100644 (file)
@@ -610,10 +610,9 @@ namespace Step8
   {
     Vector<float> estimated_error_per_cell (triangulation.n_active_cells());
 
-    typename FunctionMap<dim>::type neumann_boundary;
     KellyErrorEstimator<dim>::estimate (dof_handler,
                                         QGauss<dim-1>(2),
-                                        neumann_boundary,
+                                        typename FunctionMap<dim>::type(),
                                         solution,
                                         estimated_error_per_cell);
 
index f5f8d7fcabce28fe9ac65ffbdf77fe2bd67b0285..52e3495a0daa7fb652f72db01a269218f5855534 100644 (file)
@@ -27,16 +27,40 @@ template <int spacedim> class Function;
 
 
 /**
- * Declare a data type which denotes a mapping between a boundary indicator
- * and the function denoting the boundary values on this part of the
- * boundary. This type is required in many functions where depending on the
- * boundary indicator, different functions are used. An example is boundary
- * value interpolation.
+ * This class declares a local typedef that denotes a mapping between a boundary indicator
+ * (see @ref GlossBoundaryIndicator) that is used to describe what kind of boundary
+ * condition holds on a particular piece of the boundary,
+ * and the function describing the actual function that provides the boundary
+ * values on this part of the boundary. This type is required in many functions
+ * in the library where, for example, we need to know about the functions $h_i(\mathbf x)$
+ * used in boundary conditions
+ * @f{align*}
+ *   \mathbf n \cdot \nabla u = h_i \qquad \qquad \text{on}\ \Gamma_i\subset\partial\Omega.
+ * @f}
+ * An example is the function KellyErrorEstimator::estimate() that allows us
+ * to provide a set of functions $h_i$ for all those boundary indicators $i$ for
+ * which the boundary condition is supposed to be of Neumann type. Of course,
+ * the same kind of principle can be applied to cases where we care about
+ * Dirichlet values, where one needs to provide a map from boundary indicator $i$
+ * to Dirichlet function $h_i$ if the boundary conditions are given as
+ * @f{align*}
+ *   u = h_i \qquad \qquad \text{on}\ \Gamma_i\subset\partial\Omega.
+ * @f}
+ * This is, for example, the case for the VectorTools::interpolate() functions.
+ *
+ * Tutorial programs step-6, step-7 and step-8 show examples of how to use
+ * function arguments of this type in situations where we actually have an empty
+ * map (i.e., we want to describe that <i>no</i> part of the boundary is a
+ * Neumann boundary). step-16 actually uses it in a case where one of the
+ * parts of the boundary uses a boundary indicator for which we want to use
+ * a function object.
  *
  * It seems odd at first to declare this typedef inside a class, rather than
  * declaring a typedef at global scope. The reason is that C++ does not allow
  * to define templated typedefs, where here in fact we want a typdef that
- * depends on the space dimension.
+ * depends on the space dimension. (Defining templated typedefs is something that
+ * is possible starting with the C++11 standard, but that wasn't possible within
+ * the C++98 standard in place when this programming pattern was conceived.)
  *
  * @ingroup functions
  * @author Wolfgang Bangerth, Ralf Hartmann, 2001

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.