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
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
67VTK_ABI_NAMESPACE_BEGIN
69class vtkDoubleArray;
70class vtkIntArray;
72
73class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
74{
75public:
78 void PrintSelf(ostream& os, vtkIndent indent) override;
79
91 vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
92
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
172protected:
175
177
179
180 char* FileName;
181 int fid;
182
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 {
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
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
287 vtkInformationVector* outputVector) override;
288
289 void WriteData() override;
290
291 int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
292
295
296 int IsDouble();
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);
305 virtual void CheckBlockInfoMap();
310 char* GetCellTypeName(int t);
311
315
316 void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
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
327
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
356private:
357 vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
358 void operator=(const vtkExodusIIWriter&) = delete;
359};
360
361VTK_ABI_NAMESPACE_END
362#endif
abstract superclass for arrays of numeric data
general representation of visualization data
dynamic, self-adjusting array of double
Write Exodus II files.
int WriteSideSetInformation()
std::vector< std::vector< int > > CellToElementOffset
void StringUppercase(std::string &str)
vtkIntArray * GetBlockIdArray(const char *BlockIdArrayName, vtkUnstructuredGrid *input)
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)
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
char * StrDupWithNew(const char *s)
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 ** FlattenOutVariableNames(int nScalarArrays, const std::map< std::string, VariableInfo > &variableMap)
void ExtractPointData(const char *name, int comp, vtkDataArray *buffer)
std::vector< vtkIdType * > GlobalElementIdList
virtual int CheckParameters()
vtkTypeBool WriteOutBlockIdArray
vtkGetFilePathMacro(FileName)
int WriteCoordinateNames()
static vtkExodusIIWriter * New()
void ExtractCellData(const char *name, int comp, vtkDataArray *buffer)
char * GetCellTypeName(int t)
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:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition vtkIntArray.h:35
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:32
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