xref: /honee/src/misc.c (revision b2948607630bb74332ade58fdba4825cfa6bcd1d)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Miscellaneous utility functions
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdm.h>
13926a6279SJames Wright #include <petscsf.h>
14e419654dSJeremy L Thompson #include <petscts.h>
15e419654dSJeremy L Thompson 
16a515125bSLeila Ghaffari #include "../navierstokes.h"
179f59f36eSJames Wright #include "../qfunctions/mass.h"
18a515125bSLeila Ghaffari 
192b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
20b4c37c5cSJames Wright   Ceed         ceed = user->ceed;
21*b2948607SJames Wright   CeedVector   mult_vec;
22*b2948607SJames Wright   PetscMemType m_mem_type;
23*b2948607SJames Wright   Vec          Multiplicity, Multiplicity_loc;
24*b2948607SJames Wright 
25a515125bSLeila Ghaffari   PetscFunctionBeginUser;
26b4c37c5cSJames Wright   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
278f18bb8bSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
28a515125bSLeila Ghaffari 
29b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
30a515125bSLeila Ghaffari 
31a515125bSLeila Ghaffari   // -- Get multiplicity
32*b2948607SJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
33*b2948607SJames Wright   PetscCall(VecP2C(Multiplicity_loc, &m_mem_type, mult_vec));
34b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
35*b2948607SJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, Multiplicity_loc));
36a515125bSLeila Ghaffari 
37*b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
38*b2948607SJames Wright   PetscCall(VecZeroEntries(Multiplicity));
39*b2948607SJames Wright   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari   // -- Fix multiplicity
42*b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
43*b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
44a515125bSLeila Ghaffari 
45*b2948607SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
46*b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
47b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
48d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
49a515125bSLeila Ghaffari }
50a515125bSLeila Ghaffari 
51c56e8d5bSJames Wright // Record boundary values from initial condition
52c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
53c56e8d5bSJames Wright   Vec Qbc, boundary_mask;
54c56e8d5bSJames Wright 
55c56e8d5bSJames Wright   PetscFunctionBeginUser;
56c56e8d5bSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
57c56e8d5bSJames Wright   PetscCall(VecCopy(Q_loc, Qbc));
58c56e8d5bSJames Wright   PetscCall(VecZeroEntries(Q_loc));
59c56e8d5bSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
60c56e8d5bSJames Wright   PetscCall(VecAXPY(Qbc, -1., Q_loc));
61c56e8d5bSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
62c56e8d5bSJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
63c56e8d5bSJames Wright 
64c56e8d5bSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
65c56e8d5bSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
66c56e8d5bSJames Wright   PetscCall(VecZeroEntries(boundary_mask));
67c56e8d5bSJames Wright   PetscCall(VecSet(Q, 1.0));
68c56e8d5bSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
69c56e8d5bSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
70c56e8d5bSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
71c56e8d5bSJames Wright }
72c56e8d5bSJames Wright 
73c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
742b916ea7SJeremy L Thompson                                                   Vec grad_FVM) {
759d437337SJames Wright   Vec Qbc, boundary_mask;
76a515125bSLeila Ghaffari 
7706f41313SJames Wright   PetscFunctionBeginUser;
782eb7bf1fSJames Wright   // Mask (zero) Strong BC entries
799d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
809d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
819d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
829d437337SJames Wright 
832b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
842b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
852b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
86d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
87a515125bSLeila Ghaffari }
88a515125bSLeila Ghaffari 
89e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
90e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
91e1233009SJames Wright   PetscInt   file_step_number;
92e1233009SJames Wright   PetscInt32 token;
93e7754af5SKenneth E. Jansen   PetscReal  file_time;
94e7754af5SKenneth E. Jansen 
9506f41313SJames Wright   PetscFunctionBeginUser;
96e1233009SJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
97e1233009SJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
98e1233009SJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
99e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
100e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
101e7754af5SKenneth E. Jansen     if (time) *time = file_time;
102e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
103e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
104e7754af5SKenneth E. Jansen     PetscInt length, N;
105e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
106e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
107e7754af5SKenneth E. Jansen     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
108e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
109e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
110e7754af5SKenneth E. Jansen 
111e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
112d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
113e7754af5SKenneth E. Jansen }
114e7754af5SKenneth E. Jansen 
115a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
116c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
117a515125bSLeila Ghaffari   Vec         Qref;
118a515125bSLeila Ghaffari   PetscViewer viewer;
119a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
120e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
121a515125bSLeila Ghaffari 
12206f41313SJames Wright   PetscFunctionBeginUser;
123a515125bSLeila Ghaffari   // Read reference file
1242b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
125e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
126e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
127a515125bSLeila Ghaffari 
128a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1292b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1302b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1312b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1322b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
133a515125bSLeila Ghaffari 
134a515125bSLeila Ghaffari   // Check error
135a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1362b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
137a515125bSLeila Ghaffari   }
138a515125bSLeila Ghaffari 
139a515125bSLeila Ghaffari   // Cleanup
1402b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1412b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
142d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
143a515125bSLeila Ghaffari }
144a515125bSLeila Ghaffari 
145a515125bSLeila Ghaffari // Get error for problems with exact solutions
146c56e8d5bSJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
147a515125bSLeila Ghaffari   PetscInt  loc_nodes;
148a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
149a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
150a515125bSLeila Ghaffari 
15106f41313SJames Wright   PetscFunctionBeginUser;
152a515125bSLeila Ghaffari   // Get exact solution at final time
153*b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1542b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1552b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1562b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1592b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1602b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1612b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
162a515125bSLeila Ghaffari 
163a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
1642b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1652b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
166*b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
167d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
168a515125bSLeila Ghaffari }
169a515125bSLeila Ghaffari 
170a515125bSLeila Ghaffari // Post-processing
171c56e8d5bSJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
172a515125bSLeila Ghaffari   PetscInt          steps;
173f0784ed3SJed Brown   TSConvergedReason reason;
174a515125bSLeila Ghaffari 
17506f41313SJames Wright   PetscFunctionBeginUser;
176a515125bSLeila Ghaffari   // Print relative error
1770e1e9333SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
178c56e8d5bSJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
179a515125bSLeila Ghaffari   }
180a515125bSLeila Ghaffari 
181a515125bSLeila Ghaffari   // Print final time and number of steps
1822b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
183f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1840e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
185f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
186f0784ed3SJed Brown                           steps, (double)final_time));
187a515125bSLeila Ghaffari   }
188a515125bSLeila Ghaffari 
189a515125bSLeila Ghaffari   // Output numerical values from command line
1902b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
191a515125bSLeila Ghaffari 
192a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
1930e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
194c56e8d5bSJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
195a515125bSLeila Ghaffari   }
196d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
197a515125bSLeila Ghaffari }
198a515125bSLeila Ghaffari 
199e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
200e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
201e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2029293eaa1SJed Brown 
203a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
204a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
205a515125bSLeila Ghaffari   PetscViewer viewer;
2062b916ea7SJeremy L Thompson 
20706f41313SJames Wright   PetscFunctionBeginUser;
2082b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
209e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2102b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
211d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
212a515125bSLeila Ghaffari }
213a515125bSLeila Ghaffari 
21415a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
21515a3537eSJed Brown int FreeContextPetsc(void *data) {
2162b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
21715a3537eSJed Brown   return CEED_ERROR_SUCCESS;
21815a3537eSJed Brown }
2199f59f36eSJames Wright 
2209f59f36eSJames Wright // Return mass qfunction specification for number of components N
2219f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2229f59f36eSJames Wright   PetscFunctionBeginUser;
2239f59f36eSJames Wright   switch (N) {
2249f59f36eSJames Wright     case 1:
225b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
2269f59f36eSJames Wright       break;
2279f59f36eSJames Wright     case 5:
228b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
2299f59f36eSJames Wright       break;
230c38c977aSJames Wright     case 7:
231b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
232c38c977aSJames Wright       break;
2339f59f36eSJames Wright     case 9:
234b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
2359f59f36eSJames Wright       break;
2369f59f36eSJames Wright     case 22:
237b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
2389f59f36eSJames Wright       break;
2399f59f36eSJames Wright     default:
2406f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2419f59f36eSJames Wright   }
2429f59f36eSJames Wright 
243b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
244b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
245b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
246d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2479f59f36eSJames Wright }
248e5e81594SJames Wright 
249457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
250457a5831SJames Wright   PetscFunctionBeginUser;
251d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
252457a5831SJames Wright 
253457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
254457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
255457a5831SJames Wright 
256457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
257457a5831SJames Wright 
258457a5831SJames Wright   PetscCall(PetscFree(context));
259d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
260457a5831SJames Wright }
261457a5831SJames Wright 
262676080b4SJames Wright /*
263676080b4SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
264676080b4SJames Wright  *
265676080b4SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
266676080b4SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
267676080b4SJames Wright  *
268676080b4SJames Wright  * 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.
269676080b4SJames Wright  *
270676080b4SJames Wright  * @param[in]  comm           MPI_Comm for the program
271676080b4SJames Wright  * @param[in]  path           Path to the file
272676080b4SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
273676080b4SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
274676080b4SJames Wright  * @param[out] fp File        pointer to the opened file
275676080b4SJames Wright  */
27642454adaSJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
277676080b4SJames Wright                                  FILE **fp) {
278defe8520SJames Wright   int    ndims;
279676080b4SJames Wright   char   line[char_array_len];
280676080b4SJames Wright   char **array;
281676080b4SJames Wright 
282676080b4SJames Wright   PetscFunctionBeginUser;
283676080b4SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
284676080b4SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
285676080b4SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
286defe8520SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
287676080b4SJames Wright 
288676080b4SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
289676080b4SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
290d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
291676080b4SJames Wright }
292676080b4SJames Wright 
293676080b4SJames Wright /*
294676080b4SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
295676080b4SJames Wright  *
296676080b4SJames Wright  * 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.
297676080b4SJames Wright  *
298676080b4SJames Wright  * @param[in]  comm  MPI_Comm for the program
299676080b4SJames Wright  * @param[in]  path  Path to the file
300676080b4SJames Wright  * @param[out] nrows Number of rows
301676080b4SJames Wright  */
30242454adaSJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
303676080b4SJames Wright   const PetscInt char_array_len = 512;
304676080b4SJames Wright   PetscInt       dims[2];
305676080b4SJames Wright   FILE          *fp;
306676080b4SJames Wright 
307676080b4SJames Wright   PetscFunctionBeginUser;
30842454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
309676080b4SJames Wright   *nrows = dims[0];
310676080b4SJames Wright   PetscCall(PetscFClose(comm, fp));
311d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
312676080b4SJames Wright }
31362b7942eSJames Wright 
31442454adaSJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
315defe8520SJames Wright   PetscInt       dims[2];
316defe8520SJames Wright   int            ndims;
31762b7942eSJames Wright   FILE          *fp;
31862b7942eSJames Wright   const PetscInt char_array_len = 512;
31962b7942eSJames Wright   char           line[char_array_len];
32062b7942eSJames Wright   char         **row_array;
32162b7942eSJames Wright 
32206f41313SJames Wright   PetscFunctionBeginUser;
32342454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
32462b7942eSJames Wright 
32562b7942eSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
32662b7942eSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
32762b7942eSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3285d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
329defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
33062b7942eSJames Wright 
33162b7942eSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
33262b7942eSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
33362b7942eSJames Wright     }
33462b7942eSJames Wright   }
33562b7942eSJames Wright 
33662b7942eSJames Wright   PetscCall(PetscFClose(comm, fp));
337d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
33862b7942eSJames Wright }
3397eedc94cSJames Wright 
3407eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
3417eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
3427eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
3437eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
3447eedc94cSJames Wright static PetscClassId libCEED_classid;
3457eedc94cSJames Wright 
3467eedc94cSJames Wright PetscErrorCode RegisterLogEvents() {
3477eedc94cSJames Wright   PetscFunctionBeginUser;
3487eedc94cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
3497eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
3507eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
3517eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
3527eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
353d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3547eedc94cSJames Wright }
355defe8520SJames Wright 
356defe8520SJames Wright /**
357defe8520SJames Wright   @brief Translate array of CeedInt to PetscInt.
358defe8520SJames Wright     If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`.
359defe8520SJames Wright     Caller is responsible for freeing `array_petsc` with `free()`.
360defe8520SJames Wright 
361defe8520SJames Wright   @param[in]      num_entries  Number of array entries
362defe8520SJames Wright   @param[in,out]  array_ceed   Array of CeedInts
363defe8520SJames Wright   @param[out]     array_petsc  Array of PetscInts
364defe8520SJames Wright **/
365defe8520SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
366defe8520SJames Wright   CeedInt  int_c = 0;
367defe8520SJames Wright   PetscInt int_p = 0;
368defe8520SJames Wright 
369defe8520SJames Wright   PetscFunctionBeginUser;
370defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
371defe8520SJames Wright     *array_petsc = (PetscInt *)*array_ceed;
372defe8520SJames Wright   } else {
373defe8520SJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
374defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
375defe8520SJames Wright     free(*array_ceed);
376defe8520SJames Wright   }
377defe8520SJames Wright   *array_ceed = NULL;
378defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
379defe8520SJames Wright }
380defe8520SJames Wright 
381defe8520SJames Wright /**
382defe8520SJames Wright   @brief Translate array of PetscInt to CeedInt.
383defe8520SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
384defe8520SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
385defe8520SJames Wright 
386defe8520SJames Wright   @param[in]      num_entries  Number of array entries
387defe8520SJames Wright   @param[in,out]  array_petsc  Array of PetscInts
388defe8520SJames Wright   @param[out]     array_ceed   Array of CeedInts
389defe8520SJames Wright **/
390defe8520SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
391defe8520SJames Wright   CeedInt  int_c = 0;
392defe8520SJames Wright   PetscInt int_p = 0;
393defe8520SJames Wright 
394defe8520SJames Wright   PetscFunctionBeginUser;
395defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
396defe8520SJames Wright     *array_ceed = (CeedInt *)*array_petsc;
397defe8520SJames Wright   } else {
398defe8520SJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
399defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
400defe8520SJames Wright     PetscCall(PetscFree(*array_petsc));
401defe8520SJames Wright   }
402defe8520SJames Wright   *array_petsc = NULL;
403defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
404defe8520SJames Wright }
405926a6279SJames Wright 
406926a6279SJames Wright // Print information about the given simulation run
407926a6279SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
408b4c37c5cSJames Wright   Ceed ceed = user->ceed;
409926a6279SJames Wright   PetscFunctionBeginUser;
410926a6279SJames Wright   // Header and rank
411926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
412926a6279SJames Wright   PetscMPIInt rank, comm_size;
413926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
414926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
415926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
416926a6279SJames Wright   PetscCall(PetscPrintf(comm,
417926a6279SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
418926a6279SJames Wright                         "  MPI:\n"
419926a6279SJames Wright                         "    Host Name                          : %s\n"
420926a6279SJames Wright                         "    Total ranks                        : %d\n",
421926a6279SJames Wright                         host_name, comm_size));
422926a6279SJames Wright 
423926a6279SJames Wright   // Problem specific info
4242d49c0afSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
425926a6279SJames Wright 
426926a6279SJames Wright   // libCEED
427926a6279SJames Wright   const char *used_resource;
428926a6279SJames Wright   CeedMemType mem_type_backend;
429b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
430b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
431926a6279SJames Wright   PetscCall(PetscPrintf(comm,
432926a6279SJames Wright                         "  libCEED:\n"
433926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
434926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
435926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
436926a6279SJames Wright   // PETSc
437926a6279SJames Wright   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
438926a6279SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
439926a6279SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
440926a6279SJames Wright   MatType mat_type;
441926a6279SJames Wright   VecType vec_type;
442926a6279SJames Wright   PetscCall(DMGetMatType(user->dm, &mat_type));
443926a6279SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
444926a6279SJames Wright   PetscCall(PetscPrintf(comm,
445926a6279SJames Wright                         "  PETSc:\n"
446926a6279SJames Wright                         "    Box Faces                          : %s\n"
447926a6279SJames Wright                         "    DM MatType                         : %s\n"
448926a6279SJames Wright                         "    DM VecType                         : %s\n"
449926a6279SJames Wright                         "    Time Stepping Scheme               : %s\n",
450926a6279SJames Wright                         box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
451926a6279SJames Wright   if (user->app_ctx->cont_steps) {
452926a6279SJames Wright     PetscCall(PetscPrintf(comm,
453926a6279SJames Wright                           "  Continue:\n"
454926a6279SJames Wright                           "    Filename:                          : %s\n"
455926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
456926a6279SJames Wright                           "    Time:                              : %g\n",
457926a6279SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
458926a6279SJames Wright   }
459926a6279SJames Wright   // Mesh
460926a6279SJames Wright   const PetscInt num_comp_q = 5;
461926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
462926a6279SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
463926a6279SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
464926a6279SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
465926a6279SJames Wright   PetscCall(PetscPrintf(comm,
466926a6279SJames Wright                         "  Mesh:\n"
467926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
468926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
469926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
470926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
471dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
472dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
473926a6279SJames Wright   // -- Get Partition Statistics
474926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
475926a6279SJames Wright   {
476926a6279SJames Wright     PetscInt *gather_buffer = NULL;
477dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
478926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
479926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
480926a6279SJames Wright 
481dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
482926a6279SJames Wright     if (!rank) {
483926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
484926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
485926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
486926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
487926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
488dfeb939dSJames Wright       PetscCall(PetscPrintf(
489dfeb939dSJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
490926a6279SJames Wright           part_owned_dofs[0] / num_comp_q, part_owned_dofs[1] / num_comp_q, part_owned_dofs[2] / num_comp_q, part_owned_dof_ratio));
491926a6279SJames Wright     }
492926a6279SJames Wright 
493dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
494926a6279SJames Wright     if (!rank) {
495926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
496926a6279SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
497926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
498926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
499926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
500dfeb939dSJames Wright       PetscCall(PetscPrintf(
501dfeb939dSJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
502926a6279SJames Wright           part_local_dofs[0] / num_comp_q, part_local_dofs[1] / num_comp_q, part_local_dofs[2] / num_comp_q, part_local_dof_ratio));
503926a6279SJames Wright     }
504926a6279SJames Wright 
50545abf96eSJames Wright     if (comm_size != 1) {
506dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
507926a6279SJames Wright       {
508926a6279SJames Wright         PetscSF            sf;
509dfeb939dSJames Wright         PetscInt           nrranks, niranks;
510dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
511dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
512926a6279SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
513926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
514dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
515926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
516926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
517926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
518dfeb939dSJames Wright           num_ghost_interface_ranks++;
519dfeb939dSJames Wright         }
520dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
521dfeb939dSJames Wright           if (iranks[i] == rank) continue;
522dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
523dfeb939dSJames Wright           num_owned_interface_ranks++;
524926a6279SJames Wright         }
525926a6279SJames Wright       }
526dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
527926a6279SJames Wright       if (!rank) {
528926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
529dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
530dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
531dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
532dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
533dfeb939dSJames Wright         PetscCall(PetscPrintf(
53445abf96eSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
53545abf96eSJames Wright             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
53645abf96eSJames Wright             part_shared_dof_ratio));
537dfeb939dSJames Wright       }
538dfeb939dSJames Wright 
539dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
540dfeb939dSJames Wright       if (!rank) {
541dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
542dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
543dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
544dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
545dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
546dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
547dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
548dfeb939dSJames Wright       }
549dfeb939dSJames Wright 
550dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
551dfeb939dSJames Wright       if (!rank) {
552dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
553dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
554dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
555dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
556dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
557dfeb939dSJames Wright         PetscCall(PetscPrintf(
55845abf96eSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
55945abf96eSJames Wright             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
56045abf96eSJames Wright             part_shared_dof_ratio));
561dfeb939dSJames Wright       }
562dfeb939dSJames Wright 
563dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
564dfeb939dSJames Wright       if (!rank) {
565dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
566dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
567dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
568dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
569dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
570dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
571dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
572926a6279SJames Wright       }
57345abf96eSJames Wright     }
574926a6279SJames Wright 
575926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
576926a6279SJames Wright   }
577926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
578926a6279SJames Wright }
579