Actual source code: cgnsv.c

  1: #include <petsc/private/viewercgnsimpl.h>
  2: #include <petsc/private/dmpleximpl.h>
  3: #if defined(PETSC_HDF5_HAVE_PARALLEL)
  4:   #include <pcgnslib.h>
  5: #else
  6:   #include <cgnslib.h>
  7: #endif
  8: #include <cgns_io.h>
  9: #include <ctype.h>

 11: PetscLogEvent PETSC_VIEWER_CGNS_Open, PETSC_VIEWER_CGNS_Close, PETSC_VIEWER_CGNS_ReadMeta, PETSC_VIEWER_CGNS_WriteMeta, PETSC_VIEWER_CGNS_ReadData, PETSC_VIEWER_CGNS_WriteData;

 13: PetscErrorCode PetscViewerCGNSRegisterLogEvents_Internal()
 14: {
 15:   static PetscBool is_initialized = PETSC_FALSE;

 17:   PetscFunctionBeginUser;
 18:   if (is_initialized) PetscFunctionReturn(PETSC_SUCCESS);
 19:   PetscCall(PetscLogEventRegister("CGNSOpen", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_Open));
 20:   PetscCall(PetscLogEventRegister("CGNSClose", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_Close));
 21:   PetscCall(PetscLogEventRegister("CGNSReadMeta", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_ReadMeta));
 22:   PetscCall(PetscLogEventRegister("CGNSReadData", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_ReadData));
 23:   PetscCall(PetscLogEventRegister("CGNSWriteMeta", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_WriteMeta));
 24:   PetscCall(PetscLogEventRegister("CGNSWriteData", PETSC_VIEWER_CLASSID, &PETSC_VIEWER_CGNS_WriteData));
 25:   is_initialized = PETSC_TRUE;
 26:   PetscFunctionReturn(PETSC_SUCCESS);
 27: }

 29: static PetscErrorCode PetscViewerSetFromOptions_CGNS(PetscViewer v, PetscOptionItems PetscOptionsObject)
 30: {
 31:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)v->data;

 33:   PetscFunctionBegin;
 34:   PetscOptionsHeadBegin(PetscOptionsObject, "CGNS Viewer Options");
 35:   PetscCall(PetscOptionsInt("-viewer_cgns_batch_size", "Max number of steps to store in single file when using a template cgns:name-\%d.cgns", "", cgv->batch_size, &cgv->batch_size, NULL));
 36:   PetscOptionsHeadEnd();
 37:   PetscFunctionReturn(PETSC_SUCCESS);
 38: }

 40: static PetscErrorCode PetscViewerView_CGNS(PetscViewer v, PetscViewer viewer)
 41: {
 42:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)v->data;

 44:   PetscFunctionBegin;
 45:   if (cgv->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", cgv->filename));
 46:   PetscFunctionReturn(PETSC_SUCCESS);
 47: }

 49: static PetscErrorCode PetscViewerFileClose_CGNS(PetscViewer viewer)
 50: {
 51:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

 53:   PetscFunctionBegin;
 54:   if (cgv->output_times) {
 55:     PetscCount size, width = 32, *steps;
 56:     char      *solnames;
 57:     PetscReal *times;
 58:     cgsize_t   num_times;
 59:     PetscCall(PetscSegBufferGetSize(cgv->output_times, &size));
 60:     PetscCall(PetscSegBufferExtractInPlace(cgv->output_times, &times));
 61:     PetscCall(PetscSegBufferExtractInPlace(cgv->output_steps, &steps));
 62:     num_times = size;
 63:     PetscCallCGNSWrite(cg_biter_write(cgv->file_num, cgv->base, "TimeIterValues", num_times), viewer, 0);
 64:     PetscCallCGNS(cg_goto(cgv->file_num, cgv->base, "BaseIterativeData_t", 1, NULL));
 65:     PetscCallCGNSWrite(cg_array_write("TimeValues", CGNS_ENUMV(RealDouble), 1, &num_times, times), viewer, 0);
 66:     { // Cast output_steps to long for writing into file
 67:       int *steps_int;
 68:       PetscCall(PetscMalloc1(size, &steps_int));
 69:       for (PetscCount i = 0; i < size; i++) PetscCall(PetscCIntCast(steps[i], &steps_int[i]));
 70:       PetscCallCGNSWrite(cg_array_write("IterationValues", CGNS_ENUMV(Integer), 1, &num_times, steps_int), viewer, 0);
 71:       PetscCall(PetscFree(steps_int));
 72:     }
 73:     PetscCall(PetscSegBufferDestroy(&cgv->output_times));
 74:     PetscCallCGNSWrite(cg_ziter_write(cgv->file_num, cgv->base, cgv->zone, "ZoneIterativeData"), viewer, 0);
 75:     PetscCallCGNS(cg_goto(cgv->file_num, cgv->base, "Zone_t", cgv->zone, "ZoneIterativeData_t", 1, NULL));
 76:     PetscCall(PetscMalloc(size * width + 1, &solnames));
 77:     for (PetscCount i = 0; i < size; i++) PetscCall(PetscSNPrintf(&solnames[i * width], width + 1, "FlowSolution%-20zu", (size_t)steps[i]));
 78:     PetscCall(PetscSegBufferDestroy(&cgv->output_steps));
 79:     cgsize_t shape[2] = {(cgsize_t)width, (cgsize_t)size};
 80:     PetscCallCGNSWrite(cg_array_write("FlowSolutionPointers", CGNS_ENUMV(Character), 2, shape, solnames), viewer, 0);
 81:     // The VTK reader looks for names like FlowSolution*Pointers.
 82:     for (PetscCount i = 0; i < size; i++) PetscCall(PetscSNPrintf(&solnames[i * width], width + 1, "%-32s", "CellInfo"));
 83:     PetscCallCGNSWrite(cg_array_write("FlowSolutionCellInfoPointers", CGNS_ENUMV(Character), 2, shape, solnames), viewer, 0);
 84:     PetscCall(PetscFree(solnames));

 86:     PetscCallCGNSWrite(cg_simulation_type_write(cgv->file_num, cgv->base, CGNS_ENUMV(TimeAccurate)), viewer, 0);
 87:   }
 88:   PetscCall(PetscFree(cgv->filename));
 89: #if defined(PETSC_HDF5_HAVE_PARALLEL)
 90:   if (cgv->file_num) PetscCallCGNSClose(cgp_close(cgv->file_num), viewer, 0);
 91: #else
 92:   if (cgv->file_num) PetscCallCGNSClose(cg_close(cgv->file_num), viewer, 0);
 93: #endif
 94:   cgv->file_num = 0;
 95:   PetscCall(PetscFree(cgv->node_l2g));
 96:   PetscCall(PetscFree(cgv->nodal_field));
 97:   PetscFunctionReturn(PETSC_SUCCESS);
 98: }

