From: denis.davydov Date: Fri, 16 May 2014 14:54:57 +0000 (+0000) Subject: added new parent class for shared and distributed triangulations X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=353d182ab02f179b84059941aa1625e82240d70b;p=dealii-svn.git added new parent class for shared and distributed triangulations git-svn-id: https://svn.dealii.org/branches/branch_sharedtria@32925 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/include/deal.II/distributed/shared_tria.h b/deal.II/include/deal.II/distributed/shared_tria.h index 3168a9e271..5d432b0e46 100644 --- a/deal.II/include/deal.II/distributed/shared_tria.h +++ b/deal.II/include/deal.II/distributed/shared_tria.h @@ -44,6 +44,85 @@ template class Triangulation; namespace parallel { + template + class Triangulation : public dealii::Triangulation + { + public: + + /** + * Constructor. + */ + Triangulation (const typename dealii::Triangulation::MeshSmoothing smooth_grid = dealii::Triangulation::none, + const bool check_for_distorted_cells = false) : + dealii::Triangulation(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 & + 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 - class Triangulation : public dealii::Triangulation + class Triangulation : public dealii::parallel::Triangulation { public: typedef typename dealii::Triangulation::active_cell_iterator active_cell_iterator; diff --git a/deal.II/include/deal.II/distributed/tria.h b/deal.II/include/deal.II/distributed/tria.h index b00a01b224..22f070c7d2 100644 --- a/deal.II/include/deal.II/distributed/tria.h +++ b/deal.II/include/deal.II/distributed/tria.h @@ -27,6 +27,8 @@ #include #include +#include + #include #include #include @@ -306,7 +308,7 @@ namespace parallel * @ingroup distributed */ template - class Triangulation : public dealii::Triangulation + class Triangulation : public dealii::parallel::Triangulation { public: /** diff --git a/deal.II/source/distributed/shared_tria.cc b/deal.II/source/distributed/shared_tria.cc index 28991d5ab9..834e5cbab0 100644 --- a/deal.II/source/distributed/shared_tria.cc +++ b/deal.II/source/distributed/shared_tria.cc @@ -50,8 +50,8 @@ namespace parallel template Triangulation::Triangulation (MPI_Comm mpi_communicator, - const typename dealii::Triangulation::MeshSmoothing): - dealii::Triangulation(), + const typename dealii::Triangulation::MeshSmoothing smooth_grid): + dealii::parallel::Triangulation(smooth_grid,false), mpi_communicator (Utilities::MPI:: duplicate_communicator(mpi_communicator)), my_subdomain (Utilities::MPI::this_mpi_process (this->mpi_communicator)), diff --git a/deal.II/source/distributed/tria.cc b/deal.II/source/distributed/tria.cc index 6a22f877f5..3e69b6b329 100644 --- a/deal.II/source/distributed/tria.cc +++ b/deal.II/source/distributed/tria.cc @@ -2011,7 +2011,7 @@ namespace parallel const Settings settings_) : // do not check for distorted cells - dealii::Triangulation + dealii::parallel::Triangulation (smooth_grid, false), mpi_communicator (Utilities::MPI::