VTK  9.3.1
vtkCGNSReader.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright 2013-2014 Mickael Philit
3 // SPDX-License-Identifier: BSD-3-Clause
19 #ifndef vtkCGNSReader_h
20 #define vtkCGNSReader_h
21 
22 #include "vtkIOCGNSReaderModule.h" // for export macro
24 #include "vtkNew.h" // for vtkNew.
25 
26 #include <string> // for std::string
27 
28 namespace CGNSRead
29 {
30 VTK_ABI_NAMESPACE_BEGIN
31 class vtkCGNSMetaData;
32 VTK_ABI_NAMESPACE_END
33 }
34 
35 VTK_ABI_NAMESPACE_BEGIN
39 
40 class VTKIOCGNSREADER_EXPORT vtkCGNSReader : public vtkMultiBlockDataSetAlgorithm
41 {
42 public:
43  static vtkCGNSReader* New();
45  void PrintSelf(ostream& os, vtkIndent indent) override;
46 
48  {
49  CELL_DATA = 0,
50  FACE_DATA
51  };
52 
54 
64  vtkSetClampMacro(DataLocation, int, vtkCGNSReader::CELL_DATA, vtkCGNSReader::FACE_DATA);
65  vtkGetMacro(DataLocation, int);
67 
69 
72  vtkSetStdStringFromCharMacro(FileName);
73  vtkGetCharFromStdStringMacro(FileName);
75 
79  int CanReadFile(VTK_FILEPATH const char* filename);
80 
84  vtkDataArraySelection* GetBaseSelection();
85 
89  vtkDataArraySelection* GetFamilySelection();
90 
92 
98  int GetBaseArrayStatus(const char* name);
99  void SetBaseArrayStatus(const char* name, int status);
100  void DisableAllBases();
101  void EnableAllBases();
102  int GetNumberOfBaseArrays();
103  const char* GetBaseArrayName(int index);
105 
107 
111  int GetNumberOfFamilyArrays();
112  const char* GetFamilyArrayName(int index);
113  void SetFamilyArrayStatus(const char* name, int status);
114  int GetFamilyArrayStatus(const char* name);
115  void EnableAllFamilies();
116  void DisableAllFamilies();
118 
120 
124  int GetNumberOfPointArrays();
125  const char* GetPointArrayName(int index);
126  int GetPointArrayStatus(const char* name);
127  void SetPointArrayStatus(const char* name, int status);
128  void DisableAllPointArrays();
129  void EnableAllPointArrays();
131 
133 
137  int GetNumberOfCellArrays();
138  const char* GetCellArrayName(int index);
139  int GetCellArrayStatus(const char* name);
140  void SetCellArrayStatus(const char* name, int status);
141  void DisableAllCellArrays();
142  void EnableAllCellArrays();
144 
146 
150  int GetNumberOfFaceArrays();
151  const char* GetFaceArrayName(int index);
152  int GetFaceArrayStatus(const char* name);
153  void SetFaceArrayStatus(const char* name, int status);
154  void DisableAllFaceArrays();
155  void EnableAllFaceArrays();
157 
159 
163  vtkSetMacro(DoublePrecisionMesh, int);
164  vtkGetMacro(DoublePrecisionMesh, int);
165  vtkBooleanMacro(DoublePrecisionMesh, int);
167 
169 
173  vtkSetMacro(LoadBndPatch, bool);
174  vtkGetMacro(LoadBndPatch, bool);
175  vtkBooleanMacro(LoadBndPatch, bool);
177 
179 
183  vtkSetMacro(LoadMesh, bool);
184  vtkGetMacro(LoadMesh, bool);
185  vtkBooleanMacro(LoadMesh, bool);
187 
189 
193  vtkSetMacro(Use3DVector, bool);
194  vtkGetMacro(Use3DVector, bool);
195  vtkBooleanMacro(Use3DVector, bool);
197 
199 
207  vtkSetMacro(CreateEachSolutionAsBlock, int);
208  vtkGetMacro(CreateEachSolutionAsBlock, int);
209  vtkBooleanMacro(CreateEachSolutionAsBlock, int);
211 
213 
226  vtkSetMacro(IgnoreFlowSolutionPointers, bool);
227  vtkGetMacro(IgnoreFlowSolutionPointers, bool);
228  vtkBooleanMacro(IgnoreFlowSolutionPointers, bool);
230 
232 
238  vtkSetMacro(UseUnsteadyPattern, bool);
239  vtkGetMacro(UseUnsteadyPattern, bool);
240  vtkBooleanMacro(UseUnsteadyPattern, bool);
242 
244 
249  vtkSetMacro(UnsteadySolutionStartTimestep, int);
250  vtkGetMacro(UnsteadySolutionStartTimestep, int);
252 
254 
259  vtkSetMacro(DistributeBlocks, bool);
260  vtkGetMacro(DistributeBlocks, bool);
261  vtkBooleanMacro(DistributeBlocks, bool);
263 
265 
271  void SetCacheMesh(bool enable);
272  vtkGetMacro(CacheMesh, bool);
273  vtkBooleanMacro(CacheMesh, bool);
275 
277 
283  void SetCacheConnectivity(bool enable);
284  vtkGetMacro(CacheConnectivity, bool);
285  vtkBooleanMacro(CacheConnectivity, bool);
287 
289 
294  void SetController(vtkMultiProcessController* c);
295  vtkGetObjectMacro(Controller, vtkMultiProcessController);
297 
302  void Broadcast(vtkMultiProcessController* ctrl);
303 
307  static vtkInformationStringKey* FAMILY();
308 
309 protected:
310  vtkCGNSReader();
311  ~vtkCGNSReader() override;
312 
316 
317  int GetCurvilinearZone(
318  int base, int zone, int cell_dim, int phys_dim, void* zsize, vtkMultiBlockDataSet* mbase);
319 
320  int GetUnstructuredZone(
321  int base, int zone, int cell_dim, int phys_dim, void* zsize, vtkMultiBlockDataSet* mbase);
322 
327  int ReadUserDefinedData(int zone, vtkMultiBlockDataSet* mbase);
328 
329  vtkMultiProcessController* Controller = nullptr;
330  vtkIdType ProcRank = 0;
331  vtkIdType ProcSize = 1;
332 
335 
339 
340 private:
341  vtkCGNSReader(const vtkCGNSReader&) = delete;
342  void operator=(const vtkCGNSReader&) = delete;
343 
344  std::string FileName;
345  int DataLocation = vtkCGNSReader::CELL_DATA;
346  bool LoadBndPatch = false;
347  bool LoadMesh = true;
348  int DoublePrecisionMesh = 1;
349  int CreateEachSolutionAsBlock = 0;
350  bool IgnoreFlowSolutionPointers = false;
351  bool UseUnsteadyPattern = false;
352  bool DistributeBlocks = true;
353  bool CacheMesh = false;
354  bool CacheConnectivity = false;
355  bool Use3DVector = true;
356  int UnsteadySolutionStartTimestep = 0;
357 
358  // For internal cgio calls (low level IO)
359  int cgioNum; // cgio file reference
360  double rootId; // id of root node
361  double currentZoneId; // id of node currently being read (zone)
362 
363  unsigned int NumberOfBases = 0;
364  int ActualTimeStep = 0;
365 
366  class vtkPrivate;
367  vtkPrivate* Internals;
368  friend class vtkPrivate;
369 };
370 
371 VTK_ABI_NAMESPACE_END
372 #endif // vtkCGNSReader_h
vtkNew< vtkDataArraySelection > PointDataArraySelection
Store vtkAlgorithm input/output information.
vtkNew< vtkDataArraySelection > FaceDataArraySelection
int vtkIdType
Definition: vtkType.h:315
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:28
Store on/off settings for data arrays, etc.
This file defines functions used by vtkCGNSReader and vtkCGNSReaderInternal.
Definition: cgio_helpers.h:18
vtkCGNSReader creates a multi-block dataset and reads unstructured grids and structured meshes from b...
Definition: vtkCGNSReader.h:40
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.