VTK  9.1.0
PIOAdaptor.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: ParaView
4 Module: PIOAdaptor.h
5
6 Copyright (c) Kitware, Inc.
7 All rights reserved.
8 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html 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#ifndef PIOAdaptor_h
16#define PIOAdaptor_h
17
21
22#include "PIOData.h"
23
24#include <vector>
25
27
29{
30public:
33
34 int initializeGlobal(const char* DumpDescFile);
35 int initializeDump(int timeStep);
36
37 // Time step change requires new geometry and data
40
41 int GetNumberOfTimeSteps() { return static_cast<int>(this->CycleIndex.size()); }
42 double GetSimulationTime(int step) { return this->SimulationTime[step]; }
43 double GetCycleIndex(int step) { return this->CycleIndex[step]; }
44 double GetPIOFileIndex(int step) { return this->PIOFileIndex[step]; }
45
46 int GetNumberOfVariables() { return (int)this->variableName.size(); }
47 const char* GetVariableName(int indx) { return this->variableName[indx].c_str(); }
48 int GetNumberOfDefaultVariables() { return (int)this->variableDefault.size(); }
49 const char* GetVariableDefault(int indx) { return this->variableDefault[indx].c_str(); }
50
51 // Read pio dump file AMR as hypertree grid rather than unstructured grid
52 bool GetHyperTreeGrid() { return this->useHTG; }
53 void SetHyperTreeGrid(bool val) { this->useHTG = val; }
54
55 // Read pio dump file tracer information
56 bool GetTracers() { return this->useTracer; }
57 void SetTracers(bool val) { this->useTracer = val; }
58
59 // Read pio dump file variable data as 64 bit float
60 bool GetFloat64() { return this->useFloat64; }
61 void SetFloat64(bool val) { this->useFloat64 = val; }
62
63protected:
64 // Collect the metadata
65 int parsePIOFile(const char* DumpDescFile);
66 int collectMetaData(const char* DumpDescFile);
68 std::string trimString(const std::string& str);
69
70 // Create the unstructured grid for tracers
72
73 // Create the unstructured grid for AMR
75
77 int numberOfCells, // Number of cells all levels
78 int* cell_level, // Level within AMR
79 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
80 double* cell_center[1]); // Cell center
81
83 int numberOfCells, // Number of cells all levels
84 int* cell_level, // Level within AMR
85 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
86 double* cell_center[2]); // Cell center
87
89 int numberOfCells, // Number of cells all levels
90 int* cell_level, // Level within AMR
91 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
92 double* cell_center[3]); // Cell center
93
94 // Create the hypertree grid
96
97 int count_hypertree(int64_t curIndex, int64_t* daughter);
98
100 vtkHyperTreeGridNonOrientedCursor* treeCursor, int64_t curIndex, int64_t* daughter);
101
102 // Add variable data to the unstructured grid
105 int64_t* daughter, // Indicates top level cell or not
106 double* data[], // Data for all cells
107 int numberOfCells,
108 int numberOfComponents); // Number of components in data
109
110 // Add variable data to the hypertree grid
113 double* data[], // Data for all cells
114 int numberOfComponents); // Number of components in data
115
116 // Used in parallel reader and load balancing
118 int Rank;
120
121 // Structure to access the dump file data
123
124 // Fields of interest in dump file
125 std::list<std::string> fieldsToRead;
126
127 // Time series of dumps
128 std::string descFileName; // name.pio
129 std::string dumpBaseName; // base name to use for dumps
130 std::vector<std::string> dumpDirectory; // directories holding dumps
131 std::vector<std::string> dumpFileName; // all dump files
132
133 // Time step information
134 std::vector<double> CycleIndex; // Times as cycle index
135 std::vector<double> SimulationTime; // Times as simulation time
136 std::vector<double> PIOFileIndex; // Index into dump files
137
138 // Type of block structures to create within multiblock dataset
139 bool useHTG;
143
144 // Cell variable data and initially enabled variables
145 std::vector<std::string> variableName;
146 std::vector<std::string> variableDefault;
147
148 // Record the ordering of the cells when building the hypertree grid
149 // Needed so that the data will line up correctly
150 std::vector<int> indexNodeLeaf;
151};
152
153#endif
std::vector< double > SimulationTime
Definition: PIOAdaptor.h:135
void create_amr_UG_3D(vtkMultiBlockDataSet *grid, int numberOfCells, int *cell_level, int64_t *cell_daughter, double *cell_center[3])
void build_hypertree(vtkHyperTreeGridNonOrientedCursor *treeCursor, int64_t curIndex, int64_t *daughter)
int GetNumberOfVariables()
Definition: PIOAdaptor.h:46
void create_amr_UG_1D(vtkMultiBlockDataSet *grid, int numberOfCells, int *cell_level, int64_t *cell_daughter, double *cell_center[1])
void create_amr_UG(vtkMultiBlockDataSet *grid)
bool useTracer
Definition: PIOAdaptor.h:140
std::list< std::string > fieldsToRead
Definition: PIOAdaptor.h:125
int count_hypertree(int64_t curIndex, int64_t *daughter)
double GetSimulationTime(int step)
Definition: PIOAdaptor.h:42
void add_amr_UG_scalar(vtkMultiBlockDataSet *grid, vtkStdString varName, int64_t *daughter, double *data[], int numberOfCells, int numberOfComponents)
vtkMultiProcessController * Controller
Definition: PIOAdaptor.h:117
PIO_DATA * pioData
Definition: PIOAdaptor.h:122
bool GetHyperTreeGrid()
Definition: PIOAdaptor.h:52
void create_geometry(vtkMultiBlockDataSet *grid)
void SetTracers(bool val)
Definition: PIOAdaptor.h:57
std::string descFileName
Definition: PIOAdaptor.h:128
int parsePIOFile(const char *DumpDescFile)
double GetCycleIndex(int step)
Definition: PIOAdaptor.h:43
bool hasTracers
Definition: PIOAdaptor.h:142
int initializeDump(int timeStep)
void load_variable_data_HTG(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
void create_tracer_UG(vtkMultiBlockDataSet *grid)
bool GetTracers()
Definition: PIOAdaptor.h:56
std::vector< int > indexNodeLeaf
Definition: PIOAdaptor.h:150
std::vector< std::string > dumpFileName
Definition: PIOAdaptor.h:131
std::vector< std::string > variableDefault
Definition: PIOAdaptor.h:146
std::string trimString(const std::string &str)
void add_amr_HTG_scalar(vtkMultiBlockDataSet *grid, vtkStdString varName, double *data[], int numberOfComponents)
bool GetFloat64()
Definition: PIOAdaptor.h:60
void load_variable_data(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
std::vector< double > PIOFileIndex
Definition: PIOAdaptor.h:136
int GetNumberOfTimeSteps()
Definition: PIOAdaptor.h:41
bool useHTG
Definition: PIOAdaptor.h:139
void create_amr_UG_2D(vtkMultiBlockDataSet *grid, int numberOfCells, int *cell_level, int64_t *cell_daughter, double *cell_center[2])
int GetNumberOfDefaultVariables()
Definition: PIOAdaptor.h:48
std::vector< std::string > dumpDirectory
Definition: PIOAdaptor.h:130
double GetPIOFileIndex(int step)
Definition: PIOAdaptor.h:44
void collectVariableMetaData()
int initializeGlobal(const char *DumpDescFile)
const char * GetVariableDefault(int indx)
Definition: PIOAdaptor.h:49
void SetFloat64(bool val)
Definition: PIOAdaptor.h:61
PIOAdaptor(vtkMultiProcessController *ctrl)
std::vector< double > CycleIndex
Definition: PIOAdaptor.h:134
bool useFloat64
Definition: PIOAdaptor.h:141
int collectMetaData(const char *DumpDescFile)
int TotalRank
Definition: PIOAdaptor.h:119
void SetHyperTreeGrid(bool val)
Definition: PIOAdaptor.h:53
const char * GetVariableName(int indx)
Definition: PIOAdaptor.h:47
std::string dumpBaseName
Definition: PIOAdaptor.h:129
std::vector< std::string > variableName
Definition: PIOAdaptor.h:145
void load_variable_data_UG(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
void create_amr_HTG(vtkMultiBlockDataSet *grid)
Store on/off settings for data arrays for a vtkSource.
Objects for traversal a HyperTreeGrid.
Composite dataset that organizes datasets into blocks.
Multiprocessing communication superclass.
Wrapper around std::string to keep symbols short.
Definition: vtkStdString.h:105