]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Cleanup interface of parallel::Triangulation without MPI
authorMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sun, 30 Aug 2015 15:17:35 +0000 (17:17 +0200)
committerMartin Kronbichler <kronbichler@lnm.mw.tum.de>
Sun, 30 Aug 2015 16:03:23 +0000 (18:03 +0200)
include/deal.II/distributed/tria.h
include/deal.II/distributed/tria_base.h
source/distributed/tria.cc
source/distributed/tria_base.cc

index 7b24492e16bccee0d5e9745ea0e3e905617f7c6c..00c22c20c61a25a2f01265b55bc9aa587625d5cb 100644 (file)
@@ -1079,7 +1079,7 @@ namespace parallel
      * actually be created if we don't have p4est available.
      */
     template <int dim, int spacedim = dim>
-    class Triangulation : public dealii::Triangulation<dim,spacedim>
+    class Triangulation : public dealii::parallel::Triangulation<dim,spacedim>
     {
     private:
       /**
@@ -1087,21 +1087,6 @@ namespace parallel
        */
       Triangulation ();
 
-    public:
-
-      /**
-       * Destructor.
-       */
-      virtual ~Triangulation ();
-
-      /**
-       * Return the subdomain id of those cells that are owned by the current
-       * processor. All cells in the triangulation that do not have this
-       * subdomain id are either owned by another processor or have children
-       * that only exist on other processors.
-       */
-      types::subdomain_id locally_owned_subdomain () const;
-
     };
   }
 }
index b5f1ade171d90dd9c7cfb1f4e1c3a3009f01c94c..6ef81639cdf18f35c2d4a835c3961696d6c729c2 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2008 - 2013 by the deal.II authors
+// Copyright (C) 2008 - 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -56,25 +56,19 @@ namespace parallel
     /**
      * Constructor.
      */
-#ifdef DEAL_II_WITH_MPI
     Triangulation (MPI_Comm mpi_communicator,
                    const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid = (dealii::Triangulation<dim,spacedim>::none),
                    const bool check_for_distorted_cells = false);
-#else
-    Triangulation ();
-#endif
 
     /**
      * Destructor.
      */
     virtual ~Triangulation ();
 
-#ifdef DEAL_II_WITH_MPI
     /**
      * Return MPI communicator used by this triangulation.
      */
     virtual MPI_Comm get_communicator () const;
-#endif
 
     /**
      * Implementation of the same function as in the base class.
@@ -143,14 +137,12 @@ namespace parallel
 
 
   protected:
-#ifdef DEAL_II_WITH_MPI
     /**
      * MPI communicator to be used for the triangulation. We create a unique
      * communicator for this class, which is a duplicate of the one passed
      * to the constructor.
      */
     MPI_Comm mpi_communicator;
-#endif
 
     /**
      * The subdomain id to be used for the current processor.
index c25994341066e6070ae7ce6598ed735590a1352e..9b1687724918b1087b917ca956ee268ce8a6a490 100644 (file)
@@ -4590,37 +4590,11 @@ namespace parallel
   {
     template <int dim, int spacedim>
     Triangulation<dim,spacedim>::Triangulation ()
+      :
+      dealii::parallel::Triangulation<dim,spacedim>(MPI_COMM_SELF)
     {
       Assert (false, ExcNotImplemented());
     }
-
-
-    template <int dim, int spacedim>
-    Triangulation<dim,spacedim>::~Triangulation ()
-    {
-      Assert (false, ExcNotImplemented());
-    }
-
-
-
-    template <int dim, int spacedim>
-    types::subdomain_id
-    Triangulation<dim,spacedim>::locally_owned_subdomain () const
-    {
-      Assert (false, ExcNotImplemented());
-      return 0;
-    }
-
-
-#ifdef DEAL_II_WITH_MPI
-    template <int dim, int spacedim>
-    MPI_Comm
-    Triangulation<dim,spacedim>::get_communicator () const
-    {
-      Assert (false, ExcNotImplemented());
-      return MPI_COMM_WORLD;
-    }
-#endif
   }
 }
 
index d1b613afe01f6a26595034597d726b189e3c8304..ff948468ebcb8bafeb8fac0d85a782c42d476568 100644 (file)
@@ -37,7 +37,6 @@ DEAL_II_NAMESPACE_OPEN
 namespace parallel
 {
 
-#ifdef DEAL_II_WITH_MPI
   template <int dim, int spacedim>
   Triangulation<dim,spacedim>::Triangulation (MPI_Comm mpi_communicator,
                                               const typename dealii::Triangulation<dim,spacedim>::MeshSmoothing smooth_grid,
@@ -49,6 +48,10 @@ namespace parallel
     my_subdomain (Utilities::MPI::this_mpi_process (this->mpi_communicator)),
     n_subdomains(Utilities::MPI::n_mpi_processes(mpi_communicator))
   {
+#ifndef DEAL_II_WITH_MPI
+    Assert(false, ExcMessage("You compiled deal.II without MPI support, for "
+                             "which parallel::Triangulation is not available."));
+#endif
     number_cache.n_locally_owned_active_cells.resize (n_subdomains);
   }
 
@@ -56,6 +59,9 @@ namespace parallel
   void
   Triangulation<dim,spacedim>::copy_triangulation (const dealii::Triangulation<dim, spacedim> &old_tria)
   {
+#ifndef DEAL_II_WITH_MPI
+    Assert(false, ExcNotImplemented());
+#endif
     if (const dealii::parallel::Triangulation<dim,spacedim> *
         old_tria_x = dynamic_cast<const dealii::parallel::Triangulation<dim,spacedim> *>(&old_tria))
       {
@@ -63,22 +69,7 @@ namespace parallel
       }
   }
 
-#else
-  template <int dim, int spacedim>
-  Triangulation<dim,spacedim>::Triangulation()
-  {
-    Assert (false, ExcNotImplemented());
-  }
-
-  template <int dim, int spacedim>
-  void
-  Triangulation<dim,spacedim>::copy_triangulation (const dealii::Triangulation<dim, spacedim> &)
-  {
-    Assert (false, ExcNotImplemented());
-  }
-
 
-#endif
 
   template <int dim, int spacedim>
   std::size_t
@@ -86,9 +77,7 @@ namespace parallel
   {
     std::size_t mem=
       this->dealii::Triangulation<dim,spacedim>::memory_consumption()
-#ifdef DEAL_II_WITH_MPI
       + MemoryConsumption::memory_consumption(mpi_communicator)
-#endif
       + MemoryConsumption::memory_consumption(my_subdomain)
       + MemoryConsumption::memory_consumption(number_cache.n_locally_owned_active_cells)
       + MemoryConsumption::memory_consumption(number_cache.n_global_active_cells)
@@ -141,7 +130,6 @@ namespace parallel
     return number_cache.n_locally_owned_active_cells;
   }
 
-#ifdef DEAL_II_WITH_MPI
   template <int dim, int spacedim>
   MPI_Comm
   Triangulation<dim,spacedim>::get_communicator () const
@@ -149,6 +137,7 @@ namespace parallel
     return mpi_communicator;
   }
 
+#ifdef DEAL_II_WITH_MPI
   template <int dim, int spacedim>
   void
   Triangulation<dim,spacedim>::update_number_cache ()

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.