#include <vector>
#include <map>
+#include <list>
#include <base/exceptions.h>
#include <base/smartpointer.h>
unsigned int n_boundary_dofs () const;
/**
- * Return the number of degrees of freedom
- * located on those parts of the boundary
- * which have a boundary indicator listed
- * in the given set. The reason that a
- * @p{map} rather than a @p{set} is used is the
- * same as descibed in the section on the
- * @p{make_boundary_sparsity_pattern} function.
- */
- unsigned int n_boundary_dofs (const FunctionMap &boundary_indicators) const;
+ * Return the number of degrees
+ * of freedom located on those
+ * parts of the boundary which
+ * have a boundary indicator
+ * listed in the given set. The
+ * reason that a @p{map} rather
+ * than a @p{set} is used is the
+ * same as descibed in the
+ * section on the
+ * @p{make_boundary_sparsity_pattern}
+ * function.
+ */
+ unsigned int
+ n_boundary_dofs (const FunctionMap &boundary_indicators) const;
+
+ /**
+ * Same function, but with
+ * different data type of the
+ * argument, which is here simply
+ * a list of the boundary
+ * indicators under
+ * consideration.
+ */
+ unsigned int
+ n_boundary_dofs (const list<unsigned char> &boundary_indicators) const;
/**
* Return a constant reference to the
#include <base/exceptions.h>
#include <vector>
-#include <map>
+#include <list>
class SparsityPattern;
template <typename number> class Vector;
* numbers of the trial functions
* local to the boundary.
*
- * Prior content of @p{mapping} is deleted.
+ * Prior content of @p{mapping}
+ * is deleted.
*
* This function is not
* implemented for one
vector<unsigned int> &mapping);
/**
- * Same as the previous function, except
- * that only selected parts of the
- * boundary are considered.
+ * Same as the previous function,
+ * except that only those parts
+ * of the boundary are considered
+ * for which the boundary
+ * indicator is listed in the
+ * second argument.
*
- * See the general doc of this class for
- * more information.
+ * See the general doc of this
+ * class for more information.
*/
template <int dim>
static void
- map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
- const map<unsigned char,const Function<dim>*> &boundary_indicators,
- vector<unsigned int> &mapping);
+ map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
+ const list<unsigned char> &boundary_indicators,
+ vector<unsigned int> &mapping);
/**
#if deal_II_dimension == 1
template <>
-unsigned int DoFHandler<1>::n_boundary_dofs () const {
+unsigned int DoFHandler<1>::n_boundary_dofs () const
+{
Assert (selected_fe != 0, ExcNoFESelected());
Assert (false, ExcNotImplemented());
return 0;
};
+
template <>
-unsigned int DoFHandler<1>::n_boundary_dofs (const FunctionMap &) const {
+unsigned int DoFHandler<1>::n_boundary_dofs (const FunctionMap &) const
+{
+ Assert (selected_fe != 0, ExcNoFESelected());
+ Assert (false, ExcNotImplemented());
+ return 0;
+};
+
+
+
+template <>
+unsigned int DoFHandler<1>::n_boundary_dofs (const list<unsigned char> &) const
+{
Assert (selected_fe != 0, ExcNoFESelected());
Assert (false, ExcNotImplemented());
return 0;
template <int dim>
-unsigned int DoFHandler<dim>::n_boundary_dofs () const {
+unsigned int DoFHandler<dim>::n_boundary_dofs () const
+{
Assert (selected_fe != 0, ExcNoFESelected());
set<int> boundary_dofs;
};
+
template <int dim>
-unsigned int DoFHandler<dim>::n_boundary_dofs (const FunctionMap &boundary_indicators) const {
+unsigned int
+DoFHandler<dim>::n_boundary_dofs (const FunctionMap &boundary_indicators) const
+{
Assert (selected_fe != 0, ExcNoFESelected());
Assert (boundary_indicators.find(255) == boundary_indicators.end(),
ExcInvalidBoundaryIndicator());
};
+
template <int dim>
-const Triangulation<dim> & DoFHandler<dim>::get_tria () const {
+unsigned int
+DoFHandler<dim>::n_boundary_dofs (const list<unsigned char> &boundary_indicators) const
+{
+ Assert (selected_fe != 0, ExcNoFESelected());
+ Assert (find (boundary_indicators.begin(),
+ boundary_indicators.end(),
+ 255) ==
+ boundary_indicators.end(),
+ ExcInvalidBoundaryIndicator());
+
+ set<int> boundary_dofs;
+
+ const unsigned int dofs_per_face = selected_fe->dofs_per_face;
+ vector<unsigned int> dofs_on_face(dofs_per_face);
+ active_face_iterator face = begin_active_face (),
+ endf = end_face();
+ for (; face!=endf; ++face)
+ if (find (boundary_indicators.begin(),
+ boundary_indicators.end(),
+ face->boundary_indicator()) !=
+ boundary_indicators.end())
+ {
+ face->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>
+const Triangulation<dim> & DoFHandler<dim>::get_tria () const
+{
return *tria;
};
+
template <int dim>
void DoFHandler<dim>::distribute_dofs (const FiniteElement<dim> &ff,
const unsigned int offset)
template <>
void DoFTools::map_dof_to_boundary_indices (const DoFHandler<1> &dof_handler,
- const map<unsigned char,const Function<1>*> &,
+ const list<unsigned char> &,
vector<unsigned int> &)
{
Assert (&dof_handler.get_fe() != 0, ExcNoFESelected());
template <int dim>
-void DoFTools::map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
- const map<unsigned char,const Function<dim>*> &boundary_indicators,
- vector<unsigned int> &mapping)
+void DoFTools::map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
+ const list<unsigned char> &boundary_indicators,
+ vector<unsigned int> &mapping)
{
Assert (&dof_handler.get_fe() != 0, ExcNoFESelected());
- Assert (boundary_indicators.find(255) == boundary_indicators.end(),
+ Assert (find (boundary_indicators.begin(),
+ boundary_indicators.end(),
+ 255) ==
+ boundary_indicators.end(),
ExcInvalidBoundaryIndicator());
mapping.clear ();
typename DoFHandler<dim>::active_face_iterator face = dof_handler.begin_active_face(),
endf = dof_handler.end_face();
for (; face!=endf; ++face)
- if (boundary_indicators.find(face->boundary_indicator()) !=
+ if (find (boundary_indicators.begin(),
+ boundary_indicators.end(),
+ face->boundary_indicator()) !=
boundary_indicators.end())
{
face->get_dof_indices (dofs_on_face);
template
void
DoFTools::map_dof_to_boundary_indices (const DoFHandler<deal_II_dimension> &,
- const map<unsigned char,const Function<deal_II_dimension>*> &,
+ const list<unsigned char> &,
vector<unsigned int> &);
#endif
ExcComponentMismatch());
vector<unsigned int> dof_to_boundary_mapping;
- DoFTools::map_dof_to_boundary_indices (dof, boundary_functions,
+ list<unsigned char> selected_boundary_components;
+ for (typename map<unsigned char,const Function<dim>*>::const_iterator
+ i=boundary_functions.begin();
+ i!=boundary_functions.end(); ++i)
+ selected_boundary_components.push_back (i->first);
+
+ DoFTools::map_dof_to_boundary_indices (dof, selected_boundary_components,
dof_to_boundary_mapping);
// set up sparsity structure
<h3>deal.II</h3>
<ol>
+ <li> <p>
+ New: There is now a function
+ <code class="member">DoFHandler::n_boundary_dofs</code>
+ that takes the list of selected boundary indicators as a
+ <code class="class">list</code> of values, rather than the
+ usual <code class="class">map</code> of pairs of boundary
+ indicators and function object pointers.
+ <br>
+ (WB 2000/08/25)
+ </p>
+
<li> <p>
Changed: The
<code class="member">map_dof_to_boundary_index</code>