]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added get_used_vertices to GridTools::Cache 5238/head
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 13 Oct 2017 16:01:13 +0000 (18:01 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Fri, 13 Oct 2017 16:04:25 +0000 (18:04 +0200)
include/deal.II/grid/grid_tools_cache.h
include/deal.II/grid/grid_tools_cache_update_flags.h
source/grid/grid_tools_cache.cc
tests/grid/grid_tools_cache_03.cc [new file with mode: 0644]
tests/grid/grid_tools_cache_03.output [new file with mode: 0644]

index 138d604320ec48675480d92e8ff726b1e7f44388..c0e86ce790e36241ada0929b4dd436a19c0fa22f 100644 (file)
@@ -109,15 +109,22 @@ namespace GridTools
     const std::vector<std::vector<Tensor<1,spacedim>>>
     &get_vertex_to_cell_centers_directions() const;
 
+    /**
+     * Return the cached map of used vertices, as computed by
+     * GridTools::extract_used_vertices().
+     */
+    const std::map<unsigned int, Point<spacedim> >
+    &get_used_vertices() const;
+
     /**
      * Return a reference to the stored triangulation.
      */
-    const Triangulation<dim,spacedim> & get_triangulation() const;
+    const Triangulation<dim,spacedim> &get_triangulation() const;
 
     /**
      * Return a reference to the stored mapping.
      */
-    const Mapping<dim,spacedim> & get_mapping() const;
+    const Mapping<dim,spacedim> &get_mapping() const;
 
 #ifdef DEAL_II_WITH_NANOFLANN
     /**
@@ -155,7 +162,7 @@ namespace GridTools
 
     /**
      * Store vertex to cell center directions, as generated by
-     * GridTools::vertex_to_cell_centers_directions.
+     * GridTools::vertex_to_cell_centers_directions().
      */
     mutable std::vector<std::vector<Tensor<1,spacedim>>> vertex_to_cell_centers;
 
@@ -166,6 +173,12 @@ namespace GridTools
     mutable KDTree<spacedim> vertex_kdtree;
 #endif
 
+    /**
+     * Store the used vertices of the Triangulation, as generated by
+     * GridTools::extract_used_vertices().
+     */
+    mutable std::map<unsigned int, Point<spacedim>> used_vertices;
+
     /**
      * Storage for the status of the triangulation signal.
      */
@@ -176,7 +189,8 @@ namespace GridTools
 
   // Inline functions
   template<int dim, int spacedim>
-  inline const Triangulation<dim, spacedim>& Cache::get_triangulation() const
+  inline const Triangulation<dim, spacedim> &
+  Cache<dim,spacedim>::get_triangulation() const
   {
     return *tria;
   }
@@ -184,7 +198,8 @@ namespace GridTools
 
 
   template<int dim, int spacedim>
-  inline const Mapping<dim, spacedim>& Cache::get_mapping() const
+  inline const Mapping<dim, spacedim> &
+  Cache<dim,spacedim>::get_mapping() const
   {
     return *mapping;
   }
index 39bfd67f792fd301d6e647fd330675cf233fcb0e..dbac866d3ccc0b71429306acf152ca8b017e5402 100644 (file)
@@ -56,6 +56,11 @@ namespace GridTools
      */
     update_vertex_kdtree = 0x04,
 
+    /**
+     * Update a mapping of used vertices.
+     */
+    update_used_vertices = 0x08,
+
     /**
      * Update all objects.
      */
index 3aedd473ef4fe15a9ce19d5329e1e1a75c7f7d4b..2e7623d7ee9f00f6fa23acb155124017236bf0c6 100644 (file)
@@ -83,6 +83,20 @@ namespace GridTools
 
 
 
+  template<int dim, int spacedim>
+  const std::map<unsigned int, Point<spacedim> > &
+  Cache<dim,spacedim>::get_used_vertices() const
+  {
+    if (update_flags & update_used_vertices)
+      {
+        used_vertices = GridTools::extract_used_vertices(*tria, *mapping);
+        update_flags = update_flags & ~update_used_vertices;
+      }
+    return used_vertices;
+  }
+
+
+
 #ifdef DEAL_II_WITH_NANOFLANN
   template<int dim, int spacedim>
   const KDTree<spacedim> &Cache<dim,spacedim>::get_vertex_kdtree() const
diff --git a/tests/grid/grid_tools_cache_03.cc b/tests/grid/grid_tools_cache_03.cc
new file mode 100644 (file)
index 0000000..93108dc
--- /dev/null
@@ -0,0 +1,65 @@
+// ---------------------------------------------------------------------
+//
+// 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 and find_closest_vertex using a cache
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_tools_cache.h>
+#include <deal.II/grid/grid_out.h>
+
+
+template <int dim, int spacedim>
+void test (const Point<spacedim> &p)
+{
+  deallog << "dim: " << dim << ", spacedim: "
+          << spacedim << std::endl;
+
+  Triangulation<dim,spacedim> 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();
+
+  GridTools::Cache<dim,spacedim> cache(tria);
+
+  auto m = cache.get_used_vertices();
+
+  for (auto &e: m)
+    deallog << "Vertex: " << e.first << ": " << e.second << std::endl;
+
+  auto i = GridTools::find_closest_vertex(m,p);
+  deallog << "Closest vertex to " << p
+          << ", v["<< i << "] :"
+          << m[i] << std::endl;
+};
+
+
+int main ()
+{
+  initlog();
+  test<2,2> (Point<2>(.2,.2));
+  test<2,2> (Point<2>(.6,.9));
+  return 0;
+}
+
diff --git a/tests/grid/grid_tools_cache_03.output b/tests/grid/grid_tools_cache_03.output
new file mode 100644 (file)
index 0000000..de11874
--- /dev/null
@@ -0,0 +1,23 @@
+
+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::Closest vertex to 0.200000 0.200000, v[0] :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::Closest vertex to 0.600000 0.900000, v[7] :0.500000 1.00000

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.