VTK
vtkADIOSReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkADIOSReader.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 vtkADIOSReader_h
26 #define vtkADIOSReader_h
27 
28 #include <map> // For independently time stepped array indexing
29 #include <queue> // For post read operations
30 #include <string> // For variable name index mapping
31 #include <vector> // For independently time stepped array indexing
32 
33 #include "vtkDataObjectAlgorithm.h"
34 #include "vtkMultiProcessController.h" // For the MPI controller member
35 #include "vtkSetGet.h" // For property get/set macros
36 #include "vtkSmartPointer.h" // For the object cache
37 
38 #include "ADIOSDefs.h" // For enum definitions
39 
40 #include "vtkIOADIOSModule.h" // For export macro
41 
42 namespace ADIOS
43 {
44 class VarInfo;
45 class Reader;
46 }
47 class vtkADIOSDirTree;
48 class BaseFunctor;
49 
50 class vtkCellArray;
51 class vtkDataArray;
52 class vtkDataObject;
53 class vtkDataSet;
55 class vtkFieldData;
56 class vtkImageData;
57 class vtkPolyData;
59 
60 //----------------------------------------------------------------------------
61 
62 class VTKIOADIOS_EXPORT vtkADIOSReader : public vtkDataObjectAlgorithm
63 {
64 public:
65  static vtkADIOSReader* New(void);
67  void PrintSelf(ostream& os, vtkIndent indent) override;
68 
73  int CanReadFile(const char* name);
74 
76 
79  vtkSetStringMacro(FileName);
80  vtkGetStringMacro(FileName);
82 
84 
87  vtkGetMacro(ReadMethod, int);
88  vtkSetClampMacro(ReadMethod, int,
89  static_cast<int>(ADIOS::ReadMethod_BP),
90  static_cast<int>(ADIOS::ReadMethod_FlexPath));
91  void SetReadMethodBP() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_BP)); }
92  void SetReadMethodBPAggregate() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_BP_AGGREGATE)); }
93  void SetReadMethodDataSpaces() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_DataSpaces)); }
94  void SetReadMethodDIMES() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_DIMES)); }
95  void SetReadMethodFlexPath() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_FlexPath)); }
97 
98 
100 
103  vtkSetStringMacro(ReadMethodArguments);
104  vtkGetStringMacro(ReadMethodArguments);
106 
108 
111  void SetController(vtkMultiProcessController*);
112  vtkGetObjectMacro(Controller, vtkMultiProcessController);
114 
118  virtual int ProcessRequest(vtkInformation*, vtkInformationVector**,
120 
121 protected:
122 
126  bool OpenAndReadMetadata(void);
127 
131  void WaitForReads(void);
132 
137  template<typename T>
138  T* ReadObject(const std::string& path, int blockId);
139 
141 
147  void ReadObject(const ADIOS::VarInfo* info, const vtkADIOSDirTree *subDir,
148  vtkDataArray* data, int blockId);
149  void ReadObject(const vtkADIOSDirTree *dir, vtkCellArray* data, int blockId);
150  void ReadObject(const vtkADIOSDirTree *dir, vtkFieldData* data, int blockId);
151  void ReadObject(const vtkADIOSDirTree *dir, vtkDataSetAttributes* data, int blockId);
152  void ReadObject(const vtkADIOSDirTree *dir, vtkDataSet* data, int blockId);
153  void ReadObject(const vtkADIOSDirTree *dir, vtkImageData* data, int blockId);
154  void ReadObject(const vtkADIOSDirTree *dir, vtkPolyData* data, int blockId);
155  void ReadObject(const vtkADIOSDirTree *dir, vtkUnstructuredGrid* data, int blockId);
157 
158  char *FileName;
164 
165  // Index information for independently stepped variables
166 
167  // Map variable names to their position in the block step index
168  // [BlockId][VarName] = IndexId
169  std::vector<std::map<std::string, size_t> > BlockStepIndexIdMap;
170 
171  // [BlockId][GlobalStep][IndexId] = LocalStep
172  // Ex: The file has 30 steps, but the Variable "/Foo/Bar" in block 3 only
173  // has 2 steps, written out at global step 10 and global step 17. To
174  // lookup the local step for the variable at global time step 25:
175  //
176  // size_t idx = this->BlockStepIndexIdMap[3]["/Foo/Bar"];
177  // int localStep = this->BlockStepIndex[3][25][idx];
178  //
179  // At this point, localStep = 2, since at global step 25, local step 2 is the
180  // most recent version of "/Foo/Bar" available
181  std::vector<std::vector<std::vector<int> > > BlockStepIndex;
182 
183  // Cache the VTK objects as they are read
184  // Key = <BlockId, IndexId>, Value = <LocalStep, Object>
185  std::map<std::pair<int, size_t>,
186  std::pair<int, vtkSmartPointer<vtkObject> > >
188 
189  vtkADIOSReader();
190  virtual ~vtkADIOSReader();
191 
192  /* The design of ADIOS is such that array IO is not directly performed
193  * upon request, but instead is scheduled to be performed later, at which
194  * time all IO operations are processed at once in bulk. This creates
195  * an odd situation for data management since arrays will be allocated with
196  * junk data and scheduled to be filled, but they cannot be safely assigned
197  * to a VTK object until the data contained in them is valid, e.g. through
198  * a call to vtkUnstructuredGrid::SetPoints or similar. Similarly,
199  * they cannot have their reference could safely decremented until after
200  * they have been assigned to a vtk object. To work around this, a generic
201  * action queue is created to hold a list of arbitrary functions that need
202  * to be called in a particular order after the reads have been
203  * processed. The AddPostReadOperation prototypes use a large number of
204  * template parameters in order to permit the compiler to automatically
205  * perform the correct type deduction necessary to translate between
206  * member function signatures and the objects and arguments they get called
207  * with. This allows for arbitrary functions with arbitrary return types
208  * and arbitrary argument types to be collected into a single event queue
209  */
210 
211  // A set of operations to perform after reading is complete
212  std::queue<BaseFunctor*> PostReadOperations;
213 
214  // A set of shortcuts to allow automatic parameter deduction
215 
216  template<typename TObjectFun, typename TObjectData, typename TReturn>
217  void AddPostReadOperation(TObjectData*, TReturn (TObjectFun::*)());
218 
219  template<typename TObjectFun, typename TObjectData, typename TReturn,
220  typename TArg1Fun, typename TArg1Data>
221  void AddPostReadOperation(TObjectData*,
222  TReturn (TObjectFun::*)(TArg1Fun), TArg1Data);
223 
224  template<typename TObjectFun, typename TObjectData, typename TReturn,
225  typename TArg1Fun, typename TArg1Data,
226  typename TArg2Fun, typename TArg2Data>
227  void AddPostReadOperation(TObjectData*,
228  TReturn (TObjectFun::*)(TArg1Fun, TArg2Fun),
229  TArg1Data, TArg2Data);
230 
231  template<typename TObjectFun, typename TObjectData, typename TReturn,
232  typename TArg1Fun, typename TArg1Data,
233  typename TArg2Fun, typename TArg2Data,
234  typename TArg3Fun, typename TArg3Data>
235  void AddPostReadOperation(TObjectData*,
236  TReturn (TObjectFun::*)(TArg1Fun, TArg2Fun, TArg3Fun),
237  TArg1Data, TArg2Data, TArg3Data);
238 
239  // Used to implement vtkAlgorithm
240 
241  int FillOutputPortInformation(int, vtkInformation*);
242 
243  virtual int RequestInformation(vtkInformation *request,
244  vtkInformationVector **input,
245  vtkInformationVector *output);
246  virtual int RequestUpdateExtent(vtkInformation *request,
247  vtkInformationVector **input,
248  vtkInformationVector *output);
249  virtual int RequestData(vtkInformation *request,
250  vtkInformationVector **input,
251  vtkInformationVector *output);
252 
254  std::vector<double> TimeSteps;
255  std::map<double, size_t> TimeStepsIndex;
256 
257  double RequestStep;
261 
262 private:
263  vtkADIOSReader(const vtkADIOSReader&) = delete;
264  void operator=(const vtkADIOSReader&) = delete;
265 };
266 
267 #define DECLARE_EXPLICIT(T) \
268 template<> T* vtkADIOSReader::ReadObject<T>(const std::string& path, \
269  int blockId);
273 #undef DECLARE_EXPLICIT
274 
275 #endif
std::queue< BaseFunctor * > PostReadOperations
Store vtkAlgorithm input/output information.
#define DECLARE_EXPLICIT(T)
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
std::map< std::pair< int, size_t >, std::pair< int, vtkSmartPointer< vtkObject > > > ObjectCache
void SetReadMethodBPAggregate()
Get/Set the ADIOS read method.
Read ADIOS files.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
char * ReadMethodArguments
ReadMethod
Definition: ADIOSDefs.h:49
vtkADIOSDirTree * Tree
std::vector< std::vector< std::vector< int > > > BlockStepIndex
a simple class to control print indentation
Definition: vtkIndent.h:39
topologically and geometrically regular array of data
Definition: vtkImageData.h:45
dataset represents arbitrary combinations of all possible cell types
void SetReadMethodDataSpaces()
Get/Set the ADIOS read method.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
represent and manipulate attribute data in a dataset
void SetReadMethodFlexPath()
Get/Set the ADIOS read method.
std::map< double, size_t > TimeStepsIndex
void SetReadMethodDIMES()
Get/Set the ADIOS read method.
vtkMultiProcessController * Controller
A directory tree structure holding ADIOS data.
Superclass for algorithms that produce only data object as output.
object to represent cell connectivity
Definition: vtkCellArray.h:50
void SetReadMethodBP()
Get/Set the ADIOS read method.
Store zero or more vtkInformation instances.
ADIOS::Reader * Reader
std::vector< double > TimeSteps
general representation of visualization data
Definition: vtkDataObject.h:64
represent and manipulate fields of data
Definition: vtkFieldData.h:56
Multiprocessing communication superclass.
std::vector< std::map< std::string, size_t > > BlockStepIndexIdMap