100: PetscErrorCode PetscViewerCGNSFileOpen_Internal(PetscViewer viewer, PetscInt sequence_number)
101: {
102:   PetscViewer_CGNS *cgv          = (PetscViewer_CGNS *)viewer->data;
103:   int               cg_file_mode = -1;

105:   PetscFunctionBegin;
106:   PetscCheck((cgv->filename == NULL) ^ (sequence_number < 0), PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Expect either a template filename or non-negative sequence number");
107:   if (!cgv->filename) {
108:     char filename_numbered[PETSC_MAX_PATH_LEN];
109:     // Cast sequence_number so %d can be used also when PetscInt is 64-bit. We could upgrade the format string if users
110:     // run more than 2B time steps.
111:     PetscCall(PetscSNPrintf(filename_numbered, sizeof filename_numbered, cgv->filename_template, (int)sequence_number));
112:     PetscCall(PetscStrallocpy(filename_numbered, &cgv->filename));
113:   }
114:   switch (cgv->btype) {
115:   case FILE_MODE_READ:
116:     cg_file_mode = CG_MODE_READ;
117:     break;
118:   case FILE_MODE_WRITE:
119:     cg_file_mode = CG_MODE_WRITE;
120:     break;
121:   case FILE_MODE_UNDEFINED:
122:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
123:   default:
124:     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[cgv->btype]);
125:   }
126: #if defined(PETSC_HDF5_HAVE_PARALLEL)
127:   PetscCallCGNS(cgp_mpi_comm(PetscObjectComm((PetscObject)viewer)));
128:   PetscCallCGNSOpen(cgp_open(cgv->filename, cg_file_mode, &cgv->file_num), viewer, 0);
129: #else
130:   PetscCallCGNSOpen(cg_open(filename, cg_file_mode, &cgv->file_num), viewer, 0);
131: #endif
132:   PetscFunctionReturn(PETSC_SUCCESS);
133: }

