VTK  9.3.0
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 
65 #ifndef vtkExodusIIWriter_h
66 #define vtkExodusIIWriter_h
67 
68 #include "vtkIOExodusModule.h" // For export macro
69 #include "vtkSmartPointer.h" // For vtkSmartPointer
70 #include "vtkWriter.h"
71 
72 #include <map> // STL Header
73 #include <string> // STL Header
74 #include <vector> // STL Header
75 
76 VTK_ABI_NAMESPACE_BEGIN
77 class vtkModelMetadata;
78 class vtkDoubleArray;
79 class vtkIntArray;
81 
82 class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
83 {
84 public:
86  vtkTypeMacro(vtkExodusIIWriter, vtkWriter);
87  void PrintSelf(ostream& os, vtkIndent indent) override;
88 
100  vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
101 
111 
119  vtkSetMacro(StoreDoubles, int);
120  vtkGetMacro(StoreDoubles, int);
121 
127  vtkSetMacro(GhostLevel, int);
128  vtkGetMacro(GhostLevel, int);
129 
136  vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
137  vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
138  vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
139 
146  vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
147  vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
148  vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
149 
156  vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
157  vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
158  vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
159 
165  vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
166  vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
167  vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
168 
169  vtkSetStringMacro(BlockIdArrayName);
170  vtkGetStringMacro(BlockIdArrayName);
171 
177  vtkSetMacro(IgnoreMetaDataWarning, bool);
178  vtkGetMacro(IgnoreMetaDataWarning, bool);
179  vtkBooleanMacro(IgnoreMetaDataWarning, bool);
180 
181 protected:
183  ~vtkExodusIIWriter() override;
184 
186 
188 
189  char* FileName;
190  int fid;
191 
193  int MyRank;
194 
196 
204 
209 
211  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
212  std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
213 
214  std::vector<vtkStdString> FlattenedNames;
215  std::vector<vtkStdString> NewFlattenedNames;
216 
217  std::vector<vtkIntArray*> BlockIdList;
218 
219  struct Block
220  {
222  {
223  this->Name = nullptr;
224  this->Type = 0;
225  this->NumElements = 0;
226  this->ElementStartIndex = -1;
227  this->NodesPerElement = 0;
228  this->EntityCounts = std::vector<int>();
229  this->EntityNodeOffsets = std::vector<int>();
230  this->GridIndex = 0;
231  this->OutputIndex = -1;
232  this->NumAttributes = 0;
233  this->BlockAttributes = nullptr;
234  };
235  const char* Name;
236  int Type;
240  std::vector<int> EntityCounts;
241  std::vector<int> EntityNodeOffsets;
242  size_t GridIndex;
243  // std::vector<int> CellIndex;
246  float* BlockAttributes; // Owned by metamodel or null. Don't delete.
247  };
248  std::map<int, Block> BlockInfoMap;
249  int NumCells, NumPoints, MaxId;
250 
251  std::vector<vtkIdType*> GlobalElementIdList;
252  std::vector<vtkIdType*> GlobalNodeIdList;
253 
256 
258  {
260  int InIndex;
262  std::vector<std::string> OutNames;
263  };
264  std::map<std::string, VariableInfo> GlobalVariableMap;
265  std::map<std::string, VariableInfo> BlockVariableMap;
266  std::map<std::string, VariableInfo> NodeVariableMap;
270 
271  std::vector<std::vector<int>> CellToElementOffset;
272 
273  // By BlockId, and within block ID by element variable, with variables
274  // appearing in the same order in which they appear in OutputElementArrayNames
275 
278 
279  int BlockVariableTruthValue(int blockIdx, int varIdx);
280 
281  char* StrDupWithNew(const char* s);
283 
285  vtkInformationVector* outputVector) override;
286 
288  vtkInformationVector* outputVector);
289 
290  virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
291  vtkInformationVector* outputVector);
292 
294 
295  int RequestData(vtkInformation* request, vtkInformationVector** inputVector,
296  vtkInformationVector* outputVector) override;
297 
298  void WriteData() override;
299 
300  int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
301 
304 
305  int IsDouble();
307  int CheckParametersInternal(int numberOfProcesses, int myRank);
308  virtual int CheckParameters();
309  // If writing in parallel multiple time steps exchange after each time step
310  // if we should continue the execution. Pass local continueExecution as a
311  // parameter and return the global continueExecution.
312  virtual int GlobalContinueExecuting(int localContinueExecution);
314  virtual void CheckBlockInfoMap();
319  char* GetCellTypeName(int t);
320 
324 
325  void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
327  int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
328  std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
329 
330  std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
331  std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
332 
336 
339  int WritePoints();
349  vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
350  static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
351 
352  double ExtractGlobalData(const char* name, int comp, int ts);
353  int WriteGlobalData(int timestep, vtkDataArray* buffer);
354  void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
355  int WriteCellData(int timestep, vtkDataArray* buffer);
356  void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
357  int WritePointData(int timestep, vtkDataArray* buffer);
358 
363  virtual unsigned int GetMaxNameLength();
364 
365 private:
366  vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
367  void operator=(const vtkExodusIIWriter&) = delete;
368 };
369 
370 VTK_ABI_NAMESPACE_END
371 #endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
general representation of visualization data
Definition: vtkDataObject.h:64
dynamic, self-adjusting array of double
Write Exodus II files.
int WriteSideSetInformation()
std::vector< std::vector< int > > CellToElementOffset
void StringUppercase(std::string &str)
static vtkExodusIIWriter * New()
std::map< vtkIdType, vtkIdType > * LocalNodeIdMap
~vtkExodusIIWriter() override
void SetModelMetadata(vtkModelMetadata *)
Specify the vtkModelMetadata object which contains the Exodus file model information (metadata) absen...
int WriteVariableArrayNames()
std::map< std::string, VariableInfo > BlockVariableMap
int WriteNodeSetInformation()
vtkIdType GetNodeLocalId(vtkIdType id)
int BlockVariableTruthValue(int blockIdx, int varIdx)
int CheckParametersInternal(int numberOfProcesses, int myRank)
void ConvertVariableNames(std::map< std::string, VariableInfo > &variableMap)
virtual int GlobalContinueExecuting(int localContinueExecution)
char ** FlattenOutVariableNames(int nScalarArrays, const std::map< std::string, VariableInfo > &variableMap)
static bool SameTypeOfCells(vtkIntArray *cellToBlockId, vtkUnstructuredGrid *input)
int GetElementType(vtkIdType id)
int CreateBlockVariableMetadata(vtkModelMetadata *em)
int WriteBlockInformation()
vtkTypeBool WriteAllTimeSteps
std::vector< vtkStdString > NewFlattenedNames
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > NewFlattenedInput
vtkModelMetadata * ModelMetadata
int WritePointData(int timestep, vtkDataArray *buffer)
vtkIdType GetElementLocalId(vtkIdType id)
virtual void CheckBlockInfoMap()
virtual int RequestUpdateExtent(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
std::vector< vtkIntArray * > BlockIdList
vtkSetFilePathMacro(FileName)
Name for the output file.
int FlattenHierarchy(vtkDataObject *input, const char *name, bool &changed)
std::string CreateNameForScalarArray(const char *root, int component, int numComponents)
int WriteGlobalElementIds()
int CreateSetsMetadata(vtkModelMetadata *em)
vtkDataObject * OriginalInput
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
vtkTypeBool WriteOutGlobalNodeIdArray
std::map< std::string, VariableInfo > GlobalVariableMap
int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
int ConstructVariableInfoMaps()
std::map< std::string, VariableInfo > NodeVariableMap
double ExtractGlobalData(const char *name, int comp, int ts)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int ConstructBlockInfoMap()
virtual unsigned int GetMaxNameLength()
Get the maximum length name in the input data set.
int CreateDefaultMetadata()
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
int WriteCellData(int timestep, vtkDataArray *buffer)
int WriteInitializationParameters()
std::vector< vtkIdType * > GlobalNodeIdList
std::vector< vtkSmartPointer< vtkUnstructuredGrid > > FlattenedInput
void WriteData() override
std::map< vtkIdType, vtkIdType > * LocalElementIdMap
char * StrDupWithNew(const char *s)
char * GetCellTypeName(int t)
void ExtractPointData(const char *name, int comp, vtkDataArray *buffer)
std::vector< vtkIdType * > GlobalElementIdList
int * BlockElementVariableTruthTable
virtual int CheckParameters()
vtkTypeBool WriteOutBlockIdArray
vtkGetFilePathMacro(FileName)
int WriteCoordinateNames()
void ExtractCellData(const char *name, int comp, vtkDataArray *buffer)
vtkIntArray * GetBlockIdArray(const char *BlockIdArrayName, vtkUnstructuredGrid *input)
std::map< int, Block > BlockInfoMap
int WriteGlobalData(int timestep, vtkDataArray *buffer)
int CreateBlockIdMetadata(vtkModelMetadata *em)
vtkTypeBool WriteOutGlobalElementIdArray
vtkTypeBool ProcessRequest(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
Upstream/Downstream requests form the generalized interface through which executives invoke a algorit...
int WriteInformationRecords()
std::vector< vtkStdString > FlattenedNames
a simple class to control print indentation
Definition: vtkIndent.h:38
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:44
This class encapsulates the metadata that appear in mesh-based file formats but do not appear in vtkU...
dataset represents arbitrary combinations of all possible cell types
abstract class to write data to file(s)
Definition: vtkWriter.h:35
@ component
Definition: vtkX3D.h:175
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ name
Definition: vtkX3D.h:219
@ string
Definition: vtkX3D.h:490
std::vector< int > EntityNodeOffsets
std::vector< int > EntityCounts
std::vector< std::string > OutNames
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315