VTK  9.2.6
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 
71 #ifndef vtkExodusIIWriter_h
72 #define vtkExodusIIWriter_h
73 
74 #include "vtkIOExodusModule.h" // For export macro
75 #include "vtkSmartPointer.h" // For vtkSmartPointer
76 #include "vtkWriter.h"
77 
78 #include <map> // STL Header
79 #include <string> // STL Header
80 #include <vector> // STL Header
81 
82 class vtkModelMetadata;
83 class vtkDoubleArray;
84 class vtkIntArray;
86 
87 class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
88 {
89 public:
90  static vtkExodusIIWriter* New();
91  vtkTypeMacro(vtkExodusIIWriter, vtkWriter);
92  void PrintSelf(ostream& os, vtkIndent indent) override;
93 
104  void SetModelMetadata(vtkModelMetadata*);
105  vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
106 
114  vtkSetFilePathMacro(FileName);
115  vtkGetFilePathMacro(FileName);
116 
124  vtkSetMacro(StoreDoubles, int);
125  vtkGetMacro(StoreDoubles, int);
126 
132  vtkSetMacro(GhostLevel, int);
133  vtkGetMacro(GhostLevel, int);
134 
141  vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
142  vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
143  vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
144 
151  vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
152  vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
153  vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
154 
161  vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
162  vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
163  vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
164 
170  vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
171  vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
172  vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
173 
174  vtkSetStringMacro(BlockIdArrayName);
175  vtkGetStringMacro(BlockIdArrayName);
176 
182  vtkSetMacro(IgnoreMetaDataWarning, bool);
183  vtkGetMacro(IgnoreMetaDataWarning, bool);
184  vtkBooleanMacro(IgnoreMetaDataWarning, bool);
185 
186 protected:
188  ~vtkExodusIIWriter() override;
189 
191 
193 
194  char* FileName;
195  int fid;
196 
198  int MyRank;
199 
201 
209 
214 
216  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
217  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
218 
219  std::vector<vtkStdString> FlattenedNames;
220  std::vector<vtkStdString> NewFlattenedNames;
221 
222  std::vector<vtkIntArray*> BlockIdList;
223 
224  struct Block
225  {
227  {
228  this->Name = nullptr;
229  this->Type = 0;
230  this->NumElements = 0;
231  this->ElementStartIndex = -1;
232  this->NodesPerElement = 0;
233  this->EntityCounts = std::vector<int>();
234  this->EntityNodeOffsets = std::vector<int>();
235  this->GridIndex = 0;
236  this->OutputIndex = -1;
237  this->NumAttributes = 0;
238  this->BlockAttributes = nullptr;
239  };
240  const char* Name;
241  int Type;
245  std::vector<int> EntityCounts;
246  std::vector<int> EntityNodeOffsets;
247  size_t GridIndex;
248  // std::vector<int> CellIndex;
251  float* BlockAttributes; // Owned by metamodel or null. Don't delete.
252  };
253  std::map<int, Block> BlockInfoMap;
254  int NumCells, NumPoints, MaxId;
255 
256  std::vector<vtkIdType*> GlobalElementIdList;
257  std::vector<vtkIdType*> GlobalNodeIdList;
258 
261 
263  {
265  int InIndex;
267  std::vector<std::string> OutNames;
268  };
269  std::map<std::string, VariableInfo> GlobalVariableMap;
270  std::map<std::string, VariableInfo> BlockVariableMap;
271  std::map<std::string, VariableInfo> NodeVariableMap;
275 
276  std::vector<std::vector<int>> CellToElementOffset;
277 
278  // By BlockId, and within block ID by element variable, with variables
279  // appearing in the same order in which they appear in OutputElementArrayNames
280 
283 
284  int BlockVariableTruthValue(int blockIdx, int varIdx);
285 
286  char* StrDupWithNew(const char* s);
287  void StringUppercase(std::string& str);
288 
290  vtkInformationVector* outputVector) override;
291 
292  int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
293  vtkInformationVector* outputVector);
294 
295  virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
296  vtkInformationVector* outputVector);
297 
298  int FillInputPortInformation(int port, vtkInformation* info) override;
299 
300  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
301  vtkInformationVector* outputVector) override;
302 
303  void WriteData() override;
304 
305  int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
306 
307  int CreateNewExodusFile();
308  void CloseExodusFile();
309 
310  int IsDouble();
311  void RemoveGhostCells();
312  int CheckParametersInternal(int numberOfProcesses, int myRank);
313  virtual int CheckParameters();
314  // If writing in parallel multiple time steps exchange after each time step
315  // if we should continue the execution. Pass local continueExecution as a
316  // parameter and return the global continueExecution.
317  virtual int GlobalContinueExecuting(int localContinueExecution);
318  int CheckInputArrays();
319  virtual void CheckBlockInfoMap();
320  int ConstructBlockInfoMap();
321  int ConstructVariableInfoMaps();
322  int ParseMetadata();
323  int CreateDefaultMetadata();
324  char* GetCellTypeName(int t);
325 
326  int CreateBlockIdMetadata(vtkModelMetadata* em);
327  int CreateBlockVariableMetadata(vtkModelMetadata* em);
328  int CreateSetsMetadata(vtkModelMetadata* em);
329 
330  void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
331  char** FlattenOutVariableNames(
332  int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
333  std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
334 
335  std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
336  std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
337 
338  vtkIdType GetNodeLocalId(vtkIdType id);
339  vtkIdType GetElementLocalId(vtkIdType id);
340  int GetElementType(vtkIdType id);
341 
342  int WriteInitializationParameters();
343  int WriteInformationRecords();
344  int WritePoints();
345  int WriteCoordinateNames();
346  int WriteGlobalPointIds();
347  int WriteBlockInformation();
348  int WriteGlobalElementIds();
349  int WriteVariableArrayNames();
350  int WriteNodeSetInformation();
351  int WriteSideSetInformation();
352  int WriteProperties();
353  int WriteNextTimeStep();
354  vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
355  static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
356 
357  double ExtractGlobalData(const char* name, int comp, int ts);
358  int WriteGlobalData(int timestep, vtkDataArray* buffer);
359  void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
360  int WriteCellData(int timestep, vtkDataArray* buffer);
361  void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
362  int WritePointData(int timestep, vtkDataArray* buffer);
363 
368  virtual unsigned int GetMaxNameLength();
369 
370 private:
371  vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
372  void operator=(const vtkExodusIIWriter&) = delete;
373 };
374 
375 #endif
vtkTypeBool WriteOutBlockIdArray
std::map< std::string, VariableInfo > BlockVariableMap
std::map< std::string, VariableInfo > NodeVariableMap
int * BlockElementVariableTruthTable
Store vtkAlgorithm input/output information.
std::map< vtkIdType, vtkIdType > * LocalElementIdMap
vtkTypeBool WriteAllTimeSteps
std::vector< vtkIdType * > GlobalNodeIdList
int vtkIdType
Definition: vtkType.h:332
vtkDataObject * OriginalInput
std::vector< vtkIntArray * > BlockIdList
dynamic, self-adjusting array of double
int vtkTypeBool
Definition: vtkABI.h:69
abstract class to write data to file(s)
Definition: vtkWriter.h:45
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:45
a simple class to control print indentation
Definition: vtkIndent.h:39
vtkTypeBool WriteOutGlobalNodeIdArray
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > NewFlattenedInput
std::vector< std::vector< int > > CellToElementOffset
dataset represents arbitrary combinations of all possible cell types
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
std::vector< int > EntityCounts
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.
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > FlattenedInput
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.
std::vector< vtkStdString > NewFlattenedNames
Store zero or more vtkInformation instances.
static vtkAlgorithm * New()
std::vector< vtkStdString > FlattenedNames
virtual void WriteData()=0
std::map< std::string, VariableInfo > GlobalVariableMap
general representation of visualization data
Definition: vtkDataObject.h:65
std::map< vtkIdType, vtkIdType > * LocalNodeIdMap
vtkTypeBool WriteOutGlobalElementIdArray
vtkModelMetadata * ModelMetadata
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Upstream/Downstream requests form the generalized interface through which executives invoke a algorit...