]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Change some args from maps to lists, and add a new n_boundary_dofs function to DoFHan...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 25 Aug 2000 14:36:04 +0000 (14:36 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 25 Aug 2000 14:36:04 +0000 (14:36 +0000)
git-svn-id: https://svn.dealii.org/trunk@3274 0785d39b-7218-0410-832d-ea1e28bc413d

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

index e4efb28af0b3e1c4326436a7506a193a172fb23c..1770c7d84e36e416d7de5db2c6390d14b55a3228 100644 (file)
@@ -18,6 +18,7 @@
 
 #include <vector>
 #include <map>
+#include <list>
 #include <base/exceptions.h>
 #include <base/smartpointer.h>
 
@@ -956,15 +957,31 @@ class DoFHandler  :  public Subscriptor,
     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
index 9db699788251d658e6cebf83c771bc7988e3716e..17863219757b4836d0656841e0e42dac501bc10e 100644 (file)
@@ -16,7 +16,7 @@
 
 #include <base/exceptions.h>
 #include <vector>
-#include <map>
+#include <list>
 
 class SparsityPattern;
 template <typename number> class Vector;
@@ -801,7 +801,8 @@ class DoFTools
                                      * 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
@@ -816,18 +817,21 @@ class DoFTools
                                 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);
     
                                   
                                     /**
index 74e1f23c88c562b8211f834429767f01b93cd11a..af6ae3d4f366109b3e299249a3d2c309117fefe5 100644 (file)
@@ -1127,15 +1127,28 @@ DoFHandler<dim>::last_active_hex () const {
 #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;
@@ -1145,7 +1158,8 @@ unsigned int DoFHandler<1>::n_boundary_dofs (const FunctionMap &) const {
 
 
 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;
@@ -1165,8 +1179,11 @@ unsigned int DoFHandler<dim>::n_boundary_dofs () const {
 };    
 
 
+
 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());
@@ -1189,12 +1206,47 @@ unsigned int DoFHandler<dim>::n_boundary_dofs (const FunctionMap &boundary_indic
 };    
 
 
+
 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)
index 883f07c8b27b3dd891c59511449059566bc31a74..bc5f9657426546579aecee433cd3fbda20dcf22e 100644 (file)
@@ -1489,7 +1489,7 @@ void DoFTools::map_dof_to_boundary_indices (const DoFHandler<1> &dof_handler,
 
 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());
@@ -1532,12 +1532,15 @@ void DoFTools::map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
 
 
 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 ();
@@ -1555,7 +1558,9 @@ void DoFTools::map_dof_to_boundary_indices (const DoFHandler<dim> &dof_handler,
   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);
@@ -1670,7 +1675,7 @@ DoFTools::map_dof_to_boundary_indices (const DoFHandler<deal_II_dimension> &,
 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 
index 537691b178dd7dbd5bf5b00fad0577e53ea94f18..c05789f9d9f16fa751ab97f5b4338224df613a4d 100644 (file)
@@ -719,7 +719,13 @@ VectorTools::project_boundary_values (const DoFHandler<dim>    &dof,
          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
index 0edd600ce1e2e52cb5efc38d3998d328193f1fcb..a4eee346d3dbcc7a4344a5a29823ef347ded79fe 100644 (file)
 <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>

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.