]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Update the documentation of step-7. 580/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 20 Feb 2015 16:21:19 +0000 (10:21 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 20 Feb 2015 16:21:19 +0000 (10:21 -0600)
In particular, name the Method of Manufactured Solutions in the documentation.

examples/step-7/doc/intro.dox

index 26dcd7634def57eed773f4b03ab77f8afe873667..027df31b7a20898e2a2285bd8155eadec2e301c8 100644 (file)
@@ -16,20 +16,22 @@ There has probably never been a
 non-trivial finite element program that worked right from the start. It is
 therefore necessary to find ways to verify whether a computed solution is
 correct or not. Usually, this is done by choosing the set-up of a simulation
-such that we know the exact continuous solution and evaluate the difference
+in such a way that we know the exact continuous solution and evaluate the difference
 between continuous and computed discrete solution. If this difference
 converges to zero with the right order of convergence, this is already a good
 indication of correctness, although there may be other sources of error
 persisting which have only a small contribution to the total error or are of
-higher order.
+higher order. In the context of finite element simulations, this technique
+is often called the <i>Method of Manufactured Solution</i>.
 
 In this example, we will not go into the theories of systematic software
 verification which is a very complicated problem. Rather we will demonstrate
 the tools which deal.II can offer in this respect. This is basically centered
-around the functionality of a single function, <code>integrate_difference</code>.
+around the functionality of a single function, VectorTools::integrate_difference().
 This function computes the difference between a given continuous function and
-a finite element field in various norms on each cell. At present, the
-supported norms are the following, where $u$ denotes the continuous function
+a finite element field in various norms on each cell. At the time of writing
+this tutorial program, the norms this function can compute are the following,
+where $u$ denotes the continuous function
 and $u_h$ the finite element field, and $K$ is an element of the
 triangulation:
 @f{eqnarray*}
@@ -44,17 +46,26 @@ triangulation:
   {\| u-u_h \|}_{H^1(K)} &=& \left( {\| u-u_h \|}^2_{L_2(K)} 
                                    +{| u-u_h |}^2_{H^1(K)}    \right)^{1/2}.
 @f}
-All these norms and semi-norms can also be evaluated with weighting functions,
+(All these norms and semi-norms can also be evaluated with weighting functions,
 for example in order to exclude singularities from the determination of the
-global error. The function also works for vector-valued functions.  It should
-be noted that all these quantities are evaluated using quadrature formulas;
+global error, and the function also works for vector-valued functions.) Of
+course, like with any other integral, we can only evaluate these norms using quadrature formulas;
 the choice of the right quadrature formula is therefore crucial to the
 accurate evaluation of the error. This holds in particular for the $L_\infty$
 norm, where we evaluate the maximal deviation of numerical and exact solution
 only at the quadrature points; one should then not try to use a quadrature
-rule with points only at points where super-convergence might occur.
-
-The function <code>integrate_difference</code> evaluates the desired norm on each
+rule with points only at points where super-convergence might occur, such as
+the Gauss points of the lowest-order Gauss quadrature formula for which the
+integrals in the assembly of the matrix is correct (e.g., for linear elements,
+do not use the QGauss(2) quadrature formula). In fact, this is generally good 
+advice also for the other norms: if your quadrature points are fortuitously
+chosen at locations where the error happens to be particularly small due to
+superconvergence, the computed error will look like it is much smaller than
+it really is and may even suggest a higher convergence order. Consequently,
+we will choose a different quadrature formula for the integration of these
+error norms than for the assembly of the linear system.
+
+The function VectorTools::integrate_difference() evaluates the desired norm on each
 cell $K$ of the triangulation and returns a vector which holds these
 values for each cell. From the local values, we can then obtain the global error. For
 example, if the vector $(e_i)$ contains the local $L_2$ norms, then
@@ -87,8 +98,8 @@ boundary integrals which we have to evaluate numerically when assembling the
 right hand side vector.
 
 Before we go into programming, let's have a brief look at the mathematical
-formulation. The equation which we want to solve is Helmholtz's equation
-``with the nice sign'':
+formulation. The equation that we want to solve here is the Helmholtz equation
+"with the nice sign":
 @f[
   -\Delta u + u = f,
 @f]
@@ -101,21 +112,34 @@ on some part $\Gamma_1$ of the boundary $\Gamma$, and
   {\mathbf n}\cdot \nabla u = g_2
 @f]
 on the rest $\Gamma_2 = \Gamma \backslash \Gamma_1$.
+In our particular testcase, we will use $\Gamma_1=\Gamma \cap\{\{x=1\} \cup \{y=1\}\}$.
 
-We choose the right hand side function $f$ such that the exact solution is
+Because we want to verify the convergence of our numerical solution $u_h$,
+we want a setup so that we know the exact solution $u$. This is where
+the Method of Manufactured Solutions comes in. To this end, let us 
+choose a function
 @f[
-  u(x) = \sum_{i=1}^3 \exp\left(-\frac{|x-x_i|^2}{\sigma^2}\right)
+  \bar u(x) = \sum_{i=1}^3 \exp\left(-\frac{|x-x_i|^2}{\sigma^2}\right)
 @f]
 where the centers $x_i$ of the exponentials are 
   $x_1=(-\frac 12,\frac 12)$,
   $x_2=(-\frac 12,-\frac 12)$, and
-  $x_3=(\frac 12,-\frac 12)$.
-The half width is set to $\sigma=\frac {1}{8}$.
-
-We further choose $\Gamma_1=\Gamma \cap\{\{x=1\} \cup \{y=1\}\}$, and there
-set $g_1$ such that it resembles the exact values of $u$. Likewise, we choose
-$g_2$ on the remaining portion of the boundary to be the exact normal
-derivatives of the continuous solution.
+  $x_3=(\frac 12,-\frac 12)$,
+and the half width is set to $\sigma=\frac {1}{8}$. The method of manufactured
+solution then says: choose
+@f{align*}
+  f &= -\Delta \bar u + \bar u, \\
+  g_1 &= \bar u|_{\Gamma_1}, \\
+  g_2 &= {\mathbf n}\cdot \nabla\bar u|_{\Gamma_2}.
+@f}
+With this particular choice, we infer that of course the solution of the
+original problem happens to be $u=\bar u$. In other words, by choosing
+the right hand sides of the equation and the boundary conditions in a
+particular way, we have manufactured ourselves a problem to which we
+know the solution. This allows us then to compute the error of our
+numerical solution. In the code below, we represent $\bar u$ by the
+<code>Solution</code> class, and other classes will be used to
+denote $\bar u|_{\Gamma_1}$ and ${\mathbf n}\cdot \nabla\bar u|_{\Gamma_2}$.
 
 Using the above definitions, we can state the weak formulation of the
 equation, which reads: find $u\in H^1_g=\{v\in H^1: v|_{\Gamma_1}=g_1\}$ such

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.