]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
In 3D coarsening might delete children of lines which are still needed, if 6 or more...
authorleicht <leicht@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 28 Jun 2007 14:05:18 +0000 (14:05 +0000)
committerleicht <leicht@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 28 Jun 2007 14:05:18 +0000 (14:05 +0000)
git-svn-id: https://svn.dealii.org/trunk@14822 0785d39b-7218-0410-832d-ea1e28bc413d

tests/fail/line_coarsening_3d.cc [new file with mode: 0644]
tests/fail/line_coarsening_3d/cmp/generic [new file with mode: 0644]

diff --git a/tests/fail/line_coarsening_3d.cc b/tests/fail/line_coarsening_3d.cc
new file mode 100644 (file)
index 0000000..df198b8
--- /dev/null
@@ -0,0 +1,111 @@
+//----------------------------  line_coarsening_3d.cc  ---------------------------
+//    $Id$
+//    Version: $Name$ 
+//
+//    Copyright (C) 2007 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//----------------------------  line_coarsening_3d.cc  ---------------------------
+
+
+#include "../tests.h"
+#include <base/logstream.h>
+#include <grid/tria.h>
+#include <grid/tria_accessor.h>
+#include <grid/tria_iterator.h>
+#include <grid/grid_out.h>
+
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+
+
+std::ofstream logfile("line_coarsening_3d/output");
+
+
+// create a triangulation of a cylinder, where the line along the axis is common
+// to all cells, thus n_cells are neighboring at this line. in delete_children,
+// this is not accounted for (for n_cells>5) and the child lines are deleted,
+// although they are still needed.
+void create_star_structured_cylinder (Triangulation<3> &coarse_grid,
+                                     const unsigned int n_cells)
+{
+  Assert(n_cells>1, ExcNotImplemented());
+  
+  std::vector<Point<3> > points(2*(1+2*n_cells));
+  points[0] = Point<3>();
+  points[1] = Point<3>(1,0,0);
+  for (unsigned int i=0; i<2*n_cells-1; ++i)
+    {
+      points[2+i]=Point<3>(std::cos(deal_II_numbers::PI/n_cells*(i+1)),
+                          std::sin(deal_II_numbers::PI/n_cells*(i+1)),
+                          0);
+    }
+  
+  for (unsigned int i=0; i<2*n_cells+1; ++i)
+    {
+      points[1+2*n_cells+i]=points[i]+Point<3>(0,0,-1);
+    }
+  
+  std::vector<CellData<3> > cells(n_cells);
+    
+  for (unsigned int c=0; c<n_cells; ++c)
+    {
+      cells[c].vertices[0]   = 0;
+      cells[c].vertices[1]   = 1+2*c;
+      cells[c].vertices[2]   = 2+2*c;
+      cells[c].vertices[3]   = (3+2*c)%(2*n_cells);
+      cells[c].vertices[4]   = 0+2*n_cells+1;
+      cells[c].vertices[5]   = 1+2*c+2*n_cells+1;
+      cells[c].vertices[6]   = 2+2*c+2*n_cells+1;
+      cells[c].vertices[7]   = (3+2*c)%(2*n_cells)+2*n_cells+1;
+    }
+                                   // finally generate a triangulation
+                                   // out of this
+  coarse_grid.create_triangulation_compatibility (points, cells, SubCellData());
+}
+
+
+void check ()
+{
+  const unsigned int dim=3;
+                                   // create tria
+  Triangulation<dim> tria;
+  create_star_structured_cylinder(tria,6);
+                                  // out of the six cells, refine the first and
+                                  // fourth
+  Triangulation<3>::active_cell_iterator cell=tria.begin_active();
+  cell->set_refine_flag();
+  ++cell;
+  ++cell;
+  ++cell;
+  cell->set_refine_flag();
+  tria.execute_coarsening_and_refinement ();
+                                  // write grid to file
+  GridOut().write_ucd(tria,logfile);
+                                  // coarsen the first cell again
+  cell=tria.begin_active(1);
+  for (unsigned int c=0;c<8; ++c)
+    {
+      cell->set_coarsen_flag();
+      ++cell;
+    }
+  tria.execute_coarsening_and_refinement ();
+  
+  GridOut().write_ucd(tria,logfile);
+}
+
+
+int main()
+{
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+  check();
+  return 0;
+}
+
diff --git a/tests/fail/line_coarsening_3d/cmp/generic b/tests/fail/line_coarsening_3d/cmp/generic
new file mode 100644 (file)
index 0000000..e69de29

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.