vtkimagedatawriter.impl.hh 7.02 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
#pragma once

#include <iomanip>
#include <iostream>
#include <iterator>
#include <fstream>
#include <sstream>
#include <string>

#include <dune/geometry/referenceelements.hh>
#include <dune/geometry/type.hh>

#include "utility/enum.hh"
#include "utility/filesystem.hh"
#include "utility/string.hh"

namespace Dune { namespace experimental {

template <class GV, class DC>
void VtkImageDataWriter<GV,DC>
  ::writeSerialFile (std::string const& filename) const
{
  std::ofstream out(filename, std::ios_base::ate | std::ios::binary);
  if (format_ == Vtk::ASCII) {
    if (datatype_ == Vtk::FLOAT32)
      out << std::setprecision(std::numeric_limits<float>::digits10+2);
    else
      out << std::setprecision(std::numeric_limits<double>::digits10+2);
  }

  dataCollector_.update();

  std::vector<pos_type> offsets; // pos => offset
  out << "<VTKFile type=\"ImageData\" version=\"1.0\" "
      << "byte_order=\"" << this->getEndian() << "\" header_type=\"UInt64\""
      << (format_ == Vtk::COMPRESSED ? " compressor=\"vtkZLibDataCompressor\">\n" : ">\n");
  out << "<ImageData WholeExtent=\"";
  auto const& extent = dataCollector_.extent();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << extent[2*i] << " " << extent[2*i+1];
  }
  out << "\" Origin=\"";
  auto const& origin = dataCollector_.origin();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << origin[i];
  }
  out << "\" Spacing=\"";
  auto const& spacing = dataCollector_.spacing();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << spacing[i];
  }
  out << "\">\n";
  out << "<Piece Extent=\"";
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << extent[2*i] << " " << extent[2*i+1];
  }
  out << "\">\n";

  { // Write data associated with grid points
    auto scalar = std::find_if(pointData_.begin(), pointData_.end(), [](auto const& v) { return v.ncomps() == 1; });
    auto vector = std::find_if(pointData_.begin(), pointData_.end(), [](auto const& v) { return v.ncomps() > 1; });
    out << "<PointData" << (scalar != pointData_.end() ? " Scalars=\"" + scalar->name() + "\"" : "")
                        << (vector != pointData_.end() ? " Vectors=\"" + vector->name() + "\"" : "")
                        << ">\n";
    for (auto const& v : pointData_)
      this->writeData(out, offsets, v, Super::POINT_DATA);
    out << "</PointData>\n";
  }

  { // Write data associated with grid cells
    auto scalar = std::find_if(cellData_.begin(), cellData_.end(), [](auto const& v) { return v.ncomps() == 1; });
    auto vector = std::find_if(cellData_.begin(), cellData_.end(), [](auto const& v) { return v.ncomps() > 1; });
    out << "<CellData" << (scalar != cellData_.end() ? " Scalars=\"" + scalar->name() + "\"" : "")
                       << (vector != cellData_.end() ? " Vectors=\"" + vector->name() + "\"" : "")
                       << ">\n";
    for (auto const& v : cellData_)
      this->writeData(out, offsets, v, Super::CELL_DATA);
    out << "</CellData>\n";
  }
  out << "</Piece>\n";
  out << "</ImageData>\n";

  std::vector<std::uint64_t> blocks; // size of i'th appended block
  pos_type appended_pos = 0;
  if (is_a(format_, Vtk::APPENDED)) {
    out << "<AppendedData encoding=\"raw\">\n_";
    appended_pos = out.tellp();
    for (auto const& v : pointData_) {
      if (datatype_ == Vtk::FLOAT32)
        blocks.push_back( this->template writeDataAppended<float>(out, v, Super::POINT_DATA) );
      else
        blocks.push_back( this->template writeDataAppended<double>(out, v, Super::POINT_DATA) );
    }
    for (auto const& v : cellData_) {
      if (datatype_ == Vtk::FLOAT32)
        blocks.push_back( this->template writeDataAppended<float>(out, v, Super::CELL_DATA) );
      else
        blocks.push_back( this->template writeDataAppended<double>(out, v, Super::CELL_DATA) );
    }
    out << "</AppendedData>\n";
  }

  out << "</VTKFile>";

  // fillin offset values and block sizes
  if (is_a(format_, Vtk::APPENDED)) {
    pos_type offset = 0;
    for (std::size_t i = 0; i < offsets.size(); ++i) {
      out.seekp(offsets[i]);
      out << '"' << offset << '"';
      offset += pos_type(blocks[i]);
    }
  }
}


