]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add extract_boundary_dofs.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 17 Apr 2000 14:33:12 +0000 (14:33 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 17 Apr 2000 14:33:12 +0000 (14:33 +0000)
git-svn-id: https://svn.dealii.org/trunk@2734 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/dofs/dof_tools.h
deal.II/deal.II/source/dofs/dof_tools.cc
deal.II/doc/news/2000/c-3-0.html

index d571b1114d1d2388a54138e1e6c583a41139ccbf..657f849d6a93b1b25e63d8ecef871dd99e934854 100644 (file)
@@ -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 <int dim>
     static void
-    extract_dofs(const DoFHandler<dim> &dof_handler,
-                const vector<bool>    &select,
-                vector<bool>          &selected_dofs);
+    extract_dofs (const DoFHandler<dim> &dof_handler,
+                 const vector<bool>    &component_select,
+                 vector<bool>          &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 <int dim>
     static void
-    extract_level_dofs(const unsigned int       level,
-                      const MGDoFHandler<dim> &dof,
-                      const vector<bool>      &select,
-                      vector<bool>            &selected_dofs);
+    extract_level_dofs (const unsigned int       level,
+                       const MGDoFHandler<dim> &dof,
+                       const vector<bool>      &select,
+                       vector<bool>            &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 <int dim>
+    static void
+    extract_boundary_dofs (const DoFHandler<dim> &dof_handler,
+                          const vector<bool>    &component_select,
+                          vector<bool>          &selected_dofs);
+    
                                     /**
                                      * This function can be used when
                                      * different variables shall be
index dc8c1768af4720917539dac251adaecefd2c31be..6fb8885849d034aa0c6aab54e8b462b7fc27db30 100644 (file)
@@ -711,6 +711,74 @@ DoFTools::extract_level_dofs(const unsigned int       level,
 
 
 
+#if deal_II_dimension != 1
+
+template <int dim>
+void
+DoFTools::extract_boundary_dofs (const DoFHandler<dim> &dof_handler,
+                                const vector<bool>    &component_select,
+                                vector<bool>          &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<unsigned int> face_dof_indices (dof_handler.get_fe().dofs_per_face);
+  for (DoFHandler<dim>::active_cell_iterator cell=dof_handler.begin_active();
+       cell!=dof_handler.end(); ++cell)
+    for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
+      if (cell->at_boundary(face))
+       {
+         cell->face(face)->get_dof_indices (face_dof_indices);
+         for (unsigned int i=0; i<dof_handler.get_fe().dofs_per_face; ++i)
+           if (component_select[dof_handler.get_fe().
+                               face_system_to_component_index(i).first] == true)
+             selected_dofs[face_dof_indices[i]] = true;
+       };
+};
+
+
+#else
+
+template <>
+void
+DoFTools::extract_boundary_dofs (const DoFHandler<1> &dof_handler,
+                                const vector<bool>  &component_select,
+                                vector<bool>        &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; i<dof_handler.get_fe().dofs_per_face; ++i)
+         if (component_select[dof_handler.get_fe().
+                             face_system_to_component_index(i).first] == true)
+           selected_dofs[cell->vertex_dof_index(0,i)] = true;
+                                      // check right-most vertex
+      if (cell->neighbor(1) == dof_handler.end())
+       for (unsigned int i=0; i<dof_handler.get_fe().dofs_per_face; ++i)
+         if (component_select[dof_handler.get_fe().
+                             face_system_to_component_index(i).first] == true)
+           selected_dofs[cell->vertex_dof_index(1,i)] = true;
+    };
+};
+
+
+#endif
+
+
 
 template <int dim>
 void
@@ -1195,6 +1263,13 @@ template void DoFTools::extract_level_dofs(unsigned int level,
                                           const vector<bool>& local_select,
                                           vector<bool>& selected_dofs);
 
+#if deal_II_dimension != 1
+template
+void
+DoFTools::extract_boundary_dofs (const DoFHandler<deal_II_dimension> &,
+                                const vector<bool>                  &,
+                                vector<bool>                        &);
+#endif 
 
 template
 void
index 4aa72144ea88b6b99cafb4e7a4a65aa8ce3c5ca9..965b93476285fc95974dfc2e5d1be8566891691b 100644 (file)
   estimates the norm of the gradient on each cell from finite
   difference approximations. 
 
+  <li> New: <code class="member">DoFTools::extract_boundary_dofs</code> 
+  finds all degrees of freedom which are at the boundary and belong to
+  specified components.
+
   <li> New: <code class="member">DoFTools::compute_intergrid_constraints</code> 
   allows to use different discretization grids for different variables.
 </ol>

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.