]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Update docs.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 9 Apr 2010 04:50:20 +0000 (04:50 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 9 Apr 2010 04:50:20 +0000 (04:50 +0000)
git-svn-id: https://svn.dealii.org/trunk@20965 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_system.h

index ad1848fc0f7d7e4763e1c182106e7ce20a27390c..3cc1ffa1e75d73ba19eb0e54c42485620ba700cd 100644 (file)
@@ -28,64 +28,92 @@ DEAL_II_NAMESPACE_OPEN
  * one. To the outside world, the resulting object looks just like a usual
  * finite element object, which is composed of several other finite elements
  * that are possibly of different type. The result is then a vector-valued
- * finite element. Vector valued elements are discussed in a number of
- * tutorial programs, for example step-8, @ref step_20
- * "step-20", step-21, and in particular in the @ref vector_valued
- * module.
+ * finite element. %Vector valued elements are discussed in a number of
+ * tutorial programs, for example step-8, step-20, step-21, and in particular
+ * in the @ref vector_valued module.
  *
  * <h3>FESystem, components and blocks</h3>
  *
- * An FESystem, except in the most trivial case, produces a
- * vector-valued finite element with several components. The number of
- * components corresponds to the dimension of the function in the PDE
- * system.
+ * An FESystem, except in the most trivial case, produces a vector-valued
+ * finite element with several components. The number of components
+ * corresponds to the dimension of the solution function in the PDE system,
+ * and correspondingly also to the number of equations your PDE system
+ * has. For example, the mixed Laplace system covered in step-20 has $d+1$
+ * components in $d$ space dimensions: the scalar pressure and the $d$
+ * components of the velocity vector. Similarly, the elasticity equation
+ * covered in step-8 has $d$ components in $d$ space dimensions. In general,
+ * the number of components of a FESystem element is the
+ * accumulated number of components of all base elements times their
+ * multiplicities. A bit more on
+ * components is also given in the
+ * @ref GlossComponent "glossary entry on components".
+ *
+ * While the concept of components is important from the viewpoint of a
+ * partial differential equation, the finite element side looks a bit
+ * different Since not only FESystem, but also vector-valued elements like
+ * FE_RaviartThomas, have several components. The concept needed here is a
+ * @ref GlossBlock "block". Each block encompasses the set of degrees of
+ * freedom associated with a single base element of an FESystem, where base
+ * elements with multiplicities count multiple times. These blocks are usually
+ * addressed using the information in DoFHandler::block_info(). The number of
+ * blocks of of a FESystem object is simply the sum of all multiplicities of
+ * base elements.
  *
- * Since not only FESystem, but also vector-valued elements like
- * FE_RaviartThomas, have several components, the notion of a
- * component is of less importance in solving the finite element
- * problem. The concept needed here is a @ref GlossBlock
- * "block". A block refers to the degrees of freedom generated by a
- * single base element of an FESystem, where base elements with
- * multiplicities count multiple times. These blocks are usually
- * addressed using the information in DoFHandler::block_info(). Here,
- * we have two examples for FESystem for the Taylor-Hood element for
- * the three-dimensional Stokes problem:
+ * For example, the FESystem for the Taylor-Hood element for the
+ * three-dimensional Stokes problem can be built using the code
  *
  * @code
  * FE_Q<3> u(2);
  * FE_Q<3> p(1);
- * FESystem<3> sys1(u, 3, p, 1);
+ * FESystem<3> sys1(u,3, p,1);
  * @endcode
  *
- * This example creates an FESystem @p sys1 with four components,
- * three for the velocity components and one for the pressure, and
- * also four blocks with the degrees of freedom of each of the
- * velocity components and the pressure in a separate block each.
+ * This example creates an FESystem @p sys1 with four components, three for
+ * the velocity components and one for the pressure, and also four blocks with
+ * the degrees of freedom of each of the velocity components and the pressure
+ * in a separate block each. The number of blocks is four since the first base
+ * element is repeated three times.
+ *
+ * On the other hand, a Taylor-Hood element can also be constructed using
  *
  * @code
  * FESystem<3> U(u,3);
- * FESystem<3> sys2(U,1,p,1);
+ * FESystem<3> sys2(U,1, p,1);
  * @endcode
  *
- * The FESystem @p sys2 created here has the same four components, but
- * the degrees of freedom are distributed into only two blocks. The
- * first block has all velocity degrees of freedom from @p U, while
- * the second block contains the pressure degrees of freedom. The
- * FESystem @p U is not split accroding to its base elements. Note
- * that by blocking all velocities into one system first, we mimic a
- * block structure that would be generated by using vector-valued base
- * elements, for instance like using a mixed discretization of Darcy's
- * law using
+ * The FESystem @p sys2 created here has the same four components, but the
+ * degrees of freedom are distributed into only two blocks. The first block
+ * has all velocity degrees of freedom from @p U, while the second block
+ * contains the pressure degrees of freedom. Note that while @p U itself has 3
+ * blocks, the FESystem @p sys2 does not attempt to split @p U into its base
+ * elements but considers it a block of its own. By blocking all velocities
+ * into one system first as in @p sys2, we achieve the sam block structure
+ * that would be generated if instead of using a $Q_2^3$ element for the
+ * velocities we had used vector-valued base elements, for instance like using
+ * a mixed discretization of Darcy's law using
  *
  * @code
  * FE_RaviartThomas<3> u(1);
  * FE_DGQ<3> p(1);
- * FESystem<3> sys3(u,1,p,1);
+ * FESystem<3> sys3(u,1, p,1);
  * @endcode
  *
  * This example also produces a system with four components, but only
  * two blocks.
  *
+ * In most cases, the composed element behaves as if it were a usual
+ * element. It just has more degrees of freedom than most of the "common"
+ * elements. However the underlying structure is visible in the restriction,
+ * prolongation and interface constraint matrices, which do not couple the
+ * degrees of freedom of the base elements. E.g. the continuity requirement is
+ * imposed for the shape functions of the subobjects separately; no
+ * requirement exist between shape functions of different subobjects, i.e. in
+ * the above example: on a hanging node, the respective value of the @p u
+ * velocity is only coupled to @p u at the vertices and the line on the larger
+ * cell next to this vertex, but there is no interaction with @p v and @p w of
+ * this or the other cell.
+ *
+ *
  * <h3>Internal information on numbering of degrees of freedom</h3>
  *
  * The overall numbering of degrees of freedom is as follows: for each
@@ -107,25 +135,9 @@ DEAL_II_NAMESPACE_OPEN
  * <li> Third component on the line:
  *   <tt>p2 = s8</tt>.
  * </ul>
- * Do not rely on this numbering in your application as these
- * internals might change in future. Rather use the functions
- * @p system_to_component_index and @p component_to_system_index,
- * instead.
- *
- * In the most cases, the composed element behaves as if it were a usual element
- * with more degrees of freedom. However the underlying structure is visible in
- * the restriction, prolongation and interface constraint matrices, which do not
- * couple the degrees of freedom of the subobject. E.g. the continuity requirement
- * is imposed for the shape functions of the subobjects separately; no requirement
- * exist between shape functions of different subobjects, i.e. in the above
- * example: on a hanging node, the respective value of the @p u velocity is only
- * coupled to @p u at the vertices and the line on the larger cell next to this
- * vertex, there is no interaction with @p v and @p w of this or the other cell.
- *
- * The number of components of such a system element is the
- * accumulated number of components of all base elements times their
- * multiplicity. The number of blocks of the system is simply the sum
- * of all multiplicities.
+ * That said, you should not rely on this numbering in your application as
+ * these %internals might change in future. Rather use the functions
+ * system_to_component_index() and component_to_system_index().
  *
  * For more information on the template parameter <tt>spacedim</tt>
  * see the documentation of Triangulation.

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.