1 #pragma once 2 3 #include <iosfwd> 4 #include <map> 5 #include <memory> 6 #include <optional> 7 #include <string> 8 #include <vector> 9 10 #include <dune/common/parallel/mpihelper.hh> 11 #include <dune/vtk/filewriter.hh> 12 #include <dune/vtk/function.hh> 13 #include <dune/vtk/types.hh> 14 15 namespace Dune 16 { 17 /// Interface for file writers for the Vtk XML file formats 18 /** 19 * \tparam GV Model of Dune::GridView 20 * \tparam DC Model of \ref DataCollectorInterface 21 **/ 22 template <class GV, class DC> 23 class VtkWriterInterface 24 : public Vtk::FileWriter 25 { 26 template <class> friend class TimeseriesWriter; 27 template <class> friend class PvdWriter; 28 29 public: 30 using GridView = GV; 31 using DataCollector = DC; 32 33 protected: 34 using VtkFunction = Dune::Vtk::Function<GridView>; 35 using pos_type = typename std::ostream::pos_type; 36 37 enum PositionTypes { 38 POINT_DATA, 39 CELL_DATA 40 }; 41 42 public: 43 /// \brief Constructor, passes the gridView to the DataCollector 44 /** 45 * Creates a new VtkWriterInterface for the provided GridView. Initializes a 46 * DataCollector that is used to collect point coordinates, cell connectivity and 47 * data values. 48 * 49 * This constructor assumes, that the DataCollector can be constructed from a single argument, 50 * the passed gridView. 51 * 52 * \param gridView Implementation of Dune::GridView 53 * \param format Format of the VTK file, either Vtk::FormatTypes::BINARY, Vtk::FormatTypes::ASCII, or Vtk::COMPRESSED 54 * \param datatype Data type of a single component of the point coordinates [Vtk::DataTypes::FLOAT32] 55 * \param headertype Integer type used in binary data headers [Vtk::DataTypes::UINT32] 56 **/ VtkWriterInterface(GridView const & gridView,Vtk::FormatTypes format=Vtk::FormatTypes::BINARY,Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32,Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)57 VtkWriterInterface (GridView const& gridView, 58 Vtk::FormatTypes format = Vtk::FormatTypes::BINARY, 59 Vtk::DataTypes datatype = Vtk::DataTypes::FLOAT32, 60 Vtk::DataTypes headertype = Vtk::DataTypes::UINT32) 61 : VtkWriterInterface(std::make_shared<DataCollector>(gridView), format, datatype, headertype) 62 {} 63 64 /// \brief Constructor, wraps the passed DataCollector in a non-destroying shared_ptr VtkWriterInterface(DataCollector & dataCollector,Vtk::FormatTypes format=Vtk::FormatTypes::BINARY,Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32,Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)65 VtkWriterInterface (DataCollector& dataCollector, 66 Vtk::FormatTypes format = Vtk::FormatTypes::BINARY, 67 Vtk::DataTypes datatype = Vtk::DataTypes::FLOAT32, 68 Vtk::DataTypes headertype = Vtk::DataTypes::UINT32) 69 : VtkWriterInterface(stackobject_to_shared_ptr(dataCollector), format, datatype, headertype) 70 {} 71 72 /// \brief Constructor, stores the passed DataCollector VtkWriterInterface(std::shared_ptr<DataCollector> dataCollector,Vtk::FormatTypes format=Vtk::FormatTypes::BINARY,Vtk::DataTypes datatype=Vtk::DataTypes::FLOAT32,Vtk::DataTypes headertype=Vtk::DataTypes::UINT32)73 VtkWriterInterface (std::shared_ptr<DataCollector> dataCollector, 74 Vtk::FormatTypes format = Vtk::FormatTypes::BINARY, 75 Vtk::DataTypes datatype = Vtk::DataTypes::FLOAT32, 76 Vtk::DataTypes headertype = Vtk::DataTypes::UINT32) 77 : dataCollector_(std::move(dataCollector)) 78 { 79 setFormat(format); 80 setDatatype(datatype); 81 setHeadertype(headertype); 82 } 83 84 85 /// \brief Write the attached data to the file 86 /** 87 * \param fn Filename of the VTK file. May contain a directory and any file extension. 88 * \param dir The optional parameter specifies the directory of the partition files for parallel writes. 89 * 90 * \returns File name that is actually written. 91 **/ 92 virtual std::string write (std::string const& fn, std::optional<std::string> dir = {}) const override; 93 94 /// \brief Attach point data to the writer 95 /** 96 * Attach a global function to the writer that will be evaluated at grid points 97 * (vertices and higher order points). The global function must be 98 * assignable to the function wrapper \ref Vtk::Function. Additional argument 99 * for output datatype and number of components can be passed. See \ref Vtk::Function 100 * Constructor for possible arguments. 101 * 102 * \param fct A GridFunction, LocalFunction, or Dune::VTKFunction 103 * \param args... Additional arguments, like `name`, `numComponents`, `dataType` or `Vtk::FieldInfo` 104 **/ 105 template <class Function, class... Args> addPointData(Function && fct,Args &&...args)106 VtkWriterInterface& addPointData (Function&& fct, Args&&... args) 107 { 108 pointData_.emplace_back(std::forward<Function>(fct), std::forward<Args>(args)..., 109 datatype_, Vtk::RangeTypes::AUTO); 110 return *this; 111 } 112 113 /// \brief Attach cell data to the writer 114 /** 115 * Attach a global function to the writer that will be evaluated at cell centers. 116 * The global function must be assignable to the function wrapper \ref Vtk::Function. 117 * Additional argument for output datatype and number of components can be passed. 118 * See \ref Vtk::Function Constructor for possible arguments. 119 * 120 * \param fct A GridFunction, LocalFunction, or Dune::VTKFunction 121 * \param args... Additional arguments, like `name`, `numComponents`, `dataType` or `Vtk::FieldInfo` 122 **/ 123 template <class Function, class... Args> addCellData(Function && fct,Args &&...args)124 VtkWriterInterface& addCellData (Function&& fct, Args&&... args) 125 { 126 cellData_.emplace_back(std::forward<Function>(fct), std::forward<Args>(args)..., 127 datatype_, Vtk::RangeTypes::AUTO); 128 return *this; 129 } 130 131 132 // Sets the VTK file format setFormat(Vtk::FormatTypes format)133 void setFormat (Vtk::FormatTypes format) 134 { 135 format_ = format; 136 137 if (format_ == Vtk::FormatTypes::COMPRESSED) { 138 #if HAVE_VTK_ZLIB 139 compressor_ = Vtk::CompressorTypes::ZLIB; 140 #else 141 std::cout << "Dune is compiled without compression. Falling back to BINARY VTK output!\n"; 142 format_ = Vtk::FormatTypes::BINARY; 143 #endif 144 } else { 145 compressor_ = Vtk::CompressorTypes::NONE; 146 } 147 148 } 149 150 /// Sets the global datatype used for coordinates and other global float values setDatatype(Vtk::DataTypes datatype)151 void setDatatype (Vtk::DataTypes datatype) 152 { 153 datatype_ = datatype; 154 } 155 156 /// Sets the integer type used in binary data headers setHeadertype(Vtk::DataTypes datatype)157 void setHeadertype (Vtk::DataTypes datatype) 158 { 159 headertype_ = datatype; 160 } 161 162 /// Sets the compressor type used in binary data headers, Additionally a compression 163 /// level can be passed with level = -1 means: default compression level. Level must be in [0-9] setCompressor(Vtk::CompressorTypes compressor,int level=-1)164 void setCompressor (Vtk::CompressorTypes compressor, int level = -1) 165 { 166 compressor_ = compressor; 167 compression_level = level; 168 VTK_ASSERT(level >= -1 && level <= 9); 169 170 if (compressor_ != Vtk::CompressorTypes::NONE) 171 format_ = Vtk::FormatTypes::COMPRESSED; 172 } 173 174 private: 175 /// Write a serial VTK file in Unstructured format 176 virtual void writeSerialFile (std::ofstream& out) const = 0; 177 178 /// Write a parallel VTK file `pfilename.pvtx` in XML format, 179 /// with `size` the number of pieces and serial files given by `pfilename_p[i].vtu` 180 /// for [i] in [0,...,size). 181 virtual void writeParallelFile (std::ofstream& out, std::string const& pfilename, int size) const = 0; 182 183 /// Return the file extension of the serial file (not including the dot) 184 virtual std::string fileExtension () const = 0; 185 186 /// Write points and cells in raw/compressed format to output stream 187 virtual void writeGridAppended (std::ofstream& out, std::vector<std::uint64_t>& blocks) const = 0; 188 189 protected: 190 // Write the point or cell values given by the grid function `fct` to the 191 // output stream `out`. In case of binary format, append the streampos of XML 192 // attributes "offset" to the vector `offsets`. 193 void writeData (std::ofstream& out, 194 std::vector<pos_type>& offsets, 195 VtkFunction const& fct, 196 PositionTypes type, 197 std::optional<std::size_t> timestep = {}) const; 198 199 // Write point-data and cell-data in raw/compressed format to output stream 200 void writeDataAppended (std::ofstream& out, std::vector<std::uint64_t>& blocks) const; 201 202 // Write the coordinates of the vertices to the output stream `out`. In case 203 // of binary format, appends the streampos of XML attributes "offset" to the 204 // vector `offsets`. 205 void writePoints (std::ofstream& out, 206 std::vector<pos_type>& offsets, 207 std::optional<std::size_t> timestep = {}) const; 208 209 // Write Appended section and fillin offset values to XML attributes 210 void writeAppended (std::ofstream& out, std::vector<pos_type> const& offsets) const; 211 212 // Write the `values` in blocks (possibly compressed) to the output 213 // stream `out`. Return the written block size. 214 template <class HeaderType, class FloatType> 215 std::uint64_t writeValuesAppended (std::ofstream& out, std::vector<FloatType> const& values) const; 216 217 // Write the `values` in a space and newline separated list of ascii representations. 218 // The precision is controlled by the datatype and numerical_limits::digits10. 219 template <class T> 220 void writeValuesAscii (std::ofstream& out, std::vector<T> const& values) const; 221 222 // Write the XML file header of a VTK file `<VTKFile ...>` 223 void writeHeader (std::ofstream& out, std::string const& type) const; 224 225 /// Return PointData/CellData attributes for the name of the first scalar/vector/tensor DataArray 226 std::string getNames (std::vector<VtkFunction> const& data) const; 227 228 // Returns endianness getEndian() const229 std::string getEndian () const 230 { 231 short i = 1; 232 return (reinterpret_cast<char*>(&i)[1] == 1 ? "BigEndian" : "LittleEndian"); 233 } 234 235 // provide accessor to \ref fileExtension virtual method getFileExtension() const236 std::string getFileExtension () const 237 { 238 return fileExtension(); 239 } 240 241 // Returns the VTK file format initialized in the constructor getFormat() const242 Vtk::FormatTypes getFormat () const 243 { 244 return format_; 245 } 246 247 // Returns the global datatype used for coordinates and other global float values getDatatype() const248 Vtk::DataTypes getDatatype () const 249 { 250 return datatype_; 251 } 252 253 // Return the global MPI communicator. comm() const254 auto comm () const 255 { 256 return MPIHelper::getCollectiveCommunication(); 257 } 258 259 protected: 260 std::shared_ptr<DataCollector> dataCollector_; 261 262 Vtk::FormatTypes format_; 263 Vtk::DataTypes datatype_; 264 Vtk::DataTypes headertype_; 265 Vtk::CompressorTypes compressor_ = Vtk::CompressorTypes::NONE; 266 267 // attached data 268 std::vector<VtkFunction> pointData_; 269 std::vector<VtkFunction> cellData_; 270 271 std::size_t const block_size = 1024*32; 272 int compression_level = -1; // in [0,9], -1 ... use default value 273 }; 274 275 276 template <class Writer> 277 struct IsVtkWriter 278 { 279 template <class GV, class DC> 280 static std::uint16_t test(VtkWriterInterface<GV,DC> const&); 281 static std::uint8_t test(...); // fall-back overload 282 283 static constexpr bool value = sizeof(test(std::declval<Writer>())) > sizeof(std::uint8_t); 284 }; 285 286 } // end namespace Dune 287 288 #include "vtkwriterinterface.impl.hh" 289