135: PetscErrorCode PetscViewerCGNSCheckBatch_Internal(PetscViewer viewer)
136: {
137:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
138:   PetscCount        num_steps;

140:   PetscFunctionBegin;
141:   if (!cgv->filename_template) PetscFunctionReturn(PETSC_SUCCESS); // Batches are closed when viewer is destroyed
142:   PetscCall(PetscSegBufferGetSize(cgv->output_times, &num_steps));
143:   if (num_steps >= (PetscCount)cgv->batch_size) PetscCall(PetscViewerFileClose_CGNS(viewer));
144:   PetscFunctionReturn(PETSC_SUCCESS);
145: }

147: static PetscErrorCode PetscViewerDestroy_CGNS(PetscViewer viewer)
148: {
149:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

151:   PetscFunctionBegin;
152:   PetscCall(PetscViewerFileClose_CGNS(viewer));
153:   PetscCall(PetscFree(cgv->solution_name));
154:   PetscCall(PetscFree(cgv->filename_template));
155:   PetscCall(PetscFree(cgv));
156:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
157:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
158:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
159:   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

163: static PetscErrorCode PetscViewerFileSetMode_CGNS(PetscViewer viewer, PetscFileMode type)
164: {
165:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

167:   PetscFunctionBegin;
168:   cgv->btype = type;
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: static PetscErrorCode PetscViewerFileGetMode_CGNS(PetscViewer viewer, PetscFileMode *type)
173: {
174:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

176:   PetscFunctionBegin;
177:   *type = cgv->btype;
178:   PetscFunctionReturn(PETSC_SUCCESS);
179: }

181: static PetscErrorCode PetscViewerFileSetName_CGNS(PetscViewer viewer, const char *filename)
182: {
183:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
184:   char             *has_pattern;

186:   PetscFunctionBegin;
187: #if defined(PETSC_HDF5_HAVE_PARALLEL)
188:   if (cgv->file_num) PetscCallCGNSClose(cgp_close(cgv->file_num), viewer, 0);
189: #else
190:   if (cgv->file_num) PetscCallCGNSClose(cg_close(cgv->file_num), viewer, 0);
191: #endif
192:   PetscCall(PetscFree(cgv->filename));
193:   PetscCall(PetscFree(cgv->filename_template));
194:   PetscCall(PetscStrchr(filename, '%', &has_pattern));
195:   if (has_pattern) {
196:     PetscCall(PetscStrallocpy(filename, &cgv->filename_template));
197:     // Templated file names open lazily, once we know DMGetOutputSequenceNumber()
198:   } else {
199:     PetscCall(PetscStrallocpy(filename, &cgv->filename));
200:     PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, -1));
201:   }
202:   PetscFunctionReturn(PETSC_SUCCESS);
203: }

