VTK
vtkExodusIIWriter.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkExodusIIWriter.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
68 #ifndef vtkExodusIIWriter_h
69 #define vtkExodusIIWriter_h
70 
71 #include "vtkIOExodusModule.h" // For export macro
72 #include "vtkWriter.h"
73 #include "vtkSmartPointer.h" // For vtkSmartPointer
74 
75 #include <vector> // STL Header
76 #include <map> // STL Header
77 #include <string> // STL Header
78 
79 class vtkModelMetadata;
80 class vtkDoubleArray;
81 class vtkIntArray;
83 
84 class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
85 {
86 public:
87  static vtkExodusIIWriter *New ();
89  void PrintSelf (ostream& os, vtkIndent indent) VTK_OVERRIDE;
90 
101  void SetModelMetadata (vtkModelMetadata*);
102  vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
103 
111  vtkSetStringMacro(FileName);
112  vtkGetStringMacro(FileName);
113 
121  vtkSetMacro(StoreDoubles, int);
122  vtkGetMacro(StoreDoubles, int);
123 
129  vtkSetMacro(GhostLevel, int);
130  vtkGetMacro(GhostLevel, int);
131 
138  vtkSetMacro(WriteOutBlockIdArray, int);
139  vtkGetMacro(WriteOutBlockIdArray, int);
140  vtkBooleanMacro(WriteOutBlockIdArray, int);
141 
148  vtkSetMacro(WriteOutGlobalNodeIdArray, int);
149  vtkGetMacro(WriteOutGlobalNodeIdArray, int);
150  vtkBooleanMacro(WriteOutGlobalNodeIdArray, int);
151 
158  vtkSetMacro(WriteOutGlobalElementIdArray, int);
159  vtkGetMacro(WriteOutGlobalElementIdArray, int);
160  vtkBooleanMacro(WriteOutGlobalElementIdArray, int);
161 
167  vtkSetMacro(WriteAllTimeSteps, int);
168  vtkGetMacro(WriteAllTimeSteps, int);
169  vtkBooleanMacro(WriteAllTimeSteps, int);
170 
171  vtkSetStringMacro(BlockIdArrayName);
172  vtkGetStringMacro(BlockIdArrayName);
173 
179  vtkSetMacro(IgnoreMetaDataWarning, bool);
180  vtkGetMacro(IgnoreMetaDataWarning, bool);
181  vtkBooleanMacro(IgnoreMetaDataWarning, bool);
182 
183 protected:
185  ~vtkExodusIIWriter () VTK_OVERRIDE;
186 
187  vtkModelMetadata* ModelMetadata;
188 
189  char *BlockIdArrayName;
190 
191  char *FileName;
192  int fid;
193 
194  int NumberOfProcesses;
195  int MyRank;
196 
197  int PassDoubles;
198 
199  int StoreDoubles;
200  int GhostLevel;
201  int WriteOutBlockIdArray;
202  int WriteOutGlobalNodeIdArray;
203  int WriteOutGlobalElementIdArray;
204  int WriteAllTimeSteps;
205  int NumberOfTimeSteps;
206 
207  int CurrentTimeIndex;
208  int FileTimeOffset;
209  bool TopologyChanged;
210  bool IgnoreMetaDataWarning;
211 
212  vtkDataObject *OriginalInput;
213  std::vector< vtkSmartPointer<vtkUnstructuredGrid> > FlattenedInput;
214  std::vector< vtkSmartPointer<vtkUnstructuredGrid> > NewFlattenedInput;
215 
216  std::vector< vtkStdString > FlattenedNames;
217  std::vector< vtkStdString > NewFlattenedNames;
218 
219  std::vector< vtkIntArray* > BlockIdList;
220 
221  struct Block
222  {
223  Block ()
224  {
225  this->Name = 0;
226  this->Type = 0;
227  this->NumElements = 0;
228  this->ElementStartIndex = -1;
229  this->NodesPerElement = 0;
230  this->EntityCounts = std::vector<int>();
231  this->EntityNodeOffsets = std::vector<int>();
232  this->GridIndex = 0;
233  this->OutputIndex = -1;
234  this->NumAttributes = 0;
235  this->BlockAttributes = 0;
236  };
237  const char *Name;
238  int Type;
242  std::vector<int> EntityCounts;
243  std::vector<int> EntityNodeOffsets;
244  size_t GridIndex;
245  // std::vector<int> CellIndex;
248  float *BlockAttributes; // Owned by metamodel or null. Don't delete.
249  };
250  std::map<int, Block> BlockInfoMap;
251  int NumCells, NumPoints, MaxId;
252 
253  std::vector<vtkIdType*> GlobalElementIdList;
254  std::vector<vtkIdType*> GlobalNodeIdList;
255 
258 
260  {
262  int InIndex;
264  std::vector<std::string> OutNames;
265  };
266  std::map<std::string, VariableInfo> GlobalVariableMap;
267  std::map<std::string, VariableInfo> BlockVariableMap;
268  std::map<std::string, VariableInfo> NodeVariableMap;
272 
273  std::vector< std::vector<int> > CellToElementOffset;
274 
275  // By BlockId, and within block ID by element variable, with variables
276  // appearing in the same order in which they appear in OutputElementArrayNames
277 
280 
281  int BlockVariableTruthValue(int blockIdx, int varIdx);
282 
283  char *StrDupWithNew (const char *s);
284  void StringUppercase (std::string& str);
285 
286  int ProcessRequest (vtkInformation* request,
287  vtkInformationVector** inputVector,
288  vtkInformationVector* outputVector) VTK_OVERRIDE;
289 
290  int RequestInformation (vtkInformation* request,
291  vtkInformationVector** inputVector,
292  vtkInformationVector* outputVector);
293 
294  virtual int RequestUpdateExtent (vtkInformation* request,
295  vtkInformationVector** inputVector,
296  vtkInformationVector* outputVector);
297 
298  int FillInputPortInformation (int port, vtkInformation* info) VTK_OVERRIDE;
299 
300  int RequestData (vtkInformation* request,
301  vtkInformationVector** inputVector,
302  vtkInformationVector* outputVector) VTK_OVERRIDE;
303 
304  void WriteData () VTK_OVERRIDE;
305 
306  int FlattenHierarchy (vtkDataObject* input, const char *name, bool& changed);
307 
308  int CreateNewExodusFile ();
309  void CloseExodusFile ();
310 
311  int IsDouble ();
312  void RemoveGhostCells ();
313  int CheckParametersInternal (int NumberOfProcesses, int MyRank);
314  virtual int CheckParameters ();
315  // If writing in parallel multiple time steps exchange after each time step
316  // if we should continue the execution. Pass local continueExecution as a
317  // parameter and return the global continueExecution.
318  virtual int GlobalContinueExecuting(int localContinueExecution);
319  int CheckInputArrays ();
320  virtual void CheckBlockInfoMap();
321  int ConstructBlockInfoMap ();
322  int ConstructVariableInfoMaps ();
323  int ParseMetadata ();
324  int CreateDefaultMetadata ();
325  char *GetCellTypeName (int t);
326 
327  int CreateBlockIdMetadata(vtkModelMetadata *em);
328  int CreateBlockVariableMetadata (vtkModelMetadata* em);
329  int CreateSetsMetadata (vtkModelMetadata* em);
330 
331  void ConvertVariableNames (std::map<std::string, VariableInfo>& variableMap);
332  char **FlattenOutVariableNames (
333  int nScalarArrays,
334  const std::map<std::string, VariableInfo>& variableMap);
335  std::string CreateNameForScalarArray (const char *root,
336  int component,
337  int numComponents);
338 
339  std::map<vtkIdType, vtkIdType> *LocalNodeIdMap;
340  std::map<vtkIdType, vtkIdType> *LocalElementIdMap;
341 
342  vtkIdType GetNodeLocalId(vtkIdType id);
343  vtkIdType GetElementLocalId(vtkIdType id);
344  int GetElementType(vtkIdType id);
345 
346  int WriteInitializationParameters ();
347  int WriteInformationRecords ();
348  int WritePoints ();
349  int WriteCoordinateNames ();
350  int WriteGlobalPointIds ();
351  int WriteBlockInformation ();
352  int WriteGlobalElementIds ();
353  int WriteVariableArrayNames ();
354  int WriteNodeSetInformation ();
355  int WriteSideSetInformation ();
356  int WriteProperties ();
357  int WriteNextTimeStep ();
358  vtkIntArray* GetBlockIdArray (
359  const char* BlockIdArrayName, vtkUnstructuredGrid* input);
360  static bool SameTypeOfCells (vtkIntArray* cellToBlockId,
361  vtkUnstructuredGrid* input);
362 
363  double ExtractGlobalData (const char *name, int comp, int ts);
364  int WriteGlobalData (int timestep, vtkDataArray *buffer);
365  void ExtractCellData (const char *name, int comp, vtkDataArray *buffer);
366  int WriteCellData (int timestep, vtkDataArray *buffer);
367  void ExtractPointData (const char *name, int comp, vtkDataArray *buffer);
368  int WritePointData (int timestep, vtkDataArray *buffer);
369 
370 private:
371  vtkExodusIIWriter (const vtkExodusIIWriter&) VTK_DELETE_FUNCTION;
372  void operator= (const vtkExodusIIWriter&) VTK_DELETE_FUNCTION;
373 };
374 
375 #endif
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:41
std::map< std::string, VariableInfo > BlockVariableMap
std::map< std::string, VariableInfo > NodeVariableMap
int * BlockElementVariableTruthTable
Store vtkAlgorithm input/output information.
std::vector< std::vector< int > > CellToElementOffset
int ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Upstream/Downstream requests form the generalized interface through which executives invoke a algorit...
std::vector< vtkIdType * > GlobalNodeIdList
Hold a reference to a vtkObjectBase instance.
int vtkIdType
Definition: vtkType.h:345
dynamic, self-adjusting array of double
abstract class to write data to file(s)
Definition: vtkWriter.h:42
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:39
a simple class to control print indentation
Definition: vtkIndent.h:33
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
std::vector< int > EntityCounts
vtkGetStringMacro(ExtensionsString)
Returns a string listing all available extensions.
std::vector< int > EntityNodeOffsets
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
virtual int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
vtkSetMacro(IgnoreDriverBugs, bool)
When set known driver bugs are ignored during driver feature detection.
std::vector< std::string > OutNames
This class encapsulates the metadata that appear in mesh-based file formats but do not appear in vtkU...
std::map< int, Block > BlockInfoMap
std::vector< vtkIdType * > GlobalElementIdList
Write Exodus II files.
Store zero or more vtkInformation instances.
static vtkAlgorithm * New()
virtual void WriteData()=0
vtkBooleanMacro(IgnoreDriverBugs, bool)
When set known driver bugs are ignored during driver feature detection.
std::map< std::string, VariableInfo > GlobalVariableMap
general representation of visualization data
Definition: vtkDataObject.h:58