]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
The lexicographic <-> hierarchic numbering functions now take a FiniteElementData...
authorhartmann <hartmann@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 13 Jun 2005 11:47:55 +0000 (11:47 +0000)
committerhartmann <hartmann@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 13 Jun 2005 11:47:55 +0000 (11:47 +0000)
git-svn-id: https://svn.dealii.org/trunk@10857 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_tools.h
deal.II/deal.II/source/fe/fe_tools.cc

index b8ea6de42935e1f58fbfa5a2c3628d82a72b7dd5..aec40f3c172b759a63b36c9aefa5b307c62bb972 100644 (file)
@@ -25,7 +25,7 @@ template <typename number> class Vector;
 template <int dim> class FiniteElement;
 template <int dim> class DoFHandler;
 template <int dim> class hpDoFHandler;
-template <int dim> class FE_Q;
+template <int dim> class FiniteElementData;
 class ConstraintMatrix;
 
 
@@ -566,8 +566,13 @@ class FETools
                                      * Note that since this function
                                      * uses specifics of the
                                      * continuous finite elements, it
-                                     * can only operate on objects of
-                                     * type FE_Q().
+                                     * can only operate on
+                                     * FiniteElementData<dim> objects
+                                     * inherent in FE_Q(). However,
+                                     * this function does not take a
+                                     * FE_Q object as it is also
+                                     * invoked by the FE_Q()
+                                     * constructor.
                                      *
                                      * It is assumed that the size of
                                      * the output argument already
@@ -578,8 +583,8 @@ class FETools
                                      */
     template <int dim>
     static void
-    hierarchic_to_lexicographic_numbering (const FE_Q<dim>           &fe,
-                                          std::vector<unsigned int> &h2l);
+    hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe_data,
+                                          std::vector<unsigned int>    &h2l);
 
                                     /**
                                      * This is the reverse function
@@ -592,8 +597,8 @@ class FETools
                                      */
     template <int dim>
     static void
-    lexicographic_to_hierarchic_numbering (const FE_Q<dim>           &fe,
-                                          std::vector<unsigned int> &l2h);
+    lexicographic_to_hierarchic_numbering (const FiniteElementData<dim> &fe_data,
+                                          std::vector<unsigned int>    &l2h);
 
                                     /**
                                      * Given a name in the form which
index 3d0d6198a450c1036ace6ee899546bc824d69f08..f38b1c24a66c86fb34ebeea8d1061e7a95e1ffad 100644 (file)
@@ -1073,7 +1073,7 @@ void FETools::extrapolate(const DoFHandler<dim> &dof1,
 
 template <int dim>
 void
-FETools::hierarchic_to_lexicographic_numbering (const FE_Q<dim>           &fe,
+FETools::hierarchic_to_lexicographic_numbering (const FiniteElementData<dim> &fe,
                                                std::vector<unsigned int> &h2l)
 {
   Assert (fe.n_components() == 1, ExcInvalidFE());
@@ -1233,8 +1233,8 @@ FETools::hierarchic_to_lexicographic_numbering (const FE_Q<dim>           &fe,
 
 template <int dim>
 void
-FETools::lexicographic_to_hierarchic_numbering (const FE_Q<dim>           &fe,
-                                               std::vector<unsigned int> &l2h)
+FETools::lexicographic_to_hierarchic_numbering (const FiniteElementData<dim> &fe,
+                                               std::vector<unsigned int>    &l2h)
 {
                                   // note: this function does the
                                   // reverse operation of the
@@ -2078,13 +2078,13 @@ void FETools::extrapolate<deal_II_dimension>
 template
 void
 FETools::hierarchic_to_lexicographic_numbering<deal_II_dimension>
-(const FE_Q<deal_II_dimension> &fe,
- std::vector<unsigned int>     &h2l);
+(const FiniteElementData<deal_II_dimension> &fe,
+ std::vector<unsigned int>                  &h2l);
 template
 void
 FETools::lexicographic_to_hierarchic_numbering<deal_II_dimension>
-(const FE_Q<deal_II_dimension> &fe,
- std::vector<unsigned int>     &h2l);
+(const FiniteElementData<deal_II_dimension> &fe,
+ std::vector<unsigned int>                  &h2l);
 
 template
 FiniteElement<deal_II_dimension> *

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.