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
68#ifndef vtkExodusIIWriter_h
69#define vtkExodusIIWriter_h
70
71#include "vtkIOExodusModule.h" // For export macro
72#include "vtkSmartPointer.h" // For vtkSmartPointer
73#include "vtkWriter.h"
74
75#include <map> // STL Header
76#include <string> // STL Header
77#include <vector> // STL Header
78
80class vtkDoubleArray;
81class vtkIntArray;
83
84class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
85{
86public:
89 void PrintSelf(ostream& os, vtkIndent indent) override;
90
102 vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
103
113
121 vtkSetMacro(StoreDoubles, int);
122 vtkGetMacro(StoreDoubles, int);
123
129 vtkSetMacro(GhostLevel, int);
130 vtkGetMacro(GhostLevel, int);
131
138 vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
139 vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
140 vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
141
148 vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
149 vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
150 vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
151
158 vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
159 vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
160 vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
161
167 vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
168 vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
169 vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
170
171 vtkSetStringMacro(BlockIdArrayName);
172 vtkGetStringMacro(BlockIdArrayName);
173
179 vtkSetMacro(IgnoreMetaDataWarning, bool);
180 vtkGetMacro(IgnoreMetaDataWarning, bool);
181 vtkBooleanMacro(IgnoreMetaDataWarning, bool);
182
183protected:
186
188
190
191 char* FileName;
192 int fid;
193
196
198
206
211
213 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
214 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
215
216 std::vector<vtkStdString> FlattenedNames;
217 std::vector<vtkStdString> NewFlattenedNames;
218
219 std::vector<vtkIntArray*> BlockIdList;
220
221 struct Block
222 {
224 {
225 this->Name = nullptr;
226 this->Type = 0;
227 this->NumElements = 0;
228 this->ElementStartIndex = -1;
229 this->NodesPerElement = 0;
230 this->EntityCounts = std::vector<int>();
231 this->EntityNodeOffsets = std::vector<int>();
232 this->GridIndex = 0;
233 this->OutputIndex = -1;
234 this->NumAttributes = 0;
235 this->BlockAttributes = nullptr;
236 };
237 const char* Name;
238 int Type;
242 std::vector<int> EntityCounts;
243 std::vector<int> EntityNodeOffsets;
244 size_t GridIndex;
245 // std::vector<int> CellIndex;
248 float* BlockAttributes; // Owned by metamodel or null. Don't delete.
249 };
250 std::map<int, Block> BlockInfoMap;
251 int NumCells, NumPoints, MaxId;
252
253 std::vector<vtkIdType*> GlobalElementIdList;
254 std::vector<vtkIdType*> GlobalNodeIdList;
255
258
260 {
264 std::vector<std::string> OutNames;
265 };
266 std::map<std::string, VariableInfo> GlobalVariableMap;
267 std::map<std::string, VariableInfo> BlockVariableMap;
268 std::map<std::string, VariableInfo> NodeVariableMap;
272
273 std::vector<std::vector<int>> CellToElementOffset;
274
275 // By BlockId, and within block ID by element variable, with variables
276 // appearing in the same order in which they appear in OutputElementArrayNames
277
280
281 int BlockVariableTruthValue(int blockIdx, int varIdx);
282
283 char* StrDupWithNew(const char* s);
284 void StringUppercase(std::string& str);
285
287 vtkInformationVector* outputVector) override;
288
290 vtkInformationVector* outputVector);
291
292 virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
293 vtkInformationVector* outputVector);
294
295 int FillInputPortInformation(int port, vtkInformation* info) override;
296
298 vtkInformationVector* outputVector) override;
299
300 void WriteData() override;
301
302 int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
303
306
307 int IsDouble();
309 int CheckParametersInternal(int numberOfProcesses, int myRank);
310 virtual int CheckParameters();
311 // If writing in parallel multiple time steps exchange after each time step
312 // if we should continue the execution. Pass local continueExecution as a
313 // parameter and return the global continueExecution.
314 virtual int GlobalContinueExecuting(int localContinueExecution);
316 virtual void CheckBlockInfoMap();
321 char* GetCellTypeName(int t);
322
326
327 void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
329 int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
330 std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
331
332 std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
333 std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
334
338
351 vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
352 static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
353
354 double ExtractGlobalData(const char* name, int comp, int ts);
355 int WriteGlobalData(int timestep, vtkDataArray* buffer);
356 void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
357 int WriteCellData(int timestep, vtkDataArray* buffer);
358 void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
359 int WritePointData(int timestep, vtkDataArray* buffer);
360
365 virtual unsigned int GetMaxNameLength();
366
367private:
368 vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
369 void operator=(const vtkExodusIIWriter&) = delete;
370};
371
372#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:34
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition vtkIntArray.h:40
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:46
std::vector< int > EntityNodeOffsets
std::vector< int > EntityCounts
std::vector< std::string > OutNames
int vtkTypeBool
Definition vtkABI.h:69
int vtkIdType
Definition vtkType.h:332