VTK  9.3.0
vtkXMLReader.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
24#ifndef vtkXMLReader_h
25#define vtkXMLReader_h
26
27#include "vtkAlgorithm.h"
28#include "vtkIOXMLModule.h" // For export macro
29#include "vtkSmartPointer.h" // for vtkSmartPointer.
30
31#include <string> // for std::string
32
33VTK_ABI_NAMESPACE_BEGIN
36class vtkCommand;
37class vtkDataArray;
39class vtkDataSet;
44class vtkInformation;
45class vtkStringArray;
46
47class VTKIOXML_EXPORT vtkXMLReader : public vtkAlgorithm
48{
49public:
50 vtkTypeMacro(vtkXMLReader, vtkAlgorithm);
51 void PrintSelf(ostream& os, vtkIndent indent) override;
52
54 {
57 OTHER
58 };
59
61
67
69
72 vtkSetMacro(ReadFromInputString, vtkTypeBool);
73 vtkGetMacro(ReadFromInputString, vtkTypeBool);
74 vtkBooleanMacro(ReadFromInputString, vtkTypeBool);
75 void SetInputString(const std::string& s) { this->InputString = s; }
77
85 virtual int CanReadFile(VTK_FILEPATH const char* name);
86
88
94
96
100 vtkGetObjectMacro(PointDataArraySelection, vtkDataArraySelection);
101 vtkGetObjectMacro(CellDataArraySelection, vtkDataArraySelection);
102 vtkGetObjectMacro(ColumnArraySelection, vtkDataArraySelection);
104
106
113
115
119 const char* GetTimeDataArray(int idx) const;
120 vtkGetObjectMacro(TimeDataStringArray, vtkStringArray);
122
124
130 vtkGetStringMacro(ActiveTimeDataArrayName);
131 vtkSetStringMacro(ActiveTimeDataArrayName);
133
135
139 const char* GetPointArrayName(int index);
140 const char* GetCellArrayName(int index);
141 const char* GetColumnArrayName(int index);
143
145
149 int GetPointArrayStatus(const char* name);
150 int GetCellArrayStatus(const char* name);
151 void SetPointArrayStatus(const char* name, int status);
152 void SetCellArrayStatus(const char* name, int status);
153 int GetColumnArrayStatus(const char* name);
154 void SetColumnArrayStatus(const char* name, int status);
156
157 // For the specified port, copy the information this reader sets up in
158 // SetupOutputInformation to outInfo
159 virtual void CopyOutputInformation(vtkInformation* vtkNotUsed(outInfo), int vtkNotUsed(port)) {}
160
162
165 vtkSetMacro(TimeStep, int);
166 vtkGetMacro(TimeStep, int);
168
169 vtkGetMacro(NumberOfTimeSteps, int);
171
174 vtkGetVector2Macro(TimeStepRange, int);
175 vtkSetVector2Macro(TimeStepRange, int);
177
182 vtkXMLDataParser* GetXMLParser() { return this->XMLParser; }
183
185 vtkInformationVector* outputVector) override;
186
188
193 vtkGetObjectMacro(ReaderErrorObserver, vtkCommand);
195
197
202 vtkGetObjectMacro(ParserErrorObserver, vtkCommand);
204
205protected:
207 ~vtkXMLReader() override;
208
210
215 virtual int ReadXMLInformation();
216 virtual void ReadXMLData();
218
222 virtual const char* GetDataSetName() = 0;
223
227 virtual int CanReadFileVersion(int major, int minor);
228
232 virtual void SetupEmptyOutput() = 0;
233
237 virtual void SetupOutputInformation(vtkInformation* vtkNotUsed(outInfo)) {}
238
242 virtual void SetupOutputData();
243
248 virtual int ReadPrimaryElement(vtkXMLDataElement* ePrimary);
249
254 virtual int ReadVTKFile(vtkXMLDataElement* eVTKFile);
255
261 int GetLocalDataType(vtkXMLDataElement* da, int datatype);
262
268
274
280
282
285 virtual int OpenStream();
286 virtual void CloseStream();
287 virtual int OpenVTKFile();
288 virtual void CloseVTKFile();
289 virtual int OpenVTKString();
290 virtual void CloseVTKString();
291 virtual void CreateXMLParser();
292 virtual void DestroyXMLParser();
293 void SetupCompressor(const char* type);
294 int CanReadFileVersionString(const char* version);
296
302 virtual int CanReadFileWithDataType(const char* dsname);
303
307 vtkGetMacro(FileMajorVersion, int);
308
312 vtkGetMacro(FileMinorVersion, int);
313
315
318 int IntersectExtents(int* extent1, int* extent2, int* result);
319 int Min(int a, int b);
320 int Max(int a, int b);
321 void ComputePointDimensions(int* extent, int* dimensions);
322 void ComputePointIncrements(int* extent, vtkIdType* increments);
323 void ComputeCellDimensions(int* extent, int* dimensions);
324 void ComputeCellIncrements(int* extent, vtkIdType* increments);
325 vtkIdType GetStartTuple(int* extent, vtkIdType* increments, int i, int j, int k);
327 char** CreateStringArray(int numStrings);
328 void DestroyStringArray(int numStrings, char** strings);
330
337 virtual int ReadArrayValues(vtkXMLDataElement* da, vtkIdType arrayIndex, vtkAbstractArray* array,
338 vtkIdType startIndex, vtkIdType numValues, FieldType type = OTHER);
339
348 virtual int ReadArrayTuples(vtkXMLDataElement* da, vtkIdType arrayTupleIndex,
349 vtkAbstractArray* array, vtkIdType startTupleIndex, vtkIdType numTuples,
350 FieldType type = OTHER);
351
356
357 int SetFieldDataInfo(vtkXMLDataElement* eDSA, int association, vtkIdType numTuples,
358 vtkInformationVector*(&infoVector));
359
361
367
372 vtkObject* caller, unsigned long eid, void* clientdata, void* calldata);
373
379
390
391 // The vtkXMLDataParser instance used to hide XML reading details.
393
394 // The FieldData element representation.
396
397 // The input file's name.
398 char* FileName;
399
400 // The stream used to read the input.
401 istream* Stream;
402
403 // Whether this object is reading from a string or a file.
404 // Default is 0: read from file.
406
407 // The input string.
408 std::string InputString;
409
410 // The array selections.
415
421
427
428 // The observer to modify this object when the array selections are
429 // modified.
431
432 // Whether there was an error reading the file in RequestInformation.
434
435 // Whether there was an error reading the file in RequestData.
437
438 // incrementally fine-tuned progress updates.
439 virtual void GetProgressRange(float* range);
440 virtual void SetProgressRange(const float range[2], int curStep, int numSteps);
441 virtual void SetProgressRange(const float range[2], int curStep, const float* fractions);
442 virtual void UpdateProgressDiscrete(float progress);
443 float ProgressRange[2];
444
445 virtual int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
446 vtkInformationVector* outputVector);
447 virtual int RequestDataObject(vtkInformation* vtkNotUsed(request),
448 vtkInformationVector** vtkNotUsed(inputVector), vtkInformationVector* vtkNotUsed(outputVector))
449 {
450 return 1;
451 }
452 virtual int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
453 vtkInformationVector* outputVector);
455
456 // Whether there was an error reading the XML.
458
459 // For structured data keep track of dimensions empty of cells. For
460 // unstructured data these are always zero. This is used to support
461 // 1-D and 2-D cell data.
462 int AxesEmpty[3];
463
464 // The timestep currently being read.
468 void SetNumberOfTimeSteps(int num);
469 // buffer for reading timestep from the XML file the length is of
470 // NumberOfTimeSteps and therefore is always long enough
472 // Store the range of time steps
473 int TimeStepRange[2];
474
475 // Now we need to save what was the last time read for each kind of
476 // data to avoid rereading it that is to say we need a var for
477 // e.g. PointData/CellData/Points/Cells...
478 // See SubClass for details with member vars like PointsTimeStep/PointsOffset
479
480 // Helper function useful to know if a timestep is found in an array of timestep
481 static int IsTimeStepInArray(int timestep, int* timesteps, int length);
482
485
486 // Flag for whether DataProgressCallback should actually update
487 // progress.
489
491
493
494private:
495 // The stream used to read the input if it is in a file.
496 istream* FileStream;
497 // The stream used to read the input if it is in a string.
498 std::istringstream* StringStream;
499 int TimeStepWasReadOnce;
500
501 int FileMajorVersion;
502 int FileMinorVersion;
503
504 vtkDataObject* CurrentOutput;
505 vtkInformation* CurrentOutputInformation;
506
507 vtkXMLReader(const vtkXMLReader&) = delete;
508 void operator=(const vtkXMLReader&) = delete;
509
510 vtkCommand* ReaderErrorObserver;
511 vtkCommand* ParserErrorObserver;
512};
513
514VTK_ABI_NAMESPACE_END
515#endif
Abstract superclass for all arrays.
Superclass for all sources, filters, and sinks in VTK.
supports function callbacks
superclass for callback/observer methods
Definition vtkCommand.h:384
Store on/off settings for data arrays, etc.
abstract superclass for arrays of numeric data
general representation of visualization data
represent and manipulate attribute data in a dataset
abstract class to specify dataset behavior
Definition vtkDataSet.h:62
a simple class to control print indentation
Definition vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
abstract base class for most VTK objects
Definition vtkObject.h:58
Hold a reference to a vtkObjectBase instance.
a vtkAbstractArray subclass for strings
record modification and/or execution time
Represents an XML element and those nested inside.
Used by vtkXMLReader to parse VTK XML files.
Superclass for VTK's XML format readers.
virtual int ReadArrayTuples(vtkXMLDataElement *da, vtkIdType arrayTupleIndex, vtkAbstractArray *array, vtkIdType startTupleIndex, vtkIdType numTuples, FieldType type=OTHER)
Read an Array values starting at the given tuple index and up to numTuples taking into account the nu...
vtkXMLDataElement * FieldDataElement
vtkCallbackCommand * SelectionObserver
char ** CreateStringArray(int numStrings)
Utility methods for subclasses.
virtual void SetProgressRange(const float range[2], int curStep, int numSteps)
virtual void ConvertGhostLevelsToGhostType(FieldType, vtkAbstractArray *, vtkIdType, vtkIdType)
int GetColumnArrayStatus(const char *name)
Get/Set whether the point, cell, column or time array with the given name is to be read.
virtual int CanReadFileVersion(int major, int minor)
Test if the reader can read a file with the given version number.
virtual void GetProgressRange(float *range)
virtual void DestroyXMLParser()
Internal utility methods.
int CellDataArrayIsEnabled(vtkXMLDataElement *eCDA)
Check whether the given array element is an enabled array.
istream * Stream
void SetColumnArrayStatus(const char *name, int status)
Get/Set whether the point, cell, column or time array with the given name is to be read.
virtual int ReadArrayValues(vtkXMLDataElement *da, vtkIdType arrayIndex, vtkAbstractArray *array, vtkIdType startIndex, vtkIdType numValues, FieldType type=OTHER)
Read an Array values starting at the given index and up to numValues.
void SetReaderErrorObserver(vtkCommand *)
Set/get the ErrorObserver for the internal reader This is useful for applications that want to catch ...
const char * GetCellArrayName(int index)
Get the name of the point, cell, column or time array with the given index in the input.
vtkDataArraySelection * CellDataArraySelection
virtual void SetupOutputData()
Setup the output's data with allocation.
void ComputeCellIncrements(int *extent, vtkIdType *increments)
Utility methods for subclasses.
int CreateInformationKey(vtkXMLDataElement *eInfoKey, vtkInformation *info)
Create a vtkInformationKey from its corresponding XML representation.
void SetDataArraySelections(vtkXMLDataElement *eDSA, vtkDataArraySelection *sel)
Setup the data array selections for the input's set of arrays.
void SetParserErrorObserver(vtkCommand *)
Set/get the ErrorObserver for the internal xml parser This is useful for applications that want to ca...
virtual void SetupEmptyOutput()=0
Setup the output with no data available.
int GetNumberOfTimeDataArrays() const
Getters for time data array candidates.
virtual void SqueezeOutputArrays(vtkDataObject *)
Give concrete classes an option to squeeze any output arrays at the end of RequestData.
int GetNumberOfPointArrays()
Get the number of point, cell or column arrays available in the input.
int Max(int a, int b)
Utility methods for subclasses.
virtual void CreateXMLParser()
Internal utility methods.
vtkSmartPointer< vtkDataArray > TimeDataArray
Populated in ReadXMLInformation from the field data for the array chosen using ActiveTimeDataArrayNam...
virtual void CloseStream()
Internal utility methods.
int CanReadFileVersionString(const char *version)
Internal utility methods.
void ReadFieldData()
virtual void SetupOutputInformation(vtkInformation *vtkNotUsed(outInfo))
Setup the output's information.
vtkIdType GetStartTuple(int *extent, vtkIdType *increments, int i, int j, int k)
Utility methods for subclasses.
const char * GetColumnArrayName(int index)
Get the name of the point, cell, column or time array with the given index in the input.
vtkTypeBool ReadFromInputString
void SetupCompressor(const char *type)
Internal utility methods.
vtkInformation * GetCurrentOutputInformation()
vtkDataArraySelection * ColumnArraySelection
vtkXMLDataParser * XMLParser
virtual const char * GetDataSetName()=0
Get the name of the data set being read.
vtkXMLDataParser * GetXMLParser()
Returns the internal XML parser.
const char * GetPointArrayName(int index)
Get the name of the point, cell, column or time array with the given index in the input.
void SetInputString(const std::string &s)
Enable reading from an InputString instead of the default, a file.
void ComputePointIncrements(int *extent, vtkIdType *increments)
Utility methods for subclasses.
virtual int ReadXMLInformation()
Pipeline execution methods to be defined by subclass.
vtkAbstractArray * CreateArray(vtkXMLDataElement *da)
Create a vtkAbstractArray from its corresponding XML representation.
int GetNumberOfColumnArrays()
Get the number of point, cell or column arrays available in the input.
void ReadAttributeIndices(vtkXMLDataElement *eDSA, vtkDataSetAttributes *dsa)
Utility methods for subclasses.
int PointDataArrayIsEnabled(vtkXMLDataElement *ePDA)
Check whether the given array element is an enabled array.
vtkDataSet * GetOutputAsDataSet()
Get the output as a vtkDataSet pointer.
virtual void CloseVTKFile()
Internal utility methods.
virtual int OpenVTKFile()
Internal utility methods.
vtkTimeStamp ReadMTime
void ComputeCellDimensions(int *extent, int *dimensions)
Utility methods for subclasses.
virtual int CanReadFile(VTK_FILEPATH const char *name)
Test whether the file (type) with the given name can be read by this reader.
void SetNumberOfTimeSteps(int num)
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
static void SelectionModifiedCallback(vtkObject *caller, unsigned long eid, void *clientdata, void *calldata)
Callback registered with the SelectionObserver.
static int IsTimeStepInArray(int timestep, int *timesteps, int length)
virtual int RequestDataObject(vtkInformation *vtkNotUsed(request), vtkInformationVector **vtkNotUsed(inputVector), vtkInformationVector *vtkNotUsed(outputVector))
void ComputePointDimensions(int *extent, int *dimensions)
Utility methods for subclasses.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetNumberOfCellArrays()
Get the number of point, cell or column arrays available in the input.
void MarkIdTypeArrays(vtkXMLDataElement *da)
XML files have not consistently saved out adequate meta-data in past to correctly create vtkIdTypeArr...
~vtkXMLReader() override
virtual int CanReadFileWithDataType(const char *dsname)
This method is used by CanReadFile() to check if the reader can read an XML with the primary element ...
int IntersectExtents(int *extent1, int *extent2, int *result)
Utility methods for subclasses.
char * ActiveTimeDataArrayName
Name of the field-data array used to determine the time for the dataset being read.
virtual int ReadVTKFile(vtkXMLDataElement *eVTKFile)
Read the top-level element from the file.
virtual int OpenVTKString()
Internal utility methods.
virtual int ReadPrimaryElement(vtkXMLDataElement *ePrimary)
Read the primary element from the file.
vtkGetFilePathMacro(FileName)
Get/Set the name of the input file.
int GetLocalDataType(vtkXMLDataElement *da, int datatype)
If the IdType argument is present in the provided XMLDataElement and the provided dataType has the sa...
bool ReadInformation(vtkXMLDataElement *infoRoot, vtkInformation *info)
Populates the info object with the InformationKey children in infoRoot.
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Upstream/Downstream requests form the generalized interface through which executives invoke a algorit...
int GetCellArrayStatus(const char *name)
Get/Set whether the point, cell, column or time array with the given name is to be read.
vtkSetFilePathMacro(FileName)
Get/Set the name of the input file.
vtkDataObject * GetCurrentOutput()
int GetPointArrayStatus(const char *name)
Get/Set whether the point, cell, column or time array with the given name is to be read.
const char * GetTimeDataArray(int idx) const
Getters for time data array candidates.
virtual void UpdateProgressDiscrete(float progress)
virtual void ReadXMLData()
Pipeline execution methods to be defined by subclass.
virtual int OpenStream()
Internal utility methods.
vtkStringArray * TimeDataStringArray
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
virtual void SetProgressRange(const float range[2], int curStep, const float *fractions)
void DestroyStringArray(int numStrings, char **strings)
Utility methods for subclasses.
void SetCellArrayStatus(const char *name, int status)
Get/Set whether the point, cell, column or time array with the given name is to be read.
virtual void CloseVTKString()
Internal utility methods.
void SetPointArrayStatus(const char *name, int status)
Get/Set whether the point, cell, column or time array with the given name is to be read.
std::string InputString
vtkDataArraySelection * PointDataArraySelection
virtual void CopyOutputInformation(vtkInformation *vtkNotUsed(outInfo), int vtkNotUsed(port))
vtkDataSet * GetOutputAsDataSet(int index)
Get the output as a vtkDataSet pointer.
int SetFieldDataInfo(vtkXMLDataElement *eDSA, int association, vtkIdType numTuples, vtkInformationVector *(&infoVector))
int Min(int a, int b)
Utility methods for subclasses.
int vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
#define VTK_FILEPATH