VTK  9.3.1
vtkExodusIIWriter.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3 // SPDX-License-Identifier: BSD-3-Clause
4 
56 #ifndef vtkExodusIIWriter_h
57 #define vtkExodusIIWriter_h
58 
59 #include "vtkIOExodusModule.h" // For export macro
60 #include "vtkSmartPointer.h" // For vtkSmartPointer
61 #include "vtkWriter.h"
62 
63 #include <map> // STL Header
64 #include <string> // STL Header
65 #include <vector> // STL Header
66 
67 VTK_ABI_NAMESPACE_BEGIN
68 class vtkModelMetadata;
69 class vtkDoubleArray;
70 class vtkIntArray;
72 
73 class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
74 {
75 public:
76  static vtkExodusIIWriter* New();
77  vtkTypeMacro(vtkExodusIIWriter, vtkWriter);
78  void PrintSelf(ostream& os, vtkIndent indent) override;
79 
90  void SetModelMetadata(vtkModelMetadata*);
91  vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
92 
100  vtkSetFilePathMacro(FileName);
101  vtkGetFilePathMacro(FileName);
102 
110  vtkSetMacro(StoreDoubles, int);
111  vtkGetMacro(StoreDoubles, int);
112 
118  vtkSetMacro(GhostLevel, int);
119  vtkGetMacro(GhostLevel, int);
120 
127  vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
128  vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
129  vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
130 
137  vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
138  vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
139  vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
140 
147  vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
148  vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
149  vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
150 
156  vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
157  vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
158  vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
159 
160  vtkSetStringMacro(BlockIdArrayName);
161  vtkGetStringMacro(BlockIdArrayName);
162 
168  vtkSetMacro(IgnoreMetaDataWarning, bool);
169  vtkGetMacro(IgnoreMetaDataWarning, bool);
170  vtkBooleanMacro(IgnoreMetaDataWarning, bool);
171 
172 protected:
174  ~vtkExodusIIWriter() override;
175 
177 
179 
180  char* FileName;
181  int fid;
182 
184  int MyRank;
185 
187 
195 
200 
202  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
203  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
204 
205  std::vector<vtkStdString> FlattenedNames;
206  std::vector<vtkStdString> NewFlattenedNames;
207 
208  std::vector<vtkIntArray*> BlockIdList;
209 
210  struct Block
211  {
213  {
214  this->Name = nullptr;
215  this->Type = 0;
216  this->NumElements = 0;
217  this->ElementStartIndex = -1;
218  this->NodesPerElement = 0;
219  this->EntityCounts = std::vector<int>();
220  this->EntityNodeOffsets = std::vector<int>();
221  this->GridIndex = 0;
222  this->OutputIndex = -1;
223  this->NumAttributes = 0;
224  this->BlockAttributes = nullptr;
225  };
226  const char* Name;
227  int Type;
231  std::vector<int> EntityCounts;
232  std::vector<int> EntityNodeOffsets;
233  size_t GridIndex;
234  // std::vector<int> CellIndex;
237  float* BlockAttributes; // Owned by metamodel or null. Don't delete.
238  };
239  std::map<int, Block> BlockInfoMap;
240  int NumCells, NumPoints, MaxId;
241 
242  std::vector<vtkIdType*> GlobalElementIdList;
243  std::vector<vtkIdType*> GlobalNodeIdList;
244 
247 
249  {
251  int InIndex;
253  std::vector<std::string> OutNames;
254  };
255  std::map<std::string, VariableInfo> GlobalVariableMap;
256  std::map<std::string, VariableInfo> BlockVariableMap;
257  std::map<std::string, VariableInfo> NodeVariableMap;
261 
262  std::vector<std::vector<int>> CellToElementOffset;
263 
264  // By BlockId, and within block ID by element variable, with variables
265  // appearing in the same order in which they appear in OutputElementArrayNames
266 
269 
270  int BlockVariableTruthValue(int blockIdx, int varIdx);
271 
272  char* StrDupWithNew(const char* s);
273  void StringUppercase(std::string& str);
274 
276  vtkInformationVector* outputVector) override;
277 
278  int RequestInformation(vtkInformation* request, vtkInformationVector** inputVector,
279  vtkInformationVector* outputVector);
280 
281  virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
282  vtkInformationVector* outputVector);
283 
284  int FillInputPortInformation(int port, vtkInformation* info) override;
285 
286  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
287  vtkInformationVector* outputVector) override;
288 
289  void WriteData() override;
290 
291  int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
292 
293  int CreateNewExodusFile();
294  void CloseExodusFile();
295 
296  int IsDouble();
297  void RemoveGhostCells();
298  int CheckParametersInternal(int numberOfProcesses, int myRank);
299  virtual int CheckParameters();
300  // If writing in parallel multiple time steps exchange after each time step
301  // if we should continue the execution. Pass local continueExecution as a
302  // parameter and return the global continueExecution.
303  virtual int GlobalContinueExecuting(int localContinueExecution);
304  int CheckInputArrays();
305  virtual void CheckBlockInfoMap();
306  int ConstructBlockInfoMap();
307  int ConstructVariableInfoMaps();
308  int ParseMetadata();
309  int CreateDefaultMetadata();
310  char* GetCellTypeName(int t);
311 
312  int CreateBlockIdMetadata(vtkModelMetadata* em);
313  int CreateBlockVariableMetadata(vtkModelMetadata* em);
314  int CreateSetsMetadata(vtkModelMetadata* em);
315 
316  void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
317  char** FlattenOutVariableNames(
318  int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
319  std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
320 
321  std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
322  std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
323 
324  vtkIdType GetNodeLocalId(vtkIdType id);
325  vtkIdType GetElementLocalId(vtkIdType id);
326  int GetElementType(vtkIdType id);
327 
328  int WriteInitializationParameters();
329  int WriteInformationRecords();
330  int WritePoints();
331  int WriteCoordinateNames();
332  int WriteGlobalPointIds();
333  int WriteBlockInformation();
334  int WriteGlobalElementIds();
335  int WriteVariableArrayNames();
336  int WriteNodeSetInformation();
337  int WriteSideSetInformation();
338  int WriteProperties();
339  int WriteNextTimeStep();
340  vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
341  static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
342 
343  double ExtractGlobalData(const char* name, int comp, int ts);
344  int WriteGlobalData(int timestep, vtkDataArray* buffer);
345  void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
346  int WriteCellData(int timestep, vtkDataArray* buffer);
347  void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
348  int WritePointData(int timestep, vtkDataArray* buffer);
349 
354  virtual unsigned int GetMaxNameLength();
355 
356 private:
357  vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
358  void operator=(const vtkExodusIIWriter&) = delete;
359 };
360 
361 VTK_ABI_NAMESPACE_END
362 #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:315
vtkDataObject * OriginalInput
std::vector< vtkIntArray * > BlockIdList
dynamic, self-adjusting array of double
int vtkTypeBool
Definition: vtkABI.h:64
abstract class to write data to file(s)
Definition: vtkWriter.h:34
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:34
a simple class to control print indentation
Definition: vtkIndent.h:28
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:44
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:54
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...