205: static PetscErrorCode PetscViewerFileGetName_CGNS(PetscViewer viewer, const char **filename)
206: {
207:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

209:   PetscFunctionBegin;
210:   *filename = cgv->filename;
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*MC
215:    PETSCVIEWERCGNS - A viewer for CGNS files

217:   Level: beginner

219:   Options Database Key:
220: . -viewer_cgns_batch_size SIZE - set max number of output sequence times to write per batch

222:   Note:
223:   If the filename contains an integer format character, the CGNS viewer will created a batched output sequence. For
224:   example, one could use `-ts_monitor_solution cgns:flow-%d.cgns`. This is desirable if one wants to limit file sizes or
225:   if the job might crash/be killed by a resource manager before exiting cleanly.

227: .seealso: [](sec_viewers), `PetscViewer`, `PetscViewerCreate()`, `VecView()`, `DMView()`, `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `TSSetFromOptions()`
228: M*/
229: PetscErrorCode PetscViewerCreate_CGNS(PetscViewer v)
230: {
231:   PetscViewer_CGNS *cgv;

233:   PetscFunctionBegin;
234:   PetscCall(PetscViewerCGNSRegisterLogEvents_Internal());
235:   PetscCall(PetscNew(&cgv));

237:   v->data                   = cgv;
238:   v->ops->destroy           = PetscViewerDestroy_CGNS;
239:   v->ops->setfromoptions    = PetscViewerSetFromOptions_CGNS;
240:   v->ops->view              = PetscViewerView_CGNS;
241:   cgv->btype                = FILE_MODE_UNDEFINED;
242:   cgv->filename             = NULL;
243:   cgv->batch_size           = 20;
244:   cgv->solution_index       = -1; // Default to use the "last" solution
245:   cgv->base                 = 1;
246:   cgv->zone                 = 1;
247:   cgv->previous_output_step = -1; // No previous output step

249:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_CGNS));
250:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_CGNS));
251:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_CGNS));
252:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_CGNS));
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }

256: // Find DataArray_t node under the current node (determined by `cg_goto` and friends) that matches `name`
257: // Return the index of that array and (optionally) other data about the array
258: static inline PetscErrorCode CGNS_Find_Array(MPI_Comm comm, const char name[], int *A_index, CGNS_ENUMT(DataType_t) * data_type, int *dim, cgsize_t size[], PetscBool *found)
259: {
260:   int  narrays; // number of arrays under the current node
261:   char array_name[CGIO_MAX_NAME_LENGTH + 1];
262:   CGNS_ENUMT(DataType_t) data_type_local;
263:   int       _dim;
264:   cgsize_t  _size[12];
265:   PetscBool matches_name = PETSC_FALSE;

267:   PetscFunctionBeginUser;
268:   PetscCallCGNSRead(cg_narrays(&narrays), 0, 0);
269:   for (int i = 1; i <= narrays; i++) {
270:     PetscCallCGNSRead(cg_array_info(i, array_name, &data_type_local, &_dim, _size), 0, 0);
271:     PetscCall(PetscStrcmp(name, array_name, &matches_name));
272:     if (matches_name) {
273:       *A_index = i;
274:       if (data_type) *data_type = data_type_local;
275:       if (dim) *dim = _dim;
276:       if (size) PetscArraycpy(size, _size, _dim);
277:       if (found) *found = PETSC_TRUE;
278:       PetscFunctionReturn(PETSC_SUCCESS);
279:     }
280:   }
281:   if (found) *found = PETSC_FALSE;
282:   PetscFunctionReturn(PETSC_SUCCESS);
283: }

285: /*@
286:   PetscViewerCGNSOpen - Opens a file for CGNS input/output.

288:   Collective

290:   Input Parameters:
291: + comm - MPI communicator
292: . name - name of file
293: - type - type of file
294: .vb
295:     FILE_MODE_WRITE - create new file for binary output
296:     FILE_MODE_READ - open existing file for binary input
297:     FILE_MODE_APPEND - open existing file for binary output
298: .ve

300:   Output Parameter:
301: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file

303:   Level: beginner

305: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
306:           `DMLoad()`, `PetscFileMode`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
307: @*/
308: PetscErrorCode PetscViewerCGNSOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *viewer)
309: {
310:   PetscFunctionBegin;
311:   PetscAssertPointer(name, 2);
312:   PetscAssertPointer(viewer, 4);
313:   PetscCall(PetscViewerCreate(comm, viewer));
314:   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERCGNS));
315:   PetscCall(PetscViewerFileSetMode(*viewer, type));
316:   PetscCall(PetscViewerFileSetName(*viewer, name));
317:   PetscCall(PetscViewerSetFromOptions(*viewer));
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: /*@
322:   PetscViewerCGNSSetSolutionIndex - Set index of solution

324:   Not Collective

326:   Input Parameters:
327: + viewer      - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
328: - solution_id - Index of the solution id, or `-1` for the last solution on the file

330:   Level: intermediate

332:   Notes:
333:   By default, `solution_id` is set to `-1` to mean the last solution available in the file.
334:   If the file contains a FlowSolutionPointers node, then that array is indexed to determine which FlowSolution_t node to read from.
335:   Otherwise, `solution_id` indexes the total available FlowSolution_t nodes in the file.

337:   This solution index is used by `VecLoad()` to determine which solution to load from the file

339: .seealso: `PETSCVIEWERCGNS`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionInfo()`

341: @*/
342: PetscErrorCode PetscViewerCGNSSetSolutionIndex(PetscViewer viewer, PetscInt solution_id)
343: {
344:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

346:   PetscFunctionBeginUser;
349:   PetscCheck((solution_id != 0) && (solution_id > -2), PetscObjectComm((PetscObject)viewer), PETSC_ERR_USER_INPUT, "Solution index must be either -1 or greater than 0, not %" PetscInt_FMT, solution_id);
350:   cgv->solution_index      = solution_id;
351:   cgv->solution_file_index = 0; // Reset sol_index when solution_id changes (0 is invalid)
352:   PetscFunctionReturn(PETSC_SUCCESS);
353: }

355: /*@
356:   PetscViewerCGNSGetSolutionIndex - Get index of solution

358:   Not Collective

360:   Input Parameter:
361: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file

363:   Output Parameter:
364: . solution_id - Index of the solution id

366:   Level: intermediate

368:   Notes:
369:   By default, solution_id is set to `-1` to mean the last solution available in the file

371: .seealso: `PETSCVIEWERCGNS`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionInfo()`

373: @*/
374: PetscErrorCode PetscViewerCGNSGetSolutionIndex(PetscViewer viewer, PetscInt *solution_id)
375: {
376:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;

378:   PetscFunctionBeginUser;
380:   PetscAssertPointer(solution_id, 2);
381:   *solution_id = cgv->solution_index;
382:   PetscFunctionReturn(PETSC_SUCCESS);
383: }

385: // Gets the index for the solution in this CGNS file
386: PetscErrorCode PetscViewerCGNSGetSolutionFileIndex_Internal(PetscViewer viewer, int *sol_index)
387: {
388:   PetscViewer_CGNS *cgv  = (PetscViewer_CGNS *)viewer->data;
389:   MPI_Comm          comm = PetscObjectComm((PetscObject)viewer);
390:   int               nsols, cgns_ier;
391:   char              buffer[CGIO_MAX_NAME_LENGTH + 1];
392:   CGNS_ENUMT(GridLocation_t) gridloc; // Throwaway

394:   PetscFunctionBeginUser;
395:   if (cgv->solution_file_index > 0) {
396:     *sol_index = cgv->solution_file_index;
397:     PetscFunctionReturn(PETSC_SUCCESS);
398:   }

400:   PetscCallCGNSRead(cg_nsols(cgv->file_num, cgv->base, cgv->zone, &nsols), viewer, 0);
401:   cgns_ier = cg_goto(cgv->file_num, cgv->base, "Zone_t", cgv->zone, "ZoneIterativeData_t", 1, "FlowSolutionPointers", 0, NULL);
402:   if (cgns_ier == CG_NODE_NOT_FOUND) {
403:     // If FlowSolutionPointers does not exist, then just index off of nsols (which can include non-solution data)
404:     PetscCheck(cgv->solution_index == -1 || cgv->solution_index <= nsols, comm, PETSC_ERR_ARG_OUTOFRANGE, "CGNS Solution index (%" PetscInt_FMT ") not in [1, %d]", cgv->solution_index, nsols);

406:     cgv->solution_file_index = cgv->solution_index == -1 ? nsols : cgv->solution_index;
407:   } else {
408:     // If FlowSolutionPointers exists, then solution_id should index that array of FlowSolutions
409:     char     *pointer_id_name;
410:     PetscBool matches_name = PETSC_FALSE;
411:     int       sol_id;

413:     PetscCheck(cgns_ier == CG_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS error %d %s", cgns_ier, cg_get_error());
414:     cgns_ier = cg_goto(cgv->file_num, cgv->base, "Zone_t", cgv->zone, "ZoneIterativeData_t", 1, NULL);

416:     { // Get FlowSolutionPointer name corresponding to solution_id
417:       cgsize_t  size[12];
418:       int       dim, A_index;
419:       char     *pointer_names, *pointer_id_name_ref;
420:       PetscBool found_array;

422:       PetscCall(CGNS_Find_Array(comm, "FlowSolutionPointers", &A_index, NULL, &dim, size, &found_array));
423:       PetscCheck(found_array, comm, PETSC_ERR_SUP, "Cannot find FlowSolutionPointers array under current CGNS node");
424:       PetscCheck(cgv->solution_index == -1 || cgv->solution_index <= size[1], comm, PETSC_ERR_ARG_OUTOFRANGE, "CGNS Solution index (%" PetscInt_FMT ") not in range of FlowSolutionPointers [1, %" PRIdCGSIZE "]", cgv->solution_index, size[1]);
425:       PetscCall(PetscCalloc1(size[0] * size[1] + 1, &pointer_names)); // Need the +1 for (possibly) setting \0 for the last pointer name if it's full
426:       PetscCallCGNSRead(cg_array_read_as(1, CGNS_ENUMV(Character), pointer_names), viewer, 0);
427:       cgv->solution_file_pointer_index = cgv->solution_index == -1 ? size[1] : cgv->solution_index;
428:       pointer_id_name_ref              = &pointer_names[size[0] * (cgv->solution_file_pointer_index - 1)];
429:       { // Set last non-whitespace character of the pointer name to \0 (CGNS pads with spaces)
430:         int str_idx;
431:         for (str_idx = size[0] - 1; str_idx > 0; str_idx--) {
432:           if (!isspace((unsigned char)pointer_id_name_ref[str_idx])) break;
433:         }
434:         pointer_id_name_ref[str_idx + 1] = '\0';
435:       }
436:       PetscCall(PetscStrallocpy(pointer_id_name_ref, &pointer_id_name));
437:       PetscCall(PetscFree(pointer_names));
438:     }

440:     // Find FlowSolution_t node that matches pointer_id_name
441:     for (sol_id = 1; sol_id <= nsols; sol_id++) {
442:       PetscCallCGNSRead(cg_sol_info(cgv->file_num, cgv->base, cgv->zone, sol_id, buffer, &gridloc), viewer, 0);
443:       PetscCall(PetscStrcmp(pointer_id_name, buffer, &matches_name));
444:       if (matches_name) break;
445:     }
446:     PetscCheck(matches_name, comm, PETSC_ERR_LIB, "Cannot find FlowSolution_t node %s in file", pointer_id_name);
447:     cgv->solution_file_index = sol_id;
448:     PetscCall(PetscFree(pointer_id_name));
449:   }

451:   *sol_index = cgv->solution_file_index;
452:   PetscFunctionReturn(PETSC_SUCCESS);
453: }

455: /*@
456:   PetscViewerCGNSGetSolutionTime - Gets the solution time for the FlowSolution of the viewer

458:   Collective

460:   Input Parameter:
461: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file

463:   Output Parameters:
464: + time - Solution time of the FlowSolution_t node
465: - set  - Whether the time data is in the file

467:   Level: intermediate

469:   Notes:
470:   Reads data from a DataArray named `TimeValues` under a `BaseIterativeData_t` node

472: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSGetSolutionIteration()`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionName()`
473: @*/
474: PetscErrorCode PetscViewerCGNSGetSolutionTime(PetscViewer viewer, PetscReal *time, PetscBool *set)
475: {
476:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
477:   int               cgns_ier, A_index = 0, sol_id;
478:   PetscReal        *times;
479:   cgsize_t          size[12];

481:   PetscFunctionBeginUser;
483:   PetscAssertPointer(time, 2);
484:   PetscAssertPointer(set, 3);
485:   cgns_ier = cg_goto(cgv->file_num, cgv->base, "BaseIterativeData_t", 1, NULL);
486:   if (cgns_ier == CG_NODE_NOT_FOUND) {
487:     *set = PETSC_FALSE;
488:     PetscFunctionReturn(PETSC_SUCCESS);
489:   } else PetscCallCGNS(cgns_ier);
490:   PetscCall(CGNS_Find_Array(PetscObjectComm((PetscObject)viewer), "TimeValues", &A_index, NULL, NULL, size, set));
491:   if (!*set) PetscFunctionReturn(PETSC_SUCCESS);
492:   PetscCall(PetscMalloc1(size[0], &times));
493:   PetscCallCGNSRead(cg_array_read_as(A_index, CGNS_ENUMV(RealDouble), times), viewer, 0);
494:   PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id)); // Call to set file pointer index
495:   *time = times[cgv->solution_file_pointer_index - 1];
496:   PetscCall(PetscFree(times));
497:   PetscFunctionReturn(PETSC_SUCCESS);
498: }

