Skip to content

Commit

Permalink
Add DataOutBase::CompressionLevel::ascii
Browse files Browse the repository at this point in the history
  • Loading branch information
peterrum committed Sep 18, 2022
1 parent ab55980 commit 4ed6b6f
Show file tree
Hide file tree
Showing 4 changed files with 696 additions and 75 deletions.
6 changes: 5 additions & 1 deletion include/deal.II/base/data_out_base.h
Original file line number Diff line number Diff line change
Expand Up @@ -243,7 +243,11 @@ namespace DataOutBase
* Use the default compression algorithm. This is a compromise between
* speed and file size.
*/
default_compression
default_compression,
/**
* Output as ascii if available.
*/
ascii
};


Expand Down
190 changes: 117 additions & 73 deletions source/base/data_out_base.cc
Original file line number Diff line number Diff line change
Expand Up @@ -1922,30 +1922,38 @@ namespace
VtuStream::write_point(const unsigned int, const Point<dim> &p)
{
#ifdef DEAL_II_WITH_ZLIB
// if we want to compress, then first collect all the data in an array
for (unsigned int i = 0; i < dim; ++i)
vertices.push_back(p[i]);
for (unsigned int i = dim; i < 3; ++i)
vertices.push_back(0);
#else
// write out coordinates
stream << p;
// fill with zeroes
for (unsigned int i = dim; i < 3; ++i)
stream << " 0";
stream << '\n';
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
// if we want to compress, then first collect all the data in an array
for (unsigned int i = 0; i < dim; ++i)
vertices.push_back(p[i]);
for (unsigned int i = dim; i < 3; ++i)
vertices.push_back(0);
}
else
#endif
{
// write out coordinates
stream << p;
// fill with zeroes
for (unsigned int i = dim; i < 3; ++i)
stream << " 0";
stream << '\n';
}
}


void
VtuStream::flush_points()
{
#ifdef DEAL_II_WITH_ZLIB
// compress the data we have in memory and write them to the stream. then
// release the data
*this << vertices << '\n';
vertices.clear();
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
// compress the data we have in memory and write them to the stream.
// then release the data
*this << vertices << '\n';
vertices.clear();
}
#endif
}

Expand All @@ -1959,40 +1967,46 @@ namespace
unsigned int d3)
{
#ifdef DEAL_II_WITH_ZLIB
cells.push_back(start);
if (dim >= 1)
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
cells.push_back(start + d1);
if (dim >= 2)
cells.push_back(start);
if (dim >= 1)
{
cells.push_back(start + d2 + d1);
cells.push_back(start + d2);
if (dim >= 3)
cells.push_back(start + d1);
if (dim >= 2)
{
cells.push_back(start + d3);
cells.push_back(start + d3 + d1);
cells.push_back(start + d3 + d2 + d1);
cells.push_back(start + d3 + d2);
cells.push_back(start + d2 + d1);
cells.push_back(start + d2);
if (dim >= 3)
{
cells.push_back(start + d3);
cells.push_back(start + d3 + d1);
cells.push_back(start + d3 + d2 + d1);
cells.push_back(start + d3 + d2);
}
}
}
}
#else
stream << start;
if (dim >= 1)
else
#endif
{
stream << '\t' << start + d1;
if (dim >= 2)
stream << start;
if (dim >= 1)
{
stream << '\t' << start + d2 + d1 << '\t' << start + d2;
if (dim >= 3)
stream << '\t' << start + d1;
if (dim >= 2)
{
stream << '\t' << start + d3 << '\t' << start + d3 + d1 << '\t'
<< start + d3 + d2 + d1 << '\t' << start + d3 + d2;
stream << '\t' << start + d2 + d1 << '\t' << start + d2;
if (dim >= 3)
{
stream << '\t' << start + d3 << '\t' << start + d3 + d1
<< '\t' << start + d3 + d2 + d1 << '\t'
<< start + d3 + d2;
}
}
}
stream << '\n';
}
stream << '\n';
#endif
}

void
Expand All @@ -2006,16 +2020,22 @@ namespace
static const std::array<unsigned int, 5> table = {{0, 1, 3, 2, 4}};

#ifdef DEAL_II_WITH_ZLIB
for (unsigned int i = 0; i < n_points; ++i)
cells.push_back(
start + (reference_cell == ReferenceCells::Pyramid ? table[i] : i));
#else
for (unsigned int i = 0; i < n_points; ++i)
stream << '\t'
<< start +
(reference_cell == ReferenceCells::Pyramid ? table[i] : i);
stream << '\n';
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
for (unsigned int i = 0; i < n_points; ++i)
cells.push_back(
start + (reference_cell == ReferenceCells::Pyramid ? table[i] : i));
}
else
#endif
{
for (unsigned int i = 0; i < n_points; ++i)
stream << '\t'
<< start + (reference_cell == ReferenceCells::Pyramid ?
table[i] :
i);
stream << '\n';
}
}

