]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Added set_manifold with two arguments.
authorheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 14 Jul 2014 15:54:06 +0000 (15:54 +0000)
committerheltai <heltai@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 14 Jul 2014 15:54:06 +0000 (15:54 +0000)
git-svn-id: https://svn.dealii.org/trunk@33168 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/grid/tria.h
deal.II/source/grid/tria.cc
tests/manifold/manifold_id_05.cc
tests/manifold/manifold_id_06.cc
tests/manifold/manifold_id_07.cc

index 5af6ef9b2959496fb1db200940348a856cbbf62e..e1bb6c5b1afd0feb7281b02d6cce0c3fc7fffea2 100644 (file)
@@ -1781,59 +1781,39 @@ public:
   virtual void set_mesh_smoothing (const MeshSmoothing mesh_smoothing);
 
   /**
-   * If @p dim==spacedim, assign a boundary
-   * object to a certain part of the
-   * boundary of a the triangulation. If a
-   * face with boundary number @p number is
-   * refined, this object is used to find
-   * the location of new vertices on the
-   * boundary (see the results section of step-49
-   * for a more in-depth discussion of this, with examples).
-   * It is also used for
-   * non-linear (i.e.: non-Q1)
-   * transformations of cells to the unit
-   * cell in shape function calculations.
-   *
-   * If @p dim!=spacedim the boundary object
-   * is in fact the exact manifold that the
-   * triangulation is approximating (for
-   * example a circle approximated by a
-   * polygon triangulation). As above, the
-   * refinement is made in such a way that
-   * the new points are located on the exact
-   * manifold.
-   *
-   * Numbers of boundary objects correspond
-   * to material numbers of faces at the
-   * boundary, for instance the material id
-   * in a UCD input file. They are not
-   * necessarily consecutive but must be in
-   * the range 0-(types::boundary_id-1).
-   * Material IDs on boundaries are also
-   * called boundary indicators and are
-   * accessed with accessor functions
-   * of that name.
-   *
-   * The @p boundary_object is not copied
-   * and MUST persist until the
-   * triangulation is destroyed. This is
-   * also true for triangulations generated
-   * from this one by @p
-   * copy_triangulation.
-   *
-   * It is possible to remove or replace
-   * the boundary object during the
-   * lifetime of a non-empty
-   * triangulation. Usually, this is done
-   * before the first refinement and is
-   * dangerous afterwards. Removal of a
-   * boundary object is done by
-   * <tt>set_boundary(number)</tt>,
-   * i.e. the function of same name but
-   * only one argument. This operation then
-   * replaces the boundary object given
-   * before by a straight boundary
-   * approximation.
+   * If @p dim==spacedim, assign a boundary object to a certain part
+   * of the boundary of a the triangulation. If a face with boundary
+   * number @p number is refined, this object is used to find the
+   * location of new vertices on the boundary (see the results section
+   * of step-49 for a more in-depth discussion of this, with
+   * examples).  It is also used for non-linear (i.e.: non-Q1)
+   * transformations of cells to the unit cell in shape function
+   * calculations.
+   *
+   * If @p dim!=spacedim the boundary object is in fact the exact
+   * manifold that the triangulation is approximating (for example a
+   * circle approximated by a polygon triangulation). As above, the
+   * refinement is made in such a way that the new points are located
+   * on the exact manifold.
+   *
+   * Numbers of boundary objects correspond to material numbers of
+   * faces at the boundary, for instance the material id in a UCD
+   * input file. They are not necessarily consecutive but must be in
+   * the range 0-(types::boundary_id-1).  Material IDs on boundaries
+   * are also called boundary indicators and are accessed with
+   * accessor functions of that name.
+   *
+   * The @p boundary_object is not copied and MUST persist until the
+   * triangulation is destroyed. This is also true for triangulations
+   * generated from this one by @p copy_triangulation.
+   *
+   * It is possible to remove or replace the boundary object during
+   * the lifetime of a non-empty triangulation. Usually, this is done
+   * before the first refinement and is dangerous afterwards. Removal
+   * of a boundary object is done by <tt>set_boundary(number)</tt>,
+   * i.e. the function of same name but only one argument. This
+   * operation then replaces the boundary object given before by a
+   * straight boundary approximation.
    *
    * @ingroup boundary
    *
@@ -1854,6 +1834,35 @@ public:
    * @see @ref GlossBoundaryIndicator "Glossary entry on boundary indicators"
    */
   void set_boundary (const types::manifold_id number);
