29#ifndef vtkCGNSReaderInternal_h
30#define vtkCGNSReaderInternal_h
47#include VTK_CGNS(cgnslib.h)
48#include VTK_CGNS(cgns_io.h)
59 static const bool value =
false;
71 static const bool value =
false;
122 vtkCGNSArraySelection::const_iterator iter = other.begin();
123 for (; iter != other.end(); ++iter)
125 (*this)[iter->first] = iter->second;
133 vtkCGNSArraySelection::iterator iter = this->find(
name);
134 if (iter != this->end())
145 vtkCGNSArraySelection::iterator iter = this->find(
name);
146 return (iter != this->end());
156 for (vtkCGNSArraySelection::iterator iter = this->begin(); iter != this->end(); ++iter)
161 return iter->first.c_str();
216 this->name[0] =
'\0';
226 std::vector<CGNSRead::ZoneBCInformation>
bcs;
230 this->name[0] =
'\0';
278 std::vector<CGNSRead::FamilyInformation>
family;
281 std::vector<CGNSRead::ZoneInformation>
zones;
313 bool Parse(
const char* cgnsFileName);
328 std::vector<double>&
GetTimes() {
return this->GlobalTime; }
349 std::vector<CGNSRead::BaseInformation> baseList;
352 std::vector<double> GlobalTime;
359 return (strncmp(nameOne, nameTwo, 32) == 0);
366 char* end =
name + strlen(
name) - 1;
367 while (end >=
name && isspace(*end))
372 assert(end >=
name && end <
name + 33);
379 std::vector<CGNSVector>& vectorList,
const char_33 name)
381 for (std::vector<CGNSVector>::iterator iter = vectorList.begin(); iter != vectorList.end();
384 if (strncmp(iter->name,
name, 31) == 0)
389 return vectorList.end();
395 for (std::vector<CGNSVariable>::const_iterator iter = varList.begin(); iter != varList.end();
398 if (strncmp(iter->name,
name, 32) == 0)
408 std::vector<CGNSRead::CGNSVector>& vectors,
const int physicalDim);
410int setUpRind(
const int cgioNum,
const double rindId,
int* rind);
417 const int cgioNum,
const double parentId,
const char* label,
double*
id,
const char*
name = NULL);
420 const cgsize_t* srcStart,
const cgsize_t* srcEnd,
const cgsize_t* srcStride,
421 const cgsize_t* memStart,
const cgsize_t* memEnd,
const cgsize_t* memStride,
422 const cgsize_t* memDim,
vtkIdType* localElements);
425 const cgsize_t* srcStart,
const cgsize_t* srcEnd,
const cgsize_t* srcStride,
426 const cgsize_t* memStart,
const cgsize_t* memEnd,
const cgsize_t* memStride,
427 const cgsize_t* memDim,
vtkIdType* localElementsIdx);
430 CGNS_ENUMT(ElementType_t) elemType,
bool& higherOrderWarning,
bool& cgnsOrderFlag);
436template <
typename T,
typename Y>
437int get_XYZ_mesh(
const int cgioNum,
const std::vector<double>& gridChildId,
438 const std::size_t& nCoordsArray,
const int cellDim,
const vtkIdType nPts,
439 const cgsize_t* srcStart,
const cgsize_t* srcEnd,
const cgsize_t* srcStride,
440 const cgsize_t* memStart,
const cgsize_t* memEnd,
const cgsize_t* memStride,
443 T* coords =
static_cast<T*
>(
points->GetVoidPointer(0));
444 T* currentCoord =
static_cast<T*
>(&(coords[0]));
448 bool sameType =
true;
451 memset(coords, 0, 3 * nPts *
sizeof(T));
453 for (std::size_t c = 1; c <= nCoordsArray; ++c)
456 if (cgio_get_name(cgioNum, gridChildId[c - 1], coordName) != CG_OK)
459 cgio_error_message(message);
460 std::cerr <<
"get_XYZ_mesh : cgio_get_name :" << message;
465 if (cgio_get_data_type(cgioNum, gridChildId[c - 1], dataType))
470 if (strcmp(dataType,
"R8") == 0)
473 sameType = doubleType;
475 else if (strcmp(dataType,
"R4") == 0)
478 sameType = floatType;
482 std::cerr <<
"Invalid datatype for GridCoordinates\n";
487 len = strlen(coordName) - 1;
488 switch (coordName[len])
491 currentCoord =
static_cast<T*
>(&(coords[0]));
494 currentCoord =
static_cast<T*
>(&(coords[1]));
497 currentCoord =
static_cast<T*
>(&(coords[2]));
501 coordId = gridChildId[c - 1];
504 if (sameType ==
true)
506 constexpr const char* dtNameT = detail::cgns_type_name<T>();
507 if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameT, cellDim,
508 memEnd, memStart, memEnd, memStride, (
void*)currentCoord))
511 cgio_error_message(message);
512 std::cerr <<
"cgio_read_data_type :" << message;
517 constexpr const char* dtNameY = detail::cgns_type_name<Y>();
519 const cgsize_t memNoStride[3] = { 1, 1, 1 };
522 dataArray =
new Y[nPts];
525 std::cerr <<
"Error allocating buffer array\n";
528 if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameY, cellDim,
529 memDims, memStart, memDims, memNoStride, (
void*)dataArray))
533 cgio_error_message(message);
534 std::cerr <<
"Buffer array cgio_read_data_type :" << message;
539 currentCoord[memStride[0] * ii] =
static_cast<T
>(dataArray[ii]);
bool HasArray(const char *name)
bool ArrayIsEnabled(const char *name)
void Merge(const vtkCGNSArraySelection &other)
int GetArraySetting(const char *name)
const char * GetArrayName(int index)
void AddArray(const char *name, bool status=true)
void SetArrayStatus(const char *name, bool status)
vtkCGNSReader creates a multi-block dataset and reads unstructured grids, and structured meshes from ...
Multiprocessing communication superclass.
represent and manipulate 3D points
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
void CGNS2VTKorder(const vtkIdType size, const int *cells_types, vtkIdType *elements)
int get_section_start_offset(const int cgioNum, const double cgioSectionId, const 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)
int get_section_connectivity(const int cgioNum, const double cgioSectionId, const 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)
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.).
int GetVTKElemType(CGNS_ENUMT(ElementType_t) elemType, bool &higherOrderWarning, bool &cgnsOrderFlag)
void fillVectorsFromVars(std::vector< CGNSRead::CGNSVariable > &vars, std::vector< CGNSRead::CGNSVector > &vectors, const int physicalDim)
int get_XYZ_mesh(const int cgioNum, const std::vector< double > &gridChildId, const std::size_t &nCoordsArray, const int cellDim, const 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)
int setUpRind(const int cgioNum, const double rindId, int *rind)
std::vector< CGNSVector >::iterator getVectorFromName(std::vector< CGNSVector > &vectorList, const char_33 name)
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.).
int getFirstNodeId(const int cgioNum, const double parentId, const char *label, double *id, const char *name=NULL)
Find the first node with the given label.
bool compareName(const char_33 nameOne, const char_33 nameTwo)
void CGNS2VTKorderMonoElem(const vtkIdType size, const int cell_type, vtkIdType *elements)
CGNS_ENUMT(DataType_t) dt
CGNS_ENUMT(DataType_t) dt