From 56ef36f411bd74c84a3b49a9d9b26cfa4a1ef2e1 Mon Sep 17 00:00:00 2001 From: wolf Date: Mon, 17 Apr 2000 14:33:12 +0000 Subject: [PATCH] Add extract_boundary_dofs. git-svn-id: https://svn.dealii.org/trunk@2734 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/deal.II/include/dofs/dof_tools.h | 51 ++++++++++++---- deal.II/deal.II/source/dofs/dof_tools.cc | 75 ++++++++++++++++++++++++ deal.II/doc/news/2000/c-3-0.html | 4 ++ 3 files changed, 117 insertions(+), 13 deletions(-) diff --git a/deal.II/deal.II/include/dofs/dof_tools.h b/deal.II/deal.II/include/dofs/dof_tools.h index d571b1114d..657f849d6a 100644 --- a/deal.II/deal.II/include/dofs/dof_tools.h +++ b/deal.II/deal.II/include/dofs/dof_tools.h @@ -350,28 +350,53 @@ class DoFTools * are then flagged #true#, while all * others are set to #false#. * - * The size of #select# shall equal - * the number of components in the - * finite element used by #dof#. + * The size of #component_select# + * shall equal the number of + * components in the finite + * element used by #dof#. */ template static void - extract_dofs(const DoFHandler &dof_handler, - const vector &select, - vector &selected_dofs); + extract_dofs (const DoFHandler &dof_handler, + const vector &component_select, + vector &selected_dofs); /** - * Do the same thing as #extract_dofs# - * for one level of a multi-grid DoF - * numbering. + * Do the same thing as + * #extract_dofs# for one level + * of a multi-grid DoF numbering. */ template static void - extract_level_dofs(const unsigned int level, - const MGDoFHandler &dof, - const vector &select, - vector &selected_dofs); + extract_level_dofs (const unsigned int level, + const MGDoFHandler &dof, + const vector &select, + vector &selected_dofs); + /** + * Extract all degrees of freedom + * which are at the boundary and + * belong to specified components + * of the solution. The function + * returns its results in the + * last parameter which contains + * #true# is a degree of freedom + * is at the boundary and belongs + * to one of the selected + * components, and #false# + * otherwise. + * + * The size of #component_select# + * shall equal the number of + * components in the finite + * element used by #dof#. + */ + template + static void + extract_boundary_dofs (const DoFHandler &dof_handler, + const vector &component_select, + vector &selected_dofs); + /** * This function can be used when * different variables shall be diff --git a/deal.II/deal.II/source/dofs/dof_tools.cc b/deal.II/deal.II/source/dofs/dof_tools.cc index dc8c1768af..6fb8885849 100644 --- a/deal.II/deal.II/source/dofs/dof_tools.cc +++ b/deal.II/deal.II/source/dofs/dof_tools.cc @@ -711,6 +711,74 @@ DoFTools::extract_level_dofs(const unsigned int level, +#if deal_II_dimension != 1 + +template +void +DoFTools::extract_boundary_dofs (const DoFHandler &dof_handler, + const vector &component_select, + vector &selected_dofs) +{ + Assert (component_select.size() == dof_handler.get_fe().n_components(), + ExcWrongSize (component_select.size(), + dof_handler.get_fe().n_components())); + + selected_dofs.resize (dof_handler.n_dofs(), false); + vector face_dof_indices (dof_handler.get_fe().dofs_per_face); + for (DoFHandler::active_cell_iterator cell=dof_handler.begin_active(); + cell!=dof_handler.end(); ++cell) + for (unsigned int face=0; face::faces_per_cell; ++face) + if (cell->at_boundary(face)) + { + cell->face(face)->get_dof_indices (face_dof_indices); + for (unsigned int i=0; i +void +DoFTools::extract_boundary_dofs (const DoFHandler<1> &dof_handler, + const vector &component_select, + vector &selected_dofs) +{ + Assert (component_select.size() == dof_handler.get_fe().n_components(), + ExcWrongSize (component_select.size(), + dof_handler.get_fe().n_components())); + + selected_dofs.resize (dof_handler.n_dofs(), false); + + Assert (dof_handler.get_fe().dofs_per_face == dof_handler.get_fe().dofs_per_vertex, + ExcInternalError()); + + // loop over coarse grid cells + for (DoFHandler<1>::cell_iterator cell=dof_handler.begin(0); + cell!=dof_handler.end(0); ++cell) + { + // check left-most vertex + if (cell->neighbor(0) == dof_handler.end()) + for (unsigned int i=0; ivertex_dof_index(0,i)] = true; + // check right-most vertex + if (cell->neighbor(1) == dof_handler.end()) + for (unsigned int i=0; ivertex_dof_index(1,i)] = true; + }; +}; + + +#endif + + template void @@ -1195,6 +1263,13 @@ template void DoFTools::extract_level_dofs(unsigned int level, const vector& local_select, vector& selected_dofs); +#if deal_II_dimension != 1 +template +void +DoFTools::extract_boundary_dofs (const DoFHandler &, + const vector &, + vector &); +#endif template void diff --git a/deal.II/doc/news/2000/c-3-0.html b/deal.II/doc/news/2000/c-3-0.html index 4aa72144ea..965b934762 100644 --- a/deal.II/doc/news/2000/c-3-0.html +++ b/deal.II/doc/news/2000/c-3-0.html @@ -63,6 +63,10 @@ estimates the norm of the gradient on each cell from finite difference approximations. +
  • New: DoFTools::extract_boundary_dofs + finds all degrees of freedom which are at the boundary and belong to + specified components. +
  • New: DoFTools::compute_intergrid_constraints allows to use different discretization grids for different variables. -- 2.39.5