]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
added MeshSmoothing to constructor; Verbatim copy of NumberCache related functions...
authordenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 2 May 2014 08:33:59 +0000 (08:33 +0000)
committerdenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 2 May 2014 08:33:59 +0000 (08:33 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_sharedtria@32871 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/distributed/shared_tria.h
deal.II/source/distributed/shared_tria.cc

index 2ae1949971e1a4786ea9b10e52e2cd4d5ec883ec..f19b8e71c6e5b792365edf2e1659c31e0282b2d8 100644 (file)
@@ -41,6 +41,7 @@ DEAL_II_NAMESPACE_OPEN
 
 template <int, int> class Triangulation;
 
+
 namespace parallel
 {
   namespace shared
@@ -60,7 +61,9 @@ namespace parallel
       /**
        * Constructor.
        */
-      Triangulation (MPI_Comm mpi_communicator);
+      Triangulation (MPI_Comm mpi_communicator,
+                        const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing =
+                        typename Triangulation<dim,spacedim>::MeshSmoothing(Triangulation<dim,spacedim>::none));
 
       /**
        * Destructor.
@@ -98,11 +101,59 @@ namespace parallel
                                                                          const std::vector< CellData< dim > > &cells, 
                                                                          const SubCellData &subcelldata);
       
+      /**
+       * Return the number of active cells owned by each of the MPI
+       * processes that contribute to this triangulation. The element
+       * of this vector indexed by locally_owned_subdomain() equals
+       * the result of n_locally_owned_active_cells().
+       */
+      const std::vector<unsigned int> &
+      n_locally_owned_active_cells_per_processor () const;
+
+      /**
+       * Return the number of active cells in the triangulation that
+       * are locally owned, i.e. that have a subdomain_id equal to
+       * locally_owned_subdomain().
+       *
+       */
+      unsigned int n_locally_owned_active_cells () const;
+
+      /**
+       * Return the sum over all processors of the number of active
+       * cells owned by each processor. This equals the overall number
+       * of active cells in the shared triangulation.
+       */
+      types::global_dof_index n_global_active_cells () const;
+
+      /**
+       * Returns the global maximum level. This may be bigger than n_levels.
+       */
+      virtual unsigned int n_global_levels () const;
+
       
     private:
       MPI_Comm mpi_communicator;
       types::subdomain_id my_subdomain;
       types::subdomain_id num_subdomains;
+
+      struct NumberCache
+      {
+         std::vector<unsigned int> n_locally_owned_active_cells;
+         types::global_dof_index   n_global_active_cells;
+         unsigned int              n_global_levels;
+
+         NumberCache();
+      };
+
+      NumberCache number_cache;
+
+      /**
+       * Update the number_cache variable after mesh creation or
+       * refinement.
+       */
+      void update_number_cache ();
+
+
     };
   }
 }
index e01636453e3ad050c4561490e0583567424edc22..f28b69291a049f26dd5dd96afe117ef971589853 100644 (file)
@@ -40,15 +40,24 @@ namespace parallel
 {
   namespace shared
   {
+
+       template <int dim, int spacedim>
+       Triangulation<dim,spacedim>::NumberCache::NumberCache()
+         :
+         n_global_active_cells(0),
+         n_global_levels(0)
+        {}
+
     template <int dim, int spacedim>
-    Triangulation<dim,spacedim>::Triangulation (MPI_Comm mpi_communicator):
+    Triangulation<dim,spacedim>::Triangulation (MPI_Comm mpi_communicator,
+                                                   const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing):
         dealii::Triangulation<dim,spacedim>(),
        mpi_communicator (Utilities::MPI::
                           duplicate_communicator(mpi_communicator)),
         my_subdomain (Utilities::MPI::this_mpi_process (this->mpi_communicator)),
         num_subdomains(Utilities::MPI::n_mpi_processes(mpi_communicator))
     {
-
+       number_cache.n_locally_owned_active_cells.resize (num_subdomains);
     }
 
 
@@ -58,6 +67,72 @@ namespace parallel
 
     }
 
+    template <int dim, int spacedim>
+    unsigned int
+    Triangulation<dim,spacedim>::n_locally_owned_active_cells () const
+    {
+       return number_cache.n_locally_owned_active_cells[my_subdomain];
+    }
+
+    template <int dim, int spacedim>
+    unsigned int
+    Triangulation<dim,spacedim>::n_global_levels () const
+    {
+       return number_cache.n_global_levels;
+    }
+
+    template <int dim, int spacedim>
+    types::global_dof_index
+    Triangulation<dim,spacedim>::n_global_active_cells () const
+    {
+       return number_cache.n_global_active_cells;
+    }
+
+    template <int dim, int spacedim>
+    const std::vector<unsigned int> &
+    Triangulation<dim,spacedim>::n_locally_owned_active_cells_per_processor () const
+    {
+       return number_cache.n_locally_owned_active_cells;
+    }
+
+    template <int dim, int spacedim>
+    void
+    Triangulation<dim,spacedim>::update_number_cache ()
+    {
+       Assert (number_cache.n_locally_owned_active_cells.size()
+                       ==
+                               Utilities::MPI::n_mpi_processes (mpi_communicator),
+                               ExcInternalError());
+
+       std::fill (number_cache.n_locally_owned_active_cells.begin(),
+                       number_cache.n_locally_owned_active_cells.end(),
+                       0);
+
+       if (this->n_levels() > 0)
+               for (typename Triangulation<dim,spacedim>::active_cell_iterator
+                               cell = this->begin_active();
+                               cell != this->end(); ++cell)
+                       if (cell->subdomain_id() == my_subdomain)
+                               ++number_cache.n_locally_owned_active_cells[my_subdomain];
+
+       unsigned int send_value
+       = number_cache.n_locally_owned_active_cells[my_subdomain];
+       MPI_Allgather (&send_value,
+                       1,
+                       MPI_UNSIGNED,
+                       &number_cache.n_locally_owned_active_cells[0],
+                       1,
+                       MPI_UNSIGNED,
+                       mpi_communicator);
+
+       number_cache.n_global_active_cells
+       = std::accumulate (number_cache.n_locally_owned_active_cells.begin(),
+                       number_cache.n_locally_owned_active_cells.end(),
+                       /* ensure sum is computed with correct data type:*/
+                       static_cast<types::global_dof_index>(0));
+       number_cache.n_global_levels = Utilities::MPI::max(this->n_levels(), mpi_communicator);
+    }
+
 
 
     template <int dim, int spacedim>

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.