Loading...
Searching...
No Matches
VTK.hpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25
26#ifndef __BITPIT_VTK_HPP__
27#define __BITPIT_VTK_HPP__
28
29#include <cstdint>
30#include <typeinfo>
31#include <type_traits>
32#include <vector>
33#include <array>
34#include <typeindex>
35#include <unordered_map>
36
37#include "bitpit_common.hpp"
38#include "GenericIO.hpp"
39
40namespace bitpit{
41
42class VTK;
43class VTKField;
44
49enum class VTKWriteMode {
50 DEFAULT = 0,
51 NO_INCREMENT =1,
52 NO_SERIES =2
53};
54
63enum class VTKFieldType {
64 UNDEFINED = -1,
65 SCALAR = 0,
66 VECTOR = 1,
67 TENSOR = 2,
68};
69
74enum class VTKDataType {
75 UNDEFINED,
76 Int8 ,
77 Int16 ,
78 Int32 ,
79 Int64 ,
80 UInt8 ,
81 UInt16 ,
82 UInt32 ,
83 UInt64 ,
84 Float32 ,
85 Float64
86};
87
92enum class VTKFormat {
93 UNDEFINED,
94 ASCII,
95 APPENDED
96};
97
102enum class VTKLocation {
103 UNDEFINED,
104 CELL,
105 POINT
106};
107
112enum class VTKElementType {
113 UNDEFINED = -1,
114 VERTEX = 1,
115 LINE = 3,
116 TRIANGLE = 5,
117 POLYGON = 7,
118 PIXEL = 8,
119 QUAD = 9,
120 TETRA = 10,
121 VOXEL = 11,
122 HEXAHEDRON = 12,
123 WEDGE = 13,
124 PYRAMID = 14,
125 QUADRATIC_EDGE = 21,
126 QUADRATIC_TRIANGLE = 22,
127 QUADRATIC_QUAD = 23,
128 QUADRATIC_TETRA = 24,
129 QUADRATIC_HEXAHEDRON = 25,
130 POLYHEDRON = 42
131};
132
138 POINTS = 0,
139 OFFSETS = 1,
140 TYPES = 2,
141 CONNECTIVITY = 3,
142 FACE_STREAMS = 4,
143 FACE_OFFSETS = 5
144};
145
151 X_COORDS = 0,
152 Y_COORDS = 1,
153 Z_COORDS = 2
154};
155
157
158 private:
159 static std::unordered_map<std::type_index, VTKDataType> m_types;
161 public:
162 static uint8_t sizeOfType( VTKDataType type );
163
164 template<typename T>
165 static VTKDataType registerType();
166
167 template<typename T>
168 static VTKDataType registerType(VTKDataType VTKType);
169
170 static VTKDataType whichType( const std::type_info & ) ;
171
172 template<typename T, int nesting=0, typename std::enable_if<std::is_pod<T>::value && !utils::is_iterable<T>::value>::type* = nullptr>
173 static VTKDataType whichType();
174
175 template<typename T, int nesting=0, typename std::enable_if<utils::is_iterable<T>::value>::type* = nullptr>
176 static VTKDataType whichType();
177
178};
179
181 private:
182
183 public:
185 virtual ~VTKBaseContainer( ) = default;
186
187 virtual VTKBaseContainer * clone() const = 0 ;
188
189 virtual void flushData( std::fstream &, VTKFormat) =0 ;
190 virtual void absorbData( std::fstream &, VTKFormat, uint64_t, uint8_t) =0 ;
191};
192
193template<class T>
195 private:
196 std::vector<T>* m_ptr ;
198 public:
199 VTKVectorContainer( std::vector<T> &) ;
200
201 VTKVectorContainer* clone() const override ;
202
203 void flushData( std::fstream &, VTKFormat) override ;
204 void absorbData( std::fstream &, VTKFormat, uint64_t, uint8_t) override ;
205 void resize( std::true_type, uint64_t , uint8_t) ;
206 void resize( std::false_type, uint64_t , uint8_t) ;
207};
208
210
211 private:
212
213 public:
214 virtual ~VTKBaseStreamer() = default;
215
216 virtual void flushData( std::fstream &, const std::string &, VTKFormat) ;
217 virtual void absorbData( std::fstream &, const std::string &, VTKFormat, uint64_t, uint8_t, VTKDataType) ;
218
219 template<typename T>
220 void flushValue(std::fstream &, VTKFormat, const T &value) const;
221 template<typename T>
222 void flushValue(std::fstream &, VTKFormat, const T *values, int nValues) const;
223
224};
225
227
228 private:
229 std::unordered_map<std::string,std::unique_ptr<VTKBaseContainer> > m_field ;
231 public:
234 VTKNativeStreamer( VTKNativeStreamer && ) = default;
235
237 VTKNativeStreamer & operator=(VTKNativeStreamer &&other) = default;
238
239 template<class T>
240 void addData( const std::string &, std::vector<T> & ) ;
241 void removeData( const std::string & ) ;
242 bool hasData( const std::string & ) const ;
243 void flushData( std::fstream &, const std::string &, VTKFormat) override ;
244 void absorbData( std::fstream &, const std::string &, VTKFormat, uint64_t, uint8_t, VTKDataType) override ;
245};
246
248
249 //members
250 protected:
251 std::string m_name;
256 uint64_t m_offset;
257 std::fstream::pos_type m_position;
262 //methods
263 public:
264 static unsigned getComponentCount(VTKFieldType fieldType);
265
266 virtual ~VTKField() = default;
267
268 VTKField();
269 VTKField( const std::string & );
270
271 const std::string & getName() const;
272 VTKDataType getDataType() const;
274 unsigned getComponentCount() const;
275 VTKLocation getLocation() const;
277 uint64_t getOffset() const;
278 std::fstream::pos_type getPosition() const;
279 const VTKBaseStreamer & getStreamer() const;
280 bool isEnabled() const ;
281 bool hasAllMetaData() const ;
282
283 void setName( const std::string & ) ;
284 void setDataType( VTKDataType ) ;
285 void setFieldType( VTKFieldType ) ;
286 void setLocation( VTKLocation ) ;
287 void setCodification( VTKFormat ) ;
288 void setOffset( uint64_t ) ;
289 void setPosition( std::fstream::pos_type ) ;
291 void enable() ;
292 void disable() ;
293
294 void read( std::fstream&, uint64_t, uint8_t ) const ;
295 void write( std::fstream& ) const ;
296};
297
298class VTK{
299
300 protected:
302 uint64_t m_points ;
303 uint64_t m_cells ;
304 uint16_t m_procs ;
305 uint16_t m_rank ;
307 std::string m_headerType ;
309 std::vector<VTKField> m_geometry ;
312 std::vector<VTKField> m_data ;
317 public:
318 VTK( );
319 VTK( const std::string &, const std::string & );
320 virtual ~VTK( ) = default;
321
322 void setHeaderType( const std::string & );
323 const std::string & getHeaderType( ) const;
324
325 std::string getName( ) const;
326 std::string getDirectory( ) const;
327 int getCounter( ) const;
328
329 void setNames( const std::string & , const std::string & ) ;
330 void setName( const std::string & ) ;
331 void setDirectory( const std::string & ) ;
332 void setCounter( int c_=0 ) ;
333 void incrementCounter( ) ;
334 int unsetCounter( ) ;
335 void setParallel( uint16_t , uint16_t ) ;
336
337 void setCodex( VTKFormat );
338 void setGeomCodex( VTKFormat );
339 void setDataCodex( VTKFormat );
340
341 void setGeomData( VTKField &&field ) ;
342
343 VTKField& addData( VTKField &&field ) ;
344
345 VTKField& addData( const std::string &, VTKBaseStreamer* = nullptr ) ;
346
347 template<class T>
348 VTKField& addData( const std::string &, VTKFieldType, VTKLocation, VTKBaseStreamer* = nullptr ) ;
349
350 template<class T>
351 VTKField& addData( const std::string &, std::vector<T> & ) ;
352
353 template<class T>
354 VTKField& addData( const std::string &, VTKFieldType, VTKLocation, std::vector<T> & ) ;
355
356 void removeData( const std::string & ) ;
357 void enableData( const std::string & ) ;
358 void disableData( const std::string & ) ;
359
360 bool hasData( const std::string & ) const ;
361
362 std::vector<VTKField>::const_iterator getDataBegin( ) const ;
363 std::vector<VTKField>::const_iterator getDataEnd( ) const ;
364
365 std::vector<VTKField>::const_iterator getGeomDataBegin( ) const ;
366 std::vector<VTKField>::const_iterator getGeomDataEnd( ) const ;
367
368 std::size_t getDataCount( ) const ;
369 std::size_t getGeomDataCount( ) const ;
370
371 const VTKField * findData( const std::string & name ) const ;
372 const VTKField * findGeomData( const std::string & name ) const ;
373
374 void read() ;
375 virtual void readMetaInformation() = 0 ;
376 void readData() ;
377
378 void write( VTKWriteMode writeMode, double time ) ;
379 void write( VTKWriteMode writeMode=VTKWriteMode::DEFAULT ) ;
380 void write( const std::string &, VTKWriteMode writeMode, double time ) ;
381 void write( const std::string &, VTKWriteMode writeMode=VTKWriteMode::NO_INCREMENT ) ;
382
383 void writeCollection( ) const ;
384 void writeCollection( const std::string &outputName ) const ;
385 virtual void writeCollection( const std::string &outputName, const std::string &collectionName ) const = 0;
386
387 void writeTimeSeries( double time ) const;
388 void writeTimeSeries( const std::string &outputName, double time ) const ;
389 void writeTimeSeries( const std::string &outputName, const std::string &seriesName, double time ) const ;
390
391 protected:
392 //For Writing
393 virtual void writeMetaInformation() const = 0 ;
394 void writeData() ;
395
396 void writeDataHeader( std::fstream &, bool parallel=false ) const ;
397 void writeDataArray( std::fstream &, const VTKField &) const ;
398 void writePDataArray( std::fstream &, const VTKField &) const ;
399
400 FileHandler createCollectionHandler( const std::string &collectionName) const ;
401
402 //For Reading
403 void readDataHeader( std::fstream &) ;
404 bool readDataArray( std::fstream &, VTKField &) const ;
405
406 //General Purpose
407 std::vector<std::string> getFieldNames( const std::vector<VTKField> &fields ) const;
408
409 VTKField * getData( std::size_t id );
410 const VTKField * getData( std::size_t id ) const ;
411 VTKField * getGeomData( std::size_t id );
412 const VTKField * getGeomData( std::size_t id ) const ;
413
414 VTKField * _findData( const std::string & name);
415 VTKField * _findGeomData( const std::string & name);
416
417 int _findFieldIndex( const std::string &name, const std::vector<VTKField> &fields ) const;
418
419 bool isAppendedActive() const;
420 bool isASCIIActive() const;
421
422 void calcAppendedOffsets() ;
423 virtual uint64_t calcFieldSize( const VTKField &) const = 0;
424 virtual uint64_t calcFieldEntries( const VTKField &) const = 0;
425 virtual uint8_t calcFieldComponents( const VTKField &) const = 0;
426 void checkAllFields();
427
428 virtual std::string getExtension() const = 0 ;
429 virtual std::string getCollectionExtension() const;
430
431};
432
433class VTKUnstructuredGrid : public VTK {
434
435 protected:
437
438 private:
439 VTKElementType m_type ;
440 uint64_t m_nCells ;
441 const VTKField *m_connectivity ;
443 void flushData( std::fstream &, const std::string &, VTKFormat) override ;
444
445 public:
447 void setCellCount( uint64_t) ;
448 void setConnectivityField( const VTKField *connectivity) ;
449 };
450
456 public:
457 VTKUnstructuredGrid( VTKElementType elementType = VTKElementType::UNDEFINED );
458 VTKUnstructuredGrid( const std::string &, const std::string &, VTKElementType elementType = VTKElementType::UNDEFINED ) ;
459
460 protected:
461 uint64_t readConnectivityEntries( ) ;
462 uint64_t readFaceStreamEntries( ) ;
463
465
466 std::string getExtension() const override;
467
468 public:
469 using VTK::setGeomData;
470
471 void readMetaInformation() override ;
472 void writeMetaInformation() const override ;
473
474 void writeCollection( const std::string &outputName, const std::string &collectionName ) const override ;
475
476 void setDimensions( uint64_t , uint64_t , uint64_t nconn = 0 , uint64_t nfacestream = 0 ) ;
477
478 template<class T>
479 void setGeomData( VTKUnstructuredField, std::vector<T> & ) ;
481
482 template<class T>
484
485 uint64_t calcConnectivityEntries( ) const ;
486 uint64_t calcFieldSize( const VTKField &) const override ;
487 uint64_t calcFieldEntries( const VTKField &) const override ;
488 uint8_t calcFieldComponents( const VTKField &) const override ;
489
490 private:
491 int getFieldGeomId( VTKUnstructuredField field ) const ;
492
493};
494
495class VTKRectilinearGrid : public VTK{
496
497 typedef std::array<std::array<int,2>,2> extension2D_t ;
498 typedef std::array<std::array<int,2>,3> extension3D_t ;
500 protected:
502 extension3D_t m_localIndex ;
503 extension3D_t m_globalIndex ;
504 std::vector<extension3D_t> m_procIndex ;
506 public:
508 VTKRectilinearGrid( const std::string & , const std::string & );
509 VTKRectilinearGrid( const std::string & , const std::string & , VTKFormat, int, int, int, int, int, int );
510 VTKRectilinearGrid( const std::string & , const std::string & , VTKFormat, int, int, int );
511 VTKRectilinearGrid( const std::string & , const std::string & , VTKFormat, int, int, int, int );
512 VTKRectilinearGrid( const std::string & , const std::string & , VTKFormat, int, int );
513
514 protected:
515 std::string getExtension() const override;
516
517 public:
518 using VTK::setGeomData;
519
520 void readMetaInformation() override ;
521 void writeMetaInformation() const override ;
522
523 void writeCollection( const std::string &outputName, const std::string &collectionName ) const override ;
524
525 void setDimensions( int, int, int, int, int, int ) ;
526 void setDimensions( int, int, int ) ;
527 void setDimensions( int, int, int, int ) ;
528 void setDimensions( int, int ) ;
529
530 void setGlobalDimensions( int, int, int ) ;
531 void setGlobalDimensions( int, int ) ;
532
533 template<class T>
534 void setGeomData( VTKRectilinearField, std::vector<T> & ) ;
536
537 template<class T>
539
540 void setGlobalIndex( const std::vector<extension3D_t> & ) ;
541 void setGlobalIndex( const std::vector<extension2D_t> & ) ;
542
543 uint64_t calcFieldSize( const VTKField &) const override ;
544 uint64_t calcFieldEntries( const VTKField &) const override ;
545 uint8_t calcFieldComponents( const VTKField &) const override ;
546};
547
552namespace vtk{
554
555 std::string convertDataArrayToString( const VTKField & ) ;
556 std::string convertPDataArrayToString( const VTKField & ) ;
557
558 bool convertStringToDataArray( const std::string &, VTKField &) ;
559
560 std::string convertEnumToString( VTKLocation ) ;
561 std::string convertEnumToString( VTKFormat ) ;
562 std::string convertEnumToString( VTKDataType ) ;
563
564 bool convertStringToEnum( const std::string &, VTKLocation & ) ;
565 bool convertStringToEnum( const std::string &, VTKFormat & ) ;
566 bool convertStringToEnum( const std::string &, VTKDataType &) ;
567
568 template<class T>
569 void allocate( std::vector<T> &, int) ;
570
571 template<class T>
572 void allocate( T &, int) ;
573}
574
575}
576
577#include"VTK.tpp"
578#include"VTKTypes.tpp"
579#include"VTKStreamer.tpp"
580#include"VTKUtils.tpp"
581
582
583#endif
Creates file names and checks status.
An interface class to all containers that are batively supported by VTK.
Definition VTK.hpp:180
The base class to be used to derive VTK streamers form.
Definition VTK.hpp:209
void flushValue(std::fstream &, VTKFormat, const T &value) const
virtual void flushData(std::fstream &, const std::string &, VTKFormat)
virtual void absorbData(std::fstream &, const std::string &, VTKFormat, uint64_t, uint8_t, VTKDataType)
VTKField handles geometry and data field information for the VTK format.
Definition VTK.hpp:247
VTKFieldType m_fieldType
Definition VTK.hpp:252
const std::string & getName() const
Definition VTKField.cpp:167
VTKFormat getCodification() const
Definition VTKField.cpp:208
VTKDataType m_dataType
Definition VTK.hpp:253
bool isEnabled() const
Definition VTKField.cpp:241
std::fstream::pos_type getPosition() const
Definition VTKField.cpp:225
const VTKBaseStreamer & getStreamer() const
Definition VTKField.cpp:233
uint64_t getOffset() const
Definition VTKField.cpp:216
void setDataType(VTKDataType)
Definition VTKField.cpp:97
VTKBaseStreamer * m_streamer
Definition VTK.hpp:258
VTKFieldType getFieldType() const
Definition VTKField.cpp:175
VTKLocation getLocation() const
Definition VTKField.cpp:200
VTKFormat m_codification
Definition VTK.hpp:255
unsigned getComponentCount() const
Definition VTKField.cpp:183
void setStreamer(VTKBaseStreamer &)
Definition VTKField.cpp:145
uint64_t m_offset
Definition VTK.hpp:256
void setCodification(VTKFormat)
Definition VTKField.cpp:113
void write(std::fstream &) const
Definition VTKField.cpp:268
std::string m_name
Definition VTK.hpp:251
void setFieldType(VTKFieldType)
Definition VTKField.cpp:121
bool hasAllMetaData() const
Definition VTKField.cpp:249
void read(std::fstream &, uint64_t, uint8_t) const
Definition VTKField.cpp:280
VTKDataType getDataType() const
Definition VTKField.cpp:192
VTKLocation m_location
Definition VTK.hpp:254
std::fstream::pos_type m_position
Definition VTK.hpp:257
void setOffset(uint64_t)
Definition VTKField.cpp:137
void setName(const std::string &)
Definition VTKField.cpp:89
void setPosition(std::fstream::pos_type)
Definition VTKField.cpp:129
void setLocation(VTKLocation)
Definition VTKField.cpp:105
In VTKNativeStreamer all instances of classes derived from VTKBaseConatiner are stored....
Definition VTK.hpp:226
void absorbData(std::fstream &, const std::string &, VTKFormat, uint64_t, uint8_t, VTKDataType) override
bool hasData(const std::string &) const
void flushData(std::fstream &, const std::string &, VTKFormat) override
VTKNativeStreamer & operator=(const VTKNativeStreamer &other)
void addData(const std::string &, std::vector< T > &)
void removeData(const std::string &)
VTK input output for Rectilinear Meshes.
Definition VTK.hpp:495
std::vector< extension3D_t > m_procIndex
Definition VTK.hpp:504
void setGlobalDimensions(int, int, int)
void setDimensions(int, int, int, int, int, int)
uint8_t calcFieldComponents(const VTKField &) const override
extension3D_t m_globalIndex
Definition VTK.hpp:503
void setGeomData(VTKRectilinearField, std::vector< T > &)
Definition VTK.tpp:131
std::string getExtension() const override
void writeMetaInformation() const override
extension3D_t m_localIndex
Definition VTK.hpp:502
uint64_t calcFieldSize(const VTKField &) const override
void setGlobalIndex(const std::vector< extension3D_t > &)
uint64_t calcFieldEntries(const VTKField &) const override
void readMetaInformation() override
VTK data types handling.
Definition VTK.hpp:156
static VTKDataType whichType()
Definition VTKTypes.tpp:117
static uint8_t sizeOfType(VTKDataType type)
Definition VTKTypes.cpp:67
static VTKDataType registerType()
Definition VTKTypes.tpp:36
void setConnectivityField(const VTKField *connectivity)
VTK input output for Unstructured Meshes.
Definition VTK.hpp:433
uint8_t calcFieldComponents(const VTKField &) const override
HomogeneousInfoStreamer m_homogeneousInfoStreamer
Definition VTK.hpp:454
uint64_t calcConnectivityEntries() const
std::string getExtension() const override
uint64_t m_nConnectivityEntries
Definition VTK.hpp:451
uint64_t calcFieldSize(const VTKField &) const override
uint64_t m_nFaceStreamEntries
Definition VTK.hpp:452
void setElementType(VTKElementType)
void writeMetaInformation() const override
void setDimensions(uint64_t, uint64_t, uint64_t nconn=0, uint64_t nfacestream=0)
VTKUnstructuredGrid(VTKElementType elementType=VTKElementType::UNDEFINED)
VTKElementType m_elementType
Definition VTK.hpp:453
uint64_t calcFieldEntries(const VTKField &) const override
void setGeomData(VTKUnstructuredField, std::vector< T > &)
Definition VTK.tpp:94
Implementation of VTKBaseContainer in order to support natively std::vector in VTK.
Definition VTK.hpp:194
void absorbData(std::fstream &, VTKFormat, uint64_t, uint8_t) override
VTKVectorContainer * clone() const override
VTKVectorContainer(std::vector< T > &)
void flushData(std::fstream &, VTKFormat) override
void resize(std::true_type, uint64_t, uint8_t)
A base class for VTK input output.
Definition VTK.hpp:298
const VTKField * findGeomData(const std::string &name) const
Definition VTK.cpp:506
virtual std::string getCollectionExtension() const
Definition VTK.cpp:1495
uint16_t m_rank
Definition VTK.hpp:305
uint64_t m_cells
Definition VTK.hpp:303
void setDirectory(const std::string &)
Definition VTK.cpp:129
void incrementCounter()
Definition VTK.cpp:168
void setParallel(uint16_t, uint16_t)
Definition VTK.cpp:205
std::size_t getDataCount() const
Definition VTK.cpp:454
VTKField & addData(VTKField &&field)
Definition VTK.cpp:281
std::vector< VTKField > m_geometry
Definition VTK.hpp:309
int getCounter() const
Definition VTK.cpp:194
void readDataHeader(std::fstream &)
Definition VTK.cpp:1370
void setNames(const std::string &, const std::string &)
Definition VTK.cpp:108
int _findFieldIndex(const std::string &name, const std::vector< VTKField > &fields) const
Definition VTK.cpp:588
void writeTimeSeries(double time) const
Definition VTK.cpp:815
VTKField * _findData(const std::string &name)
Definition VTK.cpp:490
FileHandler m_fh
Definition VTK.hpp:301
VTKField * _findGeomData(const std::string &name)
Definition VTK.cpp:522
const VTKField * findData(const std::string &name) const
Definition VTK.cpp:474
std::vector< VTKField >::const_iterator getDataEnd() const
Definition VTK.cpp:427
void read()
Definition VTK.cpp:1234
void disableData(const std::string &)
Definition VTK.cpp:386
void writeData()
Definition VTK.cpp:947
void readData()
Definition VTK.cpp:1245
void setHeaderType(const std::string &)
Definition VTK.cpp:83
void removeData(const std::string &)
Definition VTK.cpp:330
std::string getName() const
Definition VTK.cpp:139
const std::string & getHeaderType() const
Definition VTK.cpp:99
std::string m_headerType
Definition VTK.hpp:307
bool readDataArray(std::fstream &, VTKField &) const
Definition VTK.cpp:1468
void setCounter(int c_=0)
Definition VTK.cpp:157
void writeCollection() const
Definition VTK.cpp:789
void enableData(const std::string &)
Definition VTK.cpp:368
void setCodex(VTKFormat)
Definition VTK.cpp:228
bool isASCIIActive() const
Definition VTK.cpp:619
void setDataCodex(VTKFormat)
Definition VTK.cpp:251
void setName(const std::string &)
Definition VTK.cpp:119
void setGeomCodex(VTKFormat)
Definition VTK.cpp:239
VTKFormat m_dataCodex
Definition VTK.hpp:313
std::vector< VTKField >::const_iterator getDataBegin() const
Definition VTK.cpp:418
void writePDataArray(std::fstream &, const VTKField &) const
Definition VTK.cpp:1224
int unsetCounter()
Definition VTK.cpp:182
VTKField * getGeomData(std::size_t id)
Definition VTK.cpp:573
bool hasData(const std::string &) const
Definition VTK.cpp:352
std::vector< VTKField >::const_iterator getGeomDataBegin() const
Definition VTK.cpp:436
VTKField * getData(std::size_t id)
Definition VTK.cpp:549
uint16_t m_procs
Definition VTK.hpp:304
void setGeomData(VTKField &&field)
Definition VTK.cpp:263
std::vector< VTKField >::const_iterator getGeomDataEnd() const
Definition VTK.cpp:445
void calcAppendedOffsets()
Definition VTK.cpp:635
void writeDataHeader(std::fstream &, bool parallel=false) const
Definition VTK.cpp:1146
bool isAppendedActive() const
Definition VTK.cpp:603
std::size_t getGeomDataCount() const
Definition VTK.cpp:463
void writeDataArray(std::fstream &, const VTKField &) const
Definition VTK.cpp:1212
void checkAllFields()
Definition VTK.cpp:675
VTKNativeStreamer m_nativeStreamer
Definition VTK.hpp:315
std::vector< VTKField > m_data
Definition VTK.hpp:312
FileHandler createCollectionHandler(const std::string &collectionName) const
Definition VTK.cpp:934
std::string getDirectory() const
Definition VTK.cpp:148
void write(VTKWriteMode writeMode, double time)
Definition VTK.cpp:704
uint64_t m_points
Definition VTK.hpp:302
std::vector< std::string > getFieldNames(const std::vector< VTKField > &fields) const
Definition VTK.cpp:403
VTKFormat m_geomCodex
Definition VTK.hpp:310
VTKFieldType
Definition VTK.hpp:63
VTKFormat
Definition VTK.hpp:92
VTKElementType
Definition VTK.hpp:112
VTKWriteMode
Definition VTK.hpp:49
VTKDataType
Definition VTK.hpp:74
VTKLocation
Definition VTK.hpp:102
VTKUnstructuredField
Definition VTK.hpp:137
VTKRectilinearField
Definition VTK.hpp:150
void allocate(std::vector< T > &, int)
Definition VTKUtils.tpp:20
std::string convertDataArrayToString(const VTKField &)
Definition VTKUtils.cpp:144
std::string convertEnumToString(VTKLocation)
Definition VTKUtils.cpp:191
bool convertStringToEnum(const std::string &, VTKLocation &)
Definition VTKUtils.cpp:265
std::string convertPDataArrayToString(const VTKField &)
Definition VTKUtils.cpp:171
bool convertStringToDataArray(const std::string &, VTKField &)
Definition VTKUtils.cpp:84
uint8_t getElementNodeCount(VTKElementType)
Definition VTKUtils.cpp:35
--- layout: doxygen_footer ---