//============================================================================ // Copyright (c) Kitware, Inc. // All rights reserved. // See LICENSE.txt for details. // // This software is distributed WITHOUT ANY WARRANTY; without even // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR // PURPOSE. See the above copyright notice for more information. //============================================================================ #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include namespace { struct CallForBaseTypeFunctor { template void operator()(T t, bool& success, Functor functor, const vtkm::cont::UnknownArrayHandle& array, Args&&... args) { if (!array.IsBaseComponentType()) { return; } success = true; functor(t, array, std::forward(args)...); } }; template void CallForBaseType(Functor&& functor, const vtkm::cont::UnknownArrayHandle& array, Args&&... args) { bool success = true; vtkm::ListForEach(CallForBaseTypeFunctor{}, vtkm::TypeListScalarAll{}, success, std::forward(functor), array, std::forward(args)...); if (!success) { std::ostringstream out; out << "Unrecognized base type in array to be written out.\nArray: "; array.PrintSummary(out); throw vtkm::cont::ErrorBadValue(out.str()); } } template using ArrayHandleRectilinearCoordinates = vtkm::cont::ArrayHandleCartesianProduct, vtkm::cont::ArrayHandle, vtkm::cont::ArrayHandle>; struct OutputArrayDataFunctor { template VTKM_CONT void operator()(T, const vtkm::cont::UnknownArrayHandle& array, std::ostream& out, vtkm::io::FileType fileType) const { auto componentArray = array.ExtractArrayFromComponents(); auto portal = componentArray.ReadPortal(); switch (fileType) { case vtkm::io::FileType::ASCII: this->OutputAsciiArray(portal, out); break; case vtkm::io::FileType::BINARY: this->OutputBinaryArray(portal, out); break; } } template void OutputAsciiArray(const PortalType& portal, std::ostream& out) const { using T = typename PortalType::ValueType::ComponentType; vtkm::Id numValues = portal.GetNumberOfValues(); for (vtkm::Id valueIndex = 0; valueIndex < numValues; ++valueIndex) { auto value = portal.Get(valueIndex); for (vtkm::IdComponent cIndex = 0; cIndex < value.GetNumberOfComponents(); ++cIndex) { out << ((cIndex == 0) ? "" : " "); if (std::numeric_limits::is_integer && sizeof(T) == 1) { out << static_cast(value[cIndex]); } else { out << value[cIndex]; } } out << "\n"; } } template void OutputBinaryArray(const PortalType& portal, std::ostream& out) const { using T = typename PortalType::ValueType::ComponentType; vtkm::Id numValues = portal.GetNumberOfValues(); std::vector tuple; for (vtkm::Id valueIndex = 0; valueIndex < numValues; ++valueIndex) { auto value = portal.Get(valueIndex); tuple.resize(static_cast(value.GetNumberOfComponents())); for (vtkm::IdComponent cIndex = 0; cIndex < value.GetNumberOfComponents(); ++cIndex) { tuple[static_cast(cIndex)] = value[cIndex]; } if (vtkm::io::internal::IsLittleEndian()) { vtkm::io::internal::FlipEndianness(tuple); } out.write(reinterpret_cast(tuple.data()), static_cast(tuple.size() * sizeof(T))); } } }; void OutputArrayData(const vtkm::cont::UnknownArrayHandle& array, std::ostream& out, vtkm::io::FileType fileType) { CallForBaseType(OutputArrayDataFunctor{}, array, out, fileType); } struct GetFieldTypeNameFunctor { template void operator()(Type, const vtkm::cont::UnknownArrayHandle& array, std::string& name) const { if (array.IsBaseComponentType()) { name = vtkm::io::internal::DataTypeName::Name(); } } }; std::string GetFieldTypeName(const vtkm::cont::UnknownArrayHandle& array) { std::string name; CallForBaseType(GetFieldTypeNameFunctor{}, array, name); return name; } template void WriteDimensions(std::ostream& out, const vtkm::cont::CellSetStructured& cellSet) { auto pointDimensions = cellSet.GetPointDimensions(); using VTraits = vtkm::VecTraits; out << "DIMENSIONS "; out << VTraits::GetComponent(pointDimensions, 0) << " "; out << (DIM > 1 ? VTraits::GetComponent(pointDimensions, 1) : 1) << " "; out << (DIM > 2 ? VTraits::GetComponent(pointDimensions, 2) : 1) << "\n"; } void WritePoints(std::ostream& out, const vtkm::cont::DataSet& dataSet, vtkm::io::FileType fileType) { ///\todo: support other coordinate systems int cindex = 0; auto cdata = dataSet.GetCoordinateSystem(cindex).GetData(); std::string typeName = GetFieldTypeName(cdata); vtkm::Id npoints = cdata.GetNumberOfValues(); out << "POINTS " << npoints << " " << typeName << " " << '\n'; OutputArrayData(cdata, out, fileType); } template void WriteExplicitCellsAscii(std::ostream& out, const CellSetType& cellSet) { vtkm::Id nCells = cellSet.GetNumberOfCells(); vtkm::Id conn_length = 0; for (vtkm::Id i = 0; i < nCells; ++i) { conn_length += 1 + cellSet.GetNumberOfPointsInCell(i); } out << "CELLS " << nCells << " " << conn_length << '\n'; for (vtkm::Id i = 0; i < nCells; ++i) { vtkm::cont::ArrayHandle ids; vtkm::Id nids = cellSet.GetNumberOfPointsInCell(i); cellSet.GetIndices(i, ids); out << nids; auto IdPortal = ids.ReadPortal(); for (int j = 0; j < nids; ++j) { out << " " << IdPortal.Get(j); } out << '\n'; } out << "CELL_TYPES " << nCells << '\n'; for (vtkm::Id i = 0; i < nCells; ++i) { vtkm::Id shape = cellSet.GetCellShape(i); out << shape << '\n'; } } template void WriteExplicitCellsBinary(std::ostream& out, const CellSetType& cellSet) { vtkm::Id nCells = cellSet.GetNumberOfCells(); vtkm::Id conn_length = 0; for (vtkm::Id i = 0; i < nCells; ++i) { conn_length += 1 + cellSet.GetNumberOfPointsInCell(i); } out << "CELLS " << nCells << " " << conn_length << '\n'; std::vector buffer; buffer.reserve(static_cast(conn_length)); for (vtkm::Id i = 0; i < nCells; ++i) { vtkm::cont::ArrayHandle ids; vtkm::Id nids = cellSet.GetNumberOfPointsInCell(i); cellSet.GetIndices(i, ids); buffer.push_back(static_cast(nids)); auto IdPortal = ids.ReadPortal(); for (int j = 0; j < nids; ++j) { buffer.push_back(static_cast(IdPortal.Get(j))); } } if (vtkm::io::internal::IsLittleEndian()) { vtkm::io::internal::FlipEndianness(buffer); } VTKM_ASSERT(static_cast(buffer.size()) == conn_length); out.write(reinterpret_cast(buffer.data()), static_cast(buffer.size() * sizeof(vtkm::Int32))); out << "CELL_TYPES " << nCells << '\n'; buffer.resize(0); buffer.reserve(static_cast(nCells)); for (vtkm::Id i = 0; i < nCells; ++i) { buffer.push_back(static_cast(cellSet.GetCellShape(i))); } if (vtkm::io::internal::IsLittleEndian()) { vtkm::io::internal::FlipEndianness(buffer); } VTKM_ASSERT(static_cast(buffer.size()) == nCells); out.write(reinterpret_cast(buffer.data()), static_cast(buffer.size() * sizeof(vtkm::Int32))); } template void WriteExplicitCells(std::ostream& out, const CellSetType& cellSet, vtkm::io::FileType fileType) { switch (fileType) { case vtkm::io::FileType::ASCII: WriteExplicitCellsAscii(out, cellSet); break; case vtkm::io::FileType::BINARY: WriteExplicitCellsBinary(out, cellSet); break; } } void WritePointFields(std::ostream& out, const vtkm::cont::DataSet& dataSet, vtkm::io::FileType fileType) { bool wrote_header = false; for (vtkm::Id f = 0; f < dataSet.GetNumberOfFields(); f++) { const vtkm::cont::Field field = dataSet.GetField(f); if (field.GetAssociation() != vtkm::cont::Field::Association::Points) { continue; } if (field.GetName() == dataSet.GetCoordinateSystemName()) { // Do not write out the first coordinate system as a field. continue; } vtkm::Id npoints = field.GetNumberOfValues(); int ncomps = field.GetData().GetNumberOfComponentsFlat(); if (!wrote_header) { out << "POINT_DATA " << npoints << '\n'; wrote_header = true; } std::string typeName = GetFieldTypeName(field.GetData()); std::string name = field.GetName(); for (auto& c : name) { if (std::isspace(c)) { c = '_'; } } out << "SCALARS " << name << " " << typeName << " " << ncomps << '\n'; out << "LOOKUP_TABLE default" << '\n'; OutputArrayData(field.GetData(), out, fileType); } } void WriteCellFields(std::ostream& out, const vtkm::cont::DataSet& dataSet, vtkm::io::FileType fileType) { bool wrote_header = false; for (vtkm::Id f = 0; f < dataSet.GetNumberOfFields(); f++) { const vtkm::cont::Field field = dataSet.GetField(f); if (!field.IsCellField()) { continue; } vtkm::Id ncells = field.GetNumberOfValues(); int ncomps = field.GetData().GetNumberOfComponentsFlat(); if (!wrote_header) { out << "CELL_DATA " << ncells << '\n'; wrote_header = true; } std::string typeName = GetFieldTypeName(field.GetData()); std::string name = field.GetName(); for (auto& c : name) { if (std::isspace(c)) { c = '_'; } } out << "SCALARS " << name << " " << typeName << " " << ncomps << '\n'; out << "LOOKUP_TABLE default" << '\n'; OutputArrayData(field.GetData(), out, fileType); } } template void WriteDataSetAsUnstructured(std::ostream& out, const vtkm::cont::DataSet& dataSet, const CellSetType& cellSet, vtkm::io::FileType fileType) { out << "DATASET UNSTRUCTURED_GRID" << '\n'; WritePoints(out, dataSet, fileType); WriteExplicitCells(out, cellSet, fileType); } template void WriteDataSetAsStructuredPoints(std::ostream& out, const vtkm::cont::ArrayHandleUniformPointCoordinates& points, const vtkm::cont::CellSetStructured& cellSet) { out << "DATASET STRUCTURED_POINTS\n"; WriteDimensions(out, cellSet); auto portal = points.ReadPortal(); auto origin = portal.GetOrigin(); auto spacing = portal.GetSpacing(); out << "ORIGIN " << origin[0] << " " << origin[1] << " " << origin[2] << "\n"; out << "SPACING " << spacing[0] << " " << spacing[1] << " " << spacing[2] << "\n"; } template void WriteDataSetAsRectilinearGrid(std::ostream& out, const ArrayHandleRectilinearCoordinates& points, const vtkm::cont::CellSetStructured& cellSet, vtkm::io::FileType fileType) { out << "DATASET RECTILINEAR_GRID\n"; WriteDimensions(out, cellSet); std::string typeName = vtkm::io::internal::DataTypeName::Name(); vtkm::cont::ArrayHandle dimArray; dimArray = points.GetFirstArray(); out << "X_COORDINATES " << dimArray.GetNumberOfValues() << " " << typeName << "\n"; OutputArrayData(dimArray, out, fileType); dimArray = points.GetSecondArray(); out << "Y_COORDINATES " << dimArray.GetNumberOfValues() << " " << typeName << "\n"; OutputArrayData(dimArray, out, fileType); dimArray = points.GetThirdArray(); out << "Z_COORDINATES " << dimArray.GetNumberOfValues() << " " << typeName << "\n"; OutputArrayData(dimArray, out, fileType); } template void WriteDataSetAsStructuredGrid(std::ostream& out, const vtkm::cont::DataSet& dataSet, const vtkm::cont::CellSetStructured& cellSet, vtkm::io::FileType fileType) { out << "DATASET STRUCTURED_GRID" << '\n'; WriteDimensions(out, cellSet); WritePoints(out, dataSet, fileType); } template void WriteDataSetAsStructured(std::ostream& out, const vtkm::cont::DataSet& dataSet, const vtkm::cont::CellSetStructured& cellSet, vtkm::io::FileType fileType) { ///\todo: support rectilinear // Type of structured grid (uniform, rectilinear, curvilinear) is determined by coordinate system auto coordSystem = dataSet.GetCoordinateSystem().GetData(); if (coordSystem.IsType()) { // uniform is written as "structured points" WriteDataSetAsStructuredPoints( out, coordSystem.AsArrayHandle(), cellSet); } else if (coordSystem.IsType>()) { WriteDataSetAsRectilinearGrid( out, coordSystem.AsArrayHandle>(), cellSet, fileType); } else if (coordSystem.IsType>()) { WriteDataSetAsRectilinearGrid( out, coordSystem.AsArrayHandle>(), cellSet, fileType); } else { // Curvilinear is written as "structured grid" WriteDataSetAsStructuredGrid(out, dataSet, cellSet, fileType); } } void Write(std::ostream& out, const vtkm::cont::DataSet& dataSet, vtkm::io::FileType fileType) { // The Paraview parser cannot handle scientific notation: out << std::fixed; // This causes a big problem for the dataset writer. // Fixed point and floating point are fundamentally different. // This is a workaround, but until Paraview supports parsing floats, // this is as good as we can do. #ifdef VTKM_USE_DOUBLE_PRECISION out << std::setprecision(18); #else out << std::setprecision(10); #endif out << "# vtk DataFile Version 3.0" << '\n'; out << "vtk output" << '\n'; switch (fileType) { case vtkm::io::FileType::ASCII: out << "ASCII" << '\n'; break; case vtkm::io::FileType::BINARY: out << "BINARY" << '\n'; break; } vtkm::cont::UnknownCellSet cellSet = dataSet.GetCellSet(); if (cellSet.IsType>()) { WriteDataSetAsUnstructured( out, dataSet, cellSet.AsCellSet>(), fileType); } else if (cellSet.IsType>()) { WriteDataSetAsStructured( out, dataSet, cellSet.AsCellSet>(), fileType); } else if (cellSet.IsType>()) { WriteDataSetAsStructured( out, dataSet, cellSet.AsCellSet>(), fileType); } else if (cellSet.IsType>()) { WriteDataSetAsStructured( out, dataSet, cellSet.AsCellSet>(), fileType); } else if (cellSet.IsType>()) { // these function just like explicit cell sets WriteDataSetAsUnstructured( out, dataSet, cellSet.AsCellSet>(), fileType); } else if (cellSet.IsType()) { WriteDataSetAsUnstructured( out, dataSet, cellSet.AsCellSet(), fileType); } else { throw vtkm::cont::ErrorBadType("Could not determine type to write out."); } WritePointFields(out, dataSet, fileType); WriteCellFields(out, dataSet, fileType); } } // anonymous namespace namespace vtkm { namespace io { VTKDataSetWriter::VTKDataSetWriter(const char* fileName) : FileName(fileName) { } VTKDataSetWriter::VTKDataSetWriter(const std::string& fileName) : FileName(fileName) { } void VTKDataSetWriter::WriteDataSet(const vtkm::cont::DataSet& dataSet) const { if (dataSet.GetNumberOfCoordinateSystems() < 1) { throw vtkm::cont::ErrorBadValue( "DataSet has no coordinate system, which is not supported by VTK file format."); } try { std::ofstream fileStream(this->FileName.c_str(), std::fstream::trunc | std::fstream::binary); Write(fileStream, dataSet, this->GetFileType()); fileStream.close(); } catch (std::ofstream::failure& error) { throw vtkm::io::ErrorIO(error.what()); } } vtkm::io::FileType VTKDataSetWriter::GetFileType() const { return this->FileType; } void VTKDataSetWriter::SetFileType(vtkm::io::FileType type) { this->FileType = type; } } } // namespace vtkm::io