]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Only allocate memory for the direction flags if we really need them, i.e. for the...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 15 Dec 2010 15:44:51 +0000 (15:44 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 15 Dec 2010 15:44:51 +0000 (15:44 +0000)
git-svn-id: https://svn.dealii.org/trunk@22971 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/grid/tria_levels.h
deal.II/source/grid/tria.cc
deal.II/source/grid/tria_accessor.cc
deal.II/source/grid/tria_levels.cc

index 0977a4274b893b523948d56938003c6bd2ca435d..edb82f3d04f3a3636e29a67a502bd3195cfe0443 100644 (file)
@@ -145,16 +145,26 @@ namespace internal
                                           * index their parent has.
                                           */
         std::vector<int> parents;
-       
+
                                         /**
                                          * One bool per cell to indicate the
                                          * direction of the normal
                                          * true:  use orientation from vertex
-                                         * false: revert the orientation.
-                                         * It is used for codim==1 meshes
+                                         * false: revert the orientation. See
+                                         * @ref GlossDirectionFlag .
+                                         *
+                                         * This is only used for
+                                         * codim==1 meshes.
                                          */
         std::vector<bool> direction_flags;
-       
+
+                                        /**
+                                         * The object containing the data on lines and
+                                         * related functions
+                                         */
+       TriaObjects<TriaObject<dim> > cells;
+
+
                                          /**
                                           *  Reserve enough space to accomodate
                                           *  @p total_cells cells on this level.
@@ -172,7 +182,8 @@ namespace internal
                                           */
 
         void reserve_space (const unsigned int total_cells,
-                            const unsigned int dimension);
+                            const unsigned int dimension,
+                           const unsigned int space_dimension);
 
                                          /**
                                           *  Check the memory consistency of the
@@ -191,13 +202,6 @@ namespace internal
                                           */
         unsigned int memory_consumption () const;
 
-                                        /**
-                                         * The object containing the data on lines and
-                                         * related functions
-                                         */
-       TriaObjects<TriaObject<dim> > cells;
-
-
                                          /**
                                           *  Exception
                                           */
@@ -229,15 +233,24 @@ namespace internal
         std::vector<unsigned char> refine_flags;
         std::vector<bool> coarsen_flags;
         std::vector<std::pair<int,int> > neighbors;
-        std::vector<int> parents;
         std::vector<types::subdomain_id_t> subdomain_ids;
-        std::vector<bool> direction_flags; //not use; only needed to allow compilation
+        std::vector<int> parents;
+
+                                        // The following is not used
+                                        // since we don't support
+                                        // codim=1 meshes in 3d; only
+                                        // needed to allow
+                                        // compilation
+        std::vector<bool> direction_flags;
+
+       TriaObjectsHex cells;
+
+
         void reserve_space (const unsigned int total_cells,
-                            const unsigned int dimension);
+                            const unsigned int dimension,
+                           const unsigned int space_dimension);
         void monitor_memory (const unsigned int true_dimension) const;
         unsigned int memory_consumption () const;
