]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Discuss the concept of superconvergence better. 10703/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Mon, 13 Jul 2020 19:03:01 +0000 (13:03 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Tue, 15 Sep 2020 23:10:03 +0000 (17:10 -0600)
doc/doxygen/references.bib
examples/step-7/step-7.cc

index 1955584265fd6808b237130cca587ede68d8c31d..d98adc66aecebb3971435fb8bb37c5c9f213245f 100644 (file)
   pages =   {101--129}
 }
 
+%-------------------------------------------------------------------------------
+% Step 7
+%-------------------------------------------------------------------------------
+
+@article{Li2019,
+  doi = {10.1007/s10915-019-01102-1},
+  url = {https://doi.org/10.1007/s10915-019-01102-1},
+  year = {2019},
+  month = dec,
+  publisher = {Springer Science and Business Media {LLC}},
+  volume = {82},
+  number = {1},
+  author = {Hao Li and Xiangxiong Zhang},
+  title = {Superconvergence of $C^0-Q^k$ Finite Element Method for Elliptic Equations with Approximated Coefficients},
+  journal = {Journal of Scientific Computing}
+}
+
 %-------------------------------------------------------------------------------
 % Step 14
 %-------------------------------------------------------------------------------
index 495d3dcf6aad16a3d5fe168cb2f84ddee61c3118..44245daadc70fb37f54b74d7aa1ce4c345dc3e42 100644 (file)
@@ -868,15 +868,30 @@ namespace Step7
                                         VectorTools::H1_seminorm);
 
     // Finally, we compute the maximum norm. Of course, we can't actually
-    // compute the true maximum, but only the maximum at the quadrature
-    // points. Since this depends quite sensitively on the quadrature rule
-    // being used, and since we would like to avoid false results due to
-    // super-convergence effects at some points, we use a special quadrature
-    // rule that is obtained by iterating the trapezoidal rule by the degree of
-    // the finite element times two plus one in each space direction.
-    // Note that the constructor of the QIterated class
-    // takes a one-dimensional quadrature rule and a number that tells it how
-    // often it shall use this rule in each space direction.
+    // compute the true maximum of the error over *all* points in the domain,
+    // but only the maximum over a finite set of evaluation points that, for
+    // convenience, we will still call "quadrature points" and represent by
+    // an object of type Quadrature even though we do not actually perform any
+    // integration.
+    //
+    // There is then the question of what points precisely we want to evaluate
+    // at. It turns out that the result we get depends quite sensitively on the
+    // "quadrature" points being used. There is also the issue of
+    // superconvergence: Finite element solutions are, on some meshes and for
+    // polynomial degrees $k\ge 2$, particularly accurate at the node points as
+    // well as at Gauss-Lobatto points, much more accurate than at randomly
+    // chosen points. (See
+    // @cite Li2019 and the discussion and references in Section 1.2 for more
+    // information on this.) In other words, if we are interested in finding
+    // the largest difference $u(\mathbf x)-u_h(\mathbf x)$, then we ought to
+    // look at points $\mathbf x$ that are specifically not of this "special"
+    // kind of points and we should specifically not use
+    // `QGauss(fe->degree+1)` to define where we evaluate. Rather, we use a
+    // special quadrature rule that is obtained by iterating the trapezoidal
+    // rule by the degree of the finite element times two plus one in each space
+    // direction. Note that the constructor of the QIterated class takes a
+    // one-dimensional quadrature rule and a number that tells it how often it
+    // shall repeat this rule in each space direction.
     //
     // Using this special quadrature rule, we can then try to find the maximal
     // error on each cell. Finally, we compute the global L infinity error

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.