From: Wolfgang Bangerth Date: Wed, 4 Dec 2013 14:38:16 +0000 (+0000) Subject: Instantiate several functions in GridTools for parallel::distributed::Triangulation. X-Git-Tag: v8.1.0~101 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=d2aa7bc6a6662e95a7f4ce67f06b6b8407e9eb04;p=dealii.git Instantiate several functions in GridTools for parallel::distributed::Triangulation. git-svn-id: https://svn.dealii.org/trunk@31876 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 18464c88da..93c0ecd3a3 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -246,6 +246,12 @@ inconvenience this causes.

Specific improvements

    +
  1. Fixed: Several functions in namespace GridTools were not instantiated + for parallel::distributed::Triangulation objects. This is now fixed. +
    + (Denis Davydov, Wolfgang Bangerth, 2013/12/01) +
  2. +
  3. 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 diff --git a/deal.II/source/grid/grid_tools.cc b/deal.II/source/grid/grid_tools.cc index 1316387b53..fef813ec37 100644 --- a/deal.II/source/grid/grid_tools.cc +++ b/deal.II/source/grid/grid_tools.cc @@ -59,6 +59,13 @@ namespace GridTools return tria; } + template + const Triangulation & + get_tria(const parallel::distributed::Triangulation &tria) + { + return tria; + } + template class Container, int spacedim> const Triangulation & get_tria(const Container &container) @@ -74,6 +81,13 @@ namespace GridTools return tria; } + template + Triangulation & + get_tria(parallel::distributed::Triangulation &tria) + { + return tria; + } + template class Container, int spacedim> const Triangulation & get_tria(Container &container) @@ -862,8 +876,10 @@ next_cell: find_active_cell_around_point (const Container &container, const Point &p) { - return find_active_cell_around_point(StaticMappingQ1::mapping, - container, p).first; + return + find_active_cell_around_point + (StaticMappingQ1::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)); } diff --git a/deal.II/source/grid/grid_tools.inst.in b/deal.II/source/grid/grid_tools.inst.in index 3d1076e57a..6e827d7111 100644 --- a/deal.II/source/grid/grid_tools.inst.in +++ b/deal.II/source/grid/grid_tools.inst.in @@ -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 &, + const Point &); + + template + std::vector::active_cell_iterator> + find_cells_adjacent_to_vertex(const parallel::distributed::Triangulation &, + const unsigned int); + template + parallel::distributed::Triangulation::active_cell_iterator + find_active_cell_around_point (const parallel::distributed::Triangulation &, + const Point &p); + + template + std::pair::active_cell_iterator, Point > + find_active_cell_around_point (const Mapping &, + const parallel::distributed::Triangulation &, + const Point &); + + template + std::list::cell_iterator, parallel::distributed::Triangulation::cell_iterator> > + get_finest_common_cells (const parallel::distributed::Triangulation &mesh_1, + const parallel::distributed::Triangulation &mesh_2); + + + template + bool + have_same_coarse_mesh (const parallel::distributed::Triangulation &mesh_1, + const parallel::distributed::Triangulation &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 index 0000000000..aca50e93c9 --- /dev/null +++ b/tests/mpi/find_active_cell_around_point_01.cc @@ -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 +#include +#include +#include +#include +#include +#include +#include +#include + + +#include + + +template +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 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 p; + for (unsigned int d=0; d::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 index 0000000000..c36f6b706d --- /dev/null +++ b/tests/mpi/find_active_cell_around_point_01.mpirun=4.output @@ -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