-       TriaObjectsHex cells;
-
 
                                          /**
                                           *  Exception
index dcf008cf23b1fa00a87f00aae323c6447d653dff..5a43056e79961f96a570d097d431a6785a702aeb 100644 (file)
@@ -1433,7 +1433,7 @@ namespace internal
 
                                             // reserve enough space
            triangulation.levels.push_back (new internal::Triangulation::TriaLevel<dim>);
-           triangulation.levels[0]->reserve_space (cells.size(), dim);
+           triangulation.levels[0]->reserve_space (cells.size(), dim, spacedim);
            triangulation.levels[0]->cells.reserve_space (0,cells.size());
 
                                             // make up cells
@@ -1661,7 +1661,7 @@ namespace internal
                                             // reserve enough space
            triangulation.levels.push_back (new internal::Triangulation::TriaLevel<dim>);
            triangulation.faces = new internal::Triangulation::TriaFaces<dim>;
-           triangulation.levels[0]->reserve_space (cells.size(), dim);
+           triangulation.levels[0]->reserve_space (cells.size(), dim, spacedim);
            triangulation.faces->lines.reserve_space (0,needed_lines.size());
            triangulation.levels[0]->cells.reserve_space (0,cells.size());
 
@@ -1970,7 +1970,7 @@ namespace internal
                                             // reserve enough space
            triangulation.levels.push_back (new internal::Triangulation::TriaLevel<dim>);
            triangulation.faces = new internal::Triangulation::TriaFaces<dim>;
-           triangulation.levels[0]->reserve_space (cells.size(), dim);
+           triangulation.levels[0]->reserve_space (cells.size(), dim, spacedim);
            triangulation.faces->lines.reserve_space (0,needed_lines.size());
 
                                             // make up lines
@@ -4181,7 +4181,7 @@ namespace internal
                                                 // properties
                subcells[i]->set_material_id (cell->material_id());
                subcells[i]->set_subdomain_id (cell->subdomain_id());
-       
+
                if (i%2==0)
                  subcells[i]->set_parent (cell->index ());
              }
@@ -4203,8 +4203,8 @@ namespace internal
 
            if (dim < spacedim)
              for (unsigned int c=0; c<n_children; ++c)
-               cell->child(c)->set_direction_flag (cell->direction_flag());  
-         
+               cell->child(c)->set_direction_flag (cell->direction_flag());
+
          }
 
 
@@ -4280,7 +4280,8 @@ namespace internal
                  ->reserve_space(used_cells+
                                  GeometryInfo<1>::max_children_per_cell *
                                  flagged_cells,
-                                 1);
+                                 1,
+                                 spacedim);
                                                 // reserve space for
                                                 // 2*flagged_cells new lines on
                                                 // the next higher level
@@ -4647,7 +4648,8 @@ namespace internal
                                                 // level as well as for the
                                                 // needed_cells that will be
                                                 // created on that level
-               triangulation.levels[level+1]->reserve_space (used_cells+needed_cells, 2);
+               triangulation.levels[level+1]
+                 ->reserve_space (used_cells+needed_cells, 2, spacedim);
 
                                                 // reserve space for
                                                 // needed_cells
@@ -5117,7 +5119,8 @@ namespace internal
                                                 // level as well as for the
                                                 // 8*flagged_cells that will be
                                                 // created on that level
-               triangulation.levels[level+1]->reserve_space (used_cells+new_cells, 3);
+               triangulation.levels[level+1]
+                 ->reserve_space (used_cells+new_cells, 3, spacedim);
                                                 // reserve space for
                                                 // 8*flagged_cells
                                                 // new hexes on the next higher
@@ -9568,13 +9571,13 @@ create_triangulation (const std::vector<Point<spacedim> >    &v,
 
     We store this data using an n_faces x n_faces full matrix, which is actually
     much bigger than the minimal data required, but it makes the code more readable.
-    
+
   */
   if (dim < spacedim) {
 
     Table<2,bool> correct(GeometryInfo< dim >::faces_per_cell,
                          GeometryInfo< dim >::faces_per_cell);
-    switch (dim) 
+    switch (dim)
       {
        case 1:
        {
@@ -9599,8 +9602,8 @@ create_triangulation (const std::vector<Point<spacedim> >    &v,
        default:
              Assert (false, ExcNotImplemented());
       }
-  
-    
+
+
     std::list<active_cell_iterator> this_round, next_round;
     active_cell_iterator neighbor;
 
@@ -9608,54 +9611,56 @@ create_triangulation (const std::vector<Point<spacedim> >    &v,
     begin_active()->set_direction_flag (true);
     begin_active()->set_user_flag ();
 
-    while (this_round.size() > 0) {
-      
-      for ( typename std::list<active_cell_iterator>::iterator cell = this_round.begin();
-           cell != this_round.end(); ++cell) {
-       
-       for (unsigned int i = 0; i < GeometryInfo< dim >::faces_per_cell; ++i) 
+    while (this_round.size() > 0)
+      {
+       for ( typename std::list<active_cell_iterator>::iterator cell = this_round.begin();
+             cell != this_round.end(); ++cell)
          {
-           if ( !((*cell)->face(i)->at_boundary()) )
+           for (unsigned int i = 0; i < GeometryInfo< dim >::faces_per_cell; ++i)
              {
-               neighbor = (*cell)->neighbor(i);
-                               
-               unsigned int cf = (*cell)->face_index(i);
-               unsigned int j = 0;
-               while(neighbor->face_index(j) != cf)
-                 {++j;}
-                
-               if ( (correct(i,j) && !(*cell)->direction_flag())
-                    ||
-                    (!correct(i,j) && (*cell)->direction_flag()) )
+               if ( !((*cell)->face(i)->at_boundary()) )
                  {
-                   if (neighbor->user_flag_set() == false)
+                   neighbor = (*cell)->neighbor(i);
+
+                   unsigned int cf = (*cell)->face_index(i);
+                   unsigned int j = 0;
+                   while(neighbor->face_index(j) != cf)
+                     {++j;}
+
+                   if ( (correct(i,j) && !(*cell)->direction_flag())
+                        ||
+                        (!correct(i,j) && (*cell)->direction_flag()) )
                      {
-                       neighbor->set_direction_flag (false);
-                       neighbor->set_user_flag ();
-                       next_round.push_back (neighbor);
+                       if (neighbor->user_flag_set() == false)
+                         {
+                           neighbor->set_direction_flag (false);
+                           neighbor->set_user_flag ();
+                           next_round.push_back (neighbor);
+                         }
+                       else
+                         Assert (neighbor->direction_flag() == false,
+                                 ExcNonOrientableTriangulation());
+
                      }
-                   else 
-                     Assert (neighbor->direction_flag() == false,
-                             ExcNonOrientableTriangulation());
-                   
                  }
              }
-           
          }
-      }
-                
-//Before we quit let's check that if the triangulation is disconnected
-      if (next_round.size() == 0) {
+
+                                        // Before we quit let's check
+                                        // that if the triangulation
+                                        // is disconnected that we
+                                        // still get all cells
+      if (next_round.size() == 0)
        for (active_cell_iterator cell = begin_active();
             cell != end(); ++cell)
-         if (cell->user_flag_set() == false) {
-           next_round.push_back (cell);
-           cell->set_direction_flag (true);
-           cell->set_user_flag ();
-           break;
-         }
-      }
-      
+         if (cell->user_flag_set() == false)
+           {
+             next_round.push_back (cell);
+             cell->set_direction_flag (true);
+             cell->set_user_flag ();
+             break;
+           }
+
       this_round = next_round;
       next_round.clear();
     }
index b5eacc9375f8ada7a560810f35359c767a76e141..7be71bc7ee108247ff354c4f74ec00f1940dfb0e 100644 (file)
@@ -1265,14 +1265,11 @@ CellAccessor<dim, spacedim>::set_subdomain_id (const types::subdomain_id_t new_s
 template <int dim, int spacedim>
 bool CellAccessor<dim, spacedim>::direction_flag () const
 {
+  Assert (this->used(), TriaAccessorExceptions::ExcCellNotUsed());
   if (dim==spacedim)
     return true;
   else
-    {
-      Assert (this->used(), TriaAccessorExceptions::ExcCellNotUsed());
-      return this->tria->levels[this->present_level]->direction_flags[this->present_index];
-    }
-  
+    return this->tria->levels[this->present_level]->direction_flags[this->present_index];
 }
 
 
@@ -1281,12 +1278,14 @@ template <int dim, int spacedim>
 void
 CellAccessor<dim, spacedim>::set_direction_flag (const bool new_direction_flag) const
 {
+  Assert (this->used(), TriaAccessorExceptions::ExcCellNotUsed());
   if (dim<spacedim)
-    {
-      Assert (this->used(), TriaAccessorExceptions::ExcCellNotUsed());
-      this->tria->levels[this->present_level]->direction_flags[this->present_index]
-       = new_direction_flag;
-    }
+    this->tria->levels[this->present_level]->direction_flags[this->present_index]
+      = new_direction_flag;
+  else
+    Assert (new_direction_flag == true,
+           ExcMessage ("If dim==spacedim, direction flags are always true and "
+                       "can not be set to anything else."));
 }
 
 
index 7718c84f23b997775513d80cc5b71196931befef..d110365599c6a34d8a4d883ed0661d92741bd22f 100644 (file)
@@ -25,7 +25,8 @@ namespace internal
     template <int dim>
     void
     TriaLevel<dim>::reserve_space (const unsigned int total_cells,
-                                  const unsigned int dimension)
+                                  const unsigned int dimension,
+                                  const unsigned int space_dimension)
     {
                                        // we need space for total_cells
                                        // cells. Maybe we have more already
@@ -51,10 +52,15 @@ namespace internal
                                 total_cells - subdomain_ids.size(),
                                 0);
 
-         direction_flags.reserve (total_cells);
-          direction_flags.insert (direction_flags.end(),
-                                 total_cells - direction_flags.size(),
-                                 true);
+         if (dimension < space_dimension)
+           {
+             direction_flags.reserve (total_cells);
+             direction_flags.insert (direction_flags.end(),
+                                     total_cells - direction_flags.size(),
+                                     true);
+           }
+         else
+           direction_flags.clear ();
 
           parents.reserve ((int) (total_cells + 1) / 2);
           parents.insert (parents.end (),
@@ -112,10 +118,13 @@ namespace internal
     unsigned int
     TriaLevel<dim>::memory_consumption () const
     {
-      return (MemoryConsumption::memory_consumption (coarsen_flags) +
+      return (MemoryConsumption::memory_consumption (refine_flags) +
+             MemoryConsumption::memory_consumption (coarsen_flags) +
               MemoryConsumption::memory_consumption (neighbors) +
-              MemoryConsumption::memory_consumption (cells) +
-              MemoryConsumption::memory_consumption (refine_flags));
+             MemoryConsumption::memory_consumption (subdomain_ids) +
+             MemoryConsumption::memory_consumption (parents) +
+             MemoryConsumption::memory_consumption (direction_flags) +
+              MemoryConsumption::memory_consumption (cells));
     }
 
 // This specialization should be only temporary, until the TriaObjects
@@ -123,7 +132,8 @@ namespace internal
 
     void
     TriaLevel<3>::reserve_space (const unsigned int total_cells,
-                                const unsigned int dimension)
+                                const unsigned int dimension,
+                                const unsigned int space_dimension)
     {
                                        // we need space for total_cells
                                        // cells. Maybe we have more already
@@ -149,6 +159,16 @@ namespace internal
                                 total_cells - subdomain_ids.size(),
                                 0);
 
+         if (dimension < space_dimension)
+           {
+             direction_flags.reserve (total_cells);
+             direction_flags.insert (direction_flags.end(),
+                                     total_cells - direction_flags.size(),
+                                     true);
+           }
+         else
+           direction_flags.clear ();
+
           parents.reserve ((int) (total_cells + 1) / 2);
           parents.insert (parents.end (),
                           (total_cells + 1) / 2 - parents.size (),
@@ -188,6 +208,10 @@ namespace internal
              subdomain_ids.capacity() + DEAL_II_MIN_VECTOR_CAPACITY,
               ExcMemoryWasted ("subdomain_ids",
                                subdomain_ids.size(), subdomain_ids.capacity()));
+      Assert (direction_flags.size() <=
+             direction_flags.capacity() + DEAL_II_MIN_VECTOR_CAPACITY,
+              ExcMemoryWasted ("direction_flags",
+                               direction_flags.size(), direction_flags.capacity()));
       Assert (2*true_dimension*refine_flags.size() == neighbors.size(),
               ExcMemoryInexact (refine_flags.size(), neighbors.size()));
       Assert (2*true_dimension*coarsen_flags.size() == neighbors.size(),
@@ -198,10 +222,13 @@ namespace internal
     unsigned int
     TriaLevel<3>::memory_consumption () const
     {
-      return (MemoryConsumption::memory_consumption (coarsen_flags) +
+      return (MemoryConsumption::memory_consumption (refine_flags) +
+             MemoryConsumption::memory_consumption (coarsen_flags) +
               MemoryConsumption::memory_consumption (neighbors) +
-              MemoryConsumption::memory_consumption (cells) +
-              MemoryConsumption::memory_consumption (refine_flags));
+             MemoryConsumption::memory_consumption (subdomain_ids) +
+             MemoryConsumption::memory_consumption (parents) +
+             MemoryConsumption::memory_consumption (direction_flags) +
+              MemoryConsumption::memory_consumption (cells));
     }
   }
 }

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.