* 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,
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,
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,
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.
*/
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;
*/
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.
*/
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";
, num_cells(numbers::invalid_unsigned_int)
, dimension(numbers::invalid_unsigned_int)
, space_dimension(numbers::invalid_unsigned_int)
+ , cell_type()
{}
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)
{}
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)
{}
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)
, 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"));
+ }
+ }
+}
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 "";
}
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";
}
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";