From: Luca Heltai Date: Mon, 18 Sep 2017 15:35:17 +0000 (+0200) Subject: New GridTools::extract_used_vertices function. X-Git-Tag: v9.0.0-rc1~1018^2~4 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=e4aff50dc847bd19ff3838e419318b3b93ecaf6b;p=dealii.git New GridTools::extract_used_vertices function. --- diff --git a/doc/news/changes/minor/20170918LucaHeltai b/doc/news/changes/minor/20170918LucaHeltai new file mode 100644 index 0000000000..b3d7bec05a --- /dev/null +++ b/doc/news/changes/minor/20170918LucaHeltai @@ -0,0 +1,6 @@ +New: A GridTools::extract_used_vertices() has been added +that returns a map of index-vertex position, using an +optional Mapping object to compute the +vertex locations. +
+(Luca Heltai, 2017/09/18) diff --git a/include/deal.II/grid/grid_tools.h b/include/deal.II/grid/grid_tools.h index c03c8eb75b..4585e16b59 100644 --- a/include/deal.II/grid/grid_tools.h +++ b/include/deal.II/grid/grid_tools.h @@ -560,6 +560,40 @@ namespace GridTools */ /*@{*/ + /** + * Return a map of index:Point, containing the used vertices of the + * given `container`. The key of the returned map is the global index in the + * triangulation. The used vertices are obtained by looping over all cells, + * and querying for each cell where its vertices are through the (optional) + * `mapping` argument. + * + * The size of the returned map equals Triangulation::n_used_vertices(), + * (not Triangulation::n_vertices()). If you use the default `mapping`, the + * returned map satisfies the following equality: + * + * @code + * used_vertices = extract_used_vertices(tria); + * all_vertices = tria.get_vertices(); + * + * for(auto &id_and_v : used_vertices) + * all_vertices[id_and_v.first] == id_and_v.second; // true + * @endcode + * + * Notice that the above is not satisfied for mappings that change the + * location of vertices, like MappingQEulerian. + * + * @ref ConceptMeshType "MeshType concept". + * @param container The container to extract vertices from + * @param mapping The mapping to use to compute the points locations + * + * @author Luca Heltai, 2017. + */ + template + std::map> extract_used_vertices ( + const Triangulation &container, + const Mapping &mapping = StaticMappingQ1::mapping); + + /** * Find and return the number of the used vertex (or marked vertex) in a * given mesh that is located closest to a given point. diff --git a/source/grid/grid_tools.cc b/source/grid/grid_tools.cc index e492b2fef4..1229339a78 100644 --- a/source/grid/grid_tools.cc +++ b/source/grid/grid_tools.cc @@ -4926,6 +4926,23 @@ next_cell: + template + std::map > + extract_used_vertices(const Triangulation &container, + const Mapping &mapping) + { + std::map > result; + for (const auto &cell : container.active_cell_iterators()) + { + auto vs = mapping.get_vertices(cell); + for (unsigned int i=0; ivertex_index(i)]=vs[i]; + } + Assert(result.size() == container.n_used_vertices(), + ExcInternalError()); + return result; + } + } /* namespace GridTools */ diff --git a/source/grid/grid_tools.inst.in b/source/grid/grid_tools.inst.in index 2b3f29b8e1..2fc9db2041 100644 --- a/source/grid/grid_tools.inst.in +++ b/source/grid/grid_tools.inst.in @@ -103,7 +103,6 @@ for (X : TRIANGULATION_AND_DOFHANDLERS; deal_II_dimension : DIMENSIONS ; deal_II bool have_same_coarse_mesh (const X &mesh_1, const X &mesh_2); - \} #endif @@ -120,8 +119,12 @@ for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS std::map compute_local_to_global_vertex_index_map(const parallel::distributed::Triangulation &triangulation); + template + std::map > + extract_used_vertices(const Triangulation&mesh, + const Mapping &mapping); - \} + \} #endif } diff --git a/tests/grid/grid_tools_extract_used_vertices_01.cc b/tests/grid/grid_tools_extract_used_vertices_01.cc new file mode 100644 index 0000000000..1a7fbcbaa4 --- /dev/null +++ b/tests/grid/grid_tools_extract_used_vertices_01.cc @@ -0,0 +1,61 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2001 - 2015 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. +// +// --------------------------------------------------------------------- + +// Check extract used_vertices. + +#include "../tests.h" +#include +#include +#include +#include + + +template +void test () +{ + deallog << "dim: " << dim << ", spacedim: " + << spacedim << std::endl; + + Triangulation tria; + GridGenerator::hyper_cube(tria); + + tria.refine_global(1); + tria.begin_active()->set_refine_flag(); + tria.execute_coarsening_and_refinement(); + + for (auto cell = tria.begin_active(2); cell != tria.end(2); ++cell) + cell->set_coarsen_flag(); + + tria.execute_coarsening_and_refinement(); + + auto m = GridTools::extract_used_vertices(tria); + + for (auto &e: m) + deallog << "Vertex: " << e.first << ": " << e.second << std::endl; +}; + + +int main () +{ + initlog(); + test<1,1> (); + test<1,2> (); + test<1,3> (); + test<2,2> (); + test<2,3> (); + test<3,3> (); + return 0; +} + diff --git a/tests/grid/grid_tools_extract_used_vertices_01.output b/tests/grid/grid_tools_extract_used_vertices_01.output new file mode 100644 index 0000000000..974e8334be --- /dev/null +++ b/tests/grid/grid_tools_extract_used_vertices_01.output @@ -0,0 +1,61 @@ + +DEAL::dim: 1, spacedim: 1 +DEAL::Vertex: 0: 0.00000 +DEAL::Vertex: 1: 1.00000 +DEAL::Vertex: 2: 0.500000 +DEAL::dim: 1, spacedim: 2 +DEAL::Vertex: 0: 0.00000 0.00000 +DEAL::Vertex: 1: 1.00000 0.00000 +DEAL::Vertex: 2: 0.500000 0.00000 +DEAL::dim: 1, spacedim: 3 +DEAL::Vertex: 0: 0.00000 0.00000 0.00000 +DEAL::Vertex: 1: 1.00000 0.00000 0.00000 +DEAL::Vertex: 2: 0.500000 0.00000 0.00000 +DEAL::dim: 2, spacedim: 2 +DEAL::Vertex: 0: 0.00000 0.00000 +DEAL::Vertex: 1: 1.00000 0.00000 +DEAL::Vertex: 2: 0.00000 1.00000 +DEAL::Vertex: 3: 1.00000 1.00000 +DEAL::Vertex: 4: 0.500000 0.00000 +DEAL::Vertex: 5: 0.00000 0.500000 +DEAL::Vertex: 6: 1.00000 0.500000 +DEAL::Vertex: 7: 0.500000 1.00000 +DEAL::Vertex: 8: 0.500000 0.500000 +DEAL::dim: 2, spacedim: 3 +DEAL::Vertex: 0: 0.00000 0.00000 0.00000 +DEAL::Vertex: 1: 1.00000 0.00000 0.00000 +DEAL::Vertex: 2: 0.00000 1.00000 0.00000 +DEAL::Vertex: 3: 1.00000 1.00000 0.00000 +DEAL::Vertex: 4: 0.500000 0.00000 0.00000 +DEAL::Vertex: 5: 0.00000 0.500000 0.00000 +DEAL::Vertex: 6: 1.00000 0.500000 0.00000 +DEAL::Vertex: 7: 0.500000 1.00000 0.00000 +DEAL::Vertex: 8: 0.500000 0.500000 0.00000 +DEAL::dim: 3, spacedim: 3 +DEAL::Vertex: 0: 0.00000 0.00000 0.00000 +DEAL::Vertex: 1: 1.00000 0.00000 0.00000 +DEAL::Vertex: 2: 0.00000 1.00000 0.00000 +DEAL::Vertex: 3: 1.00000 1.00000 0.00000 +DEAL::Vertex: 4: 0.00000 0.00000 1.00000 +DEAL::Vertex: 5: 1.00000 0.00000 1.00000 +DEAL::Vertex: 6: 0.00000 1.00000 1.00000 +DEAL::Vertex: 7: 1.00000 1.00000 1.00000 +DEAL::Vertex: 8: 0.500000 0.00000 0.00000 +DEAL::Vertex: 9: 0.00000 0.500000 0.00000 +DEAL::Vertex: 10: 0.00000 0.00000 0.500000 +DEAL::Vertex: 11: 1.00000 0.500000 0.00000 +DEAL::Vertex: 12: 1.00000 0.00000 0.500000 +DEAL::Vertex: 13: 0.500000 1.00000 0.00000 +DEAL::Vertex: 14: 0.00000 1.00000 0.500000 +DEAL::Vertex: 15: 1.00000 1.00000 0.500000 +DEAL::Vertex: 16: 0.500000 0.00000 1.00000 +DEAL::Vertex: 17: 0.00000 0.500000 1.00000 +DEAL::Vertex: 18: 1.00000 0.500000 1.00000 +DEAL::Vertex: 19: 0.500000 1.00000 1.00000 +DEAL::Vertex: 20: 0.500000 0.00000 0.500000 +DEAL::Vertex: 21: 0.500000 0.500000 0.00000 +DEAL::Vertex: 22: 0.00000 0.500000 0.500000 +DEAL::Vertex: 23: 1.00000 0.500000 0.500000 +DEAL::Vertex: 24: 0.500000 1.00000 0.500000 +DEAL::Vertex: 25: 0.500000 0.500000 1.00000 +DEAL::Vertex: 26: 0.500000 0.500000 0.500000