]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Split _01 test into two pieces: one that works for a sphere surface and one that...
authorbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 16 Nov 2010 19:55:27 +0000 (19:55 +0000)
committerbangerth <bangerth@0785d39b-7218-0410-832d-ea1e28bc413d>
Tue, 16 Nov 2010 19:55:27 +0000 (19:55 +0000)
git-svn-id: https://svn.dealii.org/trunk@22768 0785d39b-7218-0410-832d-ea1e28bc413d

tests/codim_one/extract_boundary_mesh_01.cc
tests/codim_one/extract_boundary_mesh_02.cc [new file with mode: 0644]

index 87a1b69ec02b083a2ceaee28c85b6690fc6fdf8f..e7183458cee17006f877630a970240d2bce245bc 100644 (file)
  We test that the order of cells and the orientation
  of the vertices is consistent between the two meshes.
 
- It's faling when boundary object is attached and refinement
- for a ball in 3D
+ Compared to the _00 test, we here test things with the surface of a
+ 3d sphere. There is no surface object attached, so the surface is
+ topologically a sphere but not geometrically. The _02 test checks the
+ same with a boundary object attached.
 */
 
 #include "../tests.h"
@@ -78,7 +80,7 @@ template <int dim, int spacedim>
 void save_mesh(const Triangulation<dim,spacedim>& tria)
 {
   GridOut grid_out;
-  grid_out.write_ucd (tria, deallog.get_file_stream());
+  grid_out.write_gnuplot (tria, deallog.get_file_stream());
 }
 
 
@@ -98,26 +100,18 @@ int main ()
     map< Triangulation<dim-1,dim>::cell_iterator,
         Triangulation<dim,dim>::face_iterator>
       surface_to_volume_mapping;
-    const HyperBallBoundary<dim> boundary_description;
     Triangulation<dim> volume_mesh;
     GridGenerator::hyper_ball(volume_mesh);
-    volume_mesh.set_boundary (0, boundary_description);
 
     volume_mesh.refine_global (1);
 
-    // save_mesh(volume_mesh);
-
     Triangulation<dim-1,dim> boundary_mesh;
 
     GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
                                      surface_to_volume_mapping);
 
-    if (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping,2))
-      deallog << "Passed.";
-    else
-      deallog << "Failed.";
-    deallog << endl;
-
+    Assert (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping,2),
+           ExcInternalError());
     save_mesh(boundary_mesh);
   }
 
diff --git a/tests/codim_one/extract_boundary_mesh_02.cc b/tests/codim_one/extract_boundary_mesh_02.cc
new file mode 100644 (file)
index 0000000..adb298b
--- /dev/null
@@ -0,0 +1,122 @@
+//----------------------------  extract_boundary_mesh_01.cc  ---------------------------
+//    $Id: bem.cc 22693 2010-11-11 20:11:47Z kanschat $
+//    Version: $Name$
+//
+//    Copyright (C) 2010 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.
+//
+//----------------------------  extract_boundary_mesh_01.cc  ---------------------------
+
+/*
+ Code for testing the function
+ GridTools::extract_boundary_mesh (...).
+ We test that the order of cells and the orientation
+ of the vertices is consistent between the two meshes.
+
+ It's faling when boundary object is attached and refinement
+ for a ball in 3D
+*/
+
+#include "../tests.h"
+
+#include <grid/tria.h>
+#include <grid/grid_generator.h>
+#include <grid/tria_boundary_lib.h>
+#include <grid/grid_out.h>
+#include <grid/grid_tools.h>
+
+using namespace std;
+
+
+template <int s_dim, int spacedim>
+bool test_vertices_orientation(const Triangulation<s_dim,spacedim> &boundary_mesh,
+                              map< typename Triangulation<s_dim,spacedim>::cell_iterator,
+                                   typename Triangulation<s_dim+1,spacedim>::face_iterator >
+                              &surface_to_volume_mapping,
+                              const int verbosity = 1)
+{
+  typename Triangulation<s_dim,spacedim>::active_cell_iterator
+    cell = boundary_mesh.begin_active(),
+    endc = boundary_mesh.end();
+  typename Triangulation<s_dim+1,spacedim>::face_iterator face;
+
+  bool success = true;
+
+  if (verbosity>1){
+    deallog << "The last column should be 0" <<endl;
+    deallog << "Vol faces" << "\t\t" << "Surf cell" <<
+      "\t\t" << "Distance" <<endl<<endl;
+  }
+
+  for (; cell!=endc; ++cell){
+
+    face = surface_to_volume_mapping [cell];
+
+    for (unsigned int k=0; k<GeometryInfo<s_dim>::vertices_per_cell; ++k)
+      {
+       Point<spacedim> diff(face->vertex(k));
+       diff -= cell->vertex(k);
+       if (verbosity>1){
+         deallog << face->vertex(k) << "\t\t";
+         deallog << cell->vertex(k) << "\t\t\t" << diff.square() << endl;
+       }
+       if (diff.square()>0){
+         success = false;
+         break;
+       }
+      }
+    if (verbosity>1) deallog << endl;
+  }
+  return success;
+}
+
+template <int dim, int spacedim>
+void save_mesh(const Triangulation<dim,spacedim>& tria)
+{
+  GridOut grid_out;
+  grid_out.write_ucd (tria, deallog.get_file_stream());
+}
+
+
+int main ()
+{
+
+  ofstream logfile("extract_boundary_mesh_01/output");
+  deallog.attach(logfile);
+  deallog.depth_console(0);
+
+
+  { // Extract the boundary of a hyper-sphere
+
+    const int dim = 3;
+    deallog << "Testing hyper_cube in dim: " << dim << "..."<< endl;
+
+    map< Triangulation<dim-1,dim>::cell_iterator,
+        Triangulation<dim,dim>::face_iterator>
+      surface_to_volume_mapping;
+    const HyperBallBoundary<dim> boundary_description;
+    Triangulation<dim> volume_mesh;
+    GridGenerator::hyper_ball(volume_mesh);
+    volume_mesh.set_boundary (0, boundary_description);
+
+    volume_mesh.refine_global (1);
+
+    // save_mesh(volume_mesh);
+
+    Triangulation<dim-1,dim> boundary_mesh;
+
+    GridTools::extract_boundary_mesh (volume_mesh, boundary_mesh,
+                                     surface_to_volume_mapping);
+
+    Assert (test_vertices_orientation(boundary_mesh, surface_to_volume_mapping,2),
+           ExcInternalError());
+    save_mesh(boundary_mesh);
+  }
+
+
+   return 0;
+ }

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.