xref: /honee/src/honee-file.c (revision 93ca29b6e094721a136ad10375037d18cdc7fb9d)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Custom file I/O functions for HONEE
6 
7 #include <honee-file.h>
8 
9 const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
10 const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
11 const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
12 
13 static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
14   PetscFunctionBeginUser;
15   *out = -13;  // appease the overzealous GCC compiler warning Gods
16   if (file_type == PETSC_INT32) {
17     PetscInt32 val;
18     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
19     *out = val;
20   } else if (file_type == PETSC_INT64) {
21     PetscInt64 val;
22     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
23     *out = val;
24   } else {
25     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
26   }
27   PetscFunctionReturn(PETSC_SUCCESS);
28 }
29 
30 // @brief Load vector from binary file, possibly with embedded solution time and step number
31 PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
32   PetscInt      file_step_number;
33   PetscInt32    token;
34   PetscReal     file_time;
35   PetscDataType file_type = PETSC_INT32;
36 
37   PetscFunctionBeginUser;
38   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
39   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
40       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
41     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
42     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
43     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
44     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
45     if (time) *time = file_time;
46     if (step_number) *step_number = file_step_number;
47   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
48     PetscInt length, N;
49     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
50     PetscCall(VecGetSize(Q, &N));
51     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
52     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
53   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
54 
55   PetscCall(VecLoad(Q, viewer));
56   PetscFunctionReturn(PETSC_SUCCESS);
57 }
58 
59 /*
60  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
61  *
62  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
63  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
64  *
65  * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space.
66  *
67  * @param[in]  comm           MPI_Comm for the program
68  * @param[in]  path           Path to the file
69  * @param[in]  char_array_len Length of the character array that should contain each line
70  * @param[out] dims           Dimensions of the file, taken from the first line of the file
71  * @param[out] fp File        pointer to the opened file
72  */
73 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
74                                  FILE **fp) {
75   int    ndims;
76   char   line[char_array_len];
77   char **array;
78 
79   PetscFunctionBeginUser;
80   PetscCall(PetscFOpen(comm, path, "r", fp));
81   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
82   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
83   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
84 
85   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
86   PetscCall(PetscStrToArrayDestroy(ndims, array));
87   PetscFunctionReturn(PETSC_SUCCESS);
88 }
89 
90 /*
91  * @brief Get the number of rows for the PHASTA file at path.
92  *
93  * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space.
94  *
95  * @param[in]  comm  MPI_Comm for the program
96  * @param[in]  path  Path to the file
97  * @param[out] nrows Number of rows
98  */
99 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
100   const PetscInt char_array_len = 512;
101   PetscInt       dims[2];
102   FILE          *fp;
103 
104   PetscFunctionBeginUser;
105   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
106   *nrows = dims[0];
107   PetscCall(PetscFClose(comm, fp));
108   PetscFunctionReturn(PETSC_SUCCESS);
109 }
110 
111 PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
112   PetscInt       dims[2];
113   FILE          *fp;
114   const PetscInt char_array_len = 512;
115   char           line[char_array_len];
116 
117   PetscFunctionBeginUser;
118   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
119 
120   for (PetscInt i = 0; i < dims[0]; i++) {
121     int    ndims;
122     char **row_array;
123 
124     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
125     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
126     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
127                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
128 
129     for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
130     PetscCall(PetscStrToArrayDestroy(ndims, row_array));
131   }
132 
133   PetscCall(PetscFClose(comm, fp));
134   PetscFunctionReturn(PETSC_SUCCESS);
135 }
136