VTK  9.3.0
PIOAdaptor.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Kitware, Inc.
3// SPDX-License-Identifier: BSD-3-Clause
4#ifndef PIOAdaptor_h
5#define PIOAdaptor_h
6
10
11#include "PIOData.h"
12
13#include <vector>
14
15VTK_ABI_NAMESPACE_BEGIN
17
18// class to hold information about chunk/material variables
20{
21public:
22 std::string prefix;
23 std::string var; // actual variable
24 std::string baseVar; // variable used to derive actual variable
25 std::string material_name; // full name of the material
27};
28
30{
31public:
34
35 int initializeGlobal(const char* DumpDescFile);
36 int initializeDump(int timeStep);
37
38 // Time step change requires new geometry and data
41
42 int GetNumberOfTimeSteps() { return static_cast<int>(this->CycleIndex.size()); }
43 double GetSimulationTime(int step) { return this->SimulationTime[step]; }
44 double GetCycleIndex(int step) { return this->CycleIndex[step]; }
45 double GetPIOFileIndex(int step) { return this->PIOFileIndex[step]; }
46
47 int GetNumberOfVariables() { return (int)this->variableName.size(); }
48 const char* GetVariableName(int indx) { return this->variableName[indx].c_str(); }
49 int GetNumberOfDefaultVariables() { return (int)this->variableDefault.size(); }
50 const char* GetVariableDefault(int indx) { return this->variableDefault[indx].c_str(); }
51
52 // Read pio dump file AMR as hypertree grid rather than unstructured grid
53 bool GetHyperTreeGrid() { return this->useHTG; }
54 void SetHyperTreeGrid(bool val) { this->useHTG = val; }
55
56 // Read pio dump file tracer information
57 bool GetTracers() { return this->useTracer; }
58 void SetTracers(bool val) { this->useTracer = val; }
59
60 // Read pio dump file variable data as 64 bit float
61 bool GetFloat64() { return this->useFloat64; }
62 void SetFloat64(bool val) { this->useFloat64 = val; }
63
64protected:
65 // Collect the metadata
66 int parsePIOFile(const char* DumpDescFile);
67 int collectMetaData(const char* DumpDescFile);
70 void addMaterialVariable(vtkStdString& pioFieldName, std::vector<std::string> matident);
72 std::string& prefix, std::string& baseVar, std::string& var, std::vector<std::string> matident);
73 std::string trimString(const std::string& str);
74
75 // Create the unstructured grid for tracers
77
78 // Create the unstructured grid for AMR
80
82 int numberOfCells, // Number of cells all levels
83 int* cell_level, // Level within AMR
84 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
85 double* cell_center[1]); // Cell center
86
88 int numberOfCells, // Number of cells all levels
89 int* cell_level, // Level within AMR
90 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
91 double* cell_center[2]); // Cell center
92
94 int numberOfCells, // Number of cells all levels
95 int* cell_level, // Level within AMR
96 int64_t* cell_daughter, // Daughter ID, 0 indicates no daughter
97 double* cell_center[3]); // Cell center
98
99 // Create the hypertree grid
101
102 int count_hypertree(int64_t curIndex, int64_t* daughter);
103
105 vtkHyperTreeGridNonOrientedCursor* treeCursor, int64_t curIndex, int64_t* daughter);
106
107 // Add variable data to the unstructured grid
110 int64_t* daughter, // Indicates top level cell or not
111 double* data[], // Data for all cells
112 int numberOfCells,
113 int numberOfComponents); // Number of components in data
114
115 // Add variable data to the hypertree grid
118 double* data[], // Data for all cells
119 int numberOfComponents); // Number of components in data
120
121 // Used in parallel reader and load balancing
123 int Rank;
125
126 // Structure to access the dump file data
128
129 // Time series of dumps
130 std::string descFileName; // name.pio
131 std::string dumpBaseName; // base name to use for dumps
132 std::vector<std::string> dumpDirectory; // directories holding dumps
133 std::vector<std::string> dumpFileName; // all dump files
134
135 // Time step information
136 std::vector<double> CycleIndex; // Times as cycle index
137 std::vector<double> SimulationTime; // Times as simulation time
138 std::vector<double> PIOFileIndex; // Index into dump files
139
140 // Type of block structures to create within multiblock dataset
141 bool useHTG;
145
146 // Cell variable data and initially enabled variables
147 std::vector<std::string> variableName;
148 std::vector<std::string> variableDefault;
149
150 // total number of cells in the mesh. needed when loading material variables.
151 // obtained by summing all values in pio field global_numcells
152 int64_t numCells;
153
154 // Record the ordering of the cells when building the hypertree grid
155 // Needed so that the data will line up correctly
156 std::vector<int> indexNodeLeaf;
157
158 // list of material variables
159 std::map<std::string, PIOMaterialVariable*> matVariables;
161
162 struct AdaptorImpl;
163 AdaptorImpl* Impl;
164};
165
166VTK_ABI_NAMESPACE_END
167#endif
std::vector< double > SimulationTime
Definition PIOAdaptor.h:137
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:47
int numMaterials
Definition PIOAdaptor.h:160
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:142
int count_hypertree(int64_t curIndex, int64_t *daughter)
double GetSimulationTime(int step)
Definition PIOAdaptor.h:43
void add_amr_UG_scalar(vtkMultiBlockDataSet *grid, vtkStdString varName, int64_t *daughter, double *data[], int numberOfCells, int numberOfComponents)
vtkMultiProcessController * Controller
Definition PIOAdaptor.h:122
PIO_DATA * pioData
Definition PIOAdaptor.h:127
bool GetHyperTreeGrid()
Definition PIOAdaptor.h:53
void create_geometry(vtkMultiBlockDataSet *grid)
void SetTracers(bool val)
Definition PIOAdaptor.h:58
std::string descFileName
Definition PIOAdaptor.h:130
std::map< std::string, PIOMaterialVariable * > matVariables
Definition PIOAdaptor.h:159
int parsePIOFile(const char *DumpDescFile)
double GetCycleIndex(int step)
Definition PIOAdaptor.h:44
void addMaterialVariable(vtkStdString &pioFieldName, std::vector< std::string > matident)
bool hasTracers
Definition PIOAdaptor.h:144
int initializeDump(int timeStep)
void load_variable_data_HTG(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
void create_tracer_UG(vtkMultiBlockDataSet *grid)
bool GetTracers()
Definition PIOAdaptor.h:57
std::vector< int > indexNodeLeaf
Definition PIOAdaptor.h:156
std::vector< std::string > dumpFileName
Definition PIOAdaptor.h:133
void addMaterialVariableEntries(std::string &prefix, std::string &baseVar, std::string &var, std::vector< std::string > matident)
std::vector< std::string > variableDefault
Definition PIOAdaptor.h:148
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:61
void load_variable_data(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
std::vector< double > PIOFileIndex
Definition PIOAdaptor.h:138
int GetNumberOfTimeSteps()
Definition PIOAdaptor.h:42
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:49
std::vector< std::string > dumpDirectory
Definition PIOAdaptor.h:132
double GetPIOFileIndex(int step)
Definition PIOAdaptor.h:45
void collectVariableMetaData()
int initializeGlobal(const char *DumpDescFile)
AdaptorImpl * Impl
Definition PIOAdaptor.h:163
const char * GetVariableDefault(int indx)
Definition PIOAdaptor.h:50
int64_t numCells
Definition PIOAdaptor.h:152
void SetFloat64(bool val)
Definition PIOAdaptor.h:62
PIOAdaptor(vtkMultiProcessController *ctrl)
std::vector< double > CycleIndex
Definition PIOAdaptor.h:136
void collectMaterialVariableMetaData()
bool useFloat64
Definition PIOAdaptor.h:143
int collectMetaData(const char *DumpDescFile)
void SetHyperTreeGrid(bool val)
Definition PIOAdaptor.h:54
const char * GetVariableName(int indx)
Definition PIOAdaptor.h:48
std::string dumpBaseName
Definition PIOAdaptor.h:131
std::vector< std::string > variableName
Definition PIOAdaptor.h:147
void load_variable_data_UG(vtkMultiBlockDataSet *grid, vtkDataArraySelection *cellSelection)
void create_amr_HTG(vtkMultiBlockDataSet *grid)
std::string material_name
Definition PIOAdaptor.h:25
std::string var
Definition PIOAdaptor.h:23
std::string prefix
Definition PIOAdaptor.h:22
uint32_t material_number
Definition PIOAdaptor.h:26
std::string baseVar
Definition PIOAdaptor.h:24
Store on/off settings for data arrays, etc.
Objects for traversal a HyperTreeGrid.
Composite dataset that organizes datasets into blocks.
Multiprocessing communication superclass.
Wrapper around std::string to keep symbols short.