]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce `patch_dim` and `patch_spacedim` in DataOut classes.
authorMarc Fehling <mafehling.git@gmail.com>
Tue, 24 Nov 2020 05:52:42 +0000 (22:52 -0700)
committerMarc Fehling <mafehling.git@gmail.com>
Thu, 27 May 2021 00:16:34 +0000 (18:16 -0600)
include/deal.II/base/data_out_base.h
include/deal.II/numerics/data_out_faces.h
include/deal.II/numerics/data_out_rotation.h
include/deal.II/numerics/data_out_stack.h
source/numerics/data_out_faces.cc
source/numerics/data_out_rotation.cc
source/numerics/data_out_stack.cc

index affac017a72413804ca6d7a8143a5089fe3e7565..dab3f2df5aef1abced3e29a17fab234fb168a714 100644 (file)
@@ -2456,8 +2456,8 @@ namespace DataOutBase
  * This class is thought as a base class to classes actually generating data
  * for output. It has two abstract virtual functions, get_patches() and
  * get_dataset_names() produce the data which is actually needed. These are
- * the only functions that need to be overloaded by a derived class.  In
- * additional to that, it has a function for each output format supported by
+ * the only functions that need to be overloaded by a derived class. In
+ * addition to that, it has a function for each output format supported by
  * the underlying base class which gets the output data using these two
  * virtual functions and passes them to the raw output functions.
  *
index 2b941aadd8f70cb08c52473b580da62bd7bcf4b7..840e3269eb659035085885dc70a8221f62cc4cf1 100644 (file)
@@ -100,28 +100,30 @@ namespace internal
  * applications certainly exist, for which the author is not imaginative
  * enough.
  *
- * @pre This class only makes sense if the first template argument,
- * <code>dim</code> equals the dimension of the DoFHandler type given as the
- * second template argument, i.e., if <code>dim ==
- * DoFHandlerType::dimension</code>. This redundancy is a historical relic
- * from the time where the library had only a single DoFHandler class and this
- * class consequently only a single template argument.
- *
  * @todo Reimplement this whole class using actual FEFaceValues and
  * MeshWorker.
  *
  * @ingroup output
  */
 template <int dim, int spacedim = dim>
-class DataOutFaces : public DataOut_DoFData<dim, dim - 1, spacedim, dim>
+class DataOutFaces : public DataOut_DoFData<dim, dim - 1, spacedim, spacedim>
 {
+  static_assert(dim == spacedim, "Not implemented for dim != spacedim.");
+
 public:
+  /**
+   * Dimension parameters for the patches.
+   */
+  static constexpr int patch_dim      = dim - 1;
+  static constexpr int patch_spacedim = spacedim;
+
   /**
    * Alias to the iterator type of the dof handler class under
    * consideration.
    */
   using cell_iterator =
-    typename DataOut_DoFData<dim, dim - 1, spacedim, dim>::cell_iterator;
+    typename DataOut_DoFData<dim, patch_dim, spacedim, patch_spacedim>::
+      cell_iterator;
 
   /**
    * Constructor determining whether a surface mesh (default) or the whole
@@ -167,8 +169,8 @@ public:
    * hp::MappingCollection in case of a DoFHandler with hp-capabilities.
    */
   virtual void
-  build_patches(const Mapping<dim> &mapping,
-                const unsigned int  n_subdivisions = 0);
+  build_patches(const Mapping<dim, spacedim> &mapping,
+                const unsigned int            n_subdivisions = 0);
 
   /**
    * Declare a way to describe a face which we would like to generate output
@@ -229,9 +231,9 @@ private:
    */
   void
   build_one_patch(
-    const FaceDescriptor *                                        cell_and_face,
-    internal::DataOutFacesImplementation::ParallelData<dim, dim> &data,
-    DataOutBase::Patch<dim - 1, spacedim> &                       patch);
+    const FaceDescriptor *cell_and_face,
+    internal::DataOutFacesImplementation::ParallelData<dim, spacedim> &data,
+    DataOutBase::Patch<patch_dim, patch_spacedim> &                    patch);
 };
 
 namespace Legacy
index ed8e907d6dd8c8b13762a88a5ae7a8200393078d..7416fc19dc4e8eb9bf47ed14b94eb363cb5672f6 100644 (file)
@@ -109,25 +109,28 @@ namespace internal
  * It is in the responsibility of the user to make sure that the radial
  * variable attains only non-negative values.
  *
