]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Instantiate several functions in GridTools for parallel::distributed::Triangulation.
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 4 Dec 2013 14:38:16 +0000 (14:38 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 4 Dec 2013 14:38:16 +0000 (14:38 +0000)
git-svn-id: https://svn.dealii.org/trunk@31876 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/source/grid/grid_tools.cc
deal.II/source/grid/grid_tools.inst.in
tests/mpi/find_active_cell_around_point_01.cc [new file with mode: 0644]
tests/mpi/find_active_cell_around_point_01.mpirun=4.output [new file with mode: 0644]

index 18464c88da72ad215b951e95bdb718431ca71300..93c0ecd3a338a3fa72822a8739b838b55eb5117d 100644 (file)
@@ -246,6 +246,12 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
+  <li> Fixed: Several functions in namespace GridTools were not instantiated
+  for parallel::distributed::Triangulation objects. This is now fixed.
+  <br>
+  (Denis Davydov, Wolfgang Bangerth, 2013/12/01)
+  </li>
+
   <li> Improved: The methods ConstraintMatrix::distribute_local_to_global
   now use scratch data that is private to each thread instead of allocating
   it for every cell anew. This gives better performance, in particular in
index 1316387b538c2a8c1f4282ff15d26c9014afe3ac..fef813ec37bfb448df673a9fb7eff95ac0b84ffc 100644 (file)
@@ -59,6 +59,13 @@ namespace GridTools
       return tria;
     }
 
+    template<int dim, int spacedim>
+    const Triangulation<dim, spacedim> &
+    get_tria(const parallel::distributed::Triangulation<dim, spacedim> &tria)
+    {
+      return tria;
+    }
+
     template<int dim, template<int, int> class Container, int spacedim>
     const Triangulation<dim,spacedim> &
     get_tria(const Container<dim,spacedim> &container)
@@ -74,6 +81,13 @@ namespace GridTools
       return tria;
     }
 
+    template<int dim, int spacedim>
+    Triangulation<dim, spacedim> &
+    get_tria(parallel::distributed::Triangulation<dim, spacedim> &tria)
+    {
+      return tria;
+    }
+
     template<int dim, template<int, int> class Container, int spacedim>
     const Triangulation<dim,spacedim> &
     get_tria(Container<dim,spacedim> &container)
@@ -862,8 +876,10 @@ next_cell:
   find_active_cell_around_point (const Container<dim,spacedim>  &container,
                                  const Point<spacedim> &p)
   {
-    return find_active_cell_around_point(StaticMappingQ1<dim,spacedim>::mapping,
-                                         container, p).first;
+    return
+      find_active_cell_around_point<dim,Container,spacedim>
+      (StaticMappingQ1<dim,spacedim>::mapping,
+       container, p).first;
   }
 
 
@@ -1417,8 +1433,8 @@ next_cell:
   have_same_coarse_mesh (const Container &mesh_1,
                          const Container &mesh_2)
   {
-    return have_same_coarse_mesh (mesh_1.get_tria(),
-                                  mesh_2.get_tria());
+    return have_same_coarse_mesh (get_tria(mesh_1),
+                                  get_tria(mesh_2));
   }
 
 
index 3d1076e57ad0886f8c79aad5dff171c146860918..6e827d7111ddf161d3f53ac5070eb23137409801 100644 (file)
@@ -69,6 +69,51 @@ for (X : TRIANGULATION_AND_DOFHANDLERS; deal_II_dimension : DIMENSIONS ; deal_II
 }
 
 
