VTK
/Users/kitware/Dashboards/MyTests/VTK_BLD_Release_docs/Utilities/Doxygen/dox/IO/ADIOS/vtkADIOSReader.h
Go to the documentation of this file.
00001 /*=========================================================================
00002 
00003   Program:   Visualization Toolkit
00004   Module:    vtkADIOSReader.h
00005 
00006   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
00007   All rights reserved.
00008   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
00009 
00010      This software is distributed WITHOUT ANY WARRANTY; without even
00011      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
00012      PURPOSE.  See the above copyright notice for more information.
00013 
00014 =========================================================================*/
00024 #ifndef vtkADIOSReader_h
00025 #define vtkADIOSReader_h
00026 
00027 #include <map>    // For independently time stepped array indexing
00028 #include <queue>  // For post read operations
00029 #include <string> // For variable name index mapping
00030 #include <vector> // For independently time stepped array indexing
00031 
00032 #include "vtkDataObjectAlgorithm.h"
00033 #include "vtkMultiProcessController.h" // For the MPI controller member
00034 #include "vtkSetGet.h"                 // For property get/set macros
00035 #include "vtkSmartPointer.h"           // For the object cache
00036 
00037 #include "ADIOSDefs.h"                 // For enum definitions
00038 
00039 #include "vtkIOADIOSModule.h"          // For export macro
00040 
00041 namespace ADIOS
00042 {
00043 class VarInfo;
00044 class Reader;
00045 }
00046 class vtkADIOSDirTree;
00047 class BaseFunctor;
00048 
00049 class vtkCellArray;
00050 class vtkDataArray;
00051 class vtkDataObject;
00052 class vtkDataSet;
00053 class vtkDataSetAttributes;
00054 class vtkFieldData;
00055 class vtkImageData;
00056 class vtkPolyData;
00057 class vtkUnstructuredGrid;
00058 
00059 //----------------------------------------------------------------------------
00060 
00061 class VTKIOADIOS_EXPORT vtkADIOSReader : public vtkDataObjectAlgorithm
00062 {
00063 public:
00064   static vtkADIOSReader* New(void);
00065   vtkTypeMacro(vtkADIOSReader,vtkDataObjectAlgorithm);
00066   virtual void PrintSelf(ostream& os, vtkIndent indent);
00067 
00070   int CanReadFile(const char* name);
00071 
00073 
00074   vtkSetStringMacro(FileName);
00075   vtkGetStringMacro(FileName);
00077 
00079 
00080   vtkGetMacro(ReadMethod, int);
00081   vtkSetClampMacro(ReadMethod, int,
00082                    static_cast<int>(ADIOS::ReadMethod_BP),
00083                    static_cast<int>(ADIOS::ReadMethod_FlexPath));
00084   void SetReadMethodBP()          { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_BP)); }
00085   void SetReadMethodBPAggregate() { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_BP_AGGREGATE)); }
00086   void SetReadMethodDataSpaces()  { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_DataSpaces)); }
00087   void SetReadMethodDIMES()       { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_DIMES)); }
00088   void SetReadMethodFlexPath()    { this->SetReadMethod(static_cast<int>(ADIOS::ReadMethod_FlexPath)); }
00090 
00091 
00093 
00094   vtkSetStringMacro(ReadMethodArguments);
00095   vtkGetStringMacro(ReadMethodArguments);
00097 
00099 
00100   void SetController(vtkMultiProcessController*);
00101   vtkGetObjectMacro(Controller, vtkMultiProcessController);
00103 
00105 
00106   virtual int ProcessRequest(vtkInformation*, vtkInformationVector**,
00107     vtkInformationVector*);
00109 
00110 protected:
00111 
00113   bool OpenAndReadMetadata(void);
00114 
00116   void WaitForReads(void);
00117 
00119 
00121   template<typename T>
00122   T* ReadObject(const std::string& path, int blockId);
00124 
00126 
00130   void ReadObject(const ADIOS::VarInfo* info, const vtkADIOSDirTree *subDir,
00131     vtkDataArray* data, int blockId);
00132   void ReadObject(const vtkADIOSDirTree *dir, vtkCellArray* data, int blockId);
00133   void ReadObject(const vtkADIOSDirTree *dir, vtkFieldData* data, int blockId);
00134   void ReadObject(const vtkADIOSDirTree *dir, vtkDataSetAttributes* data, int blockId);
00135   void ReadObject(const vtkADIOSDirTree *dir, vtkDataSet* data, int blockId);
00136   void ReadObject(const vtkADIOSDirTree *dir, vtkImageData* data, int blockId);
00137   void ReadObject(const vtkADIOSDirTree *dir, vtkPolyData* data, int blockId);
00138   void ReadObject(const vtkADIOSDirTree *dir, vtkUnstructuredGrid* data, int blockId);
00140 
00141   char *FileName;
00142   int ReadMethod;
00143   char *ReadMethodArguments;
00144   vtkADIOSDirTree *Tree;
00145   ADIOS::Reader *Reader;
00146   vtkMultiProcessController *Controller;
00147 
00148   // Index information for independently stepped variables
00149 
00150   // Map variable names to thier position in the block step index
00151   // [BlockId][VarName] = IndexId
00152   std::vector<std::map<std::string, size_t> > BlockStepIndexIdMap;
00153 
00154   // [BlockId][GlobalStep][IndexId] = LocalStep
00155   // Ex: The file has 30 steps, but the Variable "/Foo/Bar" in block 3 only
00156   // has 2 steps, written out at global step 10 and global step 17.  To
00157   // lookup the local step for the variable at global time step 25:
00158   //
00159   // size_t idx = this->BlockStepIndexIdMap[3]["/Foo/Bar"];
00160   // int localStep = this->BlockStepIndex[3][25][idx];
00161   //
00162   // At this point, localStep = 2, since at global step 25, local step 2 is the
00163   // most recent version of "/Foo/Bar" available
00164   std::vector<std::vector<std::vector<int> > > BlockStepIndex;
00165 
00166   // Cache the VTK objects as they are read
00167   // Key = <BlockId, IndexId>, Value = <LocalStep, Object>
00168   std::map<std::pair<int, size_t>,
00169            std::pair<int, vtkSmartPointer<vtkObject> > >
00170     ObjectCache;
00171 
00172   vtkADIOSReader();
00173   virtual ~vtkADIOSReader();
00174 
00175   /* The design of ADIOS is such that array IO is not directly performed
00176    * upon request, but instead is scheduled to be performed later, at which
00177    * time all IO operations are processed at once in bulk.  This creates
00178    * an odd situation for data management since arrays will be allocated with
00179    * junk data and scheduled to be filled, but they cannot be safely assigned
00180    * to a VTK object until the data contained in them is valid, e.g. through
00181    * a call to vtkUnstructuredGrid::SetPoints or similar.  Similary,
00182    * they cannot have thier reference cound safely decremented until after
00183    * they have been assigned to a vtk object.  To work around this, a generic
00184    * action queue is created to hold a list of arbitrary functions that need
00185    * to be called in a particular order after the reads have been
00186    * processed.  The AddPostReadOperation prototypes use a large number of
00187    * template parameters in order to permit the compiler to automatically
00188    * perform the correct type deduction necessary to translate between
00189    * member function signatures and the objects and arguments they get called
00190    * with.  This allows for arbitrary functions with arbitrary return types
00191    * and arbitrary argument types to be collected into a single event queue
00192    */
00193 
00194   // A set of operations to perform after reading is complete
00195   std::queue<BaseFunctor*> PostReadOperations;
00196 
00197   // A set of shortcuts to allow automatic parameter deduction
00198 
00199   template<typename TObjectFun, typename TObjectData, typename TReturn>
00200   void AddPostReadOperation(TObjectData*, TReturn (TObjectFun::*)());
00201 
00202   template<typename TObjectFun, typename TObjectData, typename TReturn,
00203     typename TArg1Fun, typename TArg1Data>
00204   void AddPostReadOperation(TObjectData*,
00205     TReturn (TObjectFun::*)(TArg1Fun), TArg1Data);
00206 
00207   template<typename TObjectFun, typename TObjectData, typename TReturn,
00208     typename TArg1Fun, typename TArg1Data,
00209     typename TArg2Fun, typename TArg2Data>
00210   void AddPostReadOperation(TObjectData*,
00211     TReturn (TObjectFun::*)(TArg1Fun, TArg2Fun),
00212     TArg1Data, TArg2Data);
00213 
00214   template<typename TObjectFun, typename TObjectData, typename TReturn,
00215     typename TArg1Fun, typename TArg1Data,
00216     typename TArg2Fun, typename TArg2Data,
00217     typename TArg3Fun, typename TArg3Data>
00218   void AddPostReadOperation(TObjectData*,
00219     TReturn (TObjectFun::*)(TArg1Fun, TArg2Fun, TArg3Fun),
00220     TArg1Data, TArg2Data, TArg3Data);
00221 
00222   // Used to implement vtkAlgorithm
00223 
00224   int FillOutputPortInformation(int, vtkInformation*);
00225 
00226   virtual int RequestInformation(vtkInformation *request,
00227                                  vtkInformationVector **input,
00228                                  vtkInformationVector *output);
00229   virtual int RequestUpdateExtent(vtkInformation *request,
00230                                   vtkInformationVector **input,
00231                                   vtkInformationVector *output);
00232   virtual int RequestData(vtkInformation *request,
00233                           vtkInformationVector **input,
00234                           vtkInformationVector *output);
00235 
00236   int NumberOfPieces;
00237   std::vector<double> TimeSteps;
00238   std::map<double, size_t> TimeStepsIndex;
00239 
00240   double RequestStep;
00241   int RequestStepIndex;
00242   int RequestNumberOfPieces;
00243   int RequestPiece;
00244 
00245 private:
00246   vtkADIOSReader(const vtkADIOSReader&);  // Not implemented.
00247   void operator=(const vtkADIOSReader&);  // Not implemented.
00248 };
00249 
00250 #define DECLARE_EXPLICIT(T) \
00251 template<> T* vtkADIOSReader::ReadObject<T>(const std::string& path, \
00252   int blockId);
00253 DECLARE_EXPLICIT(vtkImageData)
00254 DECLARE_EXPLICIT(vtkPolyData)
00255 DECLARE_EXPLICIT(vtkUnstructuredGrid)
00256 #undef DECLARE_EXPLICIT
00257 
00258 #endif