500: /*@
501:   PetscViewerCGNSGetSolutionIteration - Gets the solution iteration for the FlowSolution of the viewer

503:   Collective

505:   Input Parameter:
506: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file

508:   Output Parameters:
509: + iteration - Solution iteration of the FlowSolution_t node
510: - set       - Whether the time data is in the file

512:   Level: intermediate

514:   Notes:
515:   Reads data from a DataArray named `IterationValues` under a `BaseIterativeData_t` node

517: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSGetSolutionTime()`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionName()`
518: @*/
519: PetscErrorCode PetscViewerCGNSGetSolutionIteration(PetscViewer viewer, PetscInt *iteration, PetscBool *set)
520: {
521:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
522:   int               cgns_ier, A_index = 0, sol_id;
523:   int              *steps;
524:   cgsize_t          size[12];

526:   PetscFunctionBeginUser;
528:   PetscAssertPointer(iteration, 2);
529:   PetscAssertPointer(set, 3);
530:   cgns_ier = cg_goto(cgv->file_num, cgv->base, "BaseIterativeData_t", 1, NULL);
531:   if (cgns_ier == CG_NODE_NOT_FOUND) {
532:     *set = PETSC_FALSE;
533:     PetscFunctionReturn(PETSC_SUCCESS);
534:   } else PetscCallCGNS(cgns_ier);
535:   PetscCall(CGNS_Find_Array(PetscObjectComm((PetscObject)viewer), "IterationValues", &A_index, NULL, NULL, size, set));
536:   if (!*set) PetscFunctionReturn(PETSC_SUCCESS);
537:   PetscCall(PetscMalloc1(size[0], &steps));
538:   PetscCallCGNSRead(cg_array_read_as(A_index, CGNS_ENUMV(Integer), steps), viewer, 0);
539:   PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id)); // Call to set file pointer index
540:   *iteration = (PetscInt)steps[cgv->solution_file_pointer_index - 1];
541:   PetscCall(PetscFree(steps));
542:   PetscFunctionReturn(PETSC_SUCCESS);
543: }

