]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Slightly extend the infrastructure around base elements and components. Significantly...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 19 Sep 2001 14:02:42 +0000 (14:02 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 19 Sep 2001 14:02:42 +0000 (14:02 +0000)
git-svn-id: https://svn.dealii.org/trunk@5027 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe.h
deal.II/deal.II/include/fe/fe_base.h
deal.II/deal.II/include/fe/fe_dgq.h
deal.II/deal.II/include/fe/fe_q.h
deal.II/deal.II/include/fe/fe_system.h
deal.II/deal.II/source/fe/fe.cc
deal.II/deal.II/source/fe/fe_dgq.cc
deal.II/deal.II/source/fe/fe_q.cc
deal.II/deal.II/source/fe/fe_system.cc

index 34ae55012cd8e51fd67ed0abc1d4f73917a97284..ecfbc594dbcfe287cbb894e1c8f835e12ff14952 100644 (file)
@@ -69,24 +69,38 @@ class FiniteElement : public FiniteElementBase<dim>
     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)
index f3935c6aa857e9e79d00fb654c02750a84ae9a66..aa533bd5b3d97a628b8ae0ffa41c0b867525c4f8 100644 (file)
@@ -625,7 +625,7 @@ class FiniteElementBase : public Subscriptor,
                                      * 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}
@@ -881,6 +881,17 @@ class FiniteElementBase : public Subscriptor,
                                      * 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;
 
index 7f22635c3f14d297ee25ca30e15b80080d9c9896..8a89d33ae95660b411bcb559fb528a40a5270520 100644 (file)
@@ -86,6 +86,23 @@ class FE_DGQ : public FiniteElement<dim>
                                      */
     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)
index 40bfd075ff8fc5f9088543872a27c2f138729444..02e892b418462421ef57c6c66f67dc23c7b49ae5 100644 (file)
@@ -291,6 +291,23 @@ class FE_Q : public FiniteElement<dim>
                                      */
     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)
index 2b2b0b25b8b45248672bde855fede646eace65c3..b92787b70970310026d44f961266700bdd5edc4f 100644 (file)
@@ -197,28 +197,25 @@ class FESystem : public FiniteElement<dim>
     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)
@@ -626,14 +623,6 @@ template <> void FESystem<1>::initialize_unit_face_support_points ();
 
 /* ------------------------- inline functions ------------------------- */
 
-template<int dim>
-inline unsigned int
-FESystem<dim>::n_base_elements() const
-{
-  return base_elements.size();
-};
-
-
 
 template<int dim>
 inline unsigned int
@@ -646,19 +635,6 @@ FESystem<dim>::element_multiplicity (const unsigned int index) const
 
 
 
-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  ---------------------------*/
index d77a428651883d646cdcf11162a5a484c665225b..f3d88928142e5ed3e2d7a83b451a95a7c00a284e 100644 (file)
@@ -113,7 +113,7 @@ FiniteElementBase<dim>::FiniteElementBase (const FiniteElementData<dim> &fe_data
                 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,
@@ -444,17 +444,6 @@ FiniteElement<dim>::get_subface_data (const UpdateFlags        flags,
 
 
 
-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
index 6574b23e89cbb0b02252102baf8cae20b6277ab5..77b6340ab73ce3fb8b61d9977796984cb3193345 100644 (file)
@@ -563,6 +563,26 @@ FE_DGQ<dim>::fill_fe_subface_values (const Mapping<dim>                   &mappi
 }
 
 
+
+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
index ccb2cb0504ef5b5bdeee3a62bb77ec33fb9110d0..8de210c4665bfef7cca9fd51690f73d1ec4e22a6 100644 (file)
@@ -1191,6 +1191,25 @@ FE_Q<dim>::fill_fe_subface_values (const Mapping<dim>                   &mapping
 
 
 
+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
index 1e6c8eed2791bbf54a100c85f39636a0633de47e..191c766102156dd360314b87f6f63cf7bdbb1199 100644 (file)
@@ -664,6 +664,8 @@ FESystem<dim>::build_cell_table()
   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
@@ -793,17 +795,12 @@ template <int dim>
 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)
     {
@@ -1423,6 +1420,26 @@ FESystem<dim>::compute_restriction_is_additive_flags (const FiniteElement<dim> &
 
 
 
+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 

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.