From: Wolfgang Bangerth Date: Tue, 25 Aug 1998 16:18:33 +0000 (+0000) Subject: Move get_selected_fe to the .h file to allow for making it inline. X-Git-Tag: v8.0.0~22737 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=3e46b37ba884960e58dc22114376026e981c91b5;p=dealii.git Move get_selected_fe to the .h file to allow for making it inline. git-svn-id: https://svn.dealii.org/trunk@517 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/dofs/dof_handler.h b/deal.II/deal.II/include/dofs/dof_handler.h index 58634fe5bc..50233197a2 100644 --- a/deal.II/deal.II/include/dofs/dof_handler.h +++ b/deal.II/deal.II/include/dofs/dof_handler.h @@ -482,7 +482,7 @@ class DoFHandler : public DoFDimensionInfo { /** * Destructor */ - ~DoFHandler (); + virtual ~DoFHandler (); /** * Go through the triangulation and @@ -499,7 +499,7 @@ class DoFHandler : public DoFDimensionInfo { * don't need them after calling this * function, or if so store them. */ - void distribute_dofs (const FiniteElementBase &); + virtual void distribute_dofs (const FiniteElementBase &); /** * Renumber the degrees of freedom according @@ -513,9 +513,9 @@ class DoFHandler : public DoFDimensionInfo { * See the general documentation of this * class for more details. */ - void renumber_dofs (const RenumberingMethod method, - const bool use_constraints = false, - const vector &starting_points = vector()); + virtual void renumber_dofs (const RenumberingMethod method, + const bool use_constraints = false, + const vector &starting_points = vector()); /** * Make up the constraint matrix which @@ -1161,7 +1161,9 @@ class DoFHandler : public DoFDimensionInfo { /** * Return a constant reference to the - * selected finite element object. + * selected finite element object. This + * function is inline, so it should + * be reasonably fast. */ const FiniteElementBase & get_selected_fe () const; @@ -1227,8 +1229,27 @@ class DoFHandler : public DoFDimensionInfo { int, int, << "The dimension " << arg1 << " of the vector is wrong. " << "It should be " << arg2); - + protected: + + /** + * Address of the triangulation to + * work on. + */ + Triangulation * const tria; + + /** + * Store a copy of the finite element given + * latest for the distribution of dofs. In + * fact, since the FE base class itself has + * not much functionality, this object only + * stores numbers such as the number of + * dofs per line, but also restriction + * and prolongation matrices, etc. + */ + const FiniteElementBase *selected_fe; + + private: /** * Reserve enough space in the * #levels[]# objects to store the @@ -1288,12 +1309,6 @@ class DoFHandler : public DoFDimensionInfo { const cell_iterator &new_cell, dSMatrix &transfer_matrix) const; - /** - * Address of the triangulation to - * work on. - */ - Triangulation * const tria; - /** * Space to store the DoF numbers for the * different levels. Analogous to the @@ -1302,17 +1317,6 @@ class DoFHandler : public DoFDimensionInfo { */ vector*> levels; - /** - * Store a copy of the finite element given - * latest for the distribution of dofs. In - * fact, since the FE base class itself has - * not much functionality, this object only - * stores numbers such as the number of - * dofs per line, but also restriction - * and prolongation matrices, etc. - */ - const FiniteElementBase *selected_fe; - /** * Store the number of dofs created last * time. @@ -1337,6 +1341,20 @@ class DoFHandler : public DoFDimensionInfo { +/* ----------------------- Inline functions ---------------------------------- */ + + +template +inline +const FiniteElementBase & DoFHandler::get_selected_fe () const { + return *selected_fe; +}; + + + + + + /*---------------------------- dof.h ---------------------------*/ /* end of #ifndef __dof_H */ #endif diff --git a/deal.II/deal.II/source/dofs/dof_handler.cc b/deal.II/deal.II/source/dofs/dof_handler.cc index e92cb40417..1ff01dda01 100644 --- a/deal.II/deal.II/source/dofs/dof_handler.cc +++ b/deal.II/deal.II/source/dofs/dof_handler.cc @@ -759,13 +759,6 @@ const Triangulation & DoFHandler::get_tria () const { -template -const FiniteElementBase & DoFHandler::get_selected_fe () const { - return *selected_fe; -}; - - - template void DoFHandler::distribute_dofs (const FiniteElementBase &fe) { Assert (tria->n_levels() > 0, ExcInvalidTriangulation());