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.
}
+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
{
/**
*/
template<int dim, int spacedim>
std::string policy_to_string(const dealii::internal::DoFHandler::Policy::PolicyBase<dim,spacedim> &policy);
-
}
// ---------------------------------------------------------------------
//
-// 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.
//
-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);
// ---------------------------------------------------------------------
-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