xref: /honee/src/misc.c (revision 437258a37724337fca931cc1264457abb7b321be)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Miscellaneous utility functions
6 
7 #include <ceed.h>
8 #include <petscdm.h>
9 #include <petscsf.h>
10 #include <petscts.h>
11 
12 #include <navierstokes.h>
13 #include "../qfunctions/mass.h"
14 
15 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
16   Ceed         ceed = user->ceed;
17   CeedVector   mult_vec;
18   PetscMemType m_mem_type;
19   Vec          Multiplicity, Multiplicity_loc;
20 
21   PetscFunctionBeginUser;
22   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
23   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
24 
25   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
26 
27   // -- Get multiplicity
28   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
29   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
30   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
31   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
32 
33   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
34   PetscCall(VecZeroEntries(Multiplicity));
35   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
36 
37   // -- Fix multiplicity
38   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
39   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
40 
41   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
42   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
43   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
44   PetscFunctionReturn(PETSC_SUCCESS);
45 }
46 
47 // Record boundary values from initial condition
48 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
49   PetscFunctionBeginUser;
50   {  // Capture initial condition values in Qbc
51     Vec Qbc;
52 
53     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
54     PetscCall(VecCopy(Q_loc, Qbc));
55     PetscCall(VecZeroEntries(Q_loc));
56     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
57     PetscCall(VecAXPY(Qbc, -1., Q_loc));
58     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
59   }
60   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
61 
62   {  // Set boundary mask to zero out essential BCs
63     Vec boundary_mask, ones;
64 
65     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
66     PetscCall(DMGetGlobalVector(dm, &ones));
67     PetscCall(VecZeroEntries(boundary_mask));
68     PetscCall(VecSet(ones, 1.0));
69     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
70     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
71     PetscCall(DMRestoreGlobalVector(dm, &ones));
72   }
73   PetscFunctionReturn(PETSC_SUCCESS);
74 }
75 
76 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
77                                                   Vec grad_FVM) {
78   Vec Qbc, boundary_mask;
79 
80   PetscFunctionBeginUser;
81   // Mask (zero) Strong BC entries
82   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
83   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
84   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
85 
86   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
87   PetscCall(VecAXPY(Q_loc, 1., Qbc));
88   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91 
92 // Compare reference solution values with current test run for CI
93 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
94   Vec         Qref;
95   PetscViewer viewer;
96   PetscReal   error, Qrefnorm;
97   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
98 
99   PetscFunctionBeginUser;
100   // Read reference file
101   PetscCall(VecDuplicate(Q, &Qref));
102   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
103   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
104   PetscCall(HoneeLoadBinaryVec(viewer, Qref, NULL, NULL));
105 
106   // Compute error with respect to reference solution
107   PetscCall(VecAXPY(Q, -1.0, Qref));
108   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
109   PetscCall(VecScale(Q, 1. / Qrefnorm));
110   PetscCall(VecNorm(Q, NORM_MAX, &error));
111 
112   // Check error
113   if (error > app_ctx->test_tol) {
114     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
115   }
116 
117   // Cleanup
118   PetscCall(PetscViewerDestroy(&viewer));
119   PetscCall(VecDestroy(&Qref));
120   PetscFunctionReturn(PETSC_SUCCESS);
121 }
122 
123 // Get error for problems with exact solutions
124 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
125   PetscInt  loc_nodes;
126   Vec       Q_exact, Q_exact_loc;
127   PetscReal rel_error, norm_error, norm_exact;
128 
129   PetscFunctionBeginUser;
130   // Get exact solution at final time
131   PetscCall(DMGetGlobalVector(dm, &Q_exact));
132   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
133   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
134   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
135 
136   // Get |exact solution - obtained solution|
137   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
138   PetscCall(VecAXPY(Q, -1.0, Q_exact));
139   PetscCall(VecNorm(Q, NORM_1, &norm_error));
140 
141   rel_error = norm_error / norm_exact;
142   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
143   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
144   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 // Post-processing
149 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
150   PetscInt          steps;
151   TSConvergedReason reason;
152 
153   PetscFunctionBeginUser;
154   // Print relative error
155   if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
156     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
157   }
158 
159   // Print final time and number of steps
160   PetscCall(TSGetStepNumber(ts, &steps));
161   PetscCall(TSGetConvergedReason(ts, &reason));
162   if (user->app_ctx->test_type == TESTTYPE_NONE) {
163     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
164                           steps, (double)final_time));
165   }
166 
167   // Output numerical values from command line
168   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
169 
170   // Compare reference solution values with current test run for CI
171   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
172     PetscCall(RegressionTest(user->app_ctx, Q));
173   }
174   PetscFunctionReturn(PETSC_SUCCESS);
175 }
176 
177 // Free a plain data context that was allocated using PETSc; returning libCEED error codes
178 int FreeContextPetsc(void *data) {
179   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
180   return CEED_ERROR_SUCCESS;
181 }
182 
183 // Return mass qfunction specification for number of components N
184 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
185   PetscFunctionBeginUser;
186   switch (N) {
187     case 1:
188       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
189       break;
190     case 5:
191       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
192       break;
193     case 7:
194       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
195       break;
196     case 9:
197       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
198       break;
199     case 22:
200       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
201       break;
202     default:
203       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
204   }
205 
206   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
207   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
208   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
209   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
210   PetscFunctionReturn(PETSC_SUCCESS);
211 }
212 
213 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
214   PetscFunctionBeginUser;
215   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
216 
217   PetscCall(DMDestroy(&context->dm));
218   PetscCall(KSPDestroy(&context->ksp));
219 
220   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
221 
222   PetscCall(PetscFree(context));
223   PetscFunctionReturn(PETSC_SUCCESS);
224 }
225 
226 // Print information about the given simulation run
227 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
228   Ceed     ceed = user->ceed;
229   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
230 
231   PetscFunctionBeginUser;
232   // Header and rank
233   char        host_name[PETSC_MAX_PATH_LEN];
234   PetscMPIInt rank, comm_size;
235   PetscCall(PetscGetHostName(host_name, sizeof host_name));
236   PetscCallMPI(MPI_Comm_rank(comm, &rank));
237   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
238   PetscCall(PetscPrintf(comm,
239                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
240                         "  MPI:\n"
241                         "    Host Name                          : %s\n"
242                         "    Total ranks                        : %d\n",
243                         host_name, comm_size));
244 
245   // Problem specific info
246   PetscCall(problem->print_info(user, problem, user->app_ctx));
247 
248   // libCEED
249   const char *used_resource;
250   CeedMemType mem_type_backend;
251   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
252   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
253   PetscCall(PetscPrintf(comm,
254                         "  libCEED:\n"
255                         "    libCEED Backend                    : %s\n"
256                         "    libCEED Backend MemType            : %s\n",
257                         used_resource, CeedMemTypes[mem_type_backend]));
258   // PETSc
259   VecType vec_type;
260   char    box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
261   if (problem->dim == 2) box_faces_str[3] = '\0';
262   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
263   PetscCall(DMGetVecType(user->dm, &vec_type));
264   PetscCall(PetscPrintf(comm,
265                         "  PETSc:\n"
266                         "    Box Faces                          : %s\n"
267                         "    DM VecType                         : %s\n"
268                         "    Time Stepping Scheme               : %s\n",
269                         box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
270   {
271     char           pmat_type_str[PETSC_MAX_PATH_LEN];
272     MatType        amat_type, pmat_type;
273     Mat            Amat, Pmat;
274     TSIJacobianFn *ijacob_function;
275 
276     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
277     PetscCall(MatGetType(Amat, &amat_type));
278     PetscCall(MatGetType(Pmat, &pmat_type));
279 
280     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
281     if (!strcmp(pmat_type, MATCEED)) {
282       MatType pmat_coo_type;
283       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
284 
285       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
286       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
287       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
288     }
289     if (ijacob_function) {
290       PetscCall(PetscPrintf(comm,
291                             "    IJacobian A MatType                : %s\n"
292                             "    IJacobian P MatType                : %s\n",
293                             amat_type, pmat_type_str));
294     }
295   }
296   if (user->app_ctx->cont_steps) {
297     PetscCall(PetscPrintf(comm,
298                           "  Continue:\n"
299                           "    Filename:                          : %s\n"
300                           "    Step:                              : %" PetscInt_FMT "\n"
301                           "    Time:                              : %g\n",
302                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
303   }
304   // Mesh
305   const PetscInt num_comp_q = 5;
306   PetscInt       glob_dofs, owned_dofs, local_dofs;
307   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
308   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
309   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
310   PetscCall(PetscPrintf(comm,
311                         "  Mesh:\n"
312                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
313                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
314                         "    Global DoFs                        : %" PetscInt_FMT "\n"
315                         "    DoFs per node                      : %" PetscInt_FMT "\n"
316                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
317                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
318   // -- Get Partition Statistics
319   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
320   {
321     PetscInt *gather_buffer = NULL;
322     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
323     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
324     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
325 
326     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
327     if (!rank) {
328       PetscCall(PetscSortInt(comm_size, gather_buffer));
329       part_owned_dofs[0]             = gather_buffer[0];              // min
330       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
331       part_owned_dofs[2]             = gather_buffer[median_index];   // median
332       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
333       PetscCall(PetscPrintf(
334           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
335           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));
336     }
337 
338     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
339     if (!rank) {
340       PetscCall(PetscSortInt(comm_size, gather_buffer));
341       part_local_dofs[0]             = gather_buffer[0];              // min
342       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
343       part_local_dofs[2]             = gather_buffer[median_index];   // median
344       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
345       PetscCall(PetscPrintf(
346           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
347           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));
348     }
349 
350     if (comm_size != 1) {
351       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
352       {
353         PetscSF            sf;
354         PetscInt           nrranks, niranks;
355         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
356         const PetscMPIInt *rranks, *iranks;
357         PetscCall(DMGetSectionSF(user->dm, &sf));
358         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
359         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
360         for (PetscInt i = 0; i < nrranks; i++) {
361           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
362           num_remote_roots_total += roffset[i + 1] - roffset[i];
363           num_ghost_interface_ranks++;
364         }
365         for (PetscInt i = 0; i < niranks; i++) {
366           if (iranks[i] == rank) continue;
367           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
368           num_owned_interface_ranks++;
369         }
370       }
371       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
372       if (!rank) {
373         PetscCall(PetscSortInt(comm_size, gather_buffer));
374         part_boundary_dofs[0]           = gather_buffer[0];              // min
375         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
376         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
377         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
378         PetscCall(PetscPrintf(
379             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
380             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
381             part_shared_dof_ratio));
382       }
383 
384       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
385       if (!rank) {
386         PetscCall(PetscSortInt(comm_size, gather_buffer));
387         part_neighbors[0]              = gather_buffer[0];              // min
388         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
389         part_neighbors[2]              = gather_buffer[median_index];   // median
390         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
391         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
392                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
393       }
394 
395       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
396       if (!rank) {
397         PetscCall(PetscSortInt(comm_size, gather_buffer));
398         part_boundary_dofs[0]           = gather_buffer[0];              // min
399         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
400         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
401         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
402         PetscCall(PetscPrintf(
403             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
404             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
405             part_shared_dof_ratio));
406       }
407 
408       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
409       if (!rank) {
410         PetscCall(PetscSortInt(comm_size, gather_buffer));
411         part_neighbors[0]              = gather_buffer[0];              // min
412         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
413         part_neighbors[2]              = gather_buffer[median_index];   // median
414         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
415         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
416                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
417       }
418     }
419 
420     if (!rank) PetscCall(PetscFree(gather_buffer));
421   }
422   PetscFunctionReturn(PETSC_SUCCESS);
423 }
424