]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Remove specializations of the GridReordering template class. Implementation of new...
authorhartmann <hartmann@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 6 Jul 2005 16:42:03 +0000 (16:42 +0000)
committerhartmann <hartmann@0785d39b-7218-0410-832d-ea1e28bc413d>
Wed, 6 Jul 2005 16:42:03 +0000 (16:42 +0000)
git-svn-id: https://svn.dealii.org/trunk@11080 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/grid/grid_reordering.h
deal.II/deal.II/source/grid/grid_reordering.cc

index d58434e6c511b9a41d7a96dd6915074554b61ad1..6ae4fea00285839fadaa6bdc17f989676188920d 100644 (file)
  * the triangulation from this data using the
  * Triangulation::create_triangulation() function.
  *
- * @author Wolfgang Bangerth, 2000, Michael Anderson 2003
+ * @author Wolfgang Bangerth, 2000, Michael Anderson 2003, Ralf Hartmann 2005
  */
 template <int dim>
 class GridReordering
-{
-};
-
-
-
-/**
- * This is the specialization of the general template for 1d. In 1d,
- * there is actually nothing to be done.
- *
- * @author Wolfgang Bangerth, 2000
- */
-template <>
-class GridReordering<1>
-{
-  public:
-                                    /**
-                                     * Do nothing, since in 1d no
-                                     * reordering is necessary.
-                                     */
-    static void reorder_cells (const std::vector<CellData<1> > &);
-};
-
-
-
-/**
- * This specialization of the general template implements the
- * 2d-algorithm described in the documentation of the general
- * template.
- *
- * @author Michael Anderson, 2003
- */
-template <>
-class GridReordering<2>
 {
   public:
+    
                                     /**
                                      *  This is the main function,
                                      *  doing what is announced in
                                      *  the general documentation of
-                                     *  this class.
+                                     *  this class for dim=2 and 3
+                                     *  and doing nothing for dim=1.
                                      */
-    static void reorder_cells (std::vector<CellData<2> > &original_cells);
-};
+    static void reorder_cells (std::vector<CellData<dim> > &original_cells);
 
-
-
-/**
- * This specialization of the general template implements the
- * 3d-algorithm described in the documentation of the general
- * template.
- *
- * @author Michael Anderson, 2003
- */
-template <>
-class GridReordering<3>
-{
-  public:
                                     /**
-                                     *  This is the main function,
-                                     *  doing what is announced in
-                                     *  the general documentation of
-                                     *  this class.
+                                     * Grids generated by grid
+                                     * generators may have an
+                                     * orientation of cells which is
+                                     * the inverse of the orientation
+                                     * required by deal.II.
+                                     *
+                                     * In 3d this function checks
+                                     * whether all cells have
+                                     * negative or positiv
+                                     * measure/volume. In the former
+                                     * case, all cells are inverted.
+                                     * It does nothing in 1 and 2d.
+                                     *
+                                     * The invertion of cells might
+                                     * also work when only a subset
+                                     * of all cells have negative
+                                     * volume. However, grids
+                                     * consisting of a mixture of
+                                     * negative and positiv oriented
+                                     * cells are very likely to be
+                                     * broken. Therefore, an
+                                     * exception is thrown, in case
+                                     * cells are not uniformly
+                                     * oriented.
+                                     *
+                                     * Note, that this function
+                                     * should be called before
+                                     * reorder_cells().
                                      */
-    static void reorder_cells (std::vector<CellData<3> > &original_cells);
+    static void invert_all_cells_of_negative_grid(
+      const std::vector<Point<dim> > &all_vertices,
+      std::vector<CellData<dim> > &original_cells);    
 };
 
 
 
-
 #endif
index 74c51493c3a0967dfebbfb8ca045bb124e72000e..4d9bf6c38eafd02cd7024926b9ed4557d5b46e81 100644 (file)
@@ -13,6 +13,7 @@
 
 #include <grid/grid_reordering.h>
 #include <grid/grid_reordering_internal.h>
+#include <grid/grid_tools.h>
 
 #include <algorithm>
 #include <set>
 
 #if deal_II_dimension == 1
 
-void GridReordering<1>::reorder_cells (const std::vector<CellData<1> > &)
+template<>
+void
+GridReordering<1>::reorder_cells (std::vector<CellData<1> > &)
 {
                                   // there should not be much to do
                                   // in 1d...
 }
 
+
+template<>
+void
+GridReordering<1>::invert_all_cells_of_negative_grid(const std::vector<Point<1> > &,
+                                                    std::vector<CellData<1> > &)
+{
+                                  // nothing to be done in 1d
+}
+
 #endif
 
 
@@ -557,7 +569,9 @@ namespace internal
 } // namespace internal
 
 
-void GridReordering<2>::reorder_cells (std::vector<CellData<2> > &original_cells)
+template<>
+void
+GridReordering<2>::reorder_cells (std::vector<CellData<2> > &original_cells)
 {
                                    // check if grids are already
                                    // consistent. if so, do
@@ -569,6 +583,15 @@ void GridReordering<2>::reorder_cells (std::vector<CellData<2> > &original_cells
   internal::GridReordering2d::GridReordering().reorient(original_cells);
 }
 
+
+template<>
+void
+GridReordering<2>::invert_all_cells_of_negative_grid(const std::vector<Point<2> > &,
+                                                    std::vector<CellData<2> > &)
+{
+                                  // nothing to be done in 2d
+}
+
 #endif
 
 
@@ -1409,7 +1432,7 @@ namespace internal
 }
 
 
-
+template<>
 void
 GridReordering<3>::reorder_cells (std::vector<CellData<3> > &incubes)
 {
@@ -1422,6 +1445,43 @@ GridReordering<3>::reorder_cells (std::vector<CellData<3> > &incubes)
 }
 
 
+template<>
+void
+GridReordering<3>::invert_all_cells_of_negative_grid(
+  const std::vector<Point<3> > &all_vertices,
+  std::vector<CellData<3> > &cells)
+{
+  unsigned int n_negative_cells=0;
+  for (unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
+    if (GridTools::cell_measure(all_vertices, cells[cell_no].vertices) < 0)
+      {
+       ++n_negative_cells;
+       for (unsigned int i=0; i<4; ++i)
+         swap(cells[cell_no].vertices[i], cells[cell_no].vertices[i+4]);
+       
+                                        // check whether the
+                                        // resulting cell is now ok.
+                                        // if not, then the grid is
+                                        // seriously broken and
+                                        // should be sticked into the
+                                        // bin
+       AssertThrow(GridTools::cell_measure(all_vertices, cells[cell_no].vertices) > 0,
+                   ExcInternalError());
+      }
+
+                                  // We assuming that all cells of a
+                                  // grid have either positive or
+                                  // negative volumes but not both
+                                  // mixed. Although above reordering
+                                  // might work also on single cells,
+                                  // grids with both kind of cells
+                                  // are very likely to be
+                                  // broken. Check for this here.
+  AssertThrow(n_negative_cells==0 || n_negative_cells==cells.size(), ExcInternalError());
+}
+
+
+
        
 
 #endif // deal_II_dimension == 3

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.