VTK  9.3.20240831
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
75#ifndef vtkExodusIIWriter_h
76#define vtkExodusIIWriter_h
77
78#include "vtkIOExodusModule.h" // For export macro
79#include "vtkSmartPointer.h" // For vtkSmartPointer
80#include "vtkWriter.h"
81
82#include <map> // STL Header
83#include <string> // STL Header
84#include <vector> // STL Header
85
86VTK_ABI_NAMESPACE_BEGIN
88class vtkDoubleArray;
89class vtkIntArray;
91
92class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
93{
94public:
97 void PrintSelf(ostream& os, vtkIndent indent) override;
98
110 vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
111
121
129 vtkSetMacro(StoreDoubles, int);
130 vtkGetMacro(StoreDoubles, int);
131
137 vtkSetMacro(GhostLevel, int);
138 vtkGetMacro(GhostLevel, int);
139
146 vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
147 vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
148 vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
149
156 vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
157 vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
158 vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
159
166 vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
167 vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
168 vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
169
175 vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
176 vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
177 vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
178
179 vtkSetStringMacro(BlockIdArrayName);
180 vtkGetStringMacro(BlockIdArrayName);
181
187 vtkSetMacro(IgnoreMetaDataWarning, bool);
188 vtkGetMacro(IgnoreMetaDataWarning, bool);
189 vtkBooleanMacro(IgnoreMetaDataWarning, bool);
190
191protected:
194
196
198
199 char* FileName;
200 int fid;
201
204
206
214
219
221 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
222 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
223
224 std::vector<vtkStdString> FlattenedNames;
225 std::vector<vtkStdString> NewFlattenedNames;
226
227 std::vector<vtkIntArray*> BlockIdList;
228
229 struct Block
230 {
232 {
233 this->Name = nullptr;
234 this->Type = 0;
235 this->NumElements = 0;
236 this->ElementStartIndex = -1;
237 this->NodesPerElement = 0;
238 this->EntityCounts = std::vector<int>();
239 this->EntityNodeOffsets = std::vector<int>();
240 this->GridIndex = 0;
241 this->OutputIndex = -1;
242 this->NumAttributes = 0;
243 this->BlockAttributes = nullptr;
244 }
245 const char* Name;
246 int Type;
250 std::vector<int> EntityCounts;
251 std::vector<int> EntityNodeOffsets;
252 size_t GridIndex;
253 // std::vector<int> CellIndex;
256 float* BlockAttributes; // Owned by metamodel or null. Don't delete.
257 };
258 std::map<int, Block> BlockInfoMap;
259 int NumCells, NumPoints, MaxId;
260
261 std::vector<vtkIdType*> GlobalElementIdList;
262 std::vector<vtkIdType*> GlobalNodeIdList;
263
266
268 {
272 std::vector<std::string> OutNames;
273 };
274 std::map<std::string, VariableInfo> GlobalVariableMap;
275 std::map<std::string, VariableInfo> BlockVariableMap;
276 std::map<std::string, VariableInfo> NodeVariableMap;
280
281 std::vector<std::vector<int>> CellToElementOffset;
282
283 // By BlockId, and within block ID by element variable, with variables
284 // appearing in the same order in which they appear in OutputElementArrayNames
285
288
289 int BlockVariableTruthValue(int blockIdx, int varIdx);
290
291 char* StrDupWithNew(const char* s);
292 void StringUppercase(std::string& str);
293
295 vtkInformationVector* outputVector) override;
296
298 vtkInformationVector* outputVector);
299
300 virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
301 vtkInformationVector* outputVector);
302
303 int FillInputPortInformation(int port, vtkInformation* info) override;
304
306 vtkInformationVector* outputVector) override;
307
308 void WriteData() override;
309
310 int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
311
314
315 int IsDouble();
317 int CheckParametersInternal(int numberOfProcesses, int myRank);
318 virtual int CheckParameters();
319 // If writing in parallel multiple time steps exchange after each time step
320 // if we should continue the execution. Pass local continueExecution as a
321 // parameter and return the global continueExecution.
322 virtual int GlobalContinueExecuting(int localContinueExecution);
324 virtual void CheckBlockInfoMap();
329 char* GetCellTypeName(int t);
330
334
335 void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
337 int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
338 std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
339
340 std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
341 std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
342
346
359 vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
360 static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
361
362 double ExtractGlobalData(const char* name, int comp, int ts);
363 int WriteGlobalData(int timestep, vtkDataArray* buffer);
364 void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
365 int WriteCellData(int timestep, vtkDataArray* buffer);
366 void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
367 int WritePointData(int timestep, vtkDataArray* buffer);
368
373 virtual unsigned int GetMaxNameLength();
374
375private:
376 vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
377 void operator=(const vtkExodusIIWriter&) = delete;
378};
379
380VTK_ABI_NAMESPACE_END
381#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:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
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
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