]> https://gitweb.dealii.org/ - dealii.git/commitdiff
XDMF changes
authorTimo Heister <timo.heister@gmail.com>
Fri, 24 Jun 2022 15:10:53 +0000 (11:10 -0400)
committerTimo Heister <timo.heister@gmail.com>
Tue, 28 Jun 2022 15:27:52 +0000 (11:27 -0400)
include/deal.II/base/data_out_base.h
source/base/data_out_base.cc

index dc35266b27f4a4af668e6f0f755f0b921ec59109..7dc04473929325ef72d6630dd0922fb514135daa 100644 (file)
@@ -3322,6 +3322,18 @@ public:
    * cases where <code>solution_filename == mesh_filename</code>, and
    * <code>dim==spacedim</code>.
    */
+  XDMFEntry(const std::string &  filename,
+            const double         time,
+            const std::uint64_t  nodes,
+            const std::uint64_t  cells,
+            const unsigned int   dim,
+            const ReferenceCell &cell_type);
+
+  /**
+   * Deprecated constructor.
+   *
+   * @deprecated Use the constructor that additionally takes a ReferenceCell.
+   */
   XDMFEntry(const std::string & filename,
             const double        time,
             const std::uint64_t nodes,
@@ -3329,8 +3341,9 @@ public:
             const unsigned int  dim);
 
   /**
-   * Simplified constructor that calls the complete constructor for
-   * cases where <code>dim==spacedim</code>.
+   * Deprecated constructor.
+   *
+   * @deprecated Use the constructor that additionally takes a ReferenceCell.
    */
   XDMFEntry(const std::string & mesh_filename,
             const std::string & solution_filename,
@@ -3340,8 +3353,23 @@ public:
             const unsigned int  dim);
 
   /**
-   * Constructor that sets all members to provided parameters.
+   * Simplified constructor that calls the complete constructor for
+   * cases where <code>dim==spacedim</code>.
+   */
+  XDMFEntry(const std::string &  mesh_filename,
+            const std::string &  solution_filename,
+            const double         time,
+            const std::uint64_t  nodes,
+            const std::uint64_t  cells,
+            const unsigned int   dim,
+            const ReferenceCell &cell_type);
+
+  /**
+   * Deprecated constructor.
+   *
+   * @deprecated Use the constructor that additionally takes a ReferenceCell.
    */
+  DEAL_II_DEPRECATED
   XDMFEntry(const std::string & mesh_filename,
             const std::string & solution_filename,
             const double        time,
@@ -3350,6 +3378,18 @@ public:
             const unsigned int  dim,
             const unsigned int  spacedim);
 
+  /**
+   * Constructor that sets all members to provided parameters.
+   */
+  XDMFEntry(const std::string &  mesh_filename,
+            const std::string &  solution_filename,
+            const double         time,
+            const std::uint64_t  nodes,
+            const std::uint64_t  cells,
+            const unsigned int   dim,
+            const unsigned int   spacedim,
+            const ReferenceCell &cell_type);
+
   /**
    * Record an attribute and associated dimensionality.
    */
@@ -3366,24 +3406,23 @@ public:
   serialize(Archive &ar, const unsigned int /*version*/)
   {
     ar &valid &h5_sol_filename &h5_mesh_filename &entry_time &num_nodes
-      &num_cells &dimension &space_dimension &attribute_dims;
+      &num_cells &dimension &space_dimension &cell_type &attribute_dims;
   }
 
   /**
    * Get the XDMF content associated with this entry.
    * If the entry is not valid, this returns an empty string.
-   *
-   * @deprecated Use the overload taking an `unsigned int` and a
-   * `const ReferenceCell &` instead.
    */
-  DEAL_II_DEPRECATED
   std::string
   get_xdmf_content(const unsigned int indent_level) const;
 
   /**
    * Get the XDMF content associated with this entry.
    * If the entry is not valid, this returns an empty string.
+   *
+   * @deprecated Use the other function instead.
    */
+  DEAL_II_DEPRECATED
   std::string
   get_xdmf_content(const unsigned int   indent_level,
                    const ReferenceCell &reference_cell) const;
@@ -3430,6 +3469,12 @@ private:
    */
   unsigned int space_dimension;
 
