]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Better document ReferenceCell and the reference cell glossary entry. 13830/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 26 May 2022 03:22:59 +0000 (21:22 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 26 May 2022 14:04:32 +0000 (08:04 -0600)
doc/doxygen/headers/glossary.h
include/deal.II/base/geometry_info.h
include/deal.II/grid/reference_cell.h

index 2df60159e8e609ee157d8d4351bdfd775f2d8fcb..7c849fa6eb1413b75b6fbfbf4b43dfed68078924 100644 (file)
  *
  *
  * <dt class="glossary">@anchor GlossReferenceCell <b>Reference cell</b></dt>
- * <dd>The hypercube [0,1]<sup>dim</sup>, on which all parametric finite
- * element shape functions are defined. Many properties of the reference
- * cell are described by the GeometryInfo class.
+ * <dd>The finite element method is typically described by providing shape
+ * functions on a "reference cell" whose shape, along with the shape functions,
+ * is then mapped to the actual cells of the mesh. These reference cells are
+ * typically triangles or quadrilaterals for two-dimensional meshes;
+ * tetrahedra, hexahedra, wedges, or pyramids for three-dimensional meshes;
+ * and simple line segments in the one-dimensional case.
+ *
+ * Rather than hard-coding properties of the reference cell in all places
+ * where one wants to know about, say, the number of vertices of a cell,
+ * deal.II uses a single, central place to describe the properties of
+ * reference cells: The ReferenceCell class. In loops over the cells of a
+ * mesh, one typically asks for properties of these cells using the
+ * call `cell->reference_cell()`.
  * </dd>
  *
  *
index 674e501d57c4dd9a54ceada89722ab89bfdd0b98..4c778352d0b008bc70ac9ae16d78afa43dba01c6 100644 (file)
@@ -1438,8 +1438,12 @@ struct GeometryInfo<0>
  * This class provides dimension independent information to all topological
  * structures that make up the unit, or
  * @ref GlossReferenceCell "reference cell".
- * This class has been
+ * That said, this class only describes information about hypercube reference
+ * cells (i.e., lines, quadrilaterals, or hexahedra), which historically
+ * were the only kinds of cells supported by deal.II. This is no longer the
+ * case today, and consequently this class has been
  * superseded by the ReferenceCell class -- see there for more information.
+ * The rest of this class's documentation is therefore partly historical.
  *
  *
  * It is the one central point in the library where information about the
index 9d21ca8aa74ab2c88f8d1b60bd3ca0a8dfe8a544..76fdec26e449cdd0a484421e0d9c99b52ae4aa98 100644 (file)
@@ -68,10 +68,19 @@ namespace internal
 
 
 /**
- * A type that describes the kinds of reference cells that can be used.
- * This includes quadrilaterals and hexahedra (i.e., "hypercubes"),
- * triangles and tetrahedra (simplices), and the pyramids and wedges
- * necessary when using mixed 3d meshes.
+ * A type that describes the kinds of reference cells that can be
+ * used.  This includes quadrilaterals and hexahedra (i.e.,
+ * "hypercubes"), triangles and tetrahedra (simplices), and the
+ * pyramids and wedges necessary when using mixed 3d meshes. This
+ * class then describes geometric, topological, and other kinds of
+ * information about these kinds of reference cells. This includes how
+ * many vertices or faces a certain kind of reference cell has
+ * (topological information), where these vertices lie, what the
+ * cell's volume or center of mass is (geometric information), and how
+ * to output these cells in various output formats or what appropriate
+ * quadrature rules are. The documentation of this class is separated
+ * into a number of sections to group the many member functions into
+ * different categories such as those mentioned above.
  *
  * Objects of this type should not be created in user code, and as a
  * consequence the class does not have a user-accessible constructor
@@ -487,10 +496,11 @@ public:
   volume() const;
 
   /**
-   * Return the barycenter of the reference cell that corresponds
-   * to the current object. The function is not called `center()` because
-   * one can define the center of an object in a number of different ways
-   * whereas the barycenter of a reference cell $K$ is unambiguously defined as
+   * Return the barycenter (i.e., the center of mass) of the reference
+   * cell that corresponds to the current object. The function is not
+   * called `center()` because one can define the center of an object
+   * in a number of different ways whereas the barycenter of a
+   * reference cell $K$ is unambiguously defined as
    * @f[
    *   \mathbf x_K = \frac{1}{V} \int_K \mathbf x \; dx
    * @f]

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.