+
+  /**
+   * Assign a manifold object to a certain part of the the
+   * triangulation. If an object with manfifold number @p number is
+   * refined, this object is used to find the location of new vertices
+   * (see the results section of step-49 for a more in-depth
+   * discussion of this, with examples).  It is also used for
+   * non-linear (i.e.: non-Q1) transformations of cells to the unit
+   * cell in shape function calculations.
+   *
+   * The @p manifold_object is not copied and MUST persist until the
+   * triangulation is destroyed. This is also true for triangulations
+   * generated from this one by @p copy_triangulation.
+   *
+   * It is possible to remove or replace the boundary object during
+   * the lifetime of a non-empty triangulation. Usually, this is done
+   * before the first refinement and is dangerous afterwards. Removal
+   * of a manifold object is done by <tt>set_manifold(number)</tt>,
+   * i.e. the function of same name but only one argument. This
+   * operation then replaces the manifold object given before by a
+   * straight manifold approximation.
+   *
+   * @ingroup manifold
+   *
+   * @see @ref GlossManifoldIndicator "Glossary entry on manifold indicators"
+   */
+  void set_manifold (const types::manifold_id   number,
+                     const Boundary<dim,spacedim> &manifold_object);
+
     
   /**
    * Reset those parts of the triangulation with the given manifold_id
index 3ed89fc573caf7f392f325372b60cc63f386abdb..eeeef5adc971bd468b1e6da652958e81761dd956 100644 (file)
@@ -9591,11 +9591,19 @@ template <int dim, int spacedim>
 void
 Triangulation<dim, spacedim>::set_boundary (const types::manifold_id m_number,
                                             const Boundary<dim, spacedim> &boundary_object)
+{
+  set_manifold(m_number, boundary_object);
+}
+
+template <int dim, int spacedim>
+void
+Triangulation<dim, spacedim>::set_manifold (const types::manifold_id m_number,
+                                            const Boundary<dim, spacedim> &manifold_object)
 {
   Assert(m_number < numbers::invalid_manifold_id,
         ExcIndexRange(m_number,0,numbers::invalid_manifold_id));
 
-  manifold[m_number] = &boundary_object;
+  manifold[m_number] = &manifold_object;
 }
 
 
index e6a5a52b558d2182ee19ab55555dee206e474dfb..5a4ca3e063aeb7586c89160270fafa38f1d0bcc1 100644 (file)
@@ -41,7 +41,7 @@ void test(unsigned int ref=1){
        typename Triangulation<dim,spacedim>::active_cell_iterator cell;
        
        tria.begin_active()->face(0)->set_manifold_id(1);
-       tria.set_boundary(1,boundary);
+       tria.set_manifold(1,boundary);
        
        tria.refine_global(2);
 
index 8fa4e892ca284764c533b840c829520b93a1c83f..442a291da0d342952c73041d9458a89dad7548cb 100644 (file)
@@ -52,7 +52,7 @@ void test(unsigned int ref=1){
              if(cell->face(f)->center().distance(center)< radius)
                cell->face(f)->set_all_manifold_ids(1);
 
-       tria.set_boundary(1,boundary);
+       tria.set_manifold(1,boundary);
        tria.refine_global(2);
 
        GridOut gridout;
index 9e0662f45aee616e50daf2a495d6e73abb2af9df..204f66e0e25926a39ec469ea3814e35812b516c4 100644 (file)
@@ -43,7 +43,7 @@ void test(unsigned int ref=1){
        typename Triangulation<dim,spacedim>::active_cell_iterator cell;
        
        tria.begin_active()->face(0)->set_manifold_id(1);
-       tria.set_boundary(1,boundary);
+       tria.set_manifold(1,boundary);
        
        tria.refine_global(2);
 

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.