Skip to content
Snippets Groups Projects
Commit 5c48b6a8 authored by Stéphane Del Pino's avatar Stéphane Del Pino
Browse files

Begin VTKWriter improvements

- factor bunch of code
- rework node coordinates writting
parent 7b722d17
No related branches found
No related tags found
1 merge request!10Feature/data output
......@@ -14,63 +14,45 @@
class VTKWriter
{
private:
template <typename DataType> struct VTKType {};
template <> struct VTKType<int8_t>
{
inline const static std::string name{"Int8"};
};
template <> struct VTKType<uint8_t>
{
inline const static std::string name{"UInt8"};
};
template <> struct VTKType<int16_t>
{
inline const static std::string name{"Int16"};
};
template <> struct VTKType<uint16_t>
{
inline const static std::string name{"UInt16"};
};
template <> struct VTKType<int32_t>
{
inline const static std::string name{"Int32"};
};
template <> struct VTKType<uint32_t>
{
inline const static std::string name{"UInt32"};
};
template <> struct VTKType<int64_t>
{
inline const static std::string name{"Int64"};
};
template <> struct VTKType<uint64_t>
{
inline const static std::string name{"UInt64"};
};
template <> struct VTKType<float>
{
inline const static std::string name{"Float32"};
};
template <> struct VTKType<double>
{
inline const static std::string name{"Float64"};
};
const std::string m_base_filename;
unsigned int m_file_number;
double m_last_time;
const double m_time_period;
template <typename DataType> struct VTKType {};
template <typename DataType>
void _write_array(std::ofstream& os,
const std::string& name,
const Array<DataType>& item_value)
{
os << "<DataArray type=\"" << VTKType<DataType>::name
<< "\" Name=\"" << name << "\">\n";
for (typename Array<DataType>::index_type i=0;
i<item_value.size(); ++i) {
// The following '+' enforces integer output for char types
os << +item_value[i] << ' ';
}
os << "\n</DataArray>\n";
}
template <size_t N,
typename DataType>
void _write_array(std::ofstream& os,
const std::string& name,
const Array<TinyVector<N, DataType>>& item_value)
{
os << "<DataArray type=\"" << VTKType<DataType>::name
<< "\" Name=\"" << name << "\" NumberOfComponents=\"" << N << "\">\n";
for (typename Array<DataType>::index_type i=0;
i<item_value.size(); ++i) {
for (size_t j=0; j<N; ++j) {
// The following '+' enforces integer output for char types
os << +item_value[i][j] << ' ';
}
}
os << "\n</DataArray>\n";
}
template <typename DataType>
void _write_node_value(std::ofstream& os,
......@@ -80,7 +62,8 @@ class VTKWriter
os << "<DataArray type=\"" << VTKType<DataType>::name
<< "\" Name=\"" << name << "\">\n";
for (NodeId i=0; i<item_value.size(); ++i) {
os << item_value[i] << ' ';
// The following '+' enforces integer output for char types
os << +item_value[i] << ' ';
}
os << "\n</DataArray>\n";
}
......@@ -95,31 +78,28 @@ class VTKWriter
<< "\" Name=\"" << name << "\" NumberOfComponents=\"" << N << "\">\n";
for (NodeId i=0; i<item_value.size(); ++i) {
for (size_t j=0; j<N; ++j) {
os << item_value[i][j] << ' ';
// The following '+' enforces integer output for char types
os << +item_value[i][j] << ' ';
}
}
os << "\n</DataArray>\n";
}
template <typename DataType>
void _write_node_value(std::ofstream& os,
const std::string& name,
const CellValue<const DataType>& item_value) {}
template <typename DataType>
void _write_cell_value(std::ofstream& os,
const std::string& name,
const NodeValue<const DataType>& item_value) {}
void _write_node_value(std::ofstream&,
const std::string&,
const CellValue<const DataType>&) {}
template <typename DataType>
void _write_cell_value(std::ofstream& os,
const std::string& name,
const CellValue<const DataType>& item_value)
{
os << "<DataArray type=\"Float64\" Name=\"" << name << "\">\n";
os << "<DataArray type=\"" << VTKType<DataType>::name
<< "\" Name=\"" << name << "\">\n";
for (CellId i=0; i<item_value.size(); ++i) {
os << item_value[i] << ' ';
// The following '+' enforces integer output for char types
os << +item_value[i] << ' ';
}
os << "\n</DataArray>\n";
}
......@@ -130,15 +110,22 @@ class VTKWriter
const std::string& name,
const CellValue<const TinyVector<N, DataType>>& item_value)
{
os << "<DataArray type=\"Float64\" Name=\"" << name << "\" NumberOfComponents=\"" << N << "\">\n";
os << "<DataArray type=\"" << VTKType<DataType>::name
<< "\" Name=\"" << name << "\" NumberOfComponents=\"" << N << "\">\n";
for (CellId i=0; i<item_value.size(); ++i) {
for (size_t j=0; j<N; ++j) {
os << item_value[i][j] << ' ';
// The following '+' enforces integer output for char types
os << +item_value[i][j] << ' ';
}
}
os << "\n</DataArray>\n";
}
template <typename DataType>
void _write_cell_value(std::ofstream&,
const std::string&,
const NodeValue<const DataType>&) {}
public:
template <typename MeshType>
void write(const MeshType& mesh,
......@@ -175,94 +162,83 @@ class VTKWriter
}
fout << "</PointData>\n";
fout << "<Points>\n";
fout << "<DataArray Name=\"Positions\" NumberOfComponents=\"3\" type=\"Float64\" format=\"ascii\">\n";
{
using Rd = TinyVector<MeshType::dimension>;
const NodeValue<const Rd>& xr = mesh.xr();
if constexpr(MeshType::dimension == 1) {
for (NodeId r=0; r<mesh.numberOfNodes(); ++r) {
for (unsigned short i=0; i<1; ++i) {
fout << xr[r][i] << ' ';
}
fout << "0 0 "; // VTK requires 3 components
}
} else if constexpr (MeshType::dimension == 2) {
Array<TinyVector<3>> positions(mesh.numberOfNodes());
for (NodeId r=0; r<mesh.numberOfNodes(); ++r) {
for (unsigned short i=0; i<2; ++i) {
fout << xr[r][i] << ' ';
for (unsigned short i=0; i<MeshType::dimension; ++i) {
positions[r][i] = xr[r][i];
}
fout << "0 "; // VTK requires 3 components
for (unsigned short i=MeshType::dimension; i<3; ++i) {
positions[r][i] = 0;
}
} else {
for (NodeId r=0; r<mesh.numberOfNodes(); ++r) {
for (unsigned short i=0; i<3; ++i) {
fout << xr[r][i] << ' ';
}
_write_array(fout, "Positions", positions);
}
}
fout << '\n';
fout << "</DataArray>\n";
fout << "</Points>\n";
fout << "<Cells>\n";
fout << "<DataArray type=\"Int32\" Name=\"connectivity\" NumberOfComponents=\"1\" format=\"ascii\">\n";
{
const auto& cell_to_node_matrix
= mesh.connectivity().cellToNodeMatrix();
Array<int> connectivity(cell_to_node_matrix.numberOfEntries());
for (CellId j=0; j<mesh.numberOfCells(); ++j) {
const auto& cell_nodes = cell_to_node_matrix[j];
for (unsigned short r=0; r<cell_nodes.size(); ++r) {
fout << cell_nodes[r] << ' ';
for (size_t i=0; i<cell_to_node_matrix.numberOfEntries(); ++i) {
connectivity[i] = cell_to_node_matrix.entries()[i];
}
_write_array(fout, "connectivity", connectivity);
}
fout << '\n';
fout << "</DataArray>\n";
fout << "<DataArray type=\"UInt32\" Name=\"offsets\" NumberOfComponents=\"1\" format=\"ascii\">\n";
{
const auto& cell_to_node_matrix
= mesh.connectivity().cellToNodeMatrix();
Array<unsigned int> offsets(mesh.numberOfCells());
unsigned int offset=0;
for (CellId j=0; j<mesh.numberOfCells(); ++j) {
const auto& cell_nodes = cell_to_node_matrix[j];
offset += cell_nodes.size();
fout << offset << ' ';
offsets[j]=offset;
}
_write_array(fout, "offsets", offsets);
}
fout << '\n';
fout << "</DataArray>\n";
fout << "<DataArray type=\"Int8\" Name=\"types\" NumberOfComponents=\"1\" format=\"ascii\">\n";
{
const auto& cell_to_node_matrix
= mesh.connectivity().cellToNodeMatrix();
Array<int8_t> types(mesh.numberOfCells());
for (CellId j=0; j<mesh.numberOfCells(); ++j) {
const auto& cell_nodes = cell_to_node_matrix[j];
switch (cell_nodes.size()) {
case 2: {
fout << "3 ";
types[j]=3;
break;
}
case 3: {
fout << "5 ";
types[j]=5;
break;
}
case 4: {
if (mesh.meshDimension() == 3) {
fout << "10 ";
types[j]=10;
} else {
fout << "9 ";
types[j]=9;
}
break;
}
case 8: {
fout << "12 ";
types[j]=12;
break;
}
default: {
fout << "7 ";
types[j]=7;
break;
}
}
}
fout << '\n';
fout << "</DataArray>\n";
_write_array(fout, "types", types);
}
fout << "</Cells>\n";
fout << "</Piece>\n";
......@@ -271,6 +247,7 @@ class VTKWriter
m_file_number++;
}
VTKWriter(const std::string& base_filename,
const double time_period)
: m_base_filename(base_filename),
......@@ -282,4 +259,55 @@ class VTKWriter
~VTKWriter() = default;
};
template <>
struct VTKWriter::VTKType<int8_t>
{
inline const static std::string name{"Int8"};
};
template <> struct VTKWriter::VTKType<uint8_t>
{
inline const static std::string name{"UInt8"};
};
template <> struct VTKWriter::VTKType<int16_t>
{
inline const static std::string name{"Int16"};
};
template <> struct VTKWriter::VTKType<uint16_t>
{
inline const static std::string name{"UInt16"};
};
template <> struct VTKWriter::VTKType<int32_t>
{
inline const static std::string name{"Int32"};
};
template <> struct VTKWriter::VTKType<uint32_t>
{
inline const static std::string name{"UInt32"};
};
template <> struct VTKWriter::VTKType<int64_t>
{
inline const static std::string name{"Int64"};
};
template <> struct VTKWriter::VTKType<uint64_t>
{
inline const static std::string name{"UInt64"};
};
template <> struct VTKWriter::VTKType<float>
{
inline const static std::string name{"Float32"};
};
template <> struct VTKWriter::VTKType<double>
{
inline const static std::string name{"Float64"};
};
#endif // VTK_WRITER_HPP
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment