]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
added new parent class for shared and distributed triangulations
authordenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 16 May 2014 14:54:57 +0000 (14:54 +0000)
committerdenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 16 May 2014 14:54:57 +0000 (14:54 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_sharedtria@32925 0785d39b-7218-0410-832d-ea1e28bc413d

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

index 3168a9e271897eb1b3d325134414d2c9c039eaf6..5d432b0e46dd71056d83d3f015477a3c2687fc3e 100644 (file)
@@ -44,6 +44,85 @@ template <int, int> class Triangulation;
 
 namespace parallel
 {
+    template <int dim, int spacedim = dim>
+    class Triangulation : public dealii::Triangulation<dim,spacedim>
+    {
+    public:
+
+      /**
+       * Constructor.
+       */
+      Triangulation (const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid = dealii::Triangulation<dim,spacedim>::none,
+                        const bool check_for_distorted_cells = false) :
+                 dealii::Triangulation<dim,spacedim>(smooth_grid,check_for_distorted_cells) {};
+
+      /**
+       * Destructor.
+       */
+      virtual ~Triangulation () {};
+
+      /**
+       * Return locally owned subdomain id,
+       * which is equivalent to the rank of the current mpi process in
+       * communicator provided to class constructor.
+       */
+      virtual types::subdomain_id locally_owned_subdomain () const = 0;
+
+      /**
+       * Return MPI communicator used by this triangulation.
+       */
+      virtual MPI_Comm get_communicator () const = 0;
+
+      /**
+       * Coarsen and refine the mesh according to refinement and
+       * coarsening flags set.
+       *
+       * This step is equivalent to the dealii::Triangulation class
+       * with an addition of calling dealii::GridTools::partition_triangulation() at the end.
+       */
+      virtual void execute_coarsening_and_refinement () = 0;
+
+      /**
+        * Create a triangulation.
+        *
+        * This function also partitions triangulation based on the
+        * MPI communicator provided to constructor.
+        */
+      virtual void     create_triangulation (const std::vector< Point< spacedim > > &vertices,
+                                                                         const std::vector< CellData< dim > > &cells,
+                                                                         const SubCellData &subcelldata) = 0;
+
+      /**
+       * 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().
+       */
+      virtual const std::vector<unsigned int> &
+      n_locally_owned_active_cells_per_processor () const = 0;
+
+      /**
+       * 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().
+       *
+       */
+      virtual unsigned int n_locally_owned_active_cells () const = 0;
+
+      /**
+       * 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.
+       */
+      virtual types::global_dof_index n_global_active_cells () const = 0;
+
+      /**
+       * Returns the global maximum level. This may be bigger than n_levels.
+       */
+      virtual unsigned int n_global_levels () const = 0;
+    };
+
+
   namespace shared
   {
 
@@ -52,7 +131,7 @@ namespace parallel
     * partition trianguation using MPI.
     */
     template <int dim, int spacedim = dim>
-    class Triangulation : public dealii::Triangulation<dim,spacedim>
+    class Triangulation : public dealii::parallel::Triangulation<dim,spacedim>
     {
     public:
       typedef typename dealii::Triangulation<dim,spacedim>::active_cell_iterator active_cell_iterator;
index b00a01b224483e5e2cbe0825ded255754a79cc50..22f070c7d2cd3571f3d1da937e17394f772563f7 100644 (file)
@@ -27,6 +27,8 @@
 #include <deal.II/base/std_cxx1x/function.h>
 #include <deal.II/base/std_cxx1x/tuple.h>
 
+#include <deal.II/distributed/shared_tria.h>
+
 #include <set>
 #include <vector>
 #include <list>
@@ -306,7 +308,7 @@ namespace parallel
      * @ingroup distributed
      */
     template <int dim, int spacedim = dim>
-    class Triangulation : public dealii::Triangulation<dim,spacedim>
+    class Triangulation : public dealii::parallel::Triangulation<dim,spacedim>
     {
     public:
       /**
index 28991d5ab96fa78c7db2f636ea4206a604bbb0b6..834e5cbab0faa13e135ec01b3626275ace8cc92c 100644 (file)
@@ -50,8 +50,8 @@ namespace parallel
 
     template <int dim, int spacedim>
     Triangulation<dim,spacedim>::Triangulation (MPI_Comm mpi_communicator,
-                                                   const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing):
-        dealii::Triangulation<dim,spacedim>(),
+                                                   const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid):
+        dealii::parallel::Triangulation<dim,spacedim>(smooth_grid,false),
        mpi_communicator (Utilities::MPI::
                           duplicate_communicator(mpi_communicator)),
         my_subdomain (Utilities::MPI::this_mpi_process (this->mpi_communicator)),
index 6a22f877f596ce437f6191281bd9a6e42a3ef01e..3e69b6b329e1d45e8b02569216dd222bdc4e11e4 100644 (file)
@@ -2011,7 +2011,7 @@ namespace parallel
                    const Settings settings_)
       :
       // do not check for distorted cells
-      dealii::Triangulation<dim,spacedim>
+      dealii::parallel::Triangulation<dim,spacedim>
       (smooth_grid,
        false),
       mpi_communicator (Utilities::MPI::

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.