545: /*@
546:   PetscViewerCGNSGetSolutionName - Gets name of FlowSolution of the viewer

548:   Collective

550:   Input Parameter:
551: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file

553:   Output Parameter:
554: . name - Name of the FlowSolution_t node corresponding to the solution index

556:   Level: intermediate

558:   Notes:
559:   Currently assumes there is only one Zone in the CGNS file

561: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionTime()`
562: @*/
563: PetscErrorCode PetscViewerCGNSGetSolutionName(PetscViewer viewer, const char *name[])
564: {
565:   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
566:   int               sol_id;
567:   char              buffer[CGIO_MAX_NAME_LENGTH + 1];
568:   CGNS_ENUMT(GridLocation_t) gridloc; // Throwaway

570:   PetscFunctionBeginUser;
572:   PetscAssertPointer(name, 2);
573:   PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id));

575:   PetscCallCGNSRead(cg_sol_info(cgv->file_num, cgv->base, cgv->zone, sol_id, buffer, &gridloc), viewer, 0);
576:   if (cgv->solution_name) PetscCall(PetscFree(cgv->solution_name));
577:   PetscCall(PetscStrallocpy(buffer, &cgv->solution_name));
578:   *name = cgv->solution_name;
579:   PetscFunctionReturn(PETSC_SUCCESS);
580: }