VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkXMLUnstructuredDataWriter.h 00005 00006 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen 00007 All rights reserved. 00008 See Copyright.txt or http://www.kitware.com/Copyright.htm for details. 00009 00010 This software is distributed WITHOUT ANY WARRANTY; without even 00011 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 00012 PURPOSE. See the above copyright notice for more information. 00013 00014 =========================================================================*/ 00022 #ifndef vtkXMLUnstructuredDataWriter_h 00023 #define vtkXMLUnstructuredDataWriter_h 00024 00025 #include "vtkIOXMLModule.h" // For export macro 00026 #include "vtkXMLWriter.h" 00027 00028 class vtkPointSet; 00029 class vtkCellArray; 00030 class vtkCellIterator; 00031 class vtkDataArray; 00032 class vtkIdTypeArray; 00033 class vtkUnstructuredGrid; 00034 00035 class VTKIOXML_EXPORT vtkXMLUnstructuredDataWriter : public vtkXMLWriter 00036 { 00037 public: 00038 vtkTypeMacro(vtkXMLUnstructuredDataWriter,vtkXMLWriter); 00039 void PrintSelf(ostream& os, vtkIndent indent); 00040 00042 00044 vtkSetMacro(NumberOfPieces, int); 00045 vtkGetMacro(NumberOfPieces, int); 00047 00049 00051 vtkSetMacro(WritePiece, int); 00052 vtkGetMacro(WritePiece, int); 00054 00056 00057 vtkSetMacro(GhostLevel, int); 00058 vtkGetMacro(GhostLevel, int); 00060 00061 // See the vtkAlgorithm for a desciption of what these do 00062 int ProcessRequest(vtkInformation*, 00063 vtkInformationVector**, 00064 vtkInformationVector*); 00065 00066 protected: 00067 vtkXMLUnstructuredDataWriter(); 00068 ~vtkXMLUnstructuredDataWriter(); 00069 00070 vtkPointSet* GetInputAsPointSet(); 00071 virtual const char* GetDataSetName()=0; 00072 virtual void SetInputUpdateExtent(int piece, int numPieces, 00073 int ghostLevel); 00074 00075 virtual int WriteHeader(); 00076 virtual int WriteAPiece(); 00077 virtual int WriteFooter(); 00078 00079 virtual void AllocatePositionArrays(); 00080 virtual void DeletePositionArrays(); 00081 00082 virtual int WriteInlineMode(vtkIndent indent); 00083 virtual void WriteInlinePieceAttributes(); 00084 virtual void WriteInlinePiece(vtkIndent indent); 00085 00086 virtual void WriteAppendedPieceAttributes(int index); 00087 virtual void WriteAppendedPiece(int index, vtkIndent indent); 00088 virtual void WriteAppendedPieceData(int index); 00089 00090 void WriteCellsInline(const char* name, vtkCellIterator *cellIter, 00091 vtkIdType numCells, vtkIdType cellSizeEstimate, 00092 vtkIndent indent); 00093 00094 void WriteCellsInline(const char* name, vtkCellArray* cells, 00095 vtkDataArray* types, vtkIndent indent); 00096 00097 // New API with face infomration for polyhedron cell support. 00098 void WriteCellsInline(const char* name, vtkCellArray* cells, 00099 vtkDataArray* types, vtkIdTypeArray* faces, 00100 vtkIdTypeArray* faceOffsets, vtkIndent indent); 00101 00102 void WriteCellsInlineWorker(const char* name, vtkDataArray *types, 00103 vtkIndent indent); 00104 00105 void WriteCellsAppended(const char* name, vtkDataArray* types, 00106 vtkIndent indent, OffsetsManagerGroup *cellsManager); 00107 00108 void WriteCellsAppended(const char* name, vtkCellIterator *cellIter, 00109 vtkIdType numCells, vtkIndent indent, 00110 OffsetsManagerGroup *cellsManager); 00111 00112 void WriteCellsAppendedData(vtkCellArray* cells, vtkDataArray* types, 00113 int timestep, OffsetsManagerGroup *cellsManager); 00114 00115 void WriteCellsAppendedData(vtkCellIterator* cellIter, vtkIdType numCells, 00116 vtkIdType cellSizeEstimate, int timestep, 00117 OffsetsManagerGroup *cellsManager); 00118 00119 // New API with face infomration for polyhedron cell support. 00120 void WriteCellsAppendedData(vtkCellArray* cells, vtkDataArray* types, 00121 vtkIdTypeArray* faces,vtkIdTypeArray* faceOffsets, 00122 int timestep, OffsetsManagerGroup *cellsManager); 00123 00124 void WriteCellsAppendedDataWorker(vtkDataArray* types, int timestep, 00125 OffsetsManagerGroup *cellsManager); 00126 00127 void ConvertCells(vtkCellIterator* cellIter, vtkIdType numCells, 00128 vtkIdType cellSizeEstimate); 00129 00130 void ConvertCells(vtkCellArray* cells); 00131 00132 // For polyhedron support, conversion results are stored in Faces and FaceOffsets 00133 void ConvertFaces(vtkIdTypeArray* faces, vtkIdTypeArray* faceOffsets); 00134 00135 // Get the number of points/cells. Valid after Update has been 00136 // invoked on the input. 00137 virtual vtkIdType GetNumberOfInputPoints(); 00138 virtual vtkIdType GetNumberOfInputCells()=0; 00139 void CalculateDataFractions(float* fractions); 00140 void CalculateCellFractions(float* fractions, vtkIdType typesSize); 00141 00142 // Number of pieces used for streaming. 00143 int NumberOfPieces; 00144 00145 // Which piece to write, if not all. 00146 int WritePiece; 00147 00148 // The ghost level on each piece. 00149 int GhostLevel; 00150 00151 // Positions of attributes for each piece. 00152 vtkTypeInt64* NumberOfPointsPositions; 00153 00154 // For TimeStep support 00155 OffsetsManagerGroup *PointsOM; 00156 OffsetsManagerArray *PointDataOM; 00157 OffsetsManagerArray *CellDataOM; 00158 00159 // Hold the new cell representation arrays while writing a piece. 00160 vtkIdTypeArray* CellPoints; 00161 vtkIdTypeArray* CellOffsets; 00162 00163 int CurrentPiece; 00164 00165 // Hold the face arrays for polyhedron cells. 00166 vtkIdTypeArray* Faces; 00167 vtkIdTypeArray* FaceOffsets; 00168 00169 private: 00170 vtkXMLUnstructuredDataWriter(const vtkXMLUnstructuredDataWriter&); // Not implemented. 00171 void operator=(const vtkXMLUnstructuredDataWriter&); // Not implemented. 00172 }; 00173 00174 #endif