]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a FECollection::n_blocks function.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 20 Sep 2012 13:37:12 +0000 (13:37 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 20 Sep 2012 13:37:12 +0000 (13:37 +0000)
git-svn-id: https://svn.dealii.org/trunk@26564 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/hp/fe_collection.h
deal.II/source/hp/fe_collection.cc

index 91de6a82898589f0e556a95d03d209931724112d..68114248461d41e61b8585c72ffb2f0e298cf229 100644 (file)
@@ -53,6 +53,11 @@ DoFHandler, in particular removal of specializations.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> New: There is now a function hp::FECollection::n_blocks() in analogy to
+the existing function hp::FECollection::n_components().
+<br>
+(Wolfgang Bangerth, 2012/09/20)
+
 <li> Changed: step-8 now outputs data in VTK format, rather than GMV.
 GMV has long been dead.
 <br>
index cfe276996492b5690f3e3db97fbfc431056e93d7..ff0f550f773120c16631d262138f5b9e8d06d3f5 100644 (file)
@@ -124,6 +124,27 @@ namespace hp
                                         */
       unsigned int n_components () const;
 
+      /**
+       * Return the number of vector
+       * blocks of the finite elements in
+       * this collection. While this class ensures that all
+       * elements stored in it have the same number of vector
+       * components, there is no such guarantees for the number
+       * of blocks each element is made up of (an element may
+       * have fewer blocks than vector components; see
+       * @ref GlossBlock "the glossary" for more information).
+       * For example, you may have an FECollection object that stores
+       * one copy of an FESystem with <code>dim</code> FE_Q objects
+       * and one copy of an FE_RaviartThomas element. Both have
+       * <code>dim</code> vector components but while the former has
+       * <code>dim</code> blocks the latter has only one.
+       * Consequently, this function will throw an assertion
+       * if the number of blocks is not the same for all elements.
+       * If they are the same, this function returns the result
+       * of FiniteElement::n_blocks().
+       */
+      unsigned int n_blocks () const;
+
                                        /**
                                         * Return the maximal number of degrees
                                         * of freedom per vertex over all
index 033eb460139e823fa92b6c9945521453878e0f55..cb6ed6399bb5291dccd8df32c9a55d69f0238d05 100644 (file)
@@ -72,6 +72,23 @@ namespace hp
 
 
 
+  template <int dim, int spacedim>
+  unsigned int
+  FECollection<dim,spacedim>::n_blocks () const
+  {
+    Assert (finite_elements.size () > 0, ExcNoFiniteElements());
+
+    const unsigned int nb = finite_elements[0]->n_blocks ();
+    for (unsigned int i=1; i<finite_elements.size(); ++i)
+      Assert (finite_elements[i]->n_blocks() == nb,
+          ExcMessage ("Not all finite elements in this collection have "
+                      "the same number of components."));
+
+    return nb;
+  }
+
+
+
   template <int dim, int spacedim>
   std::size_t
   FECollection<dim,spacedim>::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.