]> https://gitweb.dealii.org/ - dealii.git/commitdiff
doxygen: simplex relates 14794/head
authorTimo Heister <timo.heister@gmail.com>
Mon, 13 Feb 2023 17:43:55 +0000 (12:43 -0500)
committerTimo Heister <timo.heister@gmail.com>
Mon, 13 Feb 2023 17:43:55 +0000 (12:43 -0500)
It turns out "@relates simplex" does not add the functionality (or link
to) the simplex group. Instead, it removes some functions completely
(GridGenerator functions are currently nowhere to be found).

Instead, I created a list inside simplex manually and just reference the
simplex group.

doc/doxygen/headers/simplex.h
include/deal.II/base/quadrature_lib.h
include/deal.II/fe/fe_pyramid_p.h
include/deal.II/fe/fe_simplex_p.h
include/deal.II/fe/fe_wedge_p.h
include/deal.II/fe/mapping_fe.h
include/deal.II/grid/grid_generator.h
include/deal.II/grid/grid_in.h

index 85a7279e7088b56c0d9613bc8077f1afd577c912..5beaf1d9b57078e8bd650515dfa9c0f040305eed 100644 (file)
  * Simplex and mixed meshes in deal.II are still experimental, i.e., work
  * in progress. Large parts of the library have been ported to be able to
  * operate on such kind of meshes. However, there are still many functions
- * that need to be generalized. You can get a good overview of the ported
+ * that need to be generalized.
+ *
+ * @section Important Simplex Functionality
+ *
+ * Here is an incomplete list of functionality related to simplex
+ * computations:
+ * - Mesh generation:
+ *   GridGenerator::implicit_function(),
+ *   GridGenerator::convert_hypercube_to_simplex_mesh(),
+ *   GridGenerator::subdivided_hyper_rectangle_with_simplices(),
+ *   GridGenerator::subdivided_hyper_cube_with_simplices()
+ * - Quadratures:
+ *   QGaussWedge, QGaussSimplex, QWitherdenVincentSimplex
+ * - FiniteElements:
+ *   FE_SimplexP, FE_SimplexDGP, FE_SimplexP_Bubbles
+ *   FE_PyramidP, FE_PyramidDGP, FE_WedgeP, FE_WedgeDGP
+ * - Mapping:
+ *   MappingFE
+ * - Other:
+ *   GridIn::read_vtk(), GridIn::read_msh(), GridIn::read_comsol_mphtxt()
+ *
+ *
+ *
+ * @section Examples
+ *
+ * You can get a good overview of the ported
  * functionalities by taking a look at the tests in the folder
  * "tests/simplex". In the following, we provide two very basic examples
- * to get started and provide some implementation details.
+ * to get you started and to provide some implementation details.
  *
- * @section simplex_reference_example_simplex Example: simplex mesh
+ * @subsection simplex_reference_example_simplex Example: simplex mesh
  *
  * The following code shows how to work with simplex meshes:
  *
  * @include step_3_simplex.cc
  *
- * @section simplex_reference_example_mixed Example: mixed mesh
+ * @subsection simplex_reference_example_mixed Example: mixed mesh
  *
  * The following code shows how to work with mixed meshes:
  *
index ea4d6bffeafd9dce7af6474e231a0f149ae6bb0c..3ab46408530cd10585e09530168e958440cb952d 100644 (file)
@@ -831,7 +831,8 @@ public:
  * depends on the numbering of the cell vertices. If you need rules that are
  * independent of the vertex numbering then use QWitherdenVincentSimplex.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim>
 class QGaussSimplex : public QSimplex<dim>
@@ -872,7 +873,8 @@ public:
  * @note Some rules (2d 2 odd and 3d 2 even) do not yet exist and instead a
  * higher-order rule is used in their place.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim>
 class QWitherdenVincentSimplex : public QSimplex<dim>
index 28981855f79421b2a1e49b249892083cd71e85af..9715745b3c8ba7f80fe08ac219288143bea61f6b 100644 (file)
@@ -29,7 +29,8 @@ DEAL_II_NAMESPACE_OPEN
  *
  * @note Only implemented for 3d.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_PyramidPoly : public dealii::FE_Poly<dim, spacedim>
@@ -59,7 +60,8 @@ public:
  * @note Currently, only linear polynomials (degree=1) are implemented. See
  * also the documentation of ScalarLagrangePolynomialPyramid.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_PyramidP : public FE_PyramidPoly<dim, spacedim>
