VTK  9.1.0
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
77#ifndef vtkExodusIIWriter_h
78#define vtkExodusIIWriter_h
79
80#include "vtkIOExodusModule.h" // For export macro
81#include "vtkSmartPointer.h" // For vtkSmartPointer
82#include "vtkWriter.h"
83
84#include <map> // STL Header
85#include <string> // STL Header
86#include <vector> // STL Header
87
89class vtkDoubleArray;
90class vtkIntArray;
92
93class VTKIOEXODUS_EXPORT vtkExodusIIWriter : public vtkWriter
94{
95public:
98 void PrintSelf(ostream& os, vtkIndent indent) override;
99
111 vtkGetObjectMacro(ModelMetadata, vtkModelMetadata);
112
122
130 vtkSetMacro(StoreDoubles, int);
131 vtkGetMacro(StoreDoubles, int);
132
138 vtkSetMacro(GhostLevel, int);
139 vtkGetMacro(GhostLevel, int);
140
147 vtkSetMacro(WriteOutBlockIdArray, vtkTypeBool);
148 vtkGetMacro(WriteOutBlockIdArray, vtkTypeBool);
149 vtkBooleanMacro(WriteOutBlockIdArray, vtkTypeBool);
150
157 vtkSetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
158 vtkGetMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
159 vtkBooleanMacro(WriteOutGlobalNodeIdArray, vtkTypeBool);
160
167 vtkSetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
168 vtkGetMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
169 vtkBooleanMacro(WriteOutGlobalElementIdArray, vtkTypeBool);
170
176 vtkSetMacro(WriteAllTimeSteps, vtkTypeBool);
177 vtkGetMacro(WriteAllTimeSteps, vtkTypeBool);
178 vtkBooleanMacro(WriteAllTimeSteps, vtkTypeBool);
179
180 vtkSetStringMacro(BlockIdArrayName);
181 vtkGetStringMacro(BlockIdArrayName);
182
188 vtkSetMacro(IgnoreMetaDataWarning, bool);
189 vtkGetMacro(IgnoreMetaDataWarning, bool);
190 vtkBooleanMacro(IgnoreMetaDataWarning, bool);
191
192protected:
195
197
199
200 char* FileName;
201 int fid;
202
205
207
215
220
222 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> FlattenedInput;
223 std::vector<vtkSmartPointer<vtkUnstructuredGrid>> NewFlattenedInput;
224
225 std::vector<vtkStdString> FlattenedNames;
226 std::vector<vtkStdString> NewFlattenedNames;
227
228 std::vector<vtkIntArray*> BlockIdList;
229
230 struct Block
231 {
233 {
234 this->Name = nullptr;
235 this->Type = 0;
236 this->NumElements = 0;
237 this->ElementStartIndex = -1;
238 this->NodesPerElement = 0;
239 this->EntityCounts = std::vector<int>();
240 this->EntityNodeOffsets = std::vector<int>();
241 this->GridIndex = 0;
242 this->OutputIndex = -1;
243 this->NumAttributes = 0;
244 this->BlockAttributes = nullptr;
245 };
246 const char* Name;
247 int Type;
251 std::vector<int> EntityCounts;
252 std::vector<int> EntityNodeOffsets;
253 size_t GridIndex;
254 // std::vector<int> CellIndex;
257 float* BlockAttributes; // Owned by metamodel or null. Don't delete.
258 };
259 std::map<int, Block> BlockInfoMap;
260 int NumCells, NumPoints, MaxId;
261
262 std::vector<vtkIdType*> GlobalElementIdList;
263 std::vector<vtkIdType*> GlobalNodeIdList;
264
267
269 {
273 std::vector<std::string> OutNames;
274 };
275 std::map<std::string, VariableInfo> GlobalVariableMap;
276 std::map<std::string, VariableInfo> BlockVariableMap;
277 std::map<std::string, VariableInfo> NodeVariableMap;
281
282 std::vector<std::vector<int>> CellToElementOffset;
283
284 // By BlockId, and within block ID by element variable, with variables
285 // appearing in the same order in which they appear in OutputElementArrayNames
286
289
290 int BlockVariableTruthValue(int blockIdx, int varIdx);
291
292 char* StrDupWithNew(const char* s);
294
296 vtkInformationVector* outputVector) override;
297
299 vtkInformationVector* outputVector);
300
301 virtual int RequestUpdateExtent(vtkInformation* request, vtkInformationVector** inputVector,
302 vtkInformationVector* outputVector);
303
305
307 vtkInformationVector* outputVector) override;
308
309 void WriteData() override;
310
311 int FlattenHierarchy(vtkDataObject* input, const char* name, bool& changed);
312
315
316 int IsDouble();
318 int CheckParametersInternal(int numberOfProcesses, int myRank);
319 virtual int CheckParameters();
320 // If writing in parallel multiple time steps exchange after each time step
321 // if we should continue the execution. Pass local continueExecution as a
322 // parameter and return the global continueExecution.
323 virtual int GlobalContinueExecuting(int localContinueExecution);
325 virtual void CheckBlockInfoMap();
330 char* GetCellTypeName(int t);
331
335
336 void ConvertVariableNames(std::map<std::string, VariableInfo>& variableMap);
338 int nScalarArrays, const std::map<std::string, VariableInfo>& variableMap);
339 std::string CreateNameForScalarArray(const char* root, int component, int numComponents);
340
341 std::map<vtkIdType, vtkIdType>* LocalNodeIdMap;
342 std::map<vtkIdType, vtkIdType>* LocalElementIdMap;
343
347
360 vtkIntArray* GetBlockIdArray(const char* BlockIdArrayName, vtkUnstructuredGrid* input);
361 static bool SameTypeOfCells(vtkIntArray* cellToBlockId, vtkUnstructuredGrid* input);
362
363 double ExtractGlobalData(const char* name, int comp, int ts);
364 int WriteGlobalData(int timestep, vtkDataArray* buffer);
365 void ExtractCellData(const char* name, int comp, vtkDataArray* buffer);
366 int WriteCellData(int timestep, vtkDataArray* buffer);
367 void ExtractPointData(const char* name, int comp, vtkDataArray* buffer);
368 int WritePointData(int timestep, vtkDataArray* buffer);
369
374 virtual unsigned int GetMaxNameLength();
375
376private:
377 vtkExodusIIWriter(const vtkExodusIIWriter&) = delete;
378 void operator=(const vtkExodusIIWriter&) = delete;
379};
380
381#endif
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
general representation of visualization data
Definition: vtkDataObject.h:69
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
int * BlockElementVariableTruthTable
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:43
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:49
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:43
@ component
Definition: vtkX3D.h:181
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
@ name
Definition: vtkX3D.h:225
@ string
Definition: vtkX3D.h:496
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