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
70
71#ifndef vtkExodusIIWriter_h
72#define vtkExodusIIWriter_h
73
74#include "vtkIOExodusModule.h" // For export macro
75#include "vtkSmartPointer.h" // For vtkSmartPointer
76#include "vtkWriter.h"
77
78#include <map> // STL Header
79#include <string> // STL Header
80#include <vector> // STL Header
81
83class vtkDoubleArray;
84class vtkIntArray;
86
87class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
88{
89public:
92 void PrintSelf(ostream& os, vtkIndent indent) override;
93
103
105 vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
106
113
116
123
124 vtkSetMacro(StoreDoubles, int);
125 vtkGetMacro(StoreDoubles, int);
126
131
132 vtkSetMacro(GhostLevel, int);
133 vtkGetMacro(GhostLevel, int);
134
140
144
150
154
160
164
169
173
174 vtkSetStringMacro(BlockIdArrayName);
175 vtkGetStringMacro(BlockIdArrayName);
176
181
182 vtkSetMacro(IgnoreMetaDataWarning, bool);
183 vtkGetMacro(IgnoreMetaDataWarning, bool);
184 vtkBooleanMacro(IgnoreMetaDataWarning, bool);
185
186protected:
189
191
193
194 char* FileName;
195 int fid;
196
199
201
209
214
216 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
217 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
218
219 std::vector<vtkStdString> FlattenedNames;
220 std::vector<vtkStdString> NewFlattenedNames;
221
222 std::vector<vtkIntArray*> BlockIdList;
223
224 struct Block
225 {
227 {
228 this->Name = nullptr;
229 this->Type = 0;
230 this->NumElements = 0;
231 this->ElementStartIndex = -1;
232 this->NodesPerElement = 0;
233 this->EntityCounts = std::vector<int>();
234 this->EntityNodeOffsets = std::vector<int>();
235 this->GridIndex = 0;
236 this->OutputIndex = -1;
237 this->NumAttributes = 0;
238 this->BlockAttributes = nullptr;
239 };
240 const char* Name;
241 int Type;
245 std::vector<int> EntityCounts;
246 std::vector<int> EntityNodeOffsets;
247 size_t GridIndex;
248 // std::vector<int> CellIndex;
251 float* BlockAttributes; // Owned by metamodel or null. Don't delete.
252 };
253 std::map<int, Block> BlockInfoMap;
255
256 std::vector<vtkIdType*> GlobalElementIdList;
257 std::vector<vtkIdType*> GlobalNodeIdList;
258
261
263 {
267 std::vector<std::string> OutNames;
268 };
269 std::map<std::string, VariableInfo> GlobalVariableMap;
270 std::map<std::string, VariableInfo> BlockVariableMap;
271 std::map<std::string, VariableInfo> NodeVariableMap;
275
276 std::vector<std::vector<int>> CellToElementOffset;
277
278 // By BlockId, and within block ID by element variable, with variables
279 // appearing in the same order in which they appear in OutputElementArrayNames
280
283
284 int BlockVariableTruthValue(int blockIdx, int varIdx);
285
286 char* StrDupWithNew(const char* s);
287 void StringUppercase(std::string& str);
288
290 vtkInformationVector* outputVector) override;
291
293 vtkInformationVector* outputVector);
294
295 virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
296 vtkInformationVector* outputVector);
297
298 int FillInputPortInformation(int port, vtkInformation* info) override;
299
301 vtkInformationVector* outputVector) override;
302
303 void WriteData() override;
304
305 int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
306
309
310 int IsDouble();
312 int CheckParametersInternal(int numberOfProcesses, int myRank);
313 virtual int CheckParameters();
314 // If writing in parallel multiple time steps exchange after each time step
315 // if we should continue the execution. Pass local continueExecution as a
316 // parameter and return the global continueExecution.
317 virtual int GlobalContinueExecuting(int localContinueExecution);
319 virtual void CheckBlockInfoMap();
324 char* GetCellTypeName(int t);
325
329
330 void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
332 int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
333 std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
334
335 std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
336 std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
337
341
355 static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
356
357 double ExtractGlobalData(const char* name, int comp, int ts);
358 int WriteGlobalData(int timestep, vtkDataArray* buffer);
359 void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
360 int WriteCellData(int timestep, vtkDataArray* buffer);
361 void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
362 int WritePointData(int timestep, vtkDataArray* buffer);
363
368 virtual unsigned int GetMaxNameLength();
369
370private:
371 vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
372 void operator=(const vtkExodusIIWriter&) = delete;
373};
374
375#endif
general representation of visualization data
dynamic, self-adjusting array of double
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:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition vtkIntArray.h:46
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
std::vector< int > EntityNodeOffsets
std::vector< int > EntityCounts
std::vector< std::string > OutNames
int vtkTypeBool
Definition vtkABI.h:69
#define vtkDataArray
int vtkIdType
Definition vtkType.h:332