+// now also instantiate these functions for parallel::distributed::Triangulation
+for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS)
+{
+
+#if deal_II_dimension <= deal_II_space_dimension
+  namespace GridTools \{
+
+  template
+    unsigned int
+    find_closest_vertex (const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
+                        const Point<deal_II_space_dimension> &);
+
+  template
+    std::vector<parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator>
+    find_cells_adjacent_to_vertex(const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
+                                 const unsigned int);
+  template
+    parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator
+    find_active_cell_around_point (const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
+                                  const Point<deal_II_space_dimension> &p);
+
+  template
+    std::pair<parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::active_cell_iterator, Point<deal_II_dimension> >
+    find_active_cell_around_point (const Mapping<deal_II_dimension, deal_II_space_dimension> &,
+                                  const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &,
+                                  const Point<deal_II_space_dimension> &);
+
+  template
+    std::list<std::pair<parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::cell_iterator, parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension>::cell_iterator> >
+    get_finest_common_cells (const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &mesh_1,
+                            const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &mesh_2);
+
+
+  template
+    bool
+    have_same_coarse_mesh (const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &mesh_1,
+                          const parallel::distributed::Triangulation<deal_II_dimension,deal_II_space_dimension> &mesh_2);
+
+  \}
+
+  #endif
+}
+
+
+
 for (deal_II_space_dimension : SPACE_DIMENSIONS)
 {
   
diff --git a/tests/mpi/find_active_cell_around_point_01.cc b/tests/mpi/find_active_cell_around_point_01.cc
new file mode 100644 (file)
index 0000000..aca50e9
--- /dev/null
@@ -0,0 +1,110 @@
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 2009 - 2013 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// make sure only one processor finds a locally-owned cell around a point
+
+#include "../tests.h"
+#include "coarse_grid_common.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/tensor.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/distributed/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/base/utilities.h>
+
+
+#include <fstream>
+
+
+template<int dim>
+void test()
+{
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+
+  if (true)
+    {
+      if (Utilities::MPI::this_mpi_process (MPI_COMM_WORLD) == 0)
+        deallog << "hyper_cube" << std::endl;
+
+      parallel::distributed::Triangulation<dim> tr(MPI_COMM_WORLD);
+
+      GridGenerator::hyper_cube(tr);
+      tr.refine_global(2);
+
+      // choose a point that is guaranteed to lie in the domain but not
+      // at the interface between cells
+      Point<dim> p;
+      for (unsigned int d=0; d<dim; ++d)
+       p[d] = 1./3;
+      
+      typename parallel::distributed::Triangulation<dim>::active_cell_iterator
+       cell = GridTools::find_active_cell_around_point (tr, p);
+
+      const unsigned int
+       n_locally_owned
+       = Utilities::MPI::sum (cell->is_locally_owned() ? 1 : 0,
+                              MPI_COMM_WORLD);
+      
+      const unsigned int
+       n_locally_owned_or_ghost
+       = Utilities::MPI::sum (!cell->is_artificial() ? 1 : 0,
+                              MPI_COMM_WORLD);
+      
+      if (myid == 0)
+        deallog << "Locally owned: "
+                << n_locally_owned
+                << std::endl
+               << "Locally owned or ghost: "
+                << n_locally_owned_or_ghost
+                << std::endl;
+    }
+}
+
+
+int main(int argc, char *argv[])
+{
+  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
+
+  unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
+
+
+  deallog.push(Utilities::int_to_string(myid));
+
+  if (myid == 0)
+    {
+      std::ofstream logfile("output");
+      deallog.attach(logfile);
+      deallog.depth_console(0);
+      deallog.threshold_double(1.e-10);
+
+      deallog.push("2d");
+      test<2>();
+      deallog.pop();
+      deallog.push("3d");
+      test<3>();
+      deallog.pop();
+    }
+  else
+    {
+      test<2>();
+      test<3>();
+    }
+}
diff --git a/tests/mpi/find_active_cell_around_point_01.mpirun=4.output b/tests/mpi/find_active_cell_around_point_01.mpirun=4.output
new file mode 100644 (file)
index 0000000..c36f6b7
--- /dev/null
@@ -0,0 +1,7 @@
+JobId linux-sh2g.site Wed Dec  4 08:36:59 2013
+DEAL:0:2d::hyper_cube
+DEAL:0:2d::Locally owned: 1
+DEAL:0:2d::Locally owned or ghost: 4
+DEAL:0:3d::hyper_cube
+DEAL:0:3d::Locally owned: 1
+DEAL:0:3d::Locally owned or ghost: 4

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.