+  /**
+   * The type of cell in deal.II language. We currently only support
+   * xdmf entries where all cells have the same type.
+   */
+  ReferenceCell cell_type;
+
   /**
    * The attributes associated with this entry and their dimension.
    */
index 41bb8954cf101549a56a16575a3f233839c020e1..eed58f7157ad17f0f86f2c3bb91c68b5e4051b5f 100644 (file)
@@ -8037,20 +8037,9 @@ DataOutInterface<dim, spacedim>::write_xdmf_file(
       xdmf_file
         << "    <Grid Name=\"CellTime\" GridType=\"Collection\" CollectionType=\"Temporal\">\n";
 
-      // Write out all the entries indented
-      const auto &patches = get_patches();
-      Assert(patches.size() > 0, DataOutBase::ExcNoPatches());
-
-      // We currently don't support writing mixed meshes:
-#ifdef DEBUG
-      for (const auto &patch : patches)
-        Assert(patch.reference_cell == patches[0].reference_cell,
-               ExcNotImplemented());
-#endif
-
       for (const auto &entry : entries)
         {
-          xdmf_file << entry.get_xdmf_content(3, patches[0].reference_cell);
+          xdmf_file << entry.get_xdmf_content(3);
         }
 
       xdmf_file << "    </Grid>\n";
@@ -9129,6 +9118,7 @@ XDMFEntry::XDMFEntry()
   , num_cells(numbers::invalid_unsigned_int)
   , dimension(numbers::invalid_unsigned_int)
   , space_dimension(numbers::invalid_unsigned_int)
+  , cell_type()
 {}
 
 
@@ -9138,7 +9128,16 @@ XDMFEntry::XDMFEntry(const std::string & filename,
                      const std::uint64_t nodes,
                      const std::uint64_t cells,
                      const unsigned int  dim)
-  : XDMFEntry(filename, filename, time, nodes, cells, dim, dim)
+  : XDMFEntry(filename, filename, time, nodes, cells, dim, dim, ReferenceCell())
+{}
+
+XDMFEntry::XDMFEntry(const std::string &  filename,
+                     const double         time,
+                     const std::uint64_t  nodes,
+                     const std::uint64_t  cells,
+                     const unsigned int   dim,
+                     const ReferenceCell &cell_type)
+  : XDMFEntry(filename, filename, time, nodes, cells, dim, dim, cell_type)
 {}
 
 
@@ -9149,7 +9148,33 @@ XDMFEntry::XDMFEntry(const std::string & mesh_filename,
                      const std::uint64_t nodes,
                      const std::uint64_t cells,
                      const unsigned int  dim)
-  : XDMFEntry(mesh_filename, solution_filename, time, nodes, cells, dim, dim)
+  : XDMFEntry(mesh_filename,
+              solution_filename,
+              time,
+              nodes,
+              cells,
+              dim,
+              dim,
+              ReferenceCell())
+{}
+
+
+
+XDMFEntry::XDMFEntry(const std::string &  mesh_filename,
+                     const std::string &  solution_filename,
+                     const double         time,
+                     const std::uint64_t  nodes,
+                     const std::uint64_t  cells,
+                     const unsigned int   dim,
+                     const ReferenceCell &cell_type)
+  : XDMFEntry(mesh_filename,
+              solution_filename,
+              time,
+              nodes,
+              cells,
+              dim,
+              dim,
+              cell_type)
 {}
 
 
@@ -9161,6 +9186,26 @@ XDMFEntry::XDMFEntry(const std::string & mesh_filename,
                      const std::uint64_t cells,
                      const unsigned int  dim,
                      const unsigned int  spacedim)
+  : XDMFEntry(mesh_filename,
+              solution_filename,
+              time,
+              nodes,
+              cells,
+              dim,
+              spacedim,
+              ReferenceCell())
+{}
+
+
+
+XDMFEntry::XDMFEntry(const std::string &  mesh_filename,
+                     const std::string &  solution_filename,
+                     const double         time,
+                     const std::uint64_t  nodes,
+                     const std::uint64_t  cells,
+                     const unsigned int   dim,
+                     const unsigned int   spacedim,
+                     const ReferenceCell &cell_type_)
   : valid(true)
   , h5_sol_filename(solution_filename)
   , h5_mesh_filename(mesh_filename)
@@ -9169,7 +9214,29 @@ XDMFEntry::XDMFEntry(const std::string & mesh_filename,
   , num_cells(cells)
   , dimension(dim)
   , space_dimension(spacedim)
-{}
+  , cell_type(cell_type_)
+{
+  if (cell_type == ReferenceCells::Invalid)
+    {
+      switch (dimension)
+        {
+          case 0:
+            cell_type = ReferenceCells::get_hypercube<0>();
+            break;
+          case 1:
+            cell_type = ReferenceCells::get_hypercube<1>();
+            break;
+          case 2:
+            cell_type = ReferenceCells::get_hypercube<2>();
+            break;
+          case 3:
+            cell_type = ReferenceCells::get_hypercube<3>();
+            break;
+          default:
+            AssertThrow(false, ExcMessage("Invalid dimension"));
+        }
+    }
+}
 
 
 
@@ -9200,34 +9267,20 @@ namespace
 
 
 std::string
-XDMFEntry::get_xdmf_content(const unsigned int indent_level) const
+XDMFEntry::get_xdmf_content(const unsigned int   indent_level,
+                            const ReferenceCell &reference_cell) const
 {
-  switch (dimension)
-    {
-      case 0:
-        return get_xdmf_content(indent_level,
-                                ReferenceCells::get_hypercube<0>());
-      case 1:
-        return get_xdmf_content(indent_level,
-                                ReferenceCells::get_hypercube<1>());
-      case 2:
-        return get_xdmf_content(indent_level,
-                                ReferenceCells::get_hypercube<2>());
-      case 3:
-        return get_xdmf_content(indent_level,
-                                ReferenceCells::get_hypercube<3>());
-      default:
-        Assert(false, ExcNotImplemented());
-    }
-
-  return "";
+  // We now store the type of cell in the XDMFEntry:
+  (void)reference_cell;
+  Assert(cell_type == reference_cell,
+         ExcMessage("Incorrect ReferenceCell type passed in."));
+  return get_xdmf_content(indent_level);
 }
 
 
 
 std::string
-XDMFEntry::get_xdmf_content(const unsigned int   indent_level,
-                            const ReferenceCell &reference_cell) const
+XDMFEntry::get_xdmf_content(const unsigned int indent_level) const
 {
   if (!valid)
     return "";
@@ -9260,26 +9313,26 @@ XDMFEntry::get_xdmf_content(const unsigned int   indent_level,
         }
       else if (dimension == 2)
         {
-          Assert(reference_cell == ReferenceCells::Quadrilateral ||
-                   reference_cell == ReferenceCells::Triangle,
+          Assert(cell_type == ReferenceCells::Quadrilateral ||
+                   cell_type == ReferenceCells::Triangle,
                  ExcNotImplemented());
 
-          if (reference_cell == ReferenceCells::Quadrilateral)
+          if (cell_type == ReferenceCells::Quadrilateral)
             {
               ss << "Quadrilateral";
             }
-          else // if (reference_cell == ReferenceCells::Triangle)
+          else // if (cell_type == ReferenceCells::Triangle)
             {
               ss << "Triangle";
             }
         }
       else if (dimension == 3)
         {
-          Assert(reference_cell == ReferenceCells::Hexahedron ||
-                   reference_cell == ReferenceCells::Tetrahedron,
+          Assert(cell_type == ReferenceCells::Hexahedron ||
+                   cell_type == ReferenceCells::Tetrahedron,
                  ExcNotImplemented());
 
-          if (reference_cell == ReferenceCells::Hexahedron)
+          if (cell_type == ReferenceCells::Hexahedron)
             {
               ss << "Hexahedron";
             }
@@ -9299,7 +9352,7 @@ XDMFEntry::get_xdmf_content(const unsigned int   indent_level,
         ss << "\">\n";
 
       ss << indent(indent_level + 2) << "<DataItem Dimensions=\"" << num_cells
-         << " " << reference_cell.n_vertices()
+         << " " << cell_type.n_vertices()
          << "\" NumberType=\"UInt\" Format=\"HDF\">\n";
 
       ss << indent(indent_level + 3) << h5_mesh_filename << ":/cells\n";

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.