@@ -121,7 +123,8 @@ public:
  * @note Currently, only linear polynomials (degree=1) are implemented. See
  * also the documentation of ScalarLagrangePolynomialPyramid.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_PyramidDGP : public FE_PyramidPoly<dim, spacedim>
index 38ac36cca240bddd4c08725c9fac6fabb10b6006..e3040534c3b55fbca0f35cd391014c0c175853ac 100644 (file)
@@ -30,7 +30,8 @@ DEAL_II_NAMESPACE_OPEN
  *
  * @note Only implemented for 2d and 3d.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_SimplexPoly : public dealii::FE_Poly<dim, spacedim>
@@ -121,7 +122,8 @@ protected:
  * the finite element space of continuous, piecewise polynomials of
  * degree $k$.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_SimplexP : public FE_SimplexPoly<dim, spacedim>
@@ -176,7 +178,8 @@ public:
  * element space of discontinuous, piecewise polynomials of degree
  * $k$.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_SimplexDGP : public FE_SimplexPoly<dim, spacedim>
index 772c713dbdf505c59ff819dc2c3d075a2e95fea9..3212d08485c612bc82bd6e7633f7e32a4a740aa9 100644 (file)
@@ -29,7 +29,8 @@ DEAL_II_NAMESPACE_OPEN
  *
  * @note Only implemented for 3d.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_WedgePoly : public dealii::FE_Poly<dim, spacedim>
@@ -60,7 +61,8 @@ public:
  *   (degree=2) are implemented. See also the documentation of
  *   ScalarLagrangePolynomialWedge.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_WedgeP : public FE_WedgePoly<dim, spacedim>
@@ -123,7 +125,8 @@ public:
  *   (degree=2) are implemented. See also the documentation of
  *   ScalarLagrangePolynomialWedge.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class FE_WedgeDGP : public FE_WedgePoly<dim, spacedim>
index 835900593fcbf2fb6192970436e5902892f6332b..fdb0ac6062ebb4ed1b1568f5a4fda16e561918c7 100644 (file)
@@ -53,7 +53,8 @@ DEAL_II_NAMESPACE_OPEN
  * @note Currently, only implemented for elements with tensor_degree==1 and
  *   n_components==1.
  *
- * @relates simplex
+ * Also see
+ * @ref simplex "Simplex support".
  */
 template <int dim, int spacedim = dim>
 class MappingFE : public Mapping<dim, spacedim>
index 9747f51a85024bd2c8193b1d98ad784218fa9bad..c557898e5e1f79ca8d89545cc40dd83d1678815d 100644 (file)
@@ -2268,7 +2268,8 @@ namespace GridGenerator
    *     out_tria.set_manifold(i, in_tria.get_manifold(i));
    * @endcode
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   template <int dim, int spacedim>
   void
@@ -2484,7 +2485,8 @@ namespace GridGenerator
    *
    * @note Currently, this function only works for `dim==spacedim`.
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   template <int dim, int spacedim>
   void
@@ -2507,7 +2509,8 @@ namespace GridGenerator
    * quadrilateral/hexahedral cells and subdivides these into 2/5
    * triangular/tetrahedral cells.
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   template <int dim, int spacedim>
   void
index 03254f44bd2fe163de96aeac43a0285f3c843899..b0b8fca34531a04219f835318cd0c114f0292498 100644 (file)
@@ -411,7 +411,8 @@ public:
    * The companion GridOut::write_vtk function can be used to write VTK files
    * compatible with this method.
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   void
   read_vtk(std::istream &in);
@@ -516,7 +517,8 @@ public:
    * Read grid data from an msh file. The %Gmsh formats are documented at
    * http://www.gmsh.info/.
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   void
   read_msh(std::istream &in);
@@ -575,7 +577,9 @@ public:
    * as a boundary or material id.  Physical surface numbers created in Gmsh,
    * which can be seen in the .geo file, become material IDs.
    *
-   * @relates simplex
+   *
+   * Also see
+   * @ref simplex "Simplex support".
    */
   void
   read_msh(const std::string &filename);
@@ -614,7 +618,8 @@ public:
    * @image html "comsol-mesh-boundary-lines.png"
    * @image html "comsol-mesh-boundary-triangles.png"
    *
-   * @relates simplex
+   * Also see
+   * @ref simplex "Simplex support".
    */
   void
   read_comsol_mphtxt(std::istream &in);

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.