]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move template function DoFHandler::n_boundary_dofs().
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 17 Nov 2016 03:35:50 +0000 (20:35 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 17 Nov 2016 13:46:17 +0000 (06:46 -0700)
This function is currently in a .cc file, and the .inst.in file actually tries to
instantiate it for all template arguments. Nonetheless, as #3599 shows, it is not
instantiated appropriately.

Rather than try to understand the exact cause, the easier solution is to just
move the function to the .h file. It turns out to have a lot of code duplicated
from the other function of same name, so shrink it down to the minimum and
simply defer to the non-templated member function.

include/deal.II/dofs/dof_handler.h
source/dofs/dof_handler.cc
source/dofs/dof_handler.inst.in

index bb89fba363524f571ea29e76561c7aa12a0102fd..046509edf02225a0c46dacf45d44a7b1417feb90 100644 (file)
@@ -1236,6 +1236,24 @@ DoFHandler<dim,spacedim>::block_info () const
 }
 
 
+template <int dim, int spacedim>
+template <typename number>
+types::global_dof_index
+DoFHandler<dim,spacedim>::n_boundary_dofs (const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_ids) const
+{
+  // extract the set of boundary ids and forget about the function object pointers
+  std::set<types::boundary_id> boundary_ids_only;
+  for (typename std::map<types::boundary_id, const Function<spacedim,number>*>::const_iterator
+       p = boundary_ids.begin();
+       p != boundary_ids.end(); ++p)
+    boundary_ids_only.insert (p->first);
+
+  // then just hand everything over to the other function that does the work
+  return n_boundary_dofs(boundary_ids_only);
+}
+
+
+
 namespace internal
 {
   /**
@@ -1247,7 +1265,6 @@ namespace internal
    */
   template<int dim, int spacedim>
   std::string policy_to_string(const dealii::internal::DoFHandler::Policy::PolicyBase<dim,spacedim> &policy);
-
 }
 
 
index d97af68e5b94cfe61af91fec646e4e5086a3fcd8..2031a70fec1ebd1ec34466fe5d2a9ed2ba30d9d7 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 1998 - 2015 by the deal.II authors
+// Copyright (C) 1998 - 2016 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -1009,49 +1009,14 @@ types::global_dof_index DoFHandler<dim,spacedim>::n_boundary_dofs () const
 
 
 
-template<int dim, int spacedim>
-template<typename number>
-types::global_dof_index
-DoFHandler<dim,spacedim>::n_boundary_dofs (const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_ids) const
-{
-  Assert (boundary_ids.find(numbers::internal_face_boundary_id) == boundary_ids.end(),
-          ExcInvalidBoundaryIndicator());
-
-  std::set<types::global_dof_index> boundary_dofs;
-
-  const unsigned int dofs_per_face = get_fe().dofs_per_face;
-  std::vector<types::global_dof_index> dofs_on_face(dofs_per_face);
-
-  // same as in the previous
-  // function, but with an additional
-  // check for the boundary indicator
-  active_cell_iterator cell = begin_active (),
-                       endc = end();
-  for (; cell!=endc; ++cell)
-    for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
-      if (cell->at_boundary(f)
-          &&
-          (boundary_ids.find(cell->face(f)->boundary_id()) !=
-           boundary_ids.end()))
-        {
-          cell->face(f)->get_dof_indices (dofs_on_face);
-          for (unsigned int i=0; i<dofs_per_face; ++i)
-            boundary_dofs.insert(dofs_on_face[i]);
-        }
-
-  return boundary_dofs.size();
-}
-
-
-
-template<int dim, int spacedim>
+template <int dim, int spacedim>
 types::global_dof_index
 DoFHandler<dim,spacedim>::n_boundary_dofs (const std::set<types::boundary_id> &boundary_ids) const
 {
   Assert (boundary_ids.find (numbers::internal_face_boundary_id) == boundary_ids.end(),
           ExcInvalidBoundaryIndicator());
 
-  std::set<int> boundary_dofs;
+  std::set<types::global_dof_index> boundary_dofs;
 
   const unsigned int dofs_per_face = get_fe().dofs_per_face;
   std::vector<types::global_dof_index> dofs_on_face(dofs_per_face);
index cacf64f2609998a54edd27df467378f18bb7875e..75aeba64f7b6180e1c91ea096df6504bd2cbc141 100644 (file)
 // ---------------------------------------------------------------------
 
 
-for (scalar: REAL_SCALARS; deal_II_dimension : DIMENSIONS)
-{
-    template types::global_dof_index DoFHandler<deal_II_dimension>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension,scalar>*> &boundary_ids) const;
-
-#if deal_II_dimension < 3
-    template types::global_dof_index DoFHandler<deal_II_dimension,deal_II_dimension+1>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension+1,scalar>*> &boundary_ids) const;
-#endif
-
-#if deal_II_dimension == 1
-    template types::global_dof_index DoFHandler<deal_II_dimension,deal_II_dimension+2>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension+2,scalar>*> &boundary_ids) const;
-#endif
-
-}
-
-for (scalar: COMPLEX_SCALARS; deal_II_dimension : DIMENSIONS)
-{
-    template types::global_dof_index DoFHandler<deal_II_dimension>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension,scalar>*> &boundary_ids) const;
-
-#if deal_II_dimension < 3
-    template types::global_dof_index DoFHandler<deal_II_dimension,deal_II_dimension+1>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension+1,scalar>*> &boundary_ids) const;
-#endif
-
-#if deal_II_dimension == 1
-    template types::global_dof_index DoFHandler<deal_II_dimension,deal_II_dimension+2>::n_boundary_dofs<scalar> (const std::map<types::boundary_id, const Function<deal_II_dimension+2,scalar>*> &boundary_ids) const;
-#endif
-
-}
-
 for (deal_II_dimension : DIMENSIONS)
 {
     namespace internal

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.