//---------------------------------------------------------------------------
// $Id$
-// Version: $Name$
//
-// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
+// Copyright (C) 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
* "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.
+ *
+ * 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:
+ *
+ * @code
+ * FE_Q<3> u(2);
+ * FE_Q<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.
+ *
+ * @code
+ * FESystem<3> U(u,3);
+ * 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
+ *
+ * @code
+ * FE_RaviartThomas<3> u(1);
+ * FE_DGQ<3> p(1);
+ * FESystem<3> sys3(u,1,p,1);
+ * @endcode
+ *
+ * This example also produces a system with four components, but only
+ * two blocks.
+ *
+ * <h3>Internal information on numbering of degrees of freedom</h3>
+ *
* The overall numbering of degrees of freedom is as follows: for each
* subobject (vertex, line, quad, or hex), the degrees of freedom are
* numbered such that we run over all subelements first, before
* of the sub-element @p fe.
*/
static std::vector<bool>
- compute_restriction_is_additive_flags (const FiniteElement<dim,spacedim> &fe,
- const unsigned int N);
+ compute_restriction_is_additive_flags (
+ const FiniteElement<dim,spacedim> &fe,
+ const unsigned int N);
/**
* Same as above for mixed elements
* with two different sub-elements.
*/
static std::vector<bool>
- compute_restriction_is_additive_flags (const FiniteElement<dim,spacedim> &fe1,
- const unsigned int N1,
- const FiniteElement<dim,spacedim> &fe2,
- const unsigned int N2);
+ compute_restriction_is_additive_flags (
+ const FiniteElement<dim,spacedim> &fe1,
+ const unsigned int N1,
+ const FiniteElement<dim,spacedim> &fe2,
+ const unsigned int N2);
/**
* Same as above for mixed elements
* with three different sub-elements.
*/
static std::vector<bool>
- compute_restriction_is_additive_flags (const FiniteElement<dim,spacedim> &fe1,
- const unsigned int N1,
- const FiniteElement<dim,spacedim> &fe2,
- const unsigned int N2,
- const FiniteElement<dim,spacedim> &fe3,
- const unsigned int N3);
+ compute_restriction_is_additive_flags (
+ const FiniteElement<dim,spacedim> &fe1,
+ const unsigned int N1,
+ const FiniteElement<dim,spacedim> &fe2,
+ const unsigned int N2,
+ const FiniteElement<dim,spacedim> &fe3,
+ const unsigned int N3);
/**
* with four different sub-elements
*/
static std::vector<bool>
- compute_restriction_is_additive_flags (const FiniteElement<dim,spacedim> &fe1,
- const unsigned int N1,
- const FiniteElement<dim,spacedim> &fe2,
- const unsigned int N2,
- const FiniteElement<dim,spacedim> &fe3,
- const unsigned int N3,
- const FiniteElement<dim,spacedim> &fe4,
- const unsigned int N4);
+ compute_restriction_is_additive_flags (
+ const FiniteElement<dim,spacedim> &fe1,
+ const unsigned int N1,
+ const FiniteElement<dim,spacedim> &fe2,
+ const unsigned int N2,
+ const FiniteElement<dim,spacedim> &fe3,
+ const unsigned int N3,
+ const FiniteElement<dim,spacedim> &fe4,
+ const unsigned int N4);
/**
* Compute the named flags for a
* functions.
*/
static std::vector<bool>
- compute_restriction_is_additive_flags (const std::vector<const FiniteElement<dim,spacedim>*> &fes,
- const std::vector<unsigned int> &multiplicities);
+ compute_restriction_is_additive_flags (
+ const std::vector<const FiniteElement<dim,spacedim>*> &fes,
+ const std::vector<unsigned int> &multiplicities);
/**
/*---------------------------- fe_system.h ---------------------------*/
#endif
/*---------------------------- fe_system.h ---------------------------*/
+