xref: /honee/src/misc.c (revision 990f1db09bd8cb1d2333f3dc8e2adfc9884b9d5b)
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;
21b2948607SJames Wright   CeedVector   mult_vec;
22b2948607SJames Wright   PetscMemType m_mem_type;
23b2948607SJames Wright   Vec          Multiplicity, Multiplicity_loc;
24b2948607SJames 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
32b2948607SJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
33a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
34b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
35a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
36a515125bSLeila Ghaffari 
37b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
38b2948607SJames Wright   PetscCall(VecZeroEntries(Multiplicity));
39b2948607SJames Wright   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
40a515125bSLeila Ghaffari 
41a515125bSLeila Ghaffari   // -- Fix multiplicity
42b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
43b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
44a515125bSLeila Ghaffari 
45b2948607SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
46b2948607SJames 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 
89*990f1db0SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
90*990f1db0SJed Brown   PetscFunctionBeginUser;
91*990f1db0SJed Brown   if (file_type == PETSC_INT32) {
92*990f1db0SJed Brown     PetscInt32 val;
93*990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
94*990f1db0SJed Brown     *out = val;
95*990f1db0SJed Brown   } else if (file_type == PETSC_INT64) {
96*990f1db0SJed Brown     PetscInt64 val;
97*990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
98*990f1db0SJed Brown     *out = val;
99*990f1db0SJed Brown   } else {
100*990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
101*990f1db0SJed Brown   }
102*990f1db0SJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
103*990f1db0SJed Brown }
104*990f1db0SJed Brown 
105e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
106e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
107e1233009SJames Wright   PetscInt      file_step_number;
108e1233009SJames Wright   PetscInt32    token;
109e7754af5SKenneth E. Jansen   PetscReal     file_time;
110*990f1db0SJed Brown   PetscDataType file_type = PETSC_INT32;
111e7754af5SKenneth E. Jansen 
11206f41313SJames Wright   PetscFunctionBeginUser;
113e1233009SJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
114e1233009SJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
115e1233009SJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
116*990f1db0SJed Brown     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
117*990f1db0SJed Brown     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
118*990f1db0SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
119e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
120e7754af5SKenneth E. Jansen     if (time) *time = file_time;
121e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
122e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
123e7754af5SKenneth E. Jansen     PetscInt length, N;
124*990f1db0SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
125e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
126e7754af5SKenneth 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);
127e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
128e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
129e7754af5SKenneth E. Jansen 
130e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
131d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
132e7754af5SKenneth E. Jansen }
133e7754af5SKenneth E. Jansen 
134a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
135c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
136a515125bSLeila Ghaffari   Vec         Qref;
137a515125bSLeila Ghaffari   PetscViewer viewer;
138a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
139e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
140a515125bSLeila Ghaffari 
14106f41313SJames Wright   PetscFunctionBeginUser;
142a515125bSLeila Ghaffari   // Read reference file
1432b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
144e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
145e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
146a515125bSLeila Ghaffari 
147a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1482b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1492b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1502b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1512b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
152a515125bSLeila Ghaffari 
153a515125bSLeila Ghaffari   // Check error
154a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1552b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
156a515125bSLeila Ghaffari   }
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari   // Cleanup
1592b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1602b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
161d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
162a515125bSLeila Ghaffari }
163a515125bSLeila Ghaffari 
164a515125bSLeila Ghaffari // Get error for problems with exact solutions
165c56e8d5bSJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
166a515125bSLeila Ghaffari   PetscInt  loc_nodes;
167a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
168a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
169a515125bSLeila Ghaffari 
17006f41313SJames Wright   PetscFunctionBeginUser;
171a515125bSLeila Ghaffari   // Get exact solution at final time
172b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1732b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1742b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1752b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
176a515125bSLeila Ghaffari 
177a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1782b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1792b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1802b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
181a515125bSLeila Ghaffari 
182a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
1832b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1842b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
185b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
186d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
187a515125bSLeila Ghaffari }
188a515125bSLeila Ghaffari 
189a515125bSLeila Ghaffari // Post-processing
190c56e8d5bSJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
191a515125bSLeila Ghaffari   PetscInt          steps;
192f0784ed3SJed Brown   TSConvergedReason reason;
193a515125bSLeila Ghaffari 
19406f41313SJames Wright   PetscFunctionBeginUser;
195a515125bSLeila Ghaffari   // Print relative error
1960e1e9333SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
197c56e8d5bSJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
198a515125bSLeila Ghaffari   }
199a515125bSLeila Ghaffari 
200a515125bSLeila Ghaffari   // Print final time and number of steps
2012b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
202f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
2030e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
204f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
205f0784ed3SJed Brown                           steps, (double)final_time));
206a515125bSLeila Ghaffari   }
207a515125bSLeila Ghaffari 
208a515125bSLeila Ghaffari   // Output numerical values from command line
2092b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
210a515125bSLeila Ghaffari 
211a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
2120e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
213c56e8d5bSJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
214a515125bSLeila Ghaffari   }
215d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
216a515125bSLeila Ghaffari }
217a515125bSLeila Ghaffari 
218e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
219e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
220e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2219293eaa1SJed Brown 
222a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
223a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
224a515125bSLeila Ghaffari   PetscViewer viewer;
2252b916ea7SJeremy L Thompson 
22606f41313SJames Wright   PetscFunctionBeginUser;
2272b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
228e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2292b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
230d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
231a515125bSLeila Ghaffari }
232a515125bSLeila Ghaffari 
23315a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
23415a3537eSJed Brown int FreeContextPetsc(void *data) {
2352b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
23615a3537eSJed Brown   return CEED_ERROR_SUCCESS;
23715a3537eSJed Brown }
2389f59f36eSJames Wright 
2399f59f36eSJames Wright // Return mass qfunction specification for number of components N
2409f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2419f59f36eSJames Wright   PetscFunctionBeginUser;
2429f59f36eSJames Wright   switch (N) {
2439f59f36eSJames Wright     case 1:
244b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
2459f59f36eSJames Wright       break;
2469f59f36eSJames Wright     case 5:
247b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
2489f59f36eSJames Wright       break;
249c38c977aSJames Wright     case 7:
250b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
251c38c977aSJames Wright       break;
2529f59f36eSJames Wright     case 9:
253b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
2549f59f36eSJames Wright       break;
2559f59f36eSJames Wright     case 22:
256b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
2579f59f36eSJames Wright       break;
2589f59f36eSJames Wright     default:
2596f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2609f59f36eSJames Wright   }
2619f59f36eSJames Wright 
262b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
263b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
264b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2653170c09fSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
266d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2679f59f36eSJames Wright }
268e5e81594SJames Wright 
269457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
270457a5831SJames Wright   PetscFunctionBeginUser;
271d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
272457a5831SJames Wright 
273457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
274457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
275457a5831SJames Wright 
276457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
277457a5831SJames Wright 
278457a5831SJames Wright   PetscCall(PetscFree(context));
279d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
280457a5831SJames Wright }
281457a5831SJames Wright 
282676080b4SJames Wright /*
283676080b4SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
284676080b4SJames Wright  *
285676080b4SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
286676080b4SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
287676080b4SJames Wright  *
288676080b4SJames 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.
289676080b4SJames Wright  *
290676080b4SJames Wright  * @param[in]  comm           MPI_Comm for the program
291676080b4SJames Wright  * @param[in]  path           Path to the file
292676080b4SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
293676080b4SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
294676080b4SJames Wright  * @param[out] fp File        pointer to the opened file
295676080b4SJames Wright  */
29642454adaSJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
297676080b4SJames Wright                                  FILE **fp) {
298defe8520SJames Wright   int    ndims;
299676080b4SJames Wright   char   line[char_array_len];
300676080b4SJames Wright   char **array;
301676080b4SJames Wright 
302676080b4SJames Wright   PetscFunctionBeginUser;
303676080b4SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
304676080b4SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
305676080b4SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
306defe8520SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
307676080b4SJames Wright 
308676080b4SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
309676080b4SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
310d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
311676080b4SJames Wright }
312676080b4SJames Wright 
313676080b4SJames Wright /*
314676080b4SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
315676080b4SJames Wright  *
316676080b4SJames 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.
317676080b4SJames Wright  *
318676080b4SJames Wright  * @param[in]  comm  MPI_Comm for the program
319676080b4SJames Wright  * @param[in]  path  Path to the file
320676080b4SJames Wright  * @param[out] nrows Number of rows
321676080b4SJames Wright  */
32242454adaSJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
323676080b4SJames Wright   const PetscInt char_array_len = 512;
324676080b4SJames Wright   PetscInt       dims[2];
325676080b4SJames Wright   FILE          *fp;
326676080b4SJames Wright 
327676080b4SJames Wright   PetscFunctionBeginUser;
32842454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
329676080b4SJames Wright   *nrows = dims[0];
330676080b4SJames Wright   PetscCall(PetscFClose(comm, fp));
331d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
332676080b4SJames Wright }
33362b7942eSJames Wright 
33442454adaSJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
335defe8520SJames Wright   PetscInt       dims[2];
336defe8520SJames Wright   int            ndims;
33762b7942eSJames Wright   FILE          *fp;
33862b7942eSJames Wright   const PetscInt char_array_len = 512;
33962b7942eSJames Wright   char           line[char_array_len];
34062b7942eSJames Wright   char         **row_array;
34162b7942eSJames Wright 
34206f41313SJames Wright   PetscFunctionBeginUser;
34342454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
34462b7942eSJames Wright 
34562b7942eSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
34662b7942eSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
34762b7942eSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3485d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
349defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
35062b7942eSJames Wright 
35162b7942eSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
35262b7942eSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
35362b7942eSJames Wright     }
35462b7942eSJames Wright   }
35562b7942eSJames Wright 
35662b7942eSJames Wright   PetscCall(PetscFClose(comm, fp));
357d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
35862b7942eSJames Wright }
3597eedc94cSJames Wright 
3607eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
3617eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
3627eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
3637eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
364ad2e713eSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Init;
365ad2e713eSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Meta;
366ad2e713eSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Train;
367ad2e713eSRiccardo Balin PetscLogEvent       FLUIDS_TrainDataCompute;
368855536edSJames Wright PetscLogEvent       FLUIDS_DifferentialFilter;
369855536edSJames Wright PetscLogEvent       FLUIDS_VelocityGradientProjection;
370855536edSJames Wright static PetscClassId libCEED_classid, onlineTrain_classid, misc_classid;
3717eedc94cSJames Wright 
3727eedc94cSJames Wright PetscErrorCode RegisterLogEvents() {
3737eedc94cSJames Wright   PetscFunctionBeginUser;
3747eedc94cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
3757eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
3767eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
3777eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
3787eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
379ad2e713eSRiccardo Balin 
380ad2e713eSRiccardo Balin   PetscCall(PetscClassIdRegister("onlineTrain", &onlineTrain_classid));
381ad2e713eSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Init", onlineTrain_classid, &FLUIDS_SmartRedis_Init));
382ad2e713eSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Meta", onlineTrain_classid, &FLUIDS_SmartRedis_Meta));
383ad2e713eSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Train", onlineTrain_classid, &FLUIDS_SmartRedis_Train));
384ad2e713eSRiccardo Balin   PetscCall(PetscLogEventRegister("TrainDataCompute", onlineTrain_classid, &FLUIDS_TrainDataCompute));
385855536edSJames Wright 
386855536edSJames Wright   PetscCall(PetscClassIdRegister("Miscellaneous", &misc_classid));
387855536edSJames Wright   PetscCall(PetscLogEventRegister("DiffFilter", misc_classid, &FLUIDS_DifferentialFilter));
388855536edSJames Wright   PetscCall(PetscLogEventRegister("VeloGradProj", misc_classid, &FLUIDS_VelocityGradientProjection));
389d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3907eedc94cSJames Wright }
391defe8520SJames Wright 
392926a6279SJames Wright // Print information about the given simulation run
393926a6279SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
394b4c37c5cSJames Wright   Ceed ceed = user->ceed;
395926a6279SJames Wright   PetscFunctionBeginUser;
396926a6279SJames Wright   // Header and rank
397926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
398926a6279SJames Wright   PetscMPIInt rank, comm_size;
399926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
400926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
401926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
402926a6279SJames Wright   PetscCall(PetscPrintf(comm,
403926a6279SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
404926a6279SJames Wright                         "  MPI:\n"
405926a6279SJames Wright                         "    Host Name                          : %s\n"
406926a6279SJames Wright                         "    Total ranks                        : %d\n",
407926a6279SJames Wright                         host_name, comm_size));
408926a6279SJames Wright 
409926a6279SJames Wright   // Problem specific info
4102d49c0afSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
411926a6279SJames Wright 
412926a6279SJames Wright   // libCEED
413926a6279SJames Wright   const char *used_resource;
414926a6279SJames Wright   CeedMemType mem_type_backend;
415b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
416b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
417926a6279SJames Wright   PetscCall(PetscPrintf(comm,
418926a6279SJames Wright                         "  libCEED:\n"
419926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
420926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
421926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
422926a6279SJames Wright   // PETSc
423926a6279SJames Wright   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
424926a6279SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
425926a6279SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
42665149b89SJames Wright   MatType amat_type = user->app_ctx->amat_type, pmat_type;
427926a6279SJames Wright   VecType vec_type;
42865149b89SJames Wright   PetscCall(DMGetMatType(user->dm, &pmat_type));
42965149b89SJames Wright   if (!amat_type) amat_type = pmat_type;
430926a6279SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
431926a6279SJames Wright   PetscCall(PetscPrintf(comm,
432926a6279SJames Wright                         "  PETSc:\n"
433926a6279SJames Wright                         "    Box Faces                          : %s\n"
43465149b89SJames Wright                         "    A MatType                          : %s\n"
43565149b89SJames Wright                         "    P MatType                          : %s\n"
436926a6279SJames Wright                         "    DM VecType                         : %s\n"
437926a6279SJames Wright                         "    Time Stepping Scheme               : %s\n",
43865149b89SJames Wright                         box_faces_str, amat_type, pmat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
439926a6279SJames Wright   if (user->app_ctx->cont_steps) {
440926a6279SJames Wright     PetscCall(PetscPrintf(comm,
441926a6279SJames Wright                           "  Continue:\n"
442926a6279SJames Wright                           "    Filename:                          : %s\n"
443926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
444926a6279SJames Wright                           "    Time:                              : %g\n",
445926a6279SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
446926a6279SJames Wright   }
447926a6279SJames Wright   // Mesh
448926a6279SJames Wright   const PetscInt num_comp_q = 5;
449926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
450926a6279SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
451926a6279SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
452926a6279SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
453926a6279SJames Wright   PetscCall(PetscPrintf(comm,
454926a6279SJames Wright                         "  Mesh:\n"
455926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
456926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
457926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
458926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
459dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
460dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
461926a6279SJames Wright   // -- Get Partition Statistics
462926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
463926a6279SJames Wright   {
464926a6279SJames Wright     PetscInt *gather_buffer = NULL;
465dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
466926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
467926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
468926a6279SJames Wright 
469dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
470926a6279SJames Wright     if (!rank) {
471926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
472926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
473926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
474926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
475926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
476dfeb939dSJames Wright       PetscCall(PetscPrintf(
477dfeb939dSJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
478926a6279SJames 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));
479926a6279SJames Wright     }
480926a6279SJames Wright 
481dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_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_local_dofs[0]             = gather_buffer[0];              // min
485926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
486926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
487926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
488dfeb939dSJames Wright       PetscCall(PetscPrintf(
489dfeb939dSJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
490926a6279SJames 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));
491926a6279SJames Wright     }
492926a6279SJames Wright 
49345abf96eSJames Wright     if (comm_size != 1) {
494dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
495926a6279SJames Wright       {
496926a6279SJames Wright         PetscSF            sf;
497dfeb939dSJames Wright         PetscInt           nrranks, niranks;
498dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
499dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
500926a6279SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
501926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
502dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
503926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
504926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
505926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
506dfeb939dSJames Wright           num_ghost_interface_ranks++;
507dfeb939dSJames Wright         }
508dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
509dfeb939dSJames Wright           if (iranks[i] == rank) continue;
510dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
511dfeb939dSJames Wright           num_owned_interface_ranks++;
512926a6279SJames Wright         }
513926a6279SJames Wright       }
514dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
515926a6279SJames Wright       if (!rank) {
516926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
517dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
518dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
519dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
520dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
521dfeb939dSJames Wright         PetscCall(PetscPrintf(
52245abf96eSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
52345abf96eSJames 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,
52445abf96eSJames Wright             part_shared_dof_ratio));
525dfeb939dSJames Wright       }
526dfeb939dSJames Wright 
527dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
528dfeb939dSJames Wright       if (!rank) {
529dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
530dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
531dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
532dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
533dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
534dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
535dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
536dfeb939dSJames Wright       }
537dfeb939dSJames Wright 
538dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
539dfeb939dSJames Wright       if (!rank) {
540dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
541dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
542dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
543dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
544dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
545dfeb939dSJames Wright         PetscCall(PetscPrintf(
54645abf96eSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
54745abf96eSJames 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,
54845abf96eSJames Wright             part_shared_dof_ratio));
549dfeb939dSJames Wright       }
550dfeb939dSJames Wright 
551dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
552dfeb939dSJames Wright       if (!rank) {
553dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
554dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
555dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
556dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
557dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
558dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
559dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
560926a6279SJames Wright       }
56145abf96eSJames Wright     }
562926a6279SJames Wright 
563926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
564926a6279SJames Wright   }
565926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
566926a6279SJames Wright }
567