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 
70 #ifndef vtkExodusIIWriter_h
71 #define vtkExodusIIWriter_h
72 
73 #include "vtkIOExodusModule.h" // For export macro
74 #include "vtkWriter.h"
75 #include "vtkSmartPointer.h" // For vtkSmartPointer
76 
77 #include <vector> // STL Header
78 #include <map> // STL Header
79 #include <string> // STL Header
80 
81 class vtkModelMetadata;
82 class vtkDoubleArray;
83 class vtkIntArray;
85 
87 {
88 public:
89  static vtkExodusIIWriter *New ();
91  void PrintSelf (ostream& os, vtkIndent indent);
92 
101  void SetModelMetadata (vtkModelMetadata*);
102  vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
103 
109  vtkSetStringMacro(FileName);
110  vtkGetStringMacro(FileName);
111 
117  vtkSetMacro(StoreDoubles, int);
118  vtkGetMacro(StoreDoubles, int);
119 
123  vtkSetMacro(GhostLevel, int);
124  vtkGetMacro(GhostLevel, int);
125 
131  vtkSetMacro(WriteOutBlockIdArray, int);
132  vtkGetMacro(WriteOutBlockIdArray, int);
133  vtkBooleanMacro(WriteOutBlockIdArray, int);
134 
139  vtkSetMacro(WriteOutGlobalNodeIdArray, int);
140  vtkGetMacro(WriteOutGlobalNodeIdArray, int);
141  vtkBooleanMacro(WriteOutGlobalNodeIdArray, int);
142 
147  vtkSetMacro(WriteOutGlobalElementIdArray, int);
148  vtkGetMacro(WriteOutGlobalElementIdArray, int);
149  vtkBooleanMacro(WriteOutGlobalElementIdArray, int);
150 
154  vtkSetMacro(WriteAllTimeSteps, int);
155  vtkGetMacro(WriteAllTimeSteps, int);
156  vtkBooleanMacro(WriteAllTimeSteps, int);
157 
158  vtkSetStringMacro(BlockIdArrayName);
159  vtkGetStringMacro(BlockIdArrayName);
160 
161 protected:
163  ~vtkExodusIIWriter ();
164 
166 
168 
169  char *FileName;
170  int fid;
171 
173  int MyRank;
174 
176 
184 
189 
190 //BTX
192  std::vector< vtkSmartPointer<vtkUnstructuredGrid> > FlattenedInput;
193  std::vector< vtkSmartPointer<vtkUnstructuredGrid> > NewFlattenedInput;
194 
195  std::vector< vtkIntArray* > BlockIdList;
196 
197  struct Block
198  {
199  Block ()
200  {
201  this->Type = 0;
202  this->NumElements = 0;
203  this->ElementStartIndex = -1;
204  this->NodesPerElement = 0;
205  this->EntityCounts = std::vector<int>();
206  this->EntityNodeOffsets = std::vector<int>();
207  this->GridIndex = 0;
208  this->OutputIndex = -1;
209  this->NumAttributes = 0;
210  this->BlockAttributes = 0;
211  };
212  int Type;
216  std::vector<int> EntityCounts;
217  std::vector<int> EntityNodeOffsets;
218  size_t GridIndex;
219  // std::vector<int> CellIndex;
222  float *BlockAttributes; // Owned by metamodel or null. Don't delete.
223  };
224  std::map<int, Block> BlockInfoMap;
225  int NumCells, NumPoints, MaxId;
226 
227  std::vector<vtkIdType*> GlobalElementIdList;
228  std::vector<vtkIdType*> GlobalNodeIdList;
229 //ETX
232 
233 //BTX
235  {
237  int InIndex;
239  std::vector<std::string> OutNames;
240  };
241  std::map<std::string, VariableInfo> GlobalVariableMap;
242  std::map<std::string, VariableInfo> BlockVariableMap;
243  std::map<std::string, VariableInfo> NodeVariableMap;
247 //ETX
248 
249 //BTX
250  std::vector< std::vector<int> > CellToElementOffset;
251 //ETX
252  // By BlockId, and within block ID by element variable, with variables
253  // appearing in the same order in which they appear in OutputElementArrayNames
254 
257 
258  int BlockVariableTruthValue(int blockIdx, int varIdx);
259 
260 //BTX
261  char *StrDupWithNew (const char *s);
262  void StringUppercase (std::string& str);
263 //ETX
264 
265  int ProcessRequest (vtkInformation* request,
266  vtkInformationVector** inputVector,
267  vtkInformationVector* outputVector);
268 
269  int RequestInformation (vtkInformation* request,
270  vtkInformationVector** inputVector,
271  vtkInformationVector* outputVector);
272 
273  virtual int RequestUpdateExtent (vtkInformation* request,
274  vtkInformationVector** inputVector,
275  vtkInformationVector* outputVector);
276 
278 
279  int RequestData (vtkInformation* request,
280  vtkInformationVector** inputVector,
281  vtkInformationVector* outputVector);
282 
283  void WriteData ();
284 
285  int FlattenHierarchy (vtkDataObject* input, bool& changed);
286 
287  int CreateNewExodusFile ();
288  void CloseExodusFile ();
289 
290  int IsDouble ();
291  void RemoveGhostCells ();
292  int CheckParametersInternal (int NumberOfProcesses, int MyRank);
293  virtual int CheckParameters ();
294  // If writing in parallel multiple time steps exchange after each time step
295  // if we should continue the execution. Pass local continueExecution as a
296  // parameter and return the global continueExecution.
297  virtual int GlobalContinueExecuting(int localContinueExecution);
298  int CheckInputArrays ();
299  virtual void CheckBlockInfoMap();
300  int ConstructBlockInfoMap ();
301  int ConstructVariableInfoMaps ();
302  int ParseMetadata ();
303  int CreateDefaultMetadata ();
304  char *GetCellTypeName (int t);
305 
306  int CreateBlockIdMetadata(vtkModelMetadata *em);
307  int CreateBlockVariableMetadata (vtkModelMetadata* em);
308  int CreateSetsMetadata (vtkModelMetadata* em);
309 
310 //BTX
311  void ConvertVariableNames (std::map<std::string, VariableInfo>& variableMap);
312  char **FlattenOutVariableNames (
313  int nScalarArrays,
314  const std::map<std::string, VariableInfo>& variableMap);
315  std::string CreateNameForScalarArray (const char *root,
316  int component,
317  int numComponents);
318 
319  std::map<vtkIdType, vtkIdType> *LocalNodeIdMap;
320  std::map<vtkIdType, vtkIdType> *LocalElementIdMap;
321 //ETX
322  vtkIdType GetNodeLocalId(vtkIdType id);
323  vtkIdType GetElementLocalId(vtkIdType id);
324  int GetElementType(vtkIdType id);
325 
326  int WriteInitializationParameters ();
327  int WriteInformationRecords ();
328  int WritePoints ();
329  int WriteCoordinateNames ();
330  int WriteGlobalPointIds ();
331  int WriteBlockInformation ();
332  int WriteGlobalElementIds ();
333  int WriteVariableArrayNames ();
334  int WriteNodeSetInformation ();
335  int WriteSideSetInformation ();
336  int WriteProperties ();
337  int WriteNextTimeStep ();
338  vtkIntArray* GetBlockIdArray (
339  const char* BlockIdArrayName, vtkUnstructuredGrid* input);
340  static bool SameTypeOfCells (vtkIntArray* cellToBlockId,
341  vtkUnstructuredGrid* input);
342 
343 //BTX
344  double ExtractGlobalData (const char *name, int comp, int ts);
345  int WriteGlobalData (int timestep, vtkDataArray *buffer);
346  void ExtractCellData (const char *name, int comp, vtkDataArray *buffer);
347  int WriteCellData (int timestep, vtkDataArray *buffer);
348  void ExtractPointData (const char *name, int comp, vtkDataArray *buffer);
349  int WritePointData (int timestep, vtkDataArray *buffer);
350 //ETX
351 
352 
353 private:
354  vtkExodusIIWriter (const vtkExodusIIWriter&); // Not Implemented
355  void operator= (const vtkExodusIIWriter&); // Not Implemented
356 };
357 
358 #endif
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > NewFlattenedInput
std::map< std::string, VariableInfo > BlockVariableMap
std::map< std::string, VariableInfo > NodeVariableMap
int * BlockElementVariableTruthTable
Store vtkAlgorithm input/output information.
std::map< vtkIdType, vtkIdType > * LocalElementIdMap
vtkDoubleArray * TimeValues
std::vector< std::vector< int > > CellToElementOffset
std::vector< vtkIdType * > GlobalNodeIdList
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > FlattenedInput
int vtkIdType
Definition: vtkType.h:247
virtual int ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
vtkDataObject * OriginalInput
dynamic, self-adjusting array of double
std::vector< vtkIntArray * > BlockIdList
abstract class to write data to file(s)
Definition: vtkWriter.h:44
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:49
a simple class to control print indentation
Definition: vtkIndent.h:38
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
std::vector< int > EntityCounts
#define VTKIOEXODUS_EXPORT
std::vector< int > EntityNodeOffsets
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
virtual int FillInputPortInformation(int port, vtkInformation *info)
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
std::map< std::string, VariableInfo > GlobalVariableMap
general representation of visualization data
Definition: vtkDataObject.h:64
std::map< vtkIdType, vtkIdType > * LocalNodeIdMap
void PrintSelf(ostream &os, vtkIndent indent)
vtkModelMetadata * ModelMetadata