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, ×));
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));
155: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
156: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
157: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
158: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
159: PetscFunctionReturn(PETSC_SUCCESS);
160: }
162: static PetscErrorCode PetscViewerFileSetMode_CGNS(PetscViewer viewer, PetscFileMode type)
163: {
164: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
166: PetscFunctionBegin;
167: cgv->btype = type;
168: PetscFunctionReturn(PETSC_SUCCESS);
169: }
171: static PetscErrorCode PetscViewerFileGetMode_CGNS(PetscViewer viewer, PetscFileMode *type)
172: {
173: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
175: PetscFunctionBegin;
176: *type = cgv->btype;
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: static PetscErrorCode PetscViewerFileSetName_CGNS(PetscViewer viewer, const char *filename)
181: {
182: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
183: char *has_pattern;
185: PetscFunctionBegin;
186: #if defined(PETSC_HDF5_HAVE_PARALLEL)
187: if (cgv->file_num) PetscCallCGNSClose(cgp_close(cgv->file_num), viewer, 0);
188: #else
189: if (cgv->file_num) PetscCallCGNSClose(cg_close(cgv->file_num), viewer, 0);
190: #endif
191: PetscCall(PetscFree(cgv->filename));
192: PetscCall(PetscFree(cgv->filename_template));
193: PetscCall(PetscStrchr(filename, '%', &has_pattern));
194: if (has_pattern) {
195: PetscCall(PetscStrallocpy(filename, &cgv->filename_template));
196: // Templated file names open lazily, once we know DMGetOutputSequenceNumber()
197: } else {
198: PetscCall(PetscStrallocpy(filename, &cgv->filename));
199: PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, -1));
200: }
201: PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: static PetscErrorCode PetscViewerFileGetName_CGNS(PetscViewer viewer, const char **filename)
205: {
206: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
208: PetscFunctionBegin;
209: *filename = cgv->filename;
210: PetscFunctionReturn(PETSC_SUCCESS);
211: }
213: /*MC
214: PETSCVIEWERCGNS - A viewer for CGNS files
216: Level: beginner
218: Options Database Key:
219: . -viewer_cgns_batch_size SIZE - set max number of output sequence times to write per batch
221: Note:
222: If the filename contains an integer format character, the CGNS viewer will created a batched output sequence. For
223: example, one could use `-ts_monitor_solution cgns:flow-%d.cgns`. This is desirable if one wants to limit file sizes or
224: if the job might crash/be killed by a resource manager before exiting cleanly.
226: .seealso: [](sec_viewers), `PetscViewer`, `PetscViewerCreate()`, `VecView()`, `DMView()`, `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `TSSetFromOptions()`
227: M*/
228: PetscErrorCode PetscViewerCreate_CGNS(PetscViewer v)
229: {
230: PetscViewer_CGNS *cgv;
232: PetscFunctionBegin;
233: PetscCall(PetscViewerCGNSRegisterLogEvents_Internal());
234: PetscCall(PetscNew(&cgv));
236: v->data = cgv;
237: v->ops->destroy = PetscViewerDestroy_CGNS;
238: v->ops->setfromoptions = PetscViewerSetFromOptions_CGNS;
239: v->ops->view = PetscViewerView_CGNS;
240: cgv->btype = FILE_MODE_UNDEFINED;
241: cgv->filename = NULL;
242: cgv->batch_size = 20;
243: cgv->solution_index = -1; // Default to use the "last" solution
244: cgv->base = 1;
245: cgv->zone = 1;
247: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_CGNS));
248: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_CGNS));
249: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_CGNS));
250: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_CGNS));
251: PetscFunctionReturn(PETSC_SUCCESS);
252: }
254: // Find DataArray_t node under the current node (determined by `cg_goto` and friends) that matches `name`
255: // Return the index of that array and (optionally) other data about the array
256: 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)
257: {
258: int narrays; // number of arrays under the current node
259: char array_name[CGIO_MAX_NAME_LENGTH + 1];
260: CGNS_ENUMT(DataType_t) data_type_local;
261: int _dim;
262: cgsize_t _size[12];
263: PetscBool matches_name = PETSC_FALSE;
265: PetscFunctionBeginUser;
266: PetscCallCGNSRead(cg_narrays(&narrays), 0, 0);
267: for (int i = 1; i <= narrays; i++) {
268: PetscCallCGNSRead(cg_array_info(i, array_name, &data_type_local, &_dim, _size), 0, 0);
269: PetscCall(PetscStrcmp(name, array_name, &matches_name));
270: if (matches_name) {
271: *A_index = i;
272: if (data_type) *data_type = data_type_local;
273: if (dim) *dim = _dim;
274: if (size) PetscArraycpy(size, _size, _dim);
275: if (found) *found = PETSC_TRUE;
276: PetscFunctionReturn(PETSC_SUCCESS);
277: }
278: }
279: if (found) *found = PETSC_FALSE;
280: PetscFunctionReturn(PETSC_SUCCESS);
281: }
283: /*@
284: PetscViewerCGNSOpen - Opens a file for CGNS input/output.
286: Collective
288: Input Parameters:
289: + comm - MPI communicator
290: . name - name of file
291: - type - type of file
292: .vb
293: FILE_MODE_WRITE - create new file for binary output
294: FILE_MODE_READ - open existing file for binary input
295: FILE_MODE_APPEND - open existing file for binary output
296: .ve
298: Output Parameter:
299: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
301: Level: beginner
303: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
304: `DMLoad()`, `PetscFileMode`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
305: @*/
306: PetscErrorCode PetscViewerCGNSOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *viewer)
307: {
308: PetscFunctionBegin;
309: PetscAssertPointer(name, 2);
310: PetscAssertPointer(viewer, 4);
311: PetscCall(PetscViewerCreate(comm, viewer));
312: PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERCGNS));
313: PetscCall(PetscViewerFileSetMode(*viewer, type));
314: PetscCall(PetscViewerFileSetName(*viewer, name));
315: PetscCall(PetscViewerSetFromOptions(*viewer));
316: PetscFunctionReturn(PETSC_SUCCESS);
317: }
319: /*@
320: PetscViewerCGNSSetSolutionIndex - Set index of solution
322: Not Collective
324: Input Parameters:
325: + viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
326: - solution_id - Index of the solution id, or `-1` for the last solution on the file
328: Level: intermediate
330: Notes:
331: By default, `solution_id` is set to `-1` to mean the last solution available in the file.
332: If the file contains a FlowSolutionPointers node, then that array is indexed to determine which FlowSolution_t node to read from.
333: Otherwise, `solution_id` indexes the total available FlowSolution_t nodes in the file.
335: This solution index is used by `VecLoad()` to determine which solution to load from the file
337: .seealso: `PETSCVIEWERCGNS`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionInfo()`
339: @*/
340: PetscErrorCode PetscViewerCGNSSetSolutionIndex(PetscViewer viewer, PetscInt solution_id)
341: {
342: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
344: PetscFunctionBeginUser;
347: 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);
348: cgv->solution_index = solution_id;
349: cgv->solution_file_index = 0; // Reset sol_index when solution_id changes (0 is invalid)
350: PetscFunctionReturn(PETSC_SUCCESS);
351: }
353: /*@
354: PetscViewerCGNSGetSolutionIndex - Get index of solution
356: Not Collective
358: Input Parameter:
359: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
361: Output Parameter:
362: . solution_id - Index of the solution id
364: Level: intermediate
366: Notes:
367: By default, solution_id is set to `-1` to mean the last solution available in the file
369: .seealso: `PETSCVIEWERCGNS`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionInfo()`
371: @*/
372: PetscErrorCode PetscViewerCGNSGetSolutionIndex(PetscViewer viewer, PetscInt *solution_id)
373: {
374: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
376: PetscFunctionBeginUser;
378: PetscAssertPointer(solution_id, 2);
379: *solution_id = cgv->solution_index;
380: PetscFunctionReturn(PETSC_SUCCESS);
381: }
383: // Gets the index for the solution in this CGNS file
384: PetscErrorCode PetscViewerCGNSGetSolutionFileIndex_Internal(PetscViewer viewer, int *sol_index)
385: {
386: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
387: MPI_Comm comm = PetscObjectComm((PetscObject)viewer);
388: int nsols, cgns_ier;
389: char buffer[CGIO_MAX_NAME_LENGTH + 1];
390: CGNS_ENUMT(GridLocation_t) gridloc; // Throwaway
392: PetscFunctionBeginUser;
393: if (cgv->solution_file_index > 0) {
394: *sol_index = cgv->solution_file_index;
395: PetscFunctionReturn(PETSC_SUCCESS);
396: }
398: PetscCallCGNSRead(cg_nsols(cgv->file_num, cgv->base, cgv->zone, &nsols), viewer, 0);
399: cgns_ier = cg_goto(cgv->file_num, cgv->base, "Zone_t", cgv->zone, "ZoneIterativeData_t", 1, "FlowSolutionPointers", 0, NULL);
400: if (cgns_ier == CG_NODE_NOT_FOUND) {
401: // If FlowSolutionPointers does not exist, then just index off of nsols (which can include non-solution data)
402: 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);
404: cgv->solution_file_index = cgv->solution_index == -1 ? nsols : cgv->solution_index;
405: } else {
406: // If FlowSolutionPointers exists, then solution_id should index that array of FlowSolutions
407: char *pointer_id_name;
408: PetscBool matches_name = PETSC_FALSE;
409: int sol_id;
411: PetscCheck(cgns_ier == CG_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS error %d %s", cgns_ier, cg_get_error());
412: cgns_ier = cg_goto(cgv->file_num, cgv->base, "Zone_t", cgv->zone, "ZoneIterativeData_t", 1, NULL);
414: { // Get FlowSolutionPointer name corresponding to solution_id
415: cgsize_t size[12];
416: int dim, A_index;
417: char *pointer_names, *pointer_id_name_ref;
418: PetscBool found_array;
420: PetscCall(CGNS_Find_Array(comm, "FlowSolutionPointers", &A_index, NULL, &dim, size, &found_array));
421: PetscCheck(found_array, comm, PETSC_ERR_SUP, "Cannot find FlowSolutionPointers array under current CGNS node");
422: 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]);
423: 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
424: PetscCallCGNSRead(cg_array_read_as(1, CGNS_ENUMV(Character), pointer_names), viewer, 0);
425: cgv->solution_file_pointer_index = cgv->solution_index == -1 ? size[1] : cgv->solution_index;
426: pointer_id_name_ref = &pointer_names[size[0] * (cgv->solution_file_pointer_index - 1)];
427: { // Set last non-whitespace character of the pointer name to \0 (CGNS pads with spaces)
428: int str_idx;
429: for (str_idx = size[0] - 1; str_idx > 0; str_idx--) {
430: if (!isspace((unsigned char)pointer_id_name_ref[str_idx])) break;
431: }
432: pointer_id_name_ref[str_idx + 1] = '\0';
433: }
434: PetscCall(PetscStrallocpy(pointer_id_name_ref, &pointer_id_name));
435: PetscCall(PetscFree(pointer_names));
436: }
438: // Find FlowSolution_t node that matches pointer_id_name
439: for (sol_id = 1; sol_id <= nsols; sol_id++) {
440: PetscCallCGNSRead(cg_sol_info(cgv->file_num, cgv->base, cgv->zone, sol_id, buffer, &gridloc), viewer, 0);
441: PetscCall(PetscStrcmp(pointer_id_name, buffer, &matches_name));
442: if (matches_name) break;
443: }
444: PetscCheck(matches_name, comm, PETSC_ERR_LIB, "Cannot find FlowSolution_t node %s in file", pointer_id_name);
445: cgv->solution_file_index = sol_id;
446: PetscCall(PetscFree(pointer_id_name));
447: }
449: *sol_index = cgv->solution_file_index;
450: PetscFunctionReturn(PETSC_SUCCESS);
451: }
453: /*@
454: PetscViewerCGNSGetSolutionTime - Gets the solution time for the FlowSolution of the viewer
456: Collective
458: Input Parameter:
459: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
461: Output Parameters:
462: + time - Solution time of the FlowSolution_t node
463: - set - Whether the time data is in the file
465: Level: intermediate
467: Notes:
468: Reads data from a DataArray named `TimeValues` under a `BaseIterativeData_t` node
470: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSGetSolutionIteration()`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionName()`
471: @*/
472: PetscErrorCode PetscViewerCGNSGetSolutionTime(PetscViewer viewer, PetscReal *time, PetscBool *set)
473: {
474: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
475: int cgns_ier, A_index = 0, sol_id;
476: PetscReal *times;
477: cgsize_t size[12];
479: PetscFunctionBeginUser;
481: PetscAssertPointer(time, 2);
482: PetscAssertPointer(set, 3);
483: cgns_ier = cg_goto(cgv->file_num, cgv->base, "BaseIterativeData_t", 1, NULL);
484: if (cgns_ier == CG_NODE_NOT_FOUND) {
485: *set = PETSC_FALSE;
486: PetscFunctionReturn(PETSC_SUCCESS);
487: } else PetscCallCGNS(cgns_ier);
488: PetscCall(CGNS_Find_Array(PetscObjectComm((PetscObject)viewer), "TimeValues", &A_index, NULL, NULL, size, set));
489: if (!set) PetscFunctionReturn(PETSC_SUCCESS);
490: PetscCall(PetscMalloc1(size[0], ×));
491: PetscCallCGNSRead(cg_array_read_as(A_index, CGNS_ENUMV(RealDouble), times), viewer, 0);
492: PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id)); // Call to set file pointer index
493: *time = times[cgv->solution_file_pointer_index - 1];
494: PetscCall(PetscFree(times));
495: PetscFunctionReturn(PETSC_SUCCESS);
496: }
498: /*@
499: PetscViewerCGNSGetSolutionIteration - Gets the solution iteration for the FlowSolution of the viewer
501: Collective
503: Input Parameter:
504: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
506: Output Parameters:
507: + iteration - Solution iteration of the FlowSolution_t node
508: - set - Whether the time data is in the file
510: Level: intermediate
512: Notes:
513: Reads data from a DataArray named `IterationValues` under a `BaseIterativeData_t` node
515: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSGetSolutionTime()`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionName()`
516: @*/
517: PetscErrorCode PetscViewerCGNSGetSolutionIteration(PetscViewer viewer, PetscInt *iteration, PetscBool *set)
518: {
519: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
520: int cgns_ier, A_index = 0, sol_id;
521: int *steps;
522: cgsize_t size[12];
524: PetscFunctionBeginUser;
526: PetscAssertPointer(iteration, 2);
527: PetscAssertPointer(set, 3);
528: cgns_ier = cg_goto(cgv->file_num, cgv->base, "BaseIterativeData_t", 1, NULL);
529: if (cgns_ier == CG_NODE_NOT_FOUND) {
530: *set = PETSC_FALSE;
531: PetscFunctionReturn(PETSC_SUCCESS);
532: } else PetscCallCGNS(cgns_ier);
533: PetscCall(CGNS_Find_Array(PetscObjectComm((PetscObject)viewer), "IterationValues", &A_index, NULL, NULL, size, set));
534: if (!set) PetscFunctionReturn(PETSC_SUCCESS);
535: PetscCall(PetscMalloc1(size[0], &steps));
536: PetscCallCGNSRead(cg_array_read_as(A_index, CGNS_ENUMV(Integer), steps), viewer, 0);
537: PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id)); // Call to set file pointer index
538: *iteration = (PetscInt)steps[cgv->solution_file_pointer_index - 1];
539: PetscCall(PetscFree(steps));
540: PetscFunctionReturn(PETSC_SUCCESS);
541: }
543: /*@
544: PetscViewerCGNSGetSolutionName - Gets name of FlowSolution of the viewer
546: Collective
548: Input Parameter:
549: . viewer - `PETSCVIEWERCGNS` `PetscViewer` for CGNS input/output to use with the specified file
551: Output Parameter:
552: . name - Name of the FlowSolution_t node corresponding to the solution index
554: Level: intermediate
556: Notes:
557: Currently assumes there is only one Zone in the CGNS file
559: .seealso: `PETSCVIEWERCGNS`, `PetscViewer`, `PetscViewerCGNSSetSolutionIndex()`, `PetscViewerCGNSGetSolutionIndex()`, `PetscViewerCGNSGetSolutionTime()`
560: @*/
561: PetscErrorCode PetscViewerCGNSGetSolutionName(PetscViewer viewer, const char *name[])
562: {
563: PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
564: int sol_id;
565: char buffer[CGIO_MAX_NAME_LENGTH + 1];
566: CGNS_ENUMT(GridLocation_t) gridloc; // Throwaway
568: PetscFunctionBeginUser;
570: PetscAssertPointer(name, 2);
571: PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &sol_id));
573: PetscCallCGNSRead(cg_sol_info(cgv->file_num, cgv->base, cgv->zone, sol_id, buffer, &gridloc), viewer, 0);
574: if (cgv->solution_name) PetscCall(PetscFree(cgv->solution_name));
575: PetscCall(PetscStrallocpy(buffer, &cgv->solution_name));
576: *name = cgv->solution_name;
577: PetscFunctionReturn(PETSC_SUCCESS);
578: }