VTK
vtkADIOSWriter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkADIOSWriter.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=========================================================================*/
25#ifndef vtkADIOSWriter_h
26#define vtkADIOSWriter_h
27
28#include <map> // For independently stepped array indexing
29#include <string> // For independently stepped array indexing
30#include <vector> // For independently stepped array indexing
31
33#include "vtkMultiProcessController.h" // For the MPI controller member
34#include "vtkSetGet.h" // For property get/set macros
35
36#include "ADIOSDefs.h" // For enum definitions
37
38#include "vtkIOADIOSModule.h" // For export macro
39
40namespace ADIOS
41{
42 class Writer;
43}
44
46class vtkCellArray;
47class vtkDataArray;
48class vtkDataObject;
49class vtkDataSet;
50class vtkFieldData;
51class vtkImageData;
52class vtkPolyData;
54
55class VTKIOADIOS_EXPORT vtkADIOSWriter : public vtkDataObjectAlgorithm
56{
57public:
60 void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
61
63
65
69 vtkSetStringMacro(FileName)
71
73
76 vtkGetMacro(TransportMethod, int);
77 vtkSetClampMacro(TransportMethod, int,
78 static_cast<int>(ADIOS::TransportMethod_NULL),
79 static_cast<int>(ADIOS::TransportMethod_NetCDF4));
80 void SetTransportMethodToNULL() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_NULL)); }
81 void SetTransportMethodToPOSIX() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_POSIX)); }
82 void SetTransportMethodToMPI() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_MPI)); }
83 void SetTransportMethodToMPILustre() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_MPI_LUSTRE)); }
84 void SetTransportMethodToMPIAggregate() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_MPI_AGGREGATE)); }
85 void SetTransportMethodToVarMerge() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_VAR_MERGE)); }
86 void SetTransportMethodToDataSpaces() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_DataSpaces)); }
87 void SetTransportMethodToDIMES() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_DIMES)); }
88 void SetTransportMethodToFlexPath() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_FlexPath)); }
89 void SetTransportMethodToPHDF5() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_PHDF5)); }
90 void SetTransportMethodToNetCDF4() { this->SetTransportMethod(static_cast<int>(ADIOS::TransportMethod_NetCDF4)); }
92
94
98 vtkSetStringMacro(TransportMethodArguments)
99 vtkGetStringMacro(TransportMethodArguments)
101
103
106 vtkGetMacro(Transform, int);
107 vtkSetClampMacro(Transform, int,
108 static_cast<int>(ADIOS::Transform_NONE),
109 static_cast<int>(ADIOS::Transform_SZIP));
110 void SetTransformToNone() { this->SetTransform(static_cast<int>(ADIOS::Transform_NONE)); }
111 void SetTransformToZLib() { this->SetTransform(static_cast<int>(ADIOS::Transform_ZLIB)); }
112 void SetTransformToBZip2() { this->SetTransform(static_cast<int>(ADIOS::Transform_BZLIB2)); }
113 void SetTransformToSZip() { this->SetTransform(static_cast<int>(ADIOS::Transform_SZIP)); }
115
117
122 vtkSetMacro(WriteAllTimeSteps, bool);
123 vtkGetMacro(WriteAllTimeSteps, bool);
124 vtkBooleanMacro(WriteAllTimeSteps, bool);
126
128
132 vtkGetObjectMacro(Controller, vtkMultiProcessController);
134
140
144 void Write() { return this->Update(); }
145
146protected:
147
149
152 void Define(const std::string& path, const vtkAbstractArray* value);
153 void Define(const std::string& path, const vtkDataArray* value);
154 void Define(const std::string& path, const vtkCellArray* value);
155 void Define(const std::string& path, const vtkFieldData* value);
156 void Define(const std::string& path, const vtkDataSet* value);
157 void Define(const std::string& path, const vtkImageData* value);
158 void Define(const std::string& path, const vtkPolyData* value);
159 void Define(const std::string& path, const vtkUnstructuredGrid* value);
161
163
169 void OpenFile();
170 void CloseFile();
172
174
177 void Write(const std::string& path, const vtkAbstractArray* value);
178 void Write(const std::string& path, const vtkDataArray* value);
179 void Write(const std::string& path, const vtkCellArray* value);
180 void Write(const std::string& path, const vtkFieldData* value);
181 void Write(const std::string& path, const vtkDataSet* value);
182 void Write(const std::string& path, const vtkImageData* value);
183 void Write(const std::string& path, const vtkPolyData* value);
184 void Write(const std::string& path, const vtkUnstructuredGrid* value);
186
187 char *FileName;
191 int Rank;
195
198
199protected:
200 // Used to implement vtkAlgorithm
201
203
204 virtual int RequestInformation(vtkInformation *request,
205 vtkInformationVector **input,
206 vtkInformationVector *output);
208 vtkInformationVector **input,
209 vtkInformationVector *output);
210 virtual int RequestData(vtkInformation *request,
211 vtkInformationVector **input,
212 vtkInformationVector *output);
213
218 std::vector<double> TimeSteps;
220 int RequestExtent[6];
221
222 // Used to determine whether or not the data getting written is stale
223 bool UpdateMTimeTable(const std::string& path, const vtkObject* value);
224 std::map<std::string, unsigned long> LastUpdated;
225private:
226 bool WriteInternal();
227
228 template<typename T>
229 bool DefineAndWrite(vtkDataObject *input);
230
231 vtkADIOSWriter(const vtkADIOSWriter&) VTK_DELETE_FUNCTION;
232 void operator=(const vtkADIOSWriter&) VTK_DELETE_FUNCTION;
233};
234
235#endif
Write ADIOS files.
void SetController(vtkMultiProcessController *)
Set the MPI controller.
void Write(const std::string &path, const vtkImageData *value)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Define(const std::string &path, const vtkCellArray *value)
virtual int ProcessRequest(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
The main interface which triggers the writer to start.
ADIOS::Writer * Writer
void Write(const std::string &path, const vtkDataArray *value)
void Write(const std::string &path, const vtkFieldData *value)
void SetTransportMethodToPHDF5()
const char * GetDefaultFileExtension()
void Write(const std::string &path, const vtkUnstructuredGrid *value)
vtkMultiProcessController * Controller
void SetTransportMethodToMPI()
char * TransportMethodArguments
int FillInputPortInformation(int port, vtkInformation *info)
Fill the input port information objects for this algorithm.
void Write()
Declare data if necessary and write the current step to the output stream.
void Define(const std::string &path, const vtkImageData *value)
void SetTransportMethodToNetCDF4()
bool UpdateMTimeTable(const std::string &path, const vtkObject *value)
void SetTransformToSZip()
void Write(const std::string &path, const vtkCellArray *value)
void SetTransportMethodToVarMerge()
void Write(const std::string &path, const vtkAbstractArray *value)
Write a previously defined VTK data type.
virtual int RequestUpdateExtent(vtkInformation *request, vtkInformationVector **input, vtkInformationVector *output)
This is called by the superclass.
void SetTransportMethodToFlexPath()
static vtkADIOSWriter * New()
void Define(const std::string &path, const vtkUnstructuredGrid *value)
std::vector< double > TimeSteps
virtual int RequestData(vtkInformation *request, vtkInformationVector **input, vtkInformationVector *output)
void OpenFile()
Open a file and prepare for writing already defined variables.
void SetTransformToZLib()
void Define(const std::string &path, const vtkAbstractArray *value)
Define a VTK data type.
void Define(const std::string &path, const vtkPolyData *value)
std::map< std::string, unsigned long > LastUpdated
void SetTransportMethodToDataSpaces()
void SetTransportMethodToMPIAggregate()
void Write(const std::string &path, const vtkDataSet *value)
void Define(const std::string &path, const vtkFieldData *value)
void Write(const std::string &path, const vtkPolyData *value)
void SetTransportMethodToPOSIX()
void SetTransformToBZip2()
void SetTransportMethodToDIMES()
void Define(const std::string &path, const vtkDataSet *value)
void Define(const std::string &path, const vtkDataArray *value)
void SetTransportMethodToMPILustre()
virtual int RequestInformation(vtkInformation *request, vtkInformationVector **input, vtkInformationVector *output)
Abstract superclass for all arrays.
virtual void Update()
object to represent cell connectivity
Definition: vtkCellArray.h:51
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
Superclass for algorithms that produce only data object as output.
general representation of visualization data
Definition: vtkDataObject.h:65
abstract class to specify dataset behavior
Definition: vtkDataSet.h:63
represent and manipulate fields of data
Definition: vtkFieldData.h:57
topologically and geometrically regular array of data
Definition: vtkImageData.h:46
a simple class to control print indentation
Definition: vtkIndent.h:40
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Multiprocessing communication superclass.
abstract base class for most VTK objects
Definition: vtkObject.h:60
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
dataset represents arbitrary combinations of all possible cell types
Transform
Definition: ADIOSDefs.h:40
@ Transform_SZIP
Definition: ADIOSDefs.h:44
@ Transform_BZLIB2
Definition: ADIOSDefs.h:43
@ Transform_ZLIB
Definition: ADIOSDefs.h:42
@ Transform_NONE
Definition: ADIOSDefs.h:41
TransportMethod
Definition: ADIOSDefs.h:24
@ TransportMethod_NULL
Definition: ADIOSDefs.h:25
@ TransportMethod_NetCDF4
Definition: ADIOSDefs.h:35
@ TransportMethod_PHDF5
Definition: ADIOSDefs.h:34
@ TransportMethod_MPI
Definition: ADIOSDefs.h:27
@ TransportMethod_MPI_AGGREGATE
Definition: ADIOSDefs.h:29
@ TransportMethod_DataSpaces
Definition: ADIOSDefs.h:31
@ TransportMethod_VAR_MERGE
Definition: ADIOSDefs.h:30
@ TransportMethod_DIMES
Definition: ADIOSDefs.h:32
@ TransportMethod_POSIX
Definition: ADIOSDefs.h:26
@ TransportMethod_MPI_LUSTRE
Definition: ADIOSDefs.h:28
@ TransportMethod_FlexPath
Definition: ADIOSDefs.h:33
@ info
Definition: vtkX3D.h:376
@ value
Definition: vtkX3D.h:220
@ port
Definition: vtkX3D.h:447
@ string
Definition: vtkX3D.h:490
vtkSetMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkBooleanMacro(IgnoreDriverBugs, bool)
Updates the extensions string.
vtkGetStringMacro(ExtensionsString)
Returns a string listing all available extensions.