VTK  9.3.20240422
vtkCGNSReaderInternal.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
13#ifndef vtkCGNSReaderInternal_h
14#define vtkCGNSReaderInternal_h
15
16#include <iostream>
17#include <map>
18#include <string.h> // for inline strcmp
19#include <string>
20#include <vector>
21
22#include "vtkCGNSReader.h"
24#include "vtkIdTypeArray.h"
26#include "vtkNew.h"
27#include "vtkPoints.h"
28
29// .clang-format off
30#include "vtk_cgns.h"
31#include VTK_CGNS(cgnslib.h)
32#include VTK_CGNS(cgns_io.h)
33// .clang-format on
34
35namespace CGNSRead
36{
37namespace detail
38{
39VTK_ABI_NAMESPACE_BEGIN
40template <typename T>
42{
43 static const bool value = false;
44};
45
46template <>
47struct is_double<double>
48{
49 static const bool value = true;
50};
51
52template <typename T>
54{
55 static const bool value = false;
56};
57
58template <>
59struct is_float<float>
60{
61 static const bool value = true;
62};
63VTK_ABI_NAMESPACE_END
64}
65
66namespace detail
67{
68VTK_ABI_NAMESPACE_BEGIN
69template <typename T>
70constexpr const char* cgns_type_name() noexcept
71{
72 return "MT";
73}
74
75template <>
76constexpr const char* cgns_type_name<float>() noexcept
77{
78 return "R4";
79}
80
81template <>
82constexpr const char* cgns_type_name<double>() noexcept
83{
84 return "R8";
85}
86
87template <>
88constexpr const char* cgns_type_name<vtkTypeInt32>() noexcept
89{
90 return "I4";
91}
92
93template <>
94constexpr const char* cgns_type_name<vtkTypeInt64>() noexcept
95{
96 return "I8";
97}
98VTK_ABI_NAMESPACE_END
99}
100
101VTK_ABI_NAMESPACE_BEGIN
102typedef char char_33[33];
103
104//------------------------------------------------------------------------------
105// Cell type to cell dimension
106const std::map<CGNS_ENUMT(ElementType_t), int> CellDimensions = {
107 { CGNS_ENUMV(ElementTypeUserDefined), -1 }, { CGNS_ENUMV(ElementTypeNull), -1 },
108 { CGNS_ENUMV(NODE), 0 }, { CGNS_ENUMV(BAR_2), 1 }, { CGNS_ENUMV(BAR_3), 1 },
109 { CGNS_ENUMV(TRI_3), 2 }, { CGNS_ENUMV(TRI_6), 2 }, { CGNS_ENUMV(QUAD_4), 2 },
110 { CGNS_ENUMV(QUAD_8), 2 }, { CGNS_ENUMV(QUAD_9), 2 }, { CGNS_ENUMV(TETRA_4), 3 },
111 { CGNS_ENUMV(TETRA_10), 3 }, { CGNS_ENUMV(PYRA_5), 3 }, { CGNS_ENUMV(PYRA_14), 3 },
112 { CGNS_ENUMV(PENTA_6), 3 }, { CGNS_ENUMV(PENTA_15), 3 }, { CGNS_ENUMV(PENTA_18), 3 },
113 { CGNS_ENUMV(HEXA_8), 3 }, { CGNS_ENUMV(HEXA_20), 3 }, { CGNS_ENUMV(HEXA_27), 3 },
114 { CGNS_ENUMV(MIXED), -1 }, { CGNS_ENUMV(PYRA_13), 3 }, { CGNS_ENUMV(NGON_n), 2 },
115 { CGNS_ENUMV(NFACE_n), 3 }, { CGNS_ENUMV(BAR_4), 1 }, { CGNS_ENUMV(TRI_9), 2 },
116 { CGNS_ENUMV(TRI_10), 2 }, { CGNS_ENUMV(QUAD_12), 2 }, { CGNS_ENUMV(QUAD_16), 2 },
117 { CGNS_ENUMV(TETRA_16), 3 }, { CGNS_ENUMV(TETRA_20), 3 }, { CGNS_ENUMV(PYRA_21), 3 },
118 { CGNS_ENUMV(PYRA_29), 3 }, { CGNS_ENUMV(PYRA_30), 3 }, { CGNS_ENUMV(PENTA_24), 3 },
119 { CGNS_ENUMV(PENTA_38), 3 }, { CGNS_ENUMV(PENTA_40), 3 }, { CGNS_ENUMV(HEXA_32), 3 },
120 { CGNS_ENUMV(HEXA_56), 3 }, { CGNS_ENUMV(HEXA_64), 3 }, { CGNS_ENUMV(BAR_5), 1 },
121 { CGNS_ENUMV(TRI_12), 2 }, { CGNS_ENUMV(TRI_15), 2 }, { CGNS_ENUMV(QUAD_P4_16), 2 },
122 { CGNS_ENUMV(QUAD_25), 2 }, { CGNS_ENUMV(TETRA_22), 3 }, { CGNS_ENUMV(TETRA_34), 3 },
123 { CGNS_ENUMV(TETRA_35), 3 }, { CGNS_ENUMV(PYRA_P4_29), 3 }, { CGNS_ENUMV(PYRA_50), 3 },
124 { CGNS_ENUMV(PYRA_55), 3 }, { CGNS_ENUMV(PENTA_33), 3 }, { CGNS_ENUMV(PENTA_66), 3 },
125 { CGNS_ENUMV(PENTA_75), 3 }, { CGNS_ENUMV(HEXA_44), 3 }, { CGNS_ENUMV(HEXA_98), 3 },
126 { CGNS_ENUMV(HEXA_125), 3 }
127};
128
129//------------------------------------------------------------------------------
130class vtkCGNSArraySelection : public std::map<std::string, bool>
131{
132public:
133 void Merge(const vtkCGNSArraySelection& other)
134 {
135 vtkCGNSArraySelection::const_iterator iter = other.begin();
136 for (; iter != other.end(); ++iter)
137 {
138 (*this)[iter->first] = iter->second;
139 }
140 }
141
142 void AddArray(const char* name, bool status = true) { (*this)[name] = status; }
143
144 bool ArrayIsEnabled(const char* name)
145 {
146 vtkCGNSArraySelection::iterator iter = this->find(name);
147 if (iter != this->end())
148 {
149 return iter->second;
150 }
151
152 // don't know anything about this array, enable it by default.
153 return true;
154 }
155
156 bool HasArray(const char* name)
157 {
158 vtkCGNSArraySelection::iterator iter = this->find(name);
159 return (iter != this->end());
160 }
161
162 int GetArraySetting(const char* name) { return this->ArrayIsEnabled(name) ? 1 : 0; }
163
164 void SetArrayStatus(const char* name, bool status) { this->AddArray(name, status); }
165
166 const char* GetArrayName(int index)
167 {
168 int cc = 0;
169 for (vtkCGNSArraySelection::iterator iter = this->begin(); iter != this->end(); ++iter)
170 {
171
172 if (cc == index)
173 {
174 return iter->first.c_str();
175 }
176 cc++;
177 }
178 return nullptr;
179 }
180
181 int GetNumberOfArrays() { return static_cast<int>(this->size()); }
182};
183
184//------------------------------------------------------------------------------
185typedef struct
186{
187 int cnt; // 0 1 or 3
188 int pos; // variable position in zone
191 CGNS_ENUMT(DataType_t) dt;
193} Variable;
194
195//------------------------------------------------------------------------------
196typedef struct
197{
200 CGNS_ENUMT(DataType_t) dt;
203
204//------------------------------------------------------------------------------
205typedef struct
206{
209 int xyzIndex[3];
210} CGNSVector;
211
212//------------------------------------------------------------------------------
213typedef struct
214{
219
220//------------------------------------------------------------------------------
222{
223public:
225 std::string family;
227 : family(32, '\0')
228 {
229 this->name[0] = '\0';
230 }
231};
232
233//------------------------------------------------------------------------------
235{
236public:
238 std::string family;
239 std::vector<CGNSRead::ZoneBCInformation> bcs;
241 : family(32, '\0')
242 {
243 this->name[0] = '\0';
244 }
245};
246
247//------------------------------------------------------------------------------
249{
250public:
251 std::string name;
252 bool isBC;
253};
254
255//------------------------------------------------------------------------------
257{
258public:
260
261 int32_t cellDim;
262 int32_t physicalDim;
263 //
265
266 std::vector<int32_t> steps;
267 std::vector<double> times;
268
269 // For unsteady meshes :
270 // if useGridPointers == True:
271 // loadGridPointers for first zone
272 // and assume every zone use the same
273 // notation
274 // else :
275 // assume only one grid is stored
276 // only first grid is read
277 //
278 // For unsteady flow
279 // if useFlowPointers == True :
280 // same behavior as GridPointers
281 // else if ( nstates > 1 ) :
282 // assume flow_solution are sorted
283 // to keep VisIt like behavior
284 // else :
285 // only first solution is read
286 //
287
288 bool useGridPointers; // for unsteady mesh
289 bool useFlowPointers; // for unsteady flow
290
291 std::vector<CGNSRead::FamilyInformation> family;
292 std::map<std::string, double> referenceState;
293
294 std::vector<CGNSRead::ZoneInformation> zones;
295
297
298 // std::vector<CGNSRead::zone> zone;
302};
303
304//==============================================================================
305
307
311bool ReadBase(vtkCGNSReader* reader, const BaseInformation& baseInfo);
313 vtkCGNSReader* reader, const BaseInformation& baseInfo, const ZoneInformation& zoneInfo);
315bool ReadPatch(vtkCGNSReader* reader, const BaseInformation&, const ZoneInformation& zoneInfo,
316 const std::string& patchFamilyname);
318
319//==============================================================================
321{
322public:
327 bool Parse(const char* cgnsFileName);
328
332 int GetNumberOfBaseNodes() { return static_cast<int>(this->baseList.size()); }
333
337 const CGNSRead::BaseInformation& GetBase(int numBase) { return this->baseList[numBase]; }
338
342 std::vector<double>& GetTimes() { return this->GlobalTime; }
343
347 void PrintSelf(std::ostream& os);
348
349 void Broadcast(vtkMultiProcessController* controller, int rank);
350
352
355 vtkCGNSMetaData() = default;
356 ~vtkCGNSMetaData() = default;
358
359private:
360 vtkCGNSMetaData(const vtkCGNSMetaData&) = delete;
361 void operator=(const vtkCGNSMetaData&) = delete;
362
363 std::vector<CGNSRead::BaseInformation> baseList;
364 std::string LastReadFilename;
365 // Not very elegant :
366 std::vector<double> GlobalTime;
367};
368
369//------------------------------------------------------------------------------
370// compare name return true if name1 == name2
371inline bool compareName(const char_33 nameOne, const char_33 nameTwo)
372{
373 return (strncmp(nameOne, nameTwo, 32) == 0);
374}
375
376//------------------------------------------------------------------------------
377// remove trailing whitespaces
379{
380 char* end = name + strlen(name) - 1;
381 while (end >= name && isspace(*end))
382 {
383 --end;
384 }
385 ++end;
386 assert(end >= name && end < name + 33);
387 *end = '\0';
388}
389
390//------------------------------------------------------------------------------
391// get vector from name
392inline std::vector<CGNSVector>::iterator getVectorFromName(
393 std::vector<CGNSVector>& vectorList, const char_33 name)
394{
395 for (std::vector<CGNSVector>::iterator iter = vectorList.begin(); iter != vectorList.end();
396 ++iter)
397 {
398 if (strncmp(iter->name, name, 31) == 0)
399 {
400 return iter;
401 }
402 }
403 return vectorList.end();
404}
405
406//------------------------------------------------------------------------------
407inline bool isACGNSVariable(const std::vector<CGNSVariable>& varList, const char_33 name)
408{
409 for (std::vector<CGNSVariable>::const_iterator iter = varList.begin(); iter != varList.end();
410 ++iter)
411 {
412 if (strncmp(iter->name, name, 32) == 0)
413 {
414 return true;
415 }
416 }
417 return false;
418}
419
420//------------------------------------------------------------------------------
421void fillVectorsFromVars(std::vector<CGNSRead::CGNSVariable>& vars,
422 std::vector<CGNSRead::CGNSVector>& vectors, int physicalDim);
423//------------------------------------------------------------------------------
424int setUpRind(int cgioNum, double rindId, int* rind);
425//------------------------------------------------------------------------------
431 int cgioNum, double parentId, const char* label, double* id, const char* name = nullptr);
432//------------------------------------------------------------------------------
433int get_section_connectivity(int cgioNum, double cgioSectionId, int dim, const cgsize_t* srcStart,
434 const cgsize_t* srcEnd, const cgsize_t* srcStride, const cgsize_t* memStart,
435 const cgsize_t* memEnd, const cgsize_t* memStride, const cgsize_t* memDim,
436 vtkIdType* localElements);
437//------------------------------------------------------------------------------
438int get_section_start_offset(int cgioNum, double cgioSectionId, int dim, const cgsize_t* srcStart,
439 const cgsize_t* srcEnd, const cgsize_t* srcStride, const cgsize_t* memStart,
440 const cgsize_t* memEnd, const cgsize_t* memStride, const cgsize_t* memDim,
441 vtkIdType* localElementsIdx);
442//------------------------------------------------------------------------------
443int get_section_parent_elements(int cgioNum, double cgioSectionId, int dim,
444 const cgsize_t* srcStart, const cgsize_t* srcEnd, const cgsize_t* srcStride,
445 const cgsize_t* memStart, const cgsize_t* memEnd, const cgsize_t* memStride,
446 const cgsize_t* memDim, vtkIdType* localElementsIdx);
447//------------------------------------------------------------------------------
449 CGNS_ENUMT(ElementType_t) elemType, bool& higherOrderWarning, bool& cgnsOrderFlag);
450//------------------------------------------------------------------------------
451void CGNS2VTKorder(vtkIdType size, const int* cells_types, vtkIdType* elements);
452//------------------------------------------------------------------------------
454 vtkIdType size, int cell_type, vtkIdType numPointsPerCell, vtkIdType* elements);
455//------------------------------------------------------------------------------
456template <typename T, typename Y>
457int get_XYZ_mesh(int cgioNum, const std::vector<double>& gridChildId,
458 const std::size_t& nCoordsArray, int cellDim, vtkIdType nPts, const cgsize_t* srcStart,
459 const cgsize_t* srcEnd, const cgsize_t* srcStride, const cgsize_t* memStart,
460 const cgsize_t* memEnd, const cgsize_t* memStride, const cgsize_t* memDims, vtkPoints* points)
461{
462 T* coords = static_cast<T*>(points->GetVoidPointer(0));
463 T* currentCoord = static_cast<T*>(&(coords[0]));
464
465 CGNSRead::char_33 coordName;
466 std::size_t len;
467 bool sameType = true;
468 double coordId;
469
470 memset(coords, 0, 3 * nPts * sizeof(T));
471
472 for (std::size_t c = 1; c <= nCoordsArray; ++c)
473 {
474 // Read CoordName
475 if (cgio_get_name(cgioNum, gridChildId[c - 1], coordName) != CG_OK)
476 {
477 char message[81];
478 cgio_error_message(message);
479 std::cerr << "get_XYZ_mesh : cgio_get_name :" << message;
480 }
481
482 // Read node data type
483 CGNSRead::char_33 dataType;
484 if (cgio_get_data_type(cgioNum, gridChildId[c - 1], dataType))
485 {
486 continue;
487 }
488
489 if (strcmp(dataType, "R8") == 0)
490 {
491 const bool doubleType = detail::is_double<T>::value;
492 sameType = doubleType;
493 }
494 else if (strcmp(dataType, "R4") == 0)
495 {
496 const bool floatType = detail::is_float<T>::value;
497 sameType = floatType;
498 }
499 else
500 {
501 std::cerr << "Invalid datatype for GridCoordinates\n";
502 continue;
503 }
504
505 // Determine direction X,Y,Z
506 len = strlen(coordName) - 1;
507 switch (coordName[len])
508 {
509 case 'X':
510 currentCoord = static_cast<T*>(&(coords[0]));
511 break;
512 case 'Y':
513 currentCoord = static_cast<T*>(&(coords[1]));
514 break;
515 case 'Z':
516 currentCoord = static_cast<T*>(&(coords[2]));
517 break;
518 }
519
520 coordId = gridChildId[c - 1];
521
522 // quick transfer of data if same data types
523 if (sameType == true)
524 {
525 constexpr const char* dtNameT = detail::cgns_type_name<T>();
526 if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameT, cellDim,
527 memEnd, memStart, memEnd, memStride, (void*)currentCoord))
528 {
529 char message[81];
530 cgio_error_message(message);
531 std::cerr << "cgio_read_data_type :" << message;
532 }
533 }
534 else
535 {
536 constexpr const char* dtNameY = detail::cgns_type_name<Y>();
537 Y* dataArray = nullptr;
538 const cgsize_t memNoStride[3] = { 1, 1, 1 };
539
540 // need to read into temp array to convert data
541 dataArray = new Y[nPts];
542 if (dataArray == nullptr)
543 {
544 std::cerr << "Error allocating buffer array\n";
545 break;
546 }
547 if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameY, cellDim,
548 memDims, memStart, memDims, memNoStride, (void*)dataArray))
549 {
550 delete[] dataArray;
551 char message[81];
552 cgio_error_message(message);
553 std::cerr << "Buffer array cgio_read_data_type :" << message;
554 break;
555 }
556 for (vtkIdType ii = 0; ii < nPts; ++ii)
557 {
558 currentCoord[memStride[0] * ii] = static_cast<T>(dataArray[ii]);
559 }
560 delete[] dataArray;
561 }
562 }
563 return 0;
564}
565VTK_ABI_NAMESPACE_END
566}
567
568#endif // vtkCGNSReaderInternal_h
569// VTK-HeaderTest-Exclude: vtkCGNSReaderInternal.h
std::vector< CGNSRead::FamilyInformation > family
vtkCGNSArraySelection PointDataArraySelection
vtkCGNSArraySelection CellDataArraySelection
std::map< std::string, double > referenceState
std::vector< CGNSRead::ZoneInformation > zones
vtkCGNSArraySelection FaceDataArraySelection
std::vector< CGNSRead::ZoneBCInformation > bcs
void Merge(const vtkCGNSArraySelection &other)
void AddArray(const char *name, bool status=true)
void SetArrayStatus(const char *name, bool status)
~vtkCGNSMetaData()=default
Constructor/Destructor.
int GetNumberOfBaseNodes()
return number of base nodes
vtkCGNSMetaData()=default
Constructor/Destructor.
const CGNSRead::BaseInformation & GetBase(int numBase)
return const reference to a base information
bool Parse(const char *cgnsFileName)
quick parsing of cgns file to get interesting information from a VTK point of view
void Broadcast(vtkMultiProcessController *controller, int rank)
std::vector< double > & GetTimes()
return reference to GlobalTime
void PrintSelf(std::ostream &os)
print object debugging purpose
vtkCGNSReader creates a multi-block dataset and reads unstructured grids and structured meshes from b...
Multiprocessing communication superclass.
represent and manipulate 3D points
Definition vtkPoints.h:139
constexpr const char * cgns_type_name< float >() noexcept
constexpr const char * cgns_type_name< vtkTypeInt64 >() noexcept
constexpr const char * cgns_type_name() noexcept
constexpr const char * cgns_type_name< double >() noexcept
constexpr const char * cgns_type_name< vtkTypeInt32 >() noexcept
This file defines functions used by vtkCGNSReader and vtkCGNSReaderInternal.
void ReorderMonoCellPointsCGNS2VTK(vtkIdType size, int cell_type, vtkIdType numPointsPerCell, vtkIdType *elements)
int get_section_parent_elements(int cgioNum, double cgioSectionId, int dim, const cgsize_t *srcStart, const cgsize_t *srcEnd, const cgsize_t *srcStride, const cgsize_t *memStart, const cgsize_t *memEnd, const cgsize_t *memStride, const cgsize_t *memDim, vtkIdType *localElementsIdx)
void CGNS2VTKorder(vtkIdType size, const int *cells_types, vtkIdType *elements)
int getFirstNodeId(int cgioNum, double parentId, const char *label, double *id, const char *name=nullptr)
Find the first node with the given label.
int get_XYZ_mesh(int cgioNum, const std::vector< double > &gridChildId, const std::size_t &nCoordsArray, int cellDim, vtkIdType nPts, const cgsize_t *srcStart, const cgsize_t *srcEnd, const cgsize_t *srcStride, const cgsize_t *memStart, const cgsize_t *memEnd, const cgsize_t *memStride, const cgsize_t *memDims, vtkPoints *points)
void removeTrailingWhiteSpaces(char_33 name)
bool isACGNSVariable(const std::vector< CGNSVariable > &varList, const char_33 name)
bool ReadGridForZone(vtkCGNSReader *reader, const BaseInformation &baseInfo, const ZoneInformation &zoneInfo)
Helpers to encapsulate all logic to read various nodes (zones, bc patches etc.).
bool ReadBase(vtkCGNSReader *reader, const BaseInformation &baseInfo)
Helpers to encapsulate all logic to read various nodes (zones, bc patches etc.).
const std::map< CGNS_ENUMT(ElementType_t), int > CellDimensions
int get_section_connectivity(int cgioNum, double cgioSectionId, int dim, const cgsize_t *srcStart, const cgsize_t *srcEnd, const cgsize_t *srcStride, const cgsize_t *memStart, const cgsize_t *memEnd, const cgsize_t *memStride, const cgsize_t *memDim, vtkIdType *localElements)
int GetVTKElemType(CGNS_ENUMT(ElementType_t) elemType, bool &higherOrderWarning, bool &cgnsOrderFlag)
void fillVectorsFromVars(std::vector< CGNSRead::CGNSVariable > &vars, std::vector< CGNSRead::CGNSVector > &vectors, int physicalDim)
int get_section_start_offset(int cgioNum, double cgioSectionId, int dim, const cgsize_t *srcStart, const cgsize_t *srcEnd, const cgsize_t *srcStride, const cgsize_t *memStart, const cgsize_t *memEnd, const cgsize_t *memStride, const cgsize_t *memDim, vtkIdType *localElementsIdx)
std::vector< CGNSVector >::iterator getVectorFromName(std::vector< CGNSVector > &vectorList, const char_33 name)
int setUpRind(int cgioNum, double rindId, int *rind)
bool ReadPatchesForBase(vtkCGNSReader *reader, const BaseInformation &)
Helpers to encapsulate all logic to read various nodes (zones, bc patches etc.).
bool ReadPatch(vtkCGNSReader *reader, const BaseInformation &, const ZoneInformation &zoneInfo, const std::string &patchFamilyname)
Helpers to encapsulate all logic to read various nodes (zones, bc patches etc.).
bool compareName(const char_33 nameOne, const char_33 nameTwo)
CGNS_ENUMT(DataType_t) dt
CGNS_ENUMT(DataType_t) dt
int vtkIdType
Definition vtkType.h:315
#define NODE