]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add test of FE_Q<dim>(2) with dof_to_support_patch_map
authorSpencer Patty <srobertp@gmail.com>
Wed, 10 Feb 2016 16:47:18 +0000 (10:47 -0600)
committerSpencer Patty <srobertp@gmail.com>
Wed, 10 Feb 2016 16:47:18 +0000 (10:47 -0600)
tests/grid/get_dof_to_support_patch_map_01.cc [new file with mode: 0644]
tests/grid/get_dof_to_support_patch_map_01.output [new file with mode: 0644]

diff --git a/tests/grid/get_dof_to_support_patch_map_01.cc b/tests/grid/get_dof_to_support_patch_map_01.cc
new file mode 100644 (file)
index 0000000..fee6dd9
--- /dev/null
@@ -0,0 +1,186 @@
+
+
+// ---------------------------------------------------------------------
+//
+// Copyright (C)  2015 - 2016  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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test GridTools::map_dof_to_support_patch ()
+
+
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/base/point.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_tools.h>
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <fstream>
+#include <string>
+
+using namespace dealii;
+
+
+template<int dim>
+void test()
+{
+  using namespace dealii;
+
+  Triangulation<dim> triangulation (Triangulation<dim>::limit_level_difference_at_vertices);
+
+  GridGenerator::hyper_cube(triangulation,0,1);
+
+  triangulation.refine_global (2);
+  {
+
+    std::vector<Point<dim> > refine_centers;
+    std::vector<Point<dim> > coarsen_centers;
+
+    if(dim==1)
+    {
+      refine_centers.push_back(Point<dim> (1./8.));
+
+      coarsen_centers.push_back(Point<dim> (5./8.));
+      coarsen_centers.push_back(Point<dim> (7./8.));
+    }
+    else if (dim==2)
+    {
+      refine_centers.push_back(Point<dim> (1./8., 7./8.));
+
+      coarsen_centers.push_back(Point<dim> (5./8., 5./8.));
+      coarsen_centers.push_back(Point<dim> (5./8., 7./8.));
+      coarsen_centers.push_back(Point<dim> (7./8., 5./8.));
+      coarsen_centers.push_back(Point<dim> (7./8., 7./8.));
+
+    }
+    else if (dim==3)
+    {
+      refine_centers.push_back(Point<dim> (1./8., 7./8., 1./8.));
+
+      coarsen_centers.push_back(Point<dim> (7./8., 7./8., 7./8.));
+      coarsen_centers.push_back(Point<dim> (5./8., 7./8., 7./8.));
+      coarsen_centers.push_back(Point<dim> (7./8., 5./8., 7./8.));
+      coarsen_centers.push_back(Point<dim> (5./8., 5./8., 7./8.));
+      coarsen_centers.push_back(Point<dim> (7./8., 7./8., 5./8.));
+      coarsen_centers.push_back(Point<dim> (5./8., 7./8., 5./8.));
+      coarsen_centers.push_back(Point<dim> (7./8., 5./8., 5./8.));
+      coarsen_centers.push_back(Point<dim> (5./8., 5./8., 5./8.));
+    }
+    else
+      Assert(false, ExcNotImplemented() );
+
+    unsigned int index = 0;
+    for (typename Triangulation<dim>::active_cell_iterator
+         cell = triangulation.begin_active();
+         cell != triangulation.end(); ++cell, ++index)
+
+    {
+      Point<dim> cell_bary = cell->barycenter();
+
+      // refine cells
+      for (unsigned int i=0; i<refine_centers.size(); ++i)
+      {
+        double diff = 0;
+        for (unsigned int d=0; d<dim; ++d)
+          diff += pow(cell_bary[d] - refine_centers[i][d], 2.);
+        diff = std::sqrt(diff);
+
+        if (diff < 1e-14)
+        {
+          cell->set_refine_flag ();
+          break;
+        }
+      }
+      // coarsen cells
+      for (unsigned int i=0; i<coarsen_centers.size(); ++i)
+      {
+        double diff = 0;
+        for (unsigned int d=0; d<dim; ++d)
+          diff += pow(cell_bary[d] - coarsen_centers[i][d], 2.);
+        diff = std::sqrt(diff);
+
+        if (diff < 1e-14)
+        {
+          cell->set_coarsen_flag ();
+          break;
+        }
+      }
+    }
+    triangulation.execute_coarsening_and_refinement ();
+
+  }
+
+
+  DoFHandler<dim> dof_handler(triangulation);
+  unsigned int iFEDeg = 2;
+  FE_Q<dim> finite_element(iFEDeg);
+  dof_handler.distribute_dofs (finite_element);
+
+  std::map< types::global_dof_index, std::vector<typename dealii::DoFHandler<dim>::active_cell_iterator> >
+  dof_to_cell_map = GridTools::get_dof_to_support_patch_map(dof_handler);
+
+  for (unsigned int i=0; i<dof_handler.n_dofs(); ++i)
+  {
+    // loop through and print out barycenter of cells in patch of certain dofs
+    // in system
+    if (i % (dim == 1 ? 1 : (dim == 2 ? 5 : 25 )) == 0)
+    {
+      deallog << "Patch around dof " << i << ": ";
+      typename std::vector<typename dealii::DoFHandler<dim>::active_cell_iterator>::iterator
+      patch_iter = dof_to_cell_map[i].begin(),
+      patch_iter_end = dof_to_cell_map[i].end();
+      for (; patch_iter != patch_iter_end ; ++patch_iter )
+      {
+        typename dealii::DoFHandler<dim>::active_cell_iterator patch_cell = *(patch_iter);
+        Point<dim> cell_bary = patch_cell->barycenter();
+        deallog << "(" ;
+        for (unsigned int d=0; d<dim-1; ++d)
+          deallog << cell_bary[d] << ", ";
+        deallog << cell_bary[dim-1] << ") ";
+      }
+      deallog << std::endl;
+    }
+
+  }
+
+
+  // clean up data
+  dof_handler.clear();
+}
+
+
+int main()
+{
+  using namespace dealii;
+
+  initlog();
+  deallog.threshold_double(1.e-10);
+
+  deallog.push("1d");
+  test<1>();
+  deallog.pop();
+
+  deallog.push("2d");
+  test<2>();
+  deallog.pop();
+
+  deallog.push("3d");
+  test<3>();
+  deallog.pop();
+}
\ No newline at end of file
diff --git a/tests/grid/get_dof_to_support_patch_map_01.output b/tests/grid/get_dof_to_support_patch_map_01.output
new file mode 100644 (file)
index 0000000..d70a58e
--- /dev/null
@@ -0,0 +1,60 @@
+
+DEAL:1d::Patch around dof 0: (0.750000) (0.375000) 
+DEAL:1d::Patch around dof 1: (0.750000) 
+DEAL:1d::Patch around dof 2: (0.750000) 
+DEAL:1d::Patch around dof 3: (0.375000) (0.187500) 
+DEAL:1d::Patch around dof 4: (0.375000) 
+DEAL:1d::Patch around dof 5: (0.0625000) 
+DEAL:1d::Patch around dof 6: (0.0625000) (0.187500) 
+DEAL:1d::Patch around dof 7: (0.0625000) 
+DEAL:1d::Patch around dof 8: (0.187500) 
+DEAL:2d::Patch around dof 0: (0.750000, 0.750000) (0.375000, 0.375000) (0.625000, 0.375000) (0.875000, 0.375000) (0.375000, 0.625000) (0.375000, 0.875000) 
+DEAL:2d::Patch around dof 5: (0.750000, 0.750000) 
+DEAL:2d::Patch around dof 10: (0.125000, 0.125000) (0.375000, 0.125000) 
+DEAL:2d::Patch around dof 15: (0.125000, 0.125000) 
+DEAL:2d::Patch around dof 20: (0.375000, 0.125000) (0.625000, 0.125000) 
+DEAL:2d::Patch around dof 25: (0.125000, 0.375000) (0.375000, 0.375000) (0.125000, 0.625000) (0.375000, 0.625000) 
+DEAL:2d::Patch around dof 30: (0.375000, 0.375000) (0.625000, 0.375000) 
+DEAL:2d::Patch around dof 35: (0.625000, 0.125000) (0.875000, 0.125000) 
+DEAL:2d::Patch around dof 40: (0.875000, 0.125000) (0.875000, 0.375000) 
+DEAL:2d::Patch around dof 45: (0.750000, 0.750000) (0.625000, 0.375000) (0.875000, 0.375000) 
+DEAL:2d::Patch around dof 50: (0.750000, 0.750000) (0.625000, 0.375000) (0.875000, 0.375000) 
+DEAL:2d::Patch around dof 55: (0.125000, 0.625000) (0.375000, 0.625000) 
+DEAL:2d::Patch around dof 60: (0.375000, 0.625000) (0.375000, 0.875000) 
+DEAL:2d::Patch around dof 65: (0.375000, 0.875000) 
+DEAL:2d::Patch around dof 70: (0.0625000, 0.812500) 
+DEAL:2d::Patch around dof 75: (0.375000, 0.875000) (0.187500, 0.812500) (0.187500, 0.937500) 
+DEAL:2d::Patch around dof 80: (0.0625000, 0.937500) 
+DEAL:2d::Patch around dof 85: (0.0625000, 0.937500) 
+DEAL:3d::Patch around dof 0: (0.750000, 0.750000, 0.750000) (0.375000, 0.375000, 0.375000) (0.625000, 0.375000, 0.375000) (0.875000, 0.375000, 0.375000) (0.375000, 0.625000, 0.375000) (0.375000, 0.875000, 0.375000) (0.625000, 0.625000, 0.375000) (0.875000, 0.625000, 0.375000) (0.625000, 0.875000, 0.375000) (0.875000, 0.875000, 0.375000) (0.375000, 0.375000, 0.625000) (0.375000, 0.375000, 0.875000) (0.625000, 0.375000, 0.625000) (0.875000, 0.375000, 0.625000) (0.625000, 0.375000, 0.875000) (0.875000, 0.375000, 0.875000) (0.375000, 0.625000, 0.625000) (0.375000, 0.875000, 0.625000) (0.375000, 0.625000, 0.875000) (0.375000, 0.875000, 0.875000) 
+DEAL:3d::Patch around dof 25: (0.750000, 0.750000, 0.750000) 
+DEAL:3d::Patch around dof 50: (0.125000, 0.125000, 0.125000) (0.125000, 0.375000, 0.125000) 
+DEAL:3d::Patch around dof 75: (0.125000, 0.375000, 0.125000) (0.375000, 0.375000, 0.125000) (0.125000, 0.375000, 0.375000) (0.375000, 0.375000, 0.375000) (0.125000, 0.625000, 0.125000) (0.375000, 0.625000, 0.125000) (0.125000, 0.625000, 0.375000) (0.375000, 0.625000, 0.375000) 
+DEAL:3d::Patch around dof 100: (0.375000, 0.375000, 0.125000) (0.375000, 0.375000, 0.375000) 
+DEAL:3d::Patch around dof 125: (0.375000, 0.125000, 0.375000) (0.625000, 0.125000, 0.375000) 
+DEAL:3d::Patch around dof 150: (0.375000, 0.375000, 0.375000) 
+DEAL:3d::Patch around dof 175: (0.875000, 0.125000, 0.125000) (0.875000, 0.375000, 0.125000) 
+DEAL:3d::Patch around dof 200: (0.875000, 0.375000, 0.125000) (0.875000, 0.375000, 0.375000) (0.875000, 0.625000, 0.125000) (0.875000, 0.625000, 0.375000) 
+DEAL:3d::Patch around dof 225: (0.875000, 0.125000, 0.375000) (0.875000, 0.125000, 0.625000) 
+DEAL:3d::Patch around dof 250: (0.125000, 0.625000, 0.125000) (0.0625000, 0.812500, 0.0625000) (0.187500, 0.812500, 0.0625000) (0.0625000, 0.812500, 0.187500) (0.187500, 0.812500, 0.187500) 
+DEAL:3d::Patch around dof 275: (0.375000, 0.625000, 0.125000) (0.625000, 0.625000, 0.125000) 
+DEAL:3d::Patch around dof 300: (0.125000, 0.625000, 0.375000) (0.125000, 0.625000, 0.625000) 
+DEAL:3d::Patch around dof 325: (0.125000, 0.875000, 0.375000) (0.125000, 0.875000, 0.625000) 
+DEAL:3d::Patch around dof 350: (0.625000, 0.625000, 0.125000) 
+DEAL:3d::Patch around dof 375: (0.625000, 0.875000, 0.125000) (0.625000, 0.875000, 0.375000) 
+DEAL:3d::Patch around dof 400: (0.875000, 0.625000, 0.375000) (0.875000, 0.875000, 0.375000) 
+DEAL:3d::Patch around dof 425: (0.125000, 0.125000, 0.625000) (0.375000, 0.125000, 0.625000) (0.125000, 0.125000, 0.875000) (0.375000, 0.125000, 0.875000) 
+DEAL:3d::Patch around dof 450: (0.125000, 0.375000, 0.625000) (0.125000, 0.375000, 0.875000) (0.125000, 0.625000, 0.625000) (0.125000, 0.625000, 0.875000) 
+DEAL:3d::Patch around dof 475: (0.125000, 0.125000, 0.875000) (0.375000, 0.125000, 0.875000) 
+DEAL:3d::Patch around dof 500: (0.125000, 0.375000, 0.875000) (0.125000, 0.625000, 0.875000) 
+DEAL:3d::Patch around dof 525: (0.625000, 0.125000, 0.625000) (0.875000, 0.125000, 0.625000) (0.625000, 0.375000, 0.625000) (0.875000, 0.375000, 0.625000) 
+DEAL:3d::Patch around dof 550: (0.625000, 0.375000, 0.625000) 
+DEAL:3d::Patch around dof 575: (0.875000, 0.125000, 0.875000) (0.875000, 0.375000, 0.875000) 
+DEAL:3d::Patch around dof 600: (0.125000, 0.625000, 0.625000) (0.125000, 0.625000, 0.875000) 
+DEAL:3d::Patch around dof 625: (0.125000, 0.875000, 0.625000) 
+DEAL:3d::Patch around dof 650: (0.750000, 0.750000, 0.750000) (0.375000, 0.625000, 0.625000) (0.375000, 0.875000, 0.625000) (0.375000, 0.625000, 0.875000) (0.375000, 0.875000, 0.875000) 
+DEAL:3d::Patch around dof 675: (0.375000, 0.875000, 0.875000) 
+DEAL:3d::Patch around dof 700: (0.0625000, 0.812500, 0.0625000) 
+DEAL:3d::Patch around dof 725: (0.0625000, 0.937500, 0.0625000) (0.187500, 0.937500, 0.0625000) 
+DEAL:3d::Patch around dof 750: (0.125000, 0.875000, 0.375000) (0.0625000, 0.812500, 0.187500) (0.187500, 0.812500, 0.187500) (0.0625000, 0.937500, 0.187500) (0.187500, 0.937500, 0.187500) 
+DEAL:3d::Patch around dof 775: (0.125000, 0.875000, 0.375000) (0.0625000, 0.812500, 0.187500) (0.187500, 0.812500, 0.187500) (0.0625000, 0.937500, 0.187500) (0.187500, 0.937500, 0.187500) 

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.