From: Tobias Leicht Date: Thu, 20 Mar 2008 10:46:18 +0000 (+0000) Subject: Extend DataOut to enable a user-decision, which cells shall be written as curved... X-Git-Tag: v8.0.0~9260 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=748868c1ba3545dafceb217cffbd9263f67e9c5a;p=dealii.git Extend DataOut to enable a user-decision, which cells shall be written as curved cells using the provided mapping: none, cells at the boundary, or all (for MappingQEulerian). git-svn-id: https://svn.dealii.org/trunk@15914 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/deal.II/include/numerics/data_out.h b/deal.II/deal.II/include/numerics/data_out.h index 213da371a8..24c80723c2 100644 --- a/deal.II/deal.II/include/numerics/data_out.h +++ b/deal.II/deal.II/include/numerics/data_out.h @@ -1114,6 +1114,18 @@ class DataOut : public DataOut_DoFData */ typedef typename DataOut_DoFData::cell_iterator cell_iterator; typedef typename DataOut_DoFData::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 * 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 */ virtual void build_patches (const Mapping &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 * 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; }; diff --git a/deal.II/deal.II/source/numerics/data_out.cc b/deal.II/deal.II/source/numerics/data_out.cc index aefb3f2322..18001c900f 100644 --- a/deal.II/deal.II/source/numerics/data_out.cc +++ b/deal.II/deal.II/source/numerics/data_out.cc @@ -647,7 +647,10 @@ void DataOut::build_some_patches (Data &data) const hp::FECollection fe_collection(this->dofs->get_fe()); const hp::MappingCollection 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; idof_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::build_some_patches (Data &data) const FEValues &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 > & q_points=fe_patch_values.get_quadrature_points(); @@ -930,18 +940,26 @@ void DataOut::build_patches (const unsigned int n_subdivisions, const unsigned int n_threads_) { build_patches (StaticMappingQ1::mapping, - n_subdivisions, n_threads_); + n_subdivisions, n_threads_, no_curved_cells); } template void DataOut::build_patches (const Mapping &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)); diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h index 4f767829a8..06d213e89e 100644 --- a/deal.II/doc/news/changes.h +++ b/deal.II/doc/news/changes.h @@ -282,7 +282,16 @@ constraints individually.

deal.II

-
    +
      + +
    1. 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. +
      + (Tobias Leicht, 2008/03/20) +

    2. +
    3. 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