#include <vector>
#include <map>
-#include <list>
+#include <set>
#include <base/exceptions.h>
#include <base/smartpointer.h>
* consideration.
*/
unsigned int
- n_boundary_dofs (const list<unsigned char> &boundary_indicators) const;
+ n_boundary_dofs (const set<unsigned char> &boundary_indicators) const;
/**
* Return a constant reference to the
#include <base/exceptions.h>
#include <vector>
-#include <list>
+#include <set>
class SparsityPattern;
template <typename number> class Vector;
* belong to specified components
* of the solution. The function
* returns its results in the
- * last parameter which contains
- * @p{true} is a degree of freedom
- * is at the boundary and belongs
- * to one of the selected
+ * last non-default-valued
+ * parameter which contains
+ * @p{true} if a degree of
+ * freedom is at the boundary and
+ * belongs to one of the selected
* components, and @p{false}
* otherwise.
*
+ * By specifying the
+ * @p{boundary_indicator}
+ * variable, you can select which
+ * boundary indicators the faces
+ * have to have on which the
+ * degrees of freedom are located
+ * that shall be extracted. If it
+ * is an empty list, then all
+ * boundary indicators are
+ * accepted.
+ *
* The size of @p{component_select}
* shall equal the number of
* components in the finite
static void
extract_boundary_dofs (const DoFHandler<dim> &dof_handler,
const vector<bool> &component_select,
- vector<bool> &selected_dofs);
+ vector<bool> &selected_dofs,
+ const set<unsigned char> &boundary_indicators = set<unsigned char>());
/**
* Select all dofs that will be
*/
template <int dim>
static void
- map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
- const list<unsigned char> &boundary_indicators,
- vector<unsigned int> &mapping);
+ map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
+ const set<unsigned char> &boundary_indicators,
+ vector<unsigned int> &mapping);
/**
template <>
-unsigned int DoFHandler<1>::n_boundary_dofs (const list<unsigned char> &) const
+unsigned int DoFHandler<1>::n_boundary_dofs (const set<unsigned char> &) const
{
Assert (selected_fe != 0, ExcNoFESelected());
Assert (false, ExcNotImplemented());
template <int dim>
unsigned int
-DoFHandler<dim>::n_boundary_dofs (const list<unsigned char> &boundary_indicators) const
+DoFHandler<dim>::n_boundary_dofs (const set<unsigned char> &boundary_indicators) const
{
Assert (selected_fe != 0, ExcNoFESelected());
- Assert (find (boundary_indicators.begin(),
- boundary_indicators.end(),
- 255) ==
- boundary_indicators.end(),
+ Assert (boundary_indicators.find (255) == boundary_indicators.end(),
ExcInvalidBoundaryIndicator());
set<int> boundary_dofs;
void
DoFTools::extract_boundary_dofs (const DoFHandler<dim> &dof_handler,
const vector<bool> &component_select,
- vector<bool> &selected_dofs)
+ vector<bool> &selected_dofs,
+ const set<unsigned char> &boundary_indicators)
{
Assert (component_select.size() == dof_handler.get_fe().n_components(),
ExcWrongSize (component_select.size(),
dof_handler.get_fe().n_components()));
+ Assert (boundary_indicators.find (255) == boundary_indicators.end(),
+ ExcInvalidBoundaryIndicator());
+ // let's see whether we have to
+ // check for certain boundary
+ // indicators or whether we can
+ // accept all
+ const bool check_boundary_indicator = (boundary_indicators.size() != 0);
+
// clear and reset array by default
// values
selected_dofs.clear ();
cell!=dof_handler.end(); ++cell)
for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
if (cell->at_boundary(face))
+ if (! check_boundary_indicator ||
+ (boundary_indicators.find (cell->face(face)->boundary_indicator())
+ != boundary_indicators.end()))
{
cell->face(face)->get_dof_indices (face_dof_indices);
for (unsigned int i=0; i<dof_handler.get_fe().dofs_per_face; ++i)
template <>
void DoFTools::map_dof_to_boundary_indices (const DoFHandler<1> &dof_handler,
- const list<unsigned char> &,
+ const set<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 list<unsigned char> &boundary_indicators,
- vector<unsigned int> &mapping)
+void DoFTools::map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
+ const set<unsigned char> &boundary_indicators,
+ vector<unsigned int> &mapping)
{
Assert (&dof_handler.get_fe() != 0, ExcNoFESelected());
- Assert (find (boundary_indicators.begin(),
- boundary_indicators.end(),
- 255) ==
- boundary_indicators.end(),
+ Assert (boundary_indicators.find (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 (find (boundary_indicators.begin(),
- boundary_indicators.end(),
- face->boundary_indicator()) !=
+ if (boundary_indicators.find (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 list<unsigned char> &,
+ const set<unsigned char> &,
vector<unsigned int> &);
#endif
ExcComponentMismatch());
vector<unsigned int> dof_to_boundary_mapping;
- list<unsigned char> selected_boundary_components;
+ set<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);
+ selected_boundary_components.insert (i->first);
DoFTools::map_dof_to_boundary_indices (dof, selected_boundary_components,
dof_to_boundary_mapping);
<h3>deal.II</h3>
<ol>
+ <li> <p>
+ Extend: <code class="member">DoFTools::extract_boundary_dofs</code>
+ now allows to also specify which boundary conditions shall be
+ considered.
+ <br>
+ (WB 2000/12/04)
+ </p>
+
+ <li> <p>
+ New: some arguments of the functions
+ <code class="member">DoFHandler::n_boundary_dofs</code>,
+ <code class="member">DoFTools::extract_boundary_dofs</code>,
+ and
+ <code class="member">DoFTools::map_dof_to_boundary_indices</code>
+ are changed from <code class="class">list</code> to
+ <code class="class">set</code>, since that resembles more
+ closely the purpose of the parameter, and makes computations
+ slightly faster.
+ <br>
+ (WB 2000/12/04)
+ </p>
+
<li> <p>
New: almost all classes that store data now have a function
<code class="member">memory_consumption</code> that returns an