template <int dim>
Expand All @@ -2025,23 +2045,31 @@ namespace
const std::vector<unsigned> &connectivity)
{
#ifdef DEAL_II_WITH_ZLIB
for (const auto &c : connectivity)
cells.push_back(start + c);
#else
for (const auto &c : connectivity)
stream << '\t' << start + c;
stream << '\n';
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
for (const auto &c : connectivity)
cells.push_back(start + c);
}
else
#endif
{
for (const auto &c : connectivity)
stream << '\t' << start + c;
stream << '\n';
}
}

void
VtuStream::flush_cells()
{
#ifdef DEAL_II_WITH_ZLIB
// compress the data we have in memory and write them to the stream. then
// release the data
*this << cells << '\n';
cells.clear();
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
// compress the data we have in memory and write them to the stream.
// then release the data
*this << cells << '\n';
cells.clear();
}
#endif
}

Expand All @@ -2051,13 +2079,18 @@ namespace
VtuStream::operator<<(const std::vector<T> &data)
{
#ifdef DEAL_II_WITH_ZLIB
// compress the data we have in memory and write them to the stream. then
// release the data
write_compressed_block(data, flags, stream);
#else
for (unsigned int i = 0; i < data.size(); ++i)
stream << data[i] << ' ';
if (flags.compression_level != DataOutBase::CompressionLevel::ascii)
{
// compress the data we have in memory and write them to the stream.
// then release the data
write_compressed_block(data, flags, stream);
}
else
#endif
{
for (unsigned int i = 0; i < data.size(); ++i)
stream << data[i] << ' ';
}

return stream;
}
Expand Down Expand Up @@ -5804,7 +5837,8 @@ namespace DataOutBase
out << "\n-->\n";
out << "<VTKFile type=\"UnstructuredGrid\" version=\"0.1\"";
#ifdef DEAL_II_WITH_ZLIB
out << " compressor=\"vtkZLibDataCompressor\"";
if (flags.compression_level != CompressionLevel::ascii)
out << " compressor=\"vtkZLibDataCompressor\"";
#endif
#ifdef DEAL_II_WORDS_BIGENDIAN
out << " byte_order=\"BigEndian\"";
Expand Down Expand Up @@ -5990,10 +6024,9 @@ namespace DataOutBase

const char *ascii_or_binary =
#ifdef DEAL_II_WITH_ZLIB
"binary";
#else
"ascii";
(flags.compression_level != CompressionLevel::ascii) ? "binary" :
#endif
"ascii";


// first count the number of cells and cells for later use
Expand Down Expand Up @@ -6087,14 +6120,19 @@ namespace DataOutBase

// this should compress well :-)
#ifdef DEAL_II_WITH_ZLIB
std::vector<uint8_t> cell_types_uint8_t(cell_types.size());
for (unsigned int i = 0; i < cell_types.size(); ++i)
cell_types_uint8_t[i] = static_cast<std::uint8_t>(cell_types[i]);
if (flags.compression_level != CompressionLevel::ascii)
{
std::vector<uint8_t> cell_types_uint8_t(cell_types.size());
for (unsigned int i = 0; i < cell_types.size(); ++i)
cell_types_uint8_t[i] = static_cast<std::uint8_t>(cell_types[i]);

vtu_out << cell_types_uint8_t;
#else
vtu_out << cell_types;
vtu_out << cell_types_uint8_t;
}
else
#endif
{
vtu_out << cell_types;
}

out << '\n';
out << " </DataArray>\n";
Expand Down Expand Up @@ -7648,6 +7686,8 @@ namespace DataOutBase
{
boost::iostreams::filtering_ostream f;

AssertThrow(compression != CompressionLevel::ascii, ExcNotImplemented());

if (compression != CompressionLevel::no_compression)
# ifdef DEAL_II_WITH_ZLIB
f.push(boost::iostreams::zlib_compressor(
Expand Down Expand Up @@ -9341,6 +9381,10 @@ DataOutReader<dim, spacedim>::read_whole_parallel_file(std::istream &in)
std::vector<char> temp_buffer(chunk_sizes[n]);
in.read(temp_buffer.data(), chunk_sizes[n]);

AssertThrow(static_cast<DataOutBase::CompressionLevel>(
header.compression) != DataOutBase::CompressionLevel::ascii,
ExcNotImplemented());

boost::iostreams::filtering_istreambuf f;
if (static_cast<DataOutBase::CompressionLevel>(header.compression) !=
DataOutBase::CompressionLevel::no_compression)
Expand Down
5 changes: 4 additions & 1 deletion tests/data_out/data_out_base_vtu_compression_levels.cc
Original file line number Diff line number Diff line change
Expand Up @@ -56,7 +56,7 @@ template <int dim, int spacedim>
void
check_all(std::ostream &log)
{
for (unsigned int i = 0; i < 4; ++i)
for (unsigned int i = 0; i < 5; ++i)
{
DataOutBase::VtkFlags flags;
switch (i)
Expand All @@ -76,6 +76,9 @@ check_all(std::ostream &log)
flags.compression_level =
DataOutBase::CompressionLevel::default_compression;
break;
case (4):
flags.compression_level = DataOutBase::CompressionLevel::ascii;
break;
default:
Assert(false, ExcInternalError());
}
Expand Down

0 comments on commit 4ed6b6f

Please sign in to comment.