template <class GV, class DC>
void VtkImageDataWriter<GV,DC>
  ::writeParallelFile (std::string const& pfilename, int size) const
{
  std::string filename = pfilename + ".p" + this->fileExtension();
  std::ofstream out(filename, std::ios_base::ate | std::ios::binary);

  out << "<VTKFile type=\"PImageData\" version=\"1.0\" "
      << "byte_order=\"" << this->getEndian() << "\" header_type=\"UInt64\""
      << (format_ == Vtk::COMPRESSED ? " compressor=\"vtkZLibDataCompressor\">\n" : ">\n");
  out << "<PImageData GhostLevel=\"0\" WholeExtent=\"";
  auto const& extent = dataCollector_.extent();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << extent[2*i] << " " << extent[2*i+1];
  }
  out << "\" Origin=\"";
  auto const& origin = dataCollector_.origin();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << origin[i];
  }
  out << "\" Spacing=\"";
  auto const& spacing = dataCollector_.spacing();
  for (int i = 0; i < 3; ++i) {
    out << (i == 0 ? "" : " ") << spacing[i];
  }
  out << "\">\n";

  { // Write data associated with grid points
    auto scalar = std::find_if(pointData_.begin(), pointData_.end(), [](auto const& v) { return v.ncomps() == 1; });
    auto vector = std::find_if(pointData_.begin(), pointData_.end(), [](auto const& v) { return v.ncomps() > 1; });
    out << "<PPointData" << (scalar != pointData_.end() ? " Scalars=\"" + scalar->name() + "\"" : "")
                         << (vector != pointData_.end() ? " Vectors=\"" + vector->name() + "\"" : "")
                         << ">\n";
    for (auto const& v : pointData_) {
      out << "<PDataArray Name=\"" << v.name() << "\" type=\"" << Vtk::Map::from_datatype[datatype_] << "\""
          << " NumberOfComponents=\"" << v.ncomps() << "\" />\n";
    }
    out << "</PPointData>\n";
  }

  { // Write data associated with grid cells
    auto scalar = std::find_if(cellData_.begin(), cellData_.end(), [](auto const& v) { return v.ncomps() == 1; });
    auto vector = std::find_if(cellData_.begin(), cellData_.end(), [](auto const& v) { return v.ncomps() > 1; });
    out << "<PCellData" << (scalar != cellData_.end() ? " Scalars=\"" + scalar->name() + "\"" : "")
                        << (vector != cellData_.end() ? " Vectors=\"" + vector->name() + "\"" : "")
                        << ">\n";
    for (auto const& v : cellData_) {
      out << "<PDataArray Name=\"" << v.name() << "\" type=\"" << Vtk::Map::from_datatype[datatype_] << "\""
          << " NumberOfComponents=\"" << v.ncomps() << "\" />\n";
    }
    out << "</PCellData>\n";
  }

  // Write piece file references
  for (int i = 0; i < size; ++i) {
    out << "<Piece Source=\"" << pfilename << "_p" << std::to_string(i) << "." << this->fileExtension() << "\" Extent=\"";
    for (int i = 0; i < 3; ++i) {
      out << (i == 0 ? "" : " ") << extent[2*i] << " " << extent[2*i+1];
    }
    out << "\" />\n";
  }

  out << "</PUnstructuredGrid>\n";
  out << "</VTKFile>";
}

}} // end namespace Dune::experimental