VTK  9.2.6
vtkCGNSReader.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCGNSReader.h
5 
6  Copyright (c) Ken Martin, Will Schrodeder, 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  =========================================================================*/
15 // Copyright 2013-2014 Mickael Philit.
16 
32 #ifndef vtkCGNSReader_h
33 #define vtkCGNSReader_h
34 
35 #include "vtkIOCGNSReaderModule.h" // for export macro
37 #include "vtkNew.h" // for vtkNew.
38 
39 #include <string> // for std::string
40 
43 
44 namespace CGNSRead
45 {
46 class vtkCGNSMetaData;
47 }
48 
50 class VTKIOCGNSREADER_EXPORT vtkCGNSReader : public vtkMultiBlockDataSetAlgorithm
51 {
52 public:
53  static vtkCGNSReader* New();
55  void PrintSelf(ostream& os, vtkIndent indent) override;
56 
58  {
59  CELL_DATA = 0,
60  FACE_DATA
61  };
62 
64 
74  vtkSetClampMacro(DataLocation, int, vtkCGNSReader::CELL_DATA, vtkCGNSReader::FACE_DATA);
75  vtkGetMacro(DataLocation, int);
77 
79 
82  vtkSetStdStringFromCharMacro(FileName);
83  vtkGetCharFromStdStringMacro(FileName);
85 
89  int CanReadFile(VTK_FILEPATH const char* filename);
90 
94  vtkDataArraySelection* GetBaseSelection();
95 
99  vtkDataArraySelection* GetFamilySelection();
100 
102 
108  int GetBaseArrayStatus(const char* name);
109  void SetBaseArrayStatus(const char* name, int status);
110  void DisableAllBases();
111  void EnableAllBases();
112  int GetNumberOfBaseArrays();
113  const char* GetBaseArrayName(int index);
115 
117 
121  int GetNumberOfFamilyArrays();
122  const char* GetFamilyArrayName(int index);
123  void SetFamilyArrayStatus(const char* name, int status);
124  int GetFamilyArrayStatus(const char* name);
125  void EnableAllFamilies();
126  void DisableAllFamilies();
128 
130 
134  int GetNumberOfPointArrays();
135  const char* GetPointArrayName(int index);
136  int GetPointArrayStatus(const char* name);
137  void SetPointArrayStatus(const char* name, int status);
138  void DisableAllPointArrays();
139  void EnableAllPointArrays();
141 
143 
147  int GetNumberOfCellArrays();
148  const char* GetCellArrayName(int index);
149  int GetCellArrayStatus(const char* name);
150  void SetCellArrayStatus(const char* name, int status);
151  void DisableAllCellArrays();
152  void EnableAllCellArrays();
154 
156 
160  int GetNumberOfFaceArrays();
161  const char* GetFaceArrayName(int index);
162  int GetFaceArrayStatus(const char* name);
163  void SetFaceArrayStatus(const char* name, int status);
164  void DisableAllFaceArrays();
165  void EnableAllFaceArrays();
167 
168  vtkSetMacro(DoublePrecisionMesh, int);
169  vtkGetMacro(DoublePrecisionMesh, int);
170  vtkBooleanMacro(DoublePrecisionMesh, int);
171 
173 
177  vtkSetMacro(LoadBndPatch, bool);
178  vtkGetMacro(LoadBndPatch, bool);
179  vtkBooleanMacro(LoadBndPatch, bool);
181 
183 
187  vtkSetMacro(LoadMesh, bool);
188  vtkGetMacro(LoadMesh, bool);
189  vtkBooleanMacro(LoadMesh, bool);
191 
193 
196  vtkSetMacro(Use3DVector, bool);
197  vtkGetMacro(Use3DVector, bool);
198  vtkBooleanMacro(Use3DVector, bool);
200 
208  vtkSetMacro(CreateEachSolutionAsBlock, int);
209  vtkGetMacro(CreateEachSolutionAsBlock, int);
210  vtkBooleanMacro(CreateEachSolutionAsBlock, int);
211 
225  vtkSetMacro(IgnoreFlowSolutionPointers, bool);
226  vtkGetMacro(IgnoreFlowSolutionPointers, bool);
227  vtkBooleanMacro(IgnoreFlowSolutionPointers, bool);
228 
235  vtkSetMacro(UseUnsteadyPattern, bool);
236  vtkGetMacro(UseUnsteadyPattern, bool);
237  vtkBooleanMacro(UseUnsteadyPattern, bool);
238 
244  vtkSetMacro(DistributeBlocks, bool);
245  vtkGetMacro(DistributeBlocks, bool);
246  vtkBooleanMacro(DistributeBlocks, bool);
247 
249 
254  void SetCacheMesh(bool enable);
255  vtkGetMacro(CacheMesh, bool);
256  vtkBooleanMacro(CacheMesh, bool);
257 
259 
264  void SetCacheConnectivity(bool enable);
265  vtkGetMacro(CacheConnectivity, bool);
266  vtkBooleanMacro(CacheConnectivity, bool);
267 
269 
274  void SetController(vtkMultiProcessController* c);
275  vtkGetObjectMacro(Controller, vtkMultiProcessController);
277 
282  void Broadcast(vtkMultiProcessController* ctrl);
283 
287  static vtkInformationStringKey* FAMILY();
288 
289 protected:
290  vtkCGNSReader();
291  ~vtkCGNSReader() override;
292 
296 
297  int GetCurvilinearZone(
298  int base, int zone, int cell_dim, int phys_dim, void* zsize, vtkMultiBlockDataSet* mbase);
299 
300  int GetUnstructuredZone(
301  int base, int zone, int cell_dim, int phys_dim, void* zsize, vtkMultiBlockDataSet* mbase);
302 
303  vtkMultiProcessController* Controller = nullptr;
304  vtkIdType ProcRank = 0;
305  vtkIdType ProcSize = 1;
306 
309 
313 
314 private:
315  vtkCGNSReader(const vtkCGNSReader&) = delete;
316  void operator=(const vtkCGNSReader&) = delete;
317 
318  std::string FileName = "";
319  int DataLocation = vtkCGNSReader::CELL_DATA;
320  bool LoadBndPatch = false; // option to set section loading for unstructured grid
321  bool LoadMesh = true; // option to enable/disable mesh loading
322  int DoublePrecisionMesh = 1; // option to set mesh loading to double precision
323  int CreateEachSolutionAsBlock = 0; // debug option to create
324  bool IgnoreFlowSolutionPointers = false;
325  bool UseUnsteadyPattern = false;
326  bool DistributeBlocks = true;
327  bool CacheMesh = false;
328  bool CacheConnectivity = false;
329  bool Use3DVector = true;
330 
331  // For internal cgio calls (low level IO)
332  int cgioNum; // cgio file reference
333  double rootId; // id of root node
334  double currentId; // id of node currently being read (zone)
335 
336  unsigned int NumberOfBases = 0;
337  int ActualTimeStep = 0;
338 
339  class vtkPrivate;
340  vtkPrivate* Internals;
341  friend class vtkPrivate;
342 };
343 
344 #endif // vtkCGNSReader_h
vtkNew< vtkDataArraySelection > PointDataArraySelection
Store vtkAlgorithm input/output information.
vtkNew< vtkDataArraySelection > FaceDataArraySelection
int vtkIdType
Definition: vtkType.h:332
Superclass for algorithms that produce only vtkMultiBlockDataSet as output.
vtkNew< vtkDataArraySelection > BaseSelection
virtual int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
vtkNew< vtkDataArraySelection > FamilySelection
Key for string values in vtkInformation.
int FillOutputPortInformation(int port, vtkInformation *info) override
Fill the output port information objects for this algorithm.
static vtkMultiBlockDataSetAlgorithm * New()
a simple class to control print indentation
Definition: vtkIndent.h:39
Store on/off settings for data arrays, etc.
vtkCGNSReader creates a multi-block dataset and reads unstructured grids and structured meshes from b...
Definition: vtkCGNSReader.h:50
vtkNew< vtkDataArraySelection > CellDataArraySelection
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
#define VTK_FILEPATH
Composite dataset that organizes datasets into blocks.
Store zero or more vtkInformation instances.
virtual int RequestInformation(vtkInformation *, vtkInformationVector **, vtkInformationVector *)
This is called by the superclass.
Multiprocessing communication superclass.