virtual ~FiniteElement ();
/**
- * Number of base elements in a mixed
- * discretization. This function returns
- * 1 for primitive elements.
+ * Number of base elements in a
+ * mixed discretization.
+ *
+ * Note that even for vector
+ * valued finite elements, the
+ * number of components needs not
+ * coincide with the number of
+ * base elements, since they may
+ * be reused. For example, if you
+ * create a @ref{FESystem} with
+ * three identical finite element
+ * classes by using the
+ * constructor that takes one
+ * finite element and a
+ * multiplicity, then the number
+ * of base elements is still one,
+ * although the number of
+ * components of the finite
+ * element is equal to the
+ * multiplicity.
*/
- virtual unsigned int n_base_elements () const;
+ virtual unsigned int n_base_elements () const = 0;
/**
* Access to base element
- * objects. By default,
- * #base_element(0)# is #this#.
- * This function is overloaded by
- * system elements to allow
- * access to the different
- * components of mixed
- * discretizations.
+ * objects. If the element is
+ * scalar, then
+ * @p{base_element(0)} is
+ * @p{this}.
*/
- virtual const FiniteElement<dim> & base_element (const unsigned int index) const;
-
+ virtual const FiniteElement<dim> & base_element (const unsigned int index) const = 0;
+
/**
* Determine an estimate for the
* memory consumption (in bytes)
* shape functions of the base
* element.
*/
- unsigned int component_to_base(unsigned int index) const;
+ unsigned int component_to_base (unsigned int index) const;
/**
* Access the @p{restriction_is_additive_flag}
* the result allows access to
* shape functions of the base
* element.
+ *
+ * This variable is set to the
+ * correct size by the
+ * constructor of this class, but
+ * needs to be initialized by
+ * derived classes, unless its
+ * size is one and the only entry
+ * is a zero, which is the case
+ * for scalar elements. In that
+ * case, the initialization by
+ * the base class is sufficient.
*/
std::vector<unsigned int> component_to_base_table;
*/
unsigned int get_degree () const;
+ /**
+ * Number of base elements in a
+ * mixed discretization. Since
+ * this is a scalar element,
+ * return one.
+ */
+ virtual unsigned int n_base_elements () const;
+
+ /**
+ * Access to base element
+ * objects. Since this element is
+ * scalar, @p{base_element(0)} is
+ * @p{this}, and all other
+ * indices throw an error.
+ */
+ virtual const FiniteElement<dim> & base_element (const unsigned int index) const;
+
/**
* Determine an estimate for the
* memory consumption (in bytes)
*/
unsigned int get_degree () const;
+ /**
+ * Number of base elements in a
+ * mixed discretization. Since
+ * this is a scalar element,
+ * return one.
+ */
+ virtual unsigned int n_base_elements () const;
+
+ /**
+ * Access to base element
+ * objects. Since this element is
+ * scalar, @p{base_element(0)} is
+ * @p{this}, and all other
+ * indices throw an error.
+ */
+ virtual const FiniteElement<dim> & base_element (const unsigned int index) const;
+
/**
* Determine an estimate for the
* memory consumption (in bytes)
virtual unsigned int n_base_elements () const;
/**
- * How often is a composing element used.
- *
+ * How often is a composing
+ * element used.
*/
unsigned int element_multiplicity (const unsigned int index) const;
/**
- * Access to a composing element.
- *
- * If you assemble your system
- * matrix, you usually will not
- * want to have an FEValues object
- * with a lot of equal entries. Ok,
- * so initialize your FEValues with
- * the @p{base_element} you get by
- * this function. In a mixed
- * discretization, you can choose
- * the different base element types
- * by index.
- *
+ * Access to a composing
+ * element. The index needs to be
+ * smaller than the number of
+ * base elements. Note that the
+ * number of base elements may in
+ * turn be smaller than the
+ * number of components of the
+ * system element, if the
+ * multiplicities are greater
+ * than one.
*/
- virtual const FiniteElement<dim> & base_element(const unsigned int index) const;
-
+ virtual const FiniteElement<dim> & base_element (const unsigned int index) const;
+
/**
* Determine an estimate for the
* memory consumption (in bytes)
/* ------------------------- inline functions ------------------------- */
-template<int dim>
-inline unsigned int
-FESystem<dim>::n_base_elements() const
-{
- return base_elements.size();
-};
-
-
template<int dim>
inline unsigned int
-template <int dim>
-inline const FiniteElement<dim> &
-FESystem<dim>::base_element (const unsigned int index) const
-{
- Assert (index < base_elements.size(),
- ExcIndexRange(index, 0, base_elements.size()));
- return *base_elements[index].first;
-};
-
-
-
-
-
/*---------------------------- fe_system.h ---------------------------*/
#endif
/*---------------------------- fe_system.h ---------------------------*/
face_system_to_component_table(dofs_per_face),
component_to_system_table(components, std::vector<unsigned>(dofs_per_cell)),
face_component_to_system_table(components, std::vector<unsigned>(dofs_per_face)),
- component_to_base_table(components),
+ component_to_base_table (components, 0),
restriction_is_additive_flags(restriction_is_additive_flags)
{
Assert(restriction_is_additive_flags.size()==fe_data.components,
-template <int dim>
-unsigned int
-FiniteElement<dim>::n_base_elements() const
-{
- // default implementation
- return 1;
-}
-
-
-
-
template <int dim>
unsigned int
FiniteElement<dim>::memory_consumption () const
}
+
+template <int dim>
+unsigned int
+FE_DGQ<dim>::n_base_elements () const
+{
+ return 1;
+};
+
+
+
+template <int dim>
+const FiniteElement<dim> &
+FE_DGQ<dim>::base_element (const unsigned int index) const
+{
+ Assert (index==0, ExcIndexRange(index, 0, 1));
+ return *this;
+};
+
+
+
template <int dim>
unsigned int
FE_DGQ<dim>::memory_consumption () const
+template <int dim>
+unsigned int
+FE_Q<dim>::n_base_elements () const
+{
+ return 1;
+};
+
+
+
+template <int dim>
+const FiniteElement<dim> &
+FE_Q<dim>::base_element (const unsigned int index) const
+{
+ Assert (index==0, ExcIndexRange(index, 0, 1));
+ return *this;
+};
+
+
+
template <int dim>
unsigned int
FE_Q<dim>::memory_consumption () const
for (unsigned int base=0 ; base < n_base_elements() ; ++base)
for (unsigned int m = 0; m < element_multiplicity(base); ++m)
component_to_base_table[total_index++] = base;
+ Assert (total_index == component_to_base_table.size(),
+ ExcInternalError());
// Initialize index table
// Multi-component base elements have
void
FESystem<dim>::build_face_table()
{
- unsigned total_index = 0;
- for (unsigned int base=0 ; base < n_base_elements() ; ++base)
- for (unsigned int m = 0; m < element_multiplicity(base); ++m)
- component_to_base_table[total_index++] = base;
-
// Initialize index table
// Multi-component base elements have
// to be thought of.
// 1. Vertices
- total_index = 0;
+ unsigned int total_index = 0;
for (unsigned int vertex_number= 0 ; vertex_number < GeometryInfo<dim>::vertices_per_face ;
++vertex_number)
{
+template <int dim>
+const FiniteElement<dim> &
+FESystem<dim>::base_element (const unsigned int index) const
+{
+ Assert (index < base_elements.size(),
+ ExcIndexRange(index, 0, base_elements.size()));
+ return *base_elements[index].first;
+};
+
+
+
+template <int dim>
+unsigned int
+FESystem<dim>::n_base_elements () const
+{
+ return base_elements.size();
+};
+
+
+
template <int dim>
unsigned int
FESystem<dim>::memory_consumption () const