]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Extend DataOut to enable a user-decision, which cells shall be written as curved...
authorTobias Leicht <tobias.leicht@dlr.de>
Thu, 20 Mar 2008 10:46:18 +0000 (10:46 +0000)
committerTobias Leicht <tobias.leicht@dlr.de>
Thu, 20 Mar 2008 10:46:18 +0000 (10:46 +0000)
git-svn-id: https://svn.dealii.org/trunk@15914 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/numerics/data_out.h
deal.II/deal.II/source/numerics/data_out.cc
deal.II/doc/news/changes.h

index 213da371a833b8b9971caa02adeb69e6ade8bb2d..24c80723c29521f7fa109d98b69a7ab47c185bba 100644 (file)
@@ -1114,6 +1114,18 @@ class DataOut : public DataOut_DoFData<DH, DH::dimension>
                                      */
     typedef typename DataOut_DoFData<DH,DH::dimension>::cell_iterator cell_iterator;
     typedef typename DataOut_DoFData<DH,DH::dimension>::active_cell_iterator active_cell_iterator;
+
+                                    /**
+                                     * Enumeration describing the region of the
+                                     * domain in which curved cells shall be
+                                     * created.
+                                     */
+    enum CurvedCellRegion
+    {
+         no_curved_cells,
+         curved_boundary,
+         curved_inner_cells
+    };
     
                                     /**
                                      * This is the central function
@@ -1153,13 +1165,24 @@ class DataOut : public DataOut_DoFData<DH, DH::dimension>
                                      * points interior of subdivided
                                      * patches which originate from
                                      * cells at the boundary of the
-                                     * domain are computed using the
+                                     * domain can be computed using the
                                      * mapping, i.e. a higher order
                                      * mapping leads to a
                                      * representation of a curved
                                      * boundary by using more
-                                     * subdivisions. The
-                                     * mapping argument can also be used
+                                     * subdivisions. Some mappings like
+                                     * MappingQEulerian result in curved cells
+                                     * in the interior of the domain. However,
+                                     * there is nor easy way to get this
+                                     * information from the Mapping. Thus the
+                                     * last argument @p curved_region take one
+                                     * of three values resulting in no curved
+                                     * cells at all, curved cells at the
+                                     * boundary (default) or curved cells in
+                                     * the whole domain.
+                                     *
+                                     * Even for non-curved cells the
+                                     * mapping argument can be used
                                      * for the Eulerian mappings (see
                                      * class MappingQ1Eulerian) where
                                      * a mapping is used not only to
@@ -1180,7 +1203,8 @@ class DataOut : public DataOut_DoFData<DH, DH::dimension>
                                      */
     virtual void build_patches (const Mapping<DH::dimension> &mapping,
                                const unsigned int n_subdivisions = 0,
-                               const unsigned int n_threads      = multithread_info.n_default_threads);
+                               const unsigned int n_threads      = multithread_info.n_default_threads,
+                               const CurvedCellRegion curved_region = curved_boundary);
     
                                     /**
                                      * Return the first cell which we
@@ -1270,6 +1294,12 @@ class DataOut : public DataOut_DoFData<DH, DH::dimension>
                                      * once and generates all patches.
                                      */
     void build_some_patches (Data &data);
+
+                                    /**
+                                     * Store in which region cells shall be
+                                     * curved, if a Mapping is provided.
+                                     */
+    CurvedCellRegion curved_cell_region;
 };
 
 
index aefb3f2322382d06b12a9952aa09b9a5f3ee3573..18001c900f110618562dde8e02341437291d1601 100644 (file)
@@ -647,7 +647,10 @@ void DataOut<dim,DH>::build_some_patches (Data &data)
   const hp::FECollection<DH::dimension>      fe_collection(this->dofs->get_fe());
   const hp::MappingCollection<DH::dimension> mapping_collection(*(data.mapping));
 
-  UpdateFlags update_flags=update_values | update_quadrature_points;
+  UpdateFlags update_flags=update_values;
+  if (curved_cell_region != no_curved_cells)
+    update_flags |= update_quadrature_points;
+  
   for (unsigned int i=0; i<this->dof_data.size(); ++i)
     if (this->dof_data[i]->postprocessor)
       update_flags |= this->dof_data[i]->postprocessor->get_needed_update_flags();
@@ -715,11 +718,18 @@ void DataOut<dim,DH>::build_some_patches (Data &data)
           const FEValues<DH::dimension> &fe_patch_values
             = x_fe_patch_values.get_present_fe_values ();
          
-                                          // if the cell is at the boundary,
+                                          // depending on the requested output
+                                          // of curved cells, if necessary
                                           // append the quadrature points to
-                                          // the last rows of the
-                                          // patch->data member
-         if (cell->at_boundary())
+                                          // the last rows of the patch->data
+                                          // member. THis is the case if we
+                                          // want to produce curved cells at
+                                          // the boundary and this cell
+                                          // actually is at the boundary, or
+                                          // else if we want to produce curved
+                                          // cells everywhere
+         if (curved_cell_region==curved_inner_cells ||
+             (curved_cell_region==curved_boundary && cell->at_boundary()))
            {
              Assert(patch->space_dim==dim, ExcInternalError());
              const std::vector<Point<dim> > & q_points=fe_patch_values.get_quadrature_points();
@@ -930,18 +940,26 @@ void DataOut<dim,DH>::build_patches (const unsigned int n_subdivisions,
                                     const unsigned int n_threads_) 
 {
   build_patches (StaticMappingQ1<DH::dimension>::mapping,
-                n_subdivisions, n_threads_);
+                n_subdivisions, n_threads_, no_curved_cells);
 }
 
 
 template <int dim, class DH>
 void DataOut<dim,DH>::build_patches (const Mapping<DH::dimension> &mapping,
                                     const unsigned int nnnn_subdivisions,
-                                    const unsigned int n_threads_) 
+                                    const unsigned int n_threads_,
+                                    const CurvedCellRegion curved_region) 
 {
   unsigned int n_subdivisions = (nnnn_subdivisions != 0)
                                ? nnnn_subdivisions
                                : this->default_subdivisions;
+                                  // store the region in which cells shall be
+                                  // curved. If only one subdivision is
+                                  // requested then there is no need to do this
+                                  // at all
+  curved_cell_region=curved_region;
+  if (n_subdivisions<2)
+    curved_cell_region=no_curved_cells;
   
   Assert (n_subdivisions >= 1,
          ExcInvalidNumberOfSubdivisions(n_subdivisions));
index 4f767829a8a78a55e94945094eb4ae57dfc0d466..06d213e89e5904ff04c4fd28df9fdacb8af5efe5 100644 (file)
@@ -282,7 +282,16 @@ constraints individually.
 <a name="deal.II"></a>
 <h3>deal.II</h3>
 
-<ol> 
+<ol>
+
+  <li> <p>Extended: DataOut has an extended interface now which enables a user
+  decision, which cells shall be written as curved cells using the provided
+  mapping: none, only cells at the boundary, or all cells. The last choice can
+  be useful when a mapping like MappingQEulerian is employed.
+  <br>
+  (Tobias Leicht, 2008/03/20)
+  </p></li>
+
   <li> <p>Changed: From now on the SolutionTransfer class does not use any user
   pointers of the underlying Triangulation object. Thus several SolutionTransfer
   instances with different DoFHandler objects can be used simultaneously. It is

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.