From 3d48a287141eb911b4888440e09c262743b4db3c Mon Sep 17 00:00:00 2001
From: Mickael PHILIT <mickey.phy@gmail.com>
Date: Wed, 4 Mar 2020 14:54:30 +0100
Subject: [PATCH] Add support for new API cgio_read_***data***_type
CGNS 4.1 removed old cgio API and now support providing memory data type.
Changes are made to keep current behavior of CGNS reading.
Data conversion is not let to CGNS with new API as it may only be stable for HDF5 files.
---
ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.cxx | 14 +--
ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.h | 3 +-
ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReader.cxx | 85 +++++++++++++------
.../CGNSReader/vtkCGNSReaderInternal.cxx | 44 +++++-----
.../CGNSReader/vtkCGNSReaderInternal.h | 53 ++++++++++--
5 files changed, 133 insertions(+), 66 deletions(-)
diff --git ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.cxx ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.cxx
index 1e8ecae8c5..44429a4dd0 100644
|
|
int readNodeStringData(int cgioNum, double nodeId, std::string& data) |
44 | 44 | |
45 | 45 | data.resize(size); |
46 | 46 | // read data |
47 | | if (cgio_read_all_data(cgioNum, nodeId, (void*)data.c_str()) != CG_OK) |
| 47 | if (cgio_read_all_data_type(cgioNum, nodeId, "C1", (void*)data.c_str()) != CG_OK) |
48 | 48 | { |
49 | 49 | return 1; |
50 | 50 | } |
… |
… |
int readNodeData<char>(int cgioNum, double nodeId, std::vector<char>& data) |
80 | 80 | data.resize(size + 1); |
81 | 81 | |
82 | 82 | // read data |
83 | | if (cgio_read_all_data(cgioNum, nodeId, &data[0]) != CG_OK) |
| 83 | if (cgio_read_all_data_type(cgioNum, nodeId, "C1", &data[0]) != CG_OK) |
84 | 84 | { |
85 | 85 | return 1; |
86 | 86 | } |
… |
… |
int readBaseIds(int cgioNum, double rootId, std::vector<double>& baseIds) |
167 | 167 | int readBaseCoreInfo(int cgioNum, double baseId, CGNSRead::BaseInformation& baseInfo) |
168 | 168 | { |
169 | 169 | CGNSRead::char_33 dataType; |
170 | | std::vector<int> mdata; |
| 170 | std::vector<int32_t> mdata; |
171 | 171 | |
172 | 172 | if (cgio_get_name(cgioNum, baseId, baseInfo.name) != CG_OK) |
173 | 173 | { |
… |
… |
int readBaseCoreInfo(int cgioNum, double baseId, CGNSRead::BaseInformation& base |
187 | 187 | return 1; |
188 | 188 | } |
189 | 189 | |
190 | | if (CGNSRead::readNodeData<int>(cgioNum, baseId, mdata) != 0) |
| 190 | if (CGNSRead::readNodeData<int32_t>(cgioNum, baseId, mdata) != 0) |
191 | 191 | { |
192 | 192 | std::cerr << "error while reading base dimension" << std::endl; |
193 | 193 | return 1; |
… |
… |
int readBaseIteration(int cgioNum, double nodeId, CGNSRead::BaseInformation& bas |
209 | 209 | bool createTimeStates = true; |
210 | 210 | bool createIterStates = true; |
211 | 211 | |
212 | | std::vector<int> ndata; |
| 212 | std::vector<int32_t> ndata; |
213 | 213 | // read node data type |
214 | 214 | if (cgio_get_data_type(cgioNum, nodeId, dataType) != CG_OK) |
215 | 215 | { |
… |
… |
int readBaseIteration(int cgioNum, double nodeId, CGNSRead::BaseInformation& bas |
222 | 222 | return 1; |
223 | 223 | } |
224 | 224 | |
225 | | if (CGNSRead::readNodeData<int>(cgioNum, nodeId, ndata) != 0) |
| 225 | if (CGNSRead::readNodeData<int32_t>(cgioNum, nodeId, ndata) != 0) |
226 | 226 | { |
227 | 227 | std::cerr << "error while reading number of state in base" << std::endl; |
228 | 228 | return 1; |
… |
… |
int readBaseIteration(int cgioNum, double nodeId, CGNSRead::BaseInformation& bas |
298 | 298 | } |
299 | 299 | |
300 | 300 | baseInfo.steps.clear(); |
301 | | CGNSRead::readNodeData<int>(cgioNum, childrenIterative[nc], baseInfo.steps); |
| 301 | CGNSRead::readNodeData<int32_t>(cgioNum, childrenIterative[nc], baseInfo.steps); |
302 | 302 | if (static_cast<int>(baseInfo.steps.size()) != nstates) |
303 | 303 | { |
304 | 304 | std::cerr << "Error reading steps node"; |
diff --git ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.h ParaViewCore/VTKExtensions/CGNSReader/cgio_helpers.h
index bef0f3fbfa..07deb54ebe 100644
|
|
inline int readNodeData(int cgioNum, double nodeId, std::vector<T>& data) |
46 | 46 | cgsize_t size = 1; |
47 | 47 | cgsize_t dimVals[12]; |
48 | 48 | int ndim; |
| 49 | constexpr const char* dtName = CGNSRead::detail::cgns_type_name<T>(); |
49 | 50 | |
50 | 51 | if (cgio_get_dimensions(cgioNum, nodeId, &ndim, dimVals) != CG_OK) |
51 | 52 | { |
… |
… |
inline int readNodeData(int cgioNum, double nodeId, std::vector<T>& data) |
65 | 66 | data.resize(size); |
66 | 67 | |
67 | 68 | // read data |
68 | | if (cgio_read_all_data(cgioNum, nodeId, &data[0]) != CG_OK) |
| 69 | if (cgio_read_all_data_type(cgioNum, nodeId, dtName, &data[0]) != CG_OK) |
69 | 70 | { |
70 | 71 | return 1; |
71 | 72 | } |
diff --git ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReader.cxx ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReader.cxx
index a1c417b810..7ca523997a 100644
|
|
int StartsWithFlowSolution(const char* s) |
526 | 526 | |
527 | 527 | return ret; |
528 | 528 | } |
| 529 | //---------------------------------------------------------------------------- |
| 530 | // Small helper |
| 531 | const char* get_data_type(const CGNS_ENUMT(DataType_t) dt) |
| 532 | { |
| 533 | const char* dataType; |
| 534 | switch (dt) |
| 535 | { |
| 536 | case CGNS_ENUMV(Integer): |
| 537 | dataType = "I4"; |
| 538 | break; |
| 539 | case CGNS_ENUMV(LongInteger): |
| 540 | dataType = "I8"; |
| 541 | break; |
| 542 | case CGNS_ENUMV(RealSingle): |
| 543 | dataType = "R4"; |
| 544 | break; |
| 545 | case CGNS_ENUMV(RealDouble): |
| 546 | dataType = "R8"; |
| 547 | break; |
| 548 | case CGNS_ENUMV(Character): |
| 549 | dataType = "C1"; |
| 550 | break; |
| 551 | default: |
| 552 | dataType = "MT"; |
| 553 | } |
| 554 | return dataType; |
| 555 | } |
529 | 556 | |
530 | 557 | //---------------------------------------------------------------------------- |
531 | 558 | vtkCGNSReader::vtkCGNSReader() |
… |
… |
int vtkCGNSReader::vtkPrivate::getGridAndSolutionNames(int base, std::string& gr |
672 | 699 | { |
673 | 700 | CGNSRead::char_33 gname; |
674 | 701 | const cgsize_t offset = static_cast<cgsize_t>(self->ActualTimeStep * 32); |
675 | | cgio_read_block_data(self->cgioNum, giterId, offset + 1, offset + 32, (void*)gname); |
| 702 | cgio_read_block_data_type( |
| 703 | self->cgioNum, giterId, offset + 1, offset + 32, "C1", (void*)gname); |
676 | 704 | gname[32] = '\0'; |
677 | 705 | // NOTE: Names or identifiers contain no spaces and capitalization |
678 | 706 | // is used to distinguish individual words making up a name. |
… |
… |
int vtkCGNSReader::vtkPrivate::getGridAndSolutionNames(int base, std::string& gr |
732 | 760 | EndsWithPointers(nodeName)) |
733 | 761 | { |
734 | 762 | CGNSRead::char_33 gname; |
735 | | cgio_read_block_data(self->cgioNum, iterChildId[cc], |
| 763 | cgio_read_block_data_type(self->cgioNum, iterChildId[cc], |
736 | 764 | (cgsize_t)(self->ActualTimeStep * 32 + 1), (cgsize_t)(self->ActualTimeStep * 32 + 32), |
737 | | (void*)gname); |
| 765 | "C1", (void*)gname); |
738 | 766 | gname[32] = '\0'; |
739 | 767 | CGNSRead::removeTrailingWhiteSpaces(gname); |
740 | 768 | std::string tmpStr = std::string(gname); |
… |
… |
int vtkCGNSReader::vtkPrivate::readSolution(const std::string& solutionNameStr, |
1197 | 1225 | continue; |
1198 | 1226 | } |
1199 | 1227 | double cgioVarId = solChildId[ff]; |
| 1228 | const char* fieldDataType = get_data_type(cgnsVars[ff].dt); |
1200 | 1229 | |
1201 | 1230 | // quick transfer of data because data types is given by cgns database |
1202 | 1231 | if (cgnsVars[ff].isComponent == false) |
1203 | 1232 | { |
1204 | | if (cgio_read_data(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, fieldSrcStride, |
1205 | | cellDim, fieldMemDims, fieldMemStart, fieldMemEnd, fieldMemStride, |
| 1233 | if (cgio_read_data_type(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, fieldSrcStride, |
| 1234 | fieldDataType, cellDim, fieldMemDims, fieldMemStart, fieldMemEnd, fieldMemStride, |
1206 | 1235 | (void*)vtkVars[ff]->GetVoidPointer(0)) != CG_OK) |
1207 | 1236 | { |
1208 | 1237 | char message[81]; |
1209 | 1238 | cgio_error_message(message); |
1210 | | vtkGenericWarningMacro(<< "cgio_read_data :" << message); |
| 1239 | vtkGenericWarningMacro(<< "cgio_read_data_type :" << message); |
1211 | 1240 | } |
1212 | 1241 | } |
1213 | 1242 | else |
1214 | 1243 | { |
1215 | | if (cgio_read_data(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, fieldSrcStride, |
1216 | | cellDim, fieldVectMemDims, fieldVectMemStart, fieldVectMemEnd, fieldVectMemStride, |
| 1244 | if (cgio_read_data_type(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, fieldSrcStride, |
| 1245 | fieldDataType, cellDim, fieldVectMemDims, fieldVectMemStart, fieldVectMemEnd, |
| 1246 | fieldVectMemStride, |
1217 | 1247 | (void*)vtkVars[ff]->GetVoidPointer(cgnsVars[ff].xyzIndex - 1)) != CG_OK) |
1218 | 1248 | { |
1219 | 1249 | char message[81]; |
1220 | 1250 | cgio_error_message(message); |
1221 | | vtkGenericWarningMacro(<< "cgio_read_data :" << message); |
| 1251 | vtkGenericWarningMacro(<< "cgio_read_data_type :" << message); |
1222 | 1252 | } |
1223 | 1253 | } |
1224 | 1254 | cgio_release_id(self->cgioNum, cgioVarId); |
… |
… |
int vtkCGNSReader::vtkPrivate::readBCData(const double nodeId, const int cellDim |
1448 | 1478 | continue; |
1449 | 1479 | } |
1450 | 1480 | double cgioVarId = varIds[ff]; |
| 1481 | const char* fieldDataType = get_data_type(cgnsVars[ff].dt); |
1451 | 1482 | |
1452 | 1483 | cgsize_t dataSize = 1; |
1453 | 1484 | cgsize_t dimVals[12]; |
… |
… |
int vtkCGNSReader::vtkPrivate::readBCData(const double nodeId, const int cellDim |
1474 | 1505 | // quick transfer of data because data types is given by cgns database |
1475 | 1506 | if (cgnsVars[ff].isComponent == false) |
1476 | 1507 | { |
1477 | | if (cgio_read_all_data( |
1478 | | self->cgioNum, cgioVarId, (void*)vtkVars[ff]->GetVoidPointer(0)) != CG_OK) |
| 1508 | if (cgio_read_all_data_type(self->cgioNum, cgioVarId, fieldDataType, |
| 1509 | (void*)vtkVars[ff]->GetVoidPointer(0)) != CG_OK) |
1479 | 1510 | { |
1480 | 1511 | char message[81]; |
1481 | 1512 | cgio_error_message(message); |
1482 | | vtkGenericWarningMacro(<< "cgio_read_data :" << message); |
| 1513 | vtkGenericWarningMacro(<< "cgio_read_all_data_type :" << message); |
1483 | 1514 | } |
1484 | 1515 | if (dataSize == 1) |
1485 | 1516 | { |
… |
… |
int vtkCGNSReader::vtkPrivate::readBCData(const double nodeId, const int cellDim |
1510 | 1541 | fieldVectMemDims[0] = fieldSrcEnd[0] * fieldVectMemStride[0]; |
1511 | 1542 | fieldVectMemEnd[0] = fieldSrcEnd[0] * fieldVectMemStride[0]; |
1512 | 1543 | |
1513 | | if (cgio_read_data(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, |
1514 | | fieldSrcStride, 1, fieldVectMemDims, fieldVectMemStart, fieldVectMemEnd, |
1515 | | fieldVectMemStride, |
| 1544 | if (cgio_read_data_type(self->cgioNum, cgioVarId, fieldSrcStart, fieldSrcEnd, |
| 1545 | fieldSrcStride, fieldDataType, 1, fieldVectMemDims, fieldVectMemStart, |
| 1546 | fieldVectMemEnd, fieldVectMemStride, |
1516 | 1547 | (void*)vtkVars[ff]->GetVoidPointer(cgnsVars[ff].xyzIndex - 1)) != CG_OK) |
1517 | 1548 | { |
1518 | 1549 | char message[81]; |
1519 | 1550 | cgio_error_message(message); |
1520 | | vtkGenericWarningMacro(<< "cgio_read_data :" << message); |
| 1551 | vtkGenericWarningMacro(<< "cgio_read_data_type :" << message); |
1521 | 1552 | } |
1522 | 1553 | if (dataSize == 1) |
1523 | 1554 | { |
… |
… |
int vtkCGNSReader::GetUnstructuredZone( |
2231 | 2262 | |
2232 | 2263 | // |
2233 | 2264 | CGNSRead::char_33 dataType; |
2234 | | std::vector<int> mdata; |
| 2265 | std::vector<vtkTypeInt32> mdata; |
2235 | 2266 | |
2236 | 2267 | if (cgio_get_name(this->cgioNum, elemIdList[sec], sectionInfoList[sec].name) != CG_OK) |
2237 | 2268 | { |
… |
… |
int vtkCGNSReader::GetUnstructuredZone( |
2246 | 2277 | vtkErrorMacro(<< "Unexpected data type for dimension data of Element\n"); |
2247 | 2278 | } |
2248 | 2279 | |
2249 | | CGNSRead::readNodeData<int>(cgioNum, elemIdList[sec], mdata); |
| 2280 | CGNSRead::readNodeData<vtkTypeInt32>(cgioNum, elemIdList[sec], mdata); |
2250 | 2281 | if (mdata.size() != 2) |
2251 | 2282 | { |
2252 | 2283 | vtkErrorMacro(<< "Unexpected data for Elements_t node\n"); |
… |
… |
int vtkCGNSReader::GetUnstructuredZone( |
2267 | 2298 | |
2268 | 2299 | if (strcmp(dataType, "I4") == 0) |
2269 | 2300 | { |
2270 | | std::vector<int> mdata2; |
2271 | | CGNSRead::readNodeData<int>(this->cgioNum, elemRangeId, mdata2); |
| 2301 | std::vector<vtkTypeInt32> mdata2; |
| 2302 | CGNSRead::readNodeData<vtkTypeInt32>(this->cgioNum, elemRangeId, mdata2); |
2272 | 2303 | if (mdata2.size() != 2) |
2273 | 2304 | { |
2274 | 2305 | vtkErrorMacro(<< "Unexpected data for ElementRange node\n"); |
… |
… |
int vtkCGNSReader::GetUnstructuredZone( |
2278 | 2309 | } |
2279 | 2310 | else if (strcmp(dataType, "I8") == 0) |
2280 | 2311 | { |
2281 | | std::vector<cglong_t> mdata2; |
2282 | | CGNSRead::readNodeData<cglong_t>(this->cgioNum, elemRangeId, mdata2); |
| 2312 | std::vector<vtkTypeInt64> mdata2; |
| 2313 | CGNSRead::readNodeData<vtkTypeInt64>(this->cgioNum, elemRangeId, mdata2); |
2283 | 2314 | if (mdata2.size() != 2) |
2284 | 2315 | { |
2285 | 2316 | vtkErrorMacro(<< "Unexpected data for ElementRange node\n"); |
… |
… |
int vtkCGNSReader::RequestData(vtkInformation* vtkNotUsed(request), |
4437 | 4468 | |
4438 | 4469 | if (strcmp(dataType, "I4") == 0) |
4439 | 4470 | { |
4440 | | std::vector<int> mdata; |
4441 | | CGNSRead::readNodeData<int>(this->cgioNum, baseChildId[zone], mdata); |
| 4471 | std::vector<vtkTypeInt32> mdata; |
| 4472 | CGNSRead::readNodeData<vtkTypeInt32>(this->cgioNum, baseChildId[zone], mdata); |
4442 | 4473 | for (std::size_t index = 0; index < mdata.size(); index++) |
4443 | 4474 | { |
4444 | 4475 | zsize[index] = static_cast<cgsize_t>(mdata[index]); |
… |
… |
int vtkCGNSReader::RequestData(vtkInformation* vtkNotUsed(request), |
4446 | 4477 | } |
4447 | 4478 | else if (strcmp(dataType, "I8") == 0) |
4448 | 4479 | { |
4449 | | std::vector<cglong_t> mdata; |
4450 | | CGNSRead::readNodeData<cglong_t>(this->cgioNum, baseChildId[zone], mdata); |
| 4480 | std::vector<vtkTypeInt64> mdata; |
| 4481 | CGNSRead::readNodeData<vtkTypeInt64>(this->cgioNum, baseChildId[zone], mdata); |
4451 | 4482 | for (std::size_t index = 0; index < mdata.size(); index++) |
4452 | 4483 | { |
4453 | 4484 | zsize[index] = static_cast<cgsize_t>(mdata[index]); |
… |
… |
int vtkCGNSReader::CanReadFile(const char* name) |
4724 | 4755 | } |
4725 | 4756 | |
4726 | 4757 | // read data |
4727 | | if (cgio_read_all_data(cgioFile, childId, &FileVersion)) |
| 4758 | if (cgio_read_all_data_type(cgioFile, childId, "R4", &FileVersion)) |
4728 | 4759 | { |
4729 | 4760 | vtkErrorMacro(<< "read CGNS version number"); |
4730 | 4761 | ierr = 0; |
diff --git ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReaderInternal.cxx ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReaderInternal.cxx
index 92c0d6ac51..a2bcf1a443 100644
|
|
int setUpRind(const int cgioNum, const double rindId, int* rind) |
36 | 36 | |
37 | 37 | if (strcmp(dataType, "I4") == 0) |
38 | 38 | { |
39 | | std::vector<int> mdata; |
40 | | CGNSRead::readNodeData<int>(cgioNum, rindId, mdata); |
| 39 | std::vector<vtkTypeInt32> mdata; |
| 40 | CGNSRead::readNodeData<vtkTypeInt32>(cgioNum, rindId, mdata); |
41 | 41 | for (std::size_t index = 0; index < mdata.size(); index++) |
42 | 42 | { |
43 | 43 | rind[index] = static_cast<int>(mdata[index]); |
… |
… |
int setUpRind(const int cgioNum, const double rindId, int* rind) |
45 | 45 | } |
46 | 46 | else if (strcmp(dataType, "I8") == 0) |
47 | 47 | { |
48 | | std::vector<cglong_t> mdata; |
49 | | CGNSRead::readNodeData<cglong_t>(cgioNum, rindId, mdata); |
| 48 | std::vector<vtkTypeInt64> mdata; |
| 49 | CGNSRead::readNodeData<vtkTypeInt64>(cgioNum, rindId, mdata); |
50 | 50 | for (std::size_t index = 0; index < mdata.size(); index++) |
51 | 51 | { |
52 | 52 | rind[index] = static_cast<int>(mdata[index]); |
… |
… |
int get_section_connectivity(const int cgioNum, const double cgioSectionId, cons |
156 | 156 | |
157 | 157 | if (sizeOfCnt == sizeof(vtkIdType)) |
158 | 158 | { |
159 | | if (cgio_read_data(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, dim, memDim, |
160 | | memStart, memEnd, memStride, (void*)localElements) != CG_OK) |
| 159 | if (cgio_read_data_type(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, dataType, dim, |
| 160 | memDim, memStart, memEnd, memStride, (void*)localElements) != CG_OK) |
161 | 161 | { |
162 | 162 | char message[81]; |
163 | 163 | cgio_error_message(message); |
164 | | std::cerr << "cgio_read_data :" << message; |
| 164 | std::cerr << "cgio_read_data_type :" << message; |
165 | 165 | return 1; |
166 | 166 | } |
167 | 167 | } |
… |
… |
int get_section_connectivity(const int cgioNum, const double cgioSectionId, cons |
181 | 181 | std::cerr << "Allocation failed for temporary connectivity array\n"; |
182 | 182 | } |
183 | 183 | |
184 | | if (cgio_read_data(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, dim, memDim, |
185 | | memStart, memEnd, memStride, (void*)data) != CG_OK) |
| 184 | if (cgio_read_data_type(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, "I4", dim, |
| 185 | memDim, memStart, memEnd, memStride, (void*)data) != CG_OK) |
186 | 186 | { |
187 | 187 | delete[] data; |
188 | 188 | char message[81]; |
189 | 189 | cgio_error_message(message); |
190 | | std::cerr << "cgio_read_data :" << message; |
| 190 | std::cerr << "cgio_read_data_type :" << message; |
191 | 191 | return 1; |
192 | 192 | } |
193 | 193 | for (cgsize_t n = 0; n < nn; n++) |
… |
… |
int get_section_connectivity(const int cgioNum, const double cgioSectionId, cons |
204 | 204 | std::cerr << "Allocation failed for temporary connectivity array\n"; |
205 | 205 | return 1; |
206 | 206 | } |
207 | | if (cgio_read_data(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, dim, memDim, |
208 | | memStart, memEnd, memStride, (void*)data) != CG_OK) |
| 207 | if (cgio_read_data_type(cgioNum, cgioElemConnectId, srcStart, srcEnd, srcStride, "I8", dim, |
| 208 | memDim, memStart, memEnd, memStride, (void*)data) != CG_OK) |
209 | 209 | { |
210 | 210 | delete[] data; |
211 | 211 | char message[81]; |
212 | 212 | cgio_error_message(message); |
213 | | std::cerr << "cgio_read_data :" << message; |
| 213 | std::cerr << "cgio_read_data_type :" << message; |
214 | 214 | return 1; |
215 | 215 | } |
216 | 216 | for (cgsize_t n = 0; n < nn; n++) |
… |
… |
int get_section_start_offset(const int cgioNum, const double cgioSectionId, cons |
258 | 258 | |
259 | 259 | if (sizeOfCnt == sizeof(vtkIdType)) |
260 | 260 | { |
261 | | if (cgio_read_data(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, dim, memDim, |
262 | | memStart, memEnd, memStride, (void*)localElementsIdx) != CG_OK) |
| 261 | if (cgio_read_data_type(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, dataType, dim, |
| 262 | memDim, memStart, memEnd, memStride, (void*)localElementsIdx) != CG_OK) |
263 | 263 | { |
264 | 264 | char message[81]; |
265 | 265 | cgio_error_message(message); |
266 | | std::cerr << "cgio_read_data :" << message; |
| 266 | std::cerr << "cgio_read_data_type :" << message; |
267 | 267 | return 1; |
268 | 268 | } |
269 | 269 | } |
… |
… |
int get_section_start_offset(const int cgioNum, const double cgioSectionId, cons |
283 | 283 | std::cerr << "Allocation failed for temporary connectivity offset array\n"; |
284 | 284 | } |
285 | 285 | |
286 | | if (cgio_read_data(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, dim, memDim, |
287 | | memStart, memEnd, memStride, (void*)data) != CG_OK) |
| 286 | if (cgio_read_data_type(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, "I4", dim, |
| 287 | memDim, memStart, memEnd, memStride, (void*)data) != CG_OK) |
288 | 288 | { |
289 | 289 | delete[] data; |
290 | 290 | char message[81]; |
291 | 291 | cgio_error_message(message); |
292 | | std::cerr << "cgio_read_data :" << message; |
| 292 | std::cerr << "cgio_read_data_type :" << message; |
293 | 293 | return 1; |
294 | 294 | } |
295 | 295 | for (cgsize_t n = 0; n < nn; n++) |
… |
… |
int get_section_start_offset(const int cgioNum, const double cgioSectionId, cons |
306 | 306 | std::cerr << "Allocation failed for temporary connectivity array\n"; |
307 | 307 | return 1; |
308 | 308 | } |
309 | | if (cgio_read_data(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, dim, memDim, |
310 | | memStart, memEnd, memStride, (void*)data) != CG_OK) |
| 309 | if (cgio_read_data_type(cgioNum, cgioElemOffsetId, srcStart, srcEnd, srcStride, "I8", dim, |
| 310 | memDim, memStart, memEnd, memStride, (void*)data) != CG_OK) |
311 | 311 | { |
312 | 312 | delete[] data; |
313 | 313 | char message[81]; |
314 | 314 | cgio_error_message(message); |
315 | | std::cerr << "cgio_read_data :" << message; |
| 315 | std::cerr << "cgio_read_data_type :" << message; |
316 | 316 | return 1; |
317 | 317 | } |
318 | 318 | for (cgsize_t n = 0; n < nn; n++) |
diff --git ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReaderInternal.h ParaViewCore/VTKExtensions/CGNSReader/vtkCGNSReaderInternal.h
index 4bcfd5c75a..3df3bae8de 100644
|
|
struct is_float<float> |
72 | 72 | }; |
73 | 73 | } |
74 | 74 | |
| 75 | namespace detail |
| 76 | { |
| 77 | template <typename T> |
| 78 | constexpr const char* cgns_type_name() noexcept |
| 79 | { |
| 80 | return "MT"; |
| 81 | } |
| 82 | |
| 83 | template <> |
| 84 | constexpr const char* cgns_type_name<float>() noexcept |
| 85 | { |
| 86 | return "R4"; |
| 87 | } |
| 88 | |
| 89 | template <> |
| 90 | constexpr const char* cgns_type_name<double>() noexcept |
| 91 | { |
| 92 | return "R8"; |
| 93 | } |
| 94 | |
| 95 | template <> |
| 96 | constexpr const char* cgns_type_name<vtkTypeInt32>() noexcept |
| 97 | { |
| 98 | return "I4"; |
| 99 | } |
| 100 | |
| 101 | template <> |
| 102 | constexpr const char* cgns_type_name<vtkTypeInt64>() noexcept |
| 103 | { |
| 104 | return "I8"; |
| 105 | } |
| 106 | } |
| 107 | |
75 | 108 | typedef char char_33[33]; |
76 | 109 | |
77 | 110 | //------------------------------------------------------------------------------ |
… |
… |
class BaseInformation |
206 | 239 | public: |
207 | 240 | char_33 name; |
208 | 241 | |
209 | | int cellDim; |
210 | | int physicalDim; |
| 242 | int32_t cellDim; |
| 243 | int32_t physicalDim; |
211 | 244 | // |
212 | 245 | int baseNumber; |
213 | 246 | |
214 | | std::vector<int> steps; |
| 247 | std::vector<int32_t> steps; |
215 | 248 | std::vector<double> times; |
216 | 249 | |
217 | 250 | // For unsteady meshes : |
… |
… |
int get_XYZ_mesh(const int cgioNum, const std::vector<double>& gridChildId, |
469 | 502 | // quick transfer of data if same data types |
470 | 503 | if (sameType == true) |
471 | 504 | { |
472 | | if (cgio_read_data(cgioNum, coordId, srcStart, srcEnd, srcStride, cellDim, memEnd, memStart, |
473 | | memEnd, memStride, (void*)currentCoord)) |
| 505 | constexpr const char* dtNameT = detail::cgns_type_name<T>(); |
| 506 | if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameT, cellDim, |
| 507 | memEnd, memStart, memEnd, memStride, (void*)currentCoord)) |
474 | 508 | { |
475 | 509 | char message[81]; |
476 | 510 | cgio_error_message(message); |
477 | | std::cerr << "cgio_read_data :" << message; |
| 511 | std::cerr << "cgio_read_data_type :" << message; |
478 | 512 | } |
479 | 513 | } |
480 | 514 | else |
481 | 515 | { |
| 516 | constexpr const char* dtNameY = detail::cgns_type_name<Y>(); |
482 | 517 | Y* dataArray = 0; |
483 | 518 | const cgsize_t memNoStride[3] = { 1, 1, 1 }; |
484 | 519 | |
… |
… |
int get_XYZ_mesh(const int cgioNum, const std::vector<double>& gridChildId, |
489 | 524 | std::cerr << "Error allocating buffer array\n"; |
490 | 525 | break; |
491 | 526 | } |
492 | | if (cgio_read_data(cgioNum, coordId, srcStart, srcEnd, srcStride, cellDim, memDims, memStart, |
493 | | memDims, memNoStride, (void*)dataArray)) |
| 527 | if (cgio_read_data_type(cgioNum, coordId, srcStart, srcEnd, srcStride, dtNameY, cellDim, |
| 528 | memDims, memStart, memDims, memNoStride, (void*)dataArray)) |
494 | 529 | { |
495 | 530 | delete[] dataArray; |
496 | 531 | char message[81]; |
497 | 532 | cgio_error_message(message); |
498 | | std::cerr << "Buffer array cgio_read_data :" << message; |
| 533 | std::cerr << "Buffer array cgio_read_data_type :" << message; |
499 | 534 | break; |
500 | 535 | } |
501 | 536 | for (vtkIdType ii = 0; ii < nPts; ++ii) |