- * @pre This class only makes sense if the first template argument,
- * <code>dim</code> equals the dimension of the DoFHandler type given as the
- * second template argument, i.e., if <code>dim ==
- * DoFHandlerType::dimension</code>. This redundancy is a historical relic
- * from the time where the library had only a single DoFHandler class and this
- * class consequently only a single template argument.
- *
  * @ingroup output
  */
 template <int dim, int spacedim = dim>
-class DataOutRotation : public DataOut_DoFData<dim, dim + 1, spacedim, dim + 1>
+class DataOutRotation
+  : public DataOut_DoFData<dim, dim + 1, spacedim, spacedim + 1>
 {
+  static_assert(dim == spacedim, "Not implemented for dim != spacedim.");
+
 public:
+  /**
+   * Dimension parameters for the patches.
+   */
+  static constexpr int patch_dim      = dim + 1;
+  static constexpr int patch_spacedim = spacedim + 1;
+
   /**
    * Typedef to the iterator type of the dof handler class under
    * consideration.
    */
   using cell_iterator =
-    typename DataOut_DoFData<dim, dim + 1, spacedim, dim + 1>::cell_iterator;
+    typename DataOut_DoFData<dim, patch_dim, spacedim, patch_spacedim>::
+      cell_iterator;
 
   /**
    * This is the central function of this class since it builds the list of
@@ -196,7 +199,7 @@ private:
   build_one_patch(
     const cell_iterator *                                                 cell,
     internal::DataOutRotationImplementation::ParallelData<dim, spacedim> &data,
-    std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> &my_patches);
+    std::vector<DataOutBase::Patch<patch_dim, patch_spacedim>> &my_patches);
 };
 
 namespace Legacy
index dbc01d1b379d7b0dc3c02ca29ea3bec88b0042ec..4a8d8c7ec46a57dea6c49bb4e3f455adb1db83e3 100644 (file)
@@ -131,9 +131,18 @@ public:
  * @ingroup output
  */
 template <int dim, int spacedim>
-class DataOutStack<dim, spacedim, void> : public DataOutInterface<dim + 1>
+class DataOutStack<dim, spacedim, void>
+  : public DataOutInterface<dim + 1, spacedim + 1>
 {
+  static_assert(dim == spacedim, "Not implemented for dim != spacedim.");
+
 public:
+  /**
+   * Dimension parameters for the patches.
+   */
+  static constexpr int patch_dim      = dim + 1;
+  static constexpr int patch_spacedim = spacedim + 1;
+
   /**
    * Data type declaring the two types of vectors which are used in this
    * class.
@@ -336,7 +345,7 @@ private:
   /**
    * List of patches of all past and present parameter value data sets.
    */
-  std::vector<dealii::DataOutBase::Patch<dim + 1, dim + 1>> patches;
+  std::vector<dealii::DataOutBase::Patch<patch_dim, patch_spacedim>> patches;
 
   /**
    * Structure holding data vectors (cell and dof data) for the present
@@ -377,7 +386,8 @@ private:
    * data in the form of Patch structures (declared in the base class
    * DataOutBase) to the actual output function.
    */
-  virtual const std::vector<dealii::DataOutBase::Patch<dim + 1, dim + 1>> &
+  virtual const std::vector<
+    dealii::DataOutBase::Patch<patch_dim, patch_spacedim>> &
   get_patches() const override;
 
 
index 1844379514980140095809a21d453a412109ed3d..fbe204cc7dc667ab20ce53f32412057ff0a51097 100644 (file)
@@ -69,8 +69,13 @@ namespace internal
     template <int dim, int spacedim>
     void
     append_patch_to_list(
-      const DataOutBase::Patch<dim - 1, spacedim> &       patch,
-      std::vector<DataOutBase::Patch<dim - 1, spacedim>> &patches)
+      const DataOutBase::Patch<DataOutFaces<dim, spacedim>::patch_dim,
+                               DataOutFaces<dim, spacedim>::patch_spacedim>
+        &patch,
+      std::vector<
+        DataOutBase::Patch<DataOutFaces<dim, spacedim>::patch_dim,
+                           DataOutFaces<dim, spacedim>::patch_spacedim>>
+        &patches)
     {
       patches.push_back(patch);
       patches.back().patch_index = patches.size() - 1;
@@ -90,9 +95,9 @@ DataOutFaces<dim, spacedim>::DataOutFaces(const bool so)
 template <int dim, int spacedim>
 void
 DataOutFaces<dim, spacedim>::build_one_patch(
-  const FaceDescriptor *                                        cell_and_face,
-  internal::DataOutFacesImplementation::ParallelData<dim, dim> &data,
-  DataOutBase::Patch<dim - 1, spacedim> &                       patch)
+  const FaceDescriptor *cell_and_face,
+  internal::DataOutFacesImplementation::ParallelData<dim, spacedim> &data,
+  DataOutBase::Patch<patch_dim, patch_spacedim> &                    patch)
 {
   Assert(cell_and_face->first->is_locally_owned(), ExcNotImplemented());
 
@@ -314,15 +319,16 @@ template <int dim, int spacedim>
 void
 DataOutFaces<dim, spacedim>::build_patches(const unsigned int n_subdivisions_)
 {
-  build_patches(StaticMappingQ1<dim>::mapping, n_subdivisions_);
+  build_patches(StaticMappingQ1<dim, spacedim>::mapping, n_subdivisions_);
 }
 
 
 
 template <int dim, int spacedim>
 void
-DataOutFaces<dim, spacedim>::build_patches(const Mapping<dim> &mapping,
-                                           const unsigned int  n_subdivisions_)
+DataOutFaces<dim, spacedim>::build_patches(
+  const Mapping<dim, spacedim> &mapping,
+  const unsigned int            n_subdivisions_)
 {
   const unsigned int n_subdivisions =
     (n_subdivisions_ != 0) ? n_subdivisions_ : this->default_subdivisions;
@@ -379,21 +385,23 @@ DataOutFaces<dim, spacedim>::build_patches(const Mapping<dim> &mapping,
     mapping,
     this->get_fes(),
     update_flags);
-  DataOutBase::Patch<dim - 1, spacedim> sample_patch;
+  DataOutBase::Patch<patch_dim, patch_spacedim> sample_patch;
   sample_patch.n_subdivisions = n_subdivisions;
   sample_patch.data.reinit(n_datasets,
-                           Utilities::fixed_power<dim - 1>(n_subdivisions + 1));
+                           Utilities::fixed_power<patch_dim>(n_subdivisions +
+                                                             1));
 
   // now build the patches in parallel
   WorkStream::run(
     all_faces.data(),
     all_faces.data() + all_faces.size(),
-    [this](const FaceDescriptor *cell_and_face,
-           internal::DataOutFacesImplementation::ParallelData<dim, dim> &data,
-           DataOutBase::Patch<dim - 1, spacedim> &patch) {
+    [this](
+      const FaceDescriptor *cell_and_face,
+      internal::DataOutFacesImplementation::ParallelData<dim, spacedim> &data,
+      DataOutBase::Patch<patch_dim, patch_spacedim> &patch) {
       this->build_one_patch(cell_and_face, data, patch);
     },
-    [this](const DataOutBase::Patch<dim - 1, spacedim> &patch) {
+    [this](const DataOutBase::Patch<patch_dim, patch_spacedim> &patch) {
       internal::DataOutFacesImplementation::append_patch_to_list<dim, spacedim>(
         patch, this->patches);
     },
index f75b1e95a0d185f41ad8c1a57c08a3afed068c2b..b3609418a6deedfc2e7ec2bf0812f58d5c35bd3e 100644 (file)
@@ -84,8 +84,14 @@ namespace internal
     template <int dim, int spacedim>
     void
     append_patch_to_list(
-      const std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> &new_patches,
-      std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> &      patches)
+      const std::vector<
+        DataOutBase::Patch<DataOutRotation<dim, spacedim>::patch_dim,
+                           DataOutRotation<dim, spacedim>::patch_spacedim>>
+        &new_patches,
+      std::vector<
+        DataOutBase::Patch<DataOutRotation<dim, spacedim>::patch_dim,
+                           DataOutRotation<dim, spacedim>::patch_spacedim>>
+        &patches)
     {
       for (unsigned int i = 0; i < new_patches.size(); ++i)
         {
@@ -103,7 +109,7 @@ void
 DataOutRotation<dim, spacedim>::build_one_patch(
   const cell_iterator *                                                 cell,
   internal::DataOutRotationImplementation::ParallelData<dim, spacedim> &data,
-  std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> &my_patches)
+  std::vector<DataOutBase::Patch<patch_dim, patch_spacedim>> &my_patches)
 {
   if (dim == 3)
     {
@@ -506,33 +512,34 @@ DataOutRotation<dim, spacedim>::build_patches(
       n_postprocessor_outputs[dataset] = 0;
 
   internal::DataOutRotationImplementation::ParallelData<dim, spacedim>
-                                                         thread_data(n_datasets,
+                                                             thread_data(n_datasets,
                 n_subdivisions,
                 n_patches_per_circle,
                 n_postprocessor_outputs,
                 StaticMappingQ1<dim, spacedim>::mapping,
                 this->get_fes(),
                 update_flags);
-  std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> new_patches(
+  std::vector<DataOutBase::Patch<patch_dim, patch_spacedim>> new_patches(
     n_patches_per_circle);
   for (unsigned int i = 0; i < new_patches.size(); ++i)
     {
       new_patches[i].n_subdivisions = n_subdivisions;
       new_patches[i].data.reinit(
-        n_datasets, Utilities::fixed_power<dim + 1>(n_subdivisions + 1));
+        n_datasets, Utilities::fixed_power<patch_dim>(n_subdivisions + 1));
     }
 
   // now build the patches in parallel
   WorkStream::run(
     all_cells.data(),
     all_cells.data() + all_cells.size(),
-    [this](const cell_iterator *cell,
-           internal::DataOutRotationImplementation::ParallelData<dim, spacedim>
-             &                                                     data,
-           std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>> &my_patches) {
+    [this](
+      const cell_iterator *cell,
+      internal::DataOutRotationImplementation::ParallelData<dim, spacedim>
+        &                                                         data,
+      std::vector<DataOutBase::Patch<patch_dim, patch_spacedim>> &my_patches) {
       this->build_one_patch(cell, data, my_patches);
     },
-    [this](const std::vector<DataOutBase::Patch<dim + 1, spacedim + 1>>
+    [this](const std::vector<DataOutBase::Patch<patch_dim, patch_spacedim>>
              &new_patches) {
       internal::DataOutRotationImplementation::append_patch_to_list<dim,
                                                                     spacedim>(
index f28044fedb8379b6636b81c734ac17424f1da656..f2011bb21f87f1a74f1e69b4cc5fa67081ca8984 100644 (file)
@@ -308,15 +308,16 @@ DataOutStack<dim, spacedim, void>::build_patches(
   // patch with n_q_points (in the plane
   // of the cells) times n_subdivisions+1 (in
   // the time direction) points
-  dealii::DataOutBase::Patch<dim + 1, dim + 1> default_patch;
+  dealii::DataOutBase::Patch<patch_dim, patch_spacedim> default_patch;
   default_patch.n_subdivisions = n_subdivisions;
   default_patch.data.reinit(n_datasets, n_q_points * (n_subdivisions + 1));
   patches.insert(patches.end(), n_patches, default_patch);
 
   // now loop over all cells and
   // actually create the patches
-  typename std::vector<dealii::DataOutBase::Patch<dim + 1, dim + 1>>::iterator
-               patch       = patches.begin() + (patches.size() - n_patches);
+  typename std::vector<
+    dealii::DataOutBase::Patch<patch_dim, patch_spacedim>>::iterator patch =
+    patches.begin() + (patches.size() - n_patches);
   unsigned int cell_number = 0;
   for (typename DoFHandler<dim, spacedim>::active_cell_iterator cell =
          dof_handler->begin_active();
@@ -451,7 +452,7 @@ template <int dim, int spacedim>
 std::size_t
 DataOutStack<dim, spacedim, void>::memory_consumption() const
 {
-  return (DataOutInterface<dim + 1>::memory_consumption() +
+  return (DataOutInterface<patch_dim, patch_spacedim>::memory_consumption() +
           MemoryConsumption::memory_consumption(parameter) +
           MemoryConsumption::memory_consumption(parameter_step) +
           MemoryConsumption::memory_consumption(dof_handler) +
@@ -463,8 +464,11 @@ DataOutStack<dim, spacedim, void>::memory_consumption() const
 
 
 template <int dim, int spacedim>
-const std::vector<dealii::DataOutBase::Patch<dim + 1, dim + 1>> &
-DataOutStack<dim, spacedim, void>::get_patches() const
+const std::vector<
+  dealii::DataOutBase::Patch<DataOutStack<dim, spacedim, void>::patch_dim,
+                             DataOutStack<dim, spacedim, void>::patch_spacedim>>
+  &
+  DataOutStack<dim, spacedim, void>::get_patches() const
 {
   return patches;
 }

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.