xref: /honee/src/misc.c (revision 039caf0d64f861d88e5a4ad7f261290c43ab5010)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Miscellaneous utility functions
10 
11 #include <ceed.h>
12 #include <petscdm.h>
13 #include <petscsf.h>
14 #include <petscts.h>
15 
16 #include <navierstokes.h>
17 #include "../qfunctions/mass.h"
18 
19 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
20   Ceed         ceed = user->ceed;
21   CeedVector   mult_vec;
22   PetscMemType m_mem_type;
23   Vec          Multiplicity, Multiplicity_loc;
24 
25   PetscFunctionBeginUser;
26   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
27   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
28 
29   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
30 
31   // -- Get multiplicity
32   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
33   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
34   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
35   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
36 
37   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
38   PetscCall(VecZeroEntries(Multiplicity));
39   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
40 
41   // -- Fix multiplicity
42   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
43   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
44 
45   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
46   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
47   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 // Record boundary values from initial condition
52 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
53   PetscFunctionBeginUser;
54   {  // Capture initial condition values in Qbc
55     Vec Qbc;
56 
57     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
58     PetscCall(VecCopy(Q_loc, Qbc));
59     PetscCall(VecZeroEntries(Q_loc));
60     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
61     PetscCall(VecAXPY(Qbc, -1., Q_loc));
62     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
63   }
64   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
65 
66   {  // Set boundary mask to zero out essential BCs
67     Vec boundary_mask, ones;
68 
69     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
70     PetscCall(DMGetGlobalVector(dm, &ones));
71     PetscCall(VecZeroEntries(boundary_mask));
72     PetscCall(VecSet(ones, 1.0));
73     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
74     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
75     PetscCall(DMRestoreGlobalVector(dm, &ones));
76   }
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
81                                                   Vec grad_FVM) {
82   Vec Qbc, boundary_mask;
83 
84   PetscFunctionBeginUser;
85   // Mask (zero) Strong BC entries
86   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
87   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
88   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
89 
90   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
91   PetscCall(VecAXPY(Q_loc, 1., Qbc));
92   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
93   PetscFunctionReturn(PETSC_SUCCESS);
94 }
95 
96 static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
97   PetscFunctionBeginUser;
98   *out = -13;  // appease the overzealous GCC compiler warning Gods
99   if (file_type == PETSC_INT32) {
100     PetscInt32 val;
101     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
102     *out = val;
103   } else if (file_type == PETSC_INT64) {
104     PetscInt64 val;
105     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
106     *out = val;
107   } else {
108     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
109   }
110   PetscFunctionReturn(PETSC_SUCCESS);
111 }
112 
113 // @brief Load vector from binary file, possibly with embedded solution time and step number
114 PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
115   PetscInt      file_step_number;
116   PetscInt32    token;
117   PetscReal     file_time;
118   PetscDataType file_type = PETSC_INT32;
119 
120   PetscFunctionBeginUser;
121   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
122   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
123       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
124     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
125     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
126     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
127     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
128     if (time) *time = file_time;
129     if (step_number) *step_number = file_step_number;
130   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
131     PetscInt length, N;
132     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
133     PetscCall(VecGetSize(Q, &N));
134     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
135     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
136   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
137 
138   PetscCall(VecLoad(Q, viewer));
139   PetscFunctionReturn(PETSC_SUCCESS);
140 }
141 
142 // Compare reference solution values with current test run for CI
143 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
144   Vec         Qref;
145   PetscViewer viewer;
146   PetscReal   error, Qrefnorm;
147   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
148 
149   PetscFunctionBeginUser;
150   // Read reference file
151   PetscCall(VecDuplicate(Q, &Qref));
152   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
153   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
154   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
155 
156   // Compute error with respect to reference solution
157   PetscCall(VecAXPY(Q, -1.0, Qref));
158   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
159   PetscCall(VecScale(Q, 1. / Qrefnorm));
160   PetscCall(VecNorm(Q, NORM_MAX, &error));
161 
162   // Check error
163   if (error > app_ctx->test_tol) {
164     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
165   }
166 
167   // Cleanup
168   PetscCall(PetscViewerDestroy(&viewer));
169   PetscCall(VecDestroy(&Qref));
170   PetscFunctionReturn(PETSC_SUCCESS);
171 }
172 
173 // Get error for problems with exact solutions
174 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
175   PetscInt  loc_nodes;
176   Vec       Q_exact, Q_exact_loc;
177   PetscReal rel_error, norm_error, norm_exact;
178 
179   PetscFunctionBeginUser;
180   // Get exact solution at final time
181   PetscCall(DMGetGlobalVector(dm, &Q_exact));
182   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
183   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
184   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
185 
186   // Get |exact solution - obtained solution|
187   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
188   PetscCall(VecAXPY(Q, -1.0, Q_exact));
189   PetscCall(VecNorm(Q, NORM_1, &norm_error));
190 
191   rel_error = norm_error / norm_exact;
192   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
193   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
194   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
195   PetscFunctionReturn(PETSC_SUCCESS);
196 }
197 
198 // Post-processing
199 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
200   PetscInt          steps;
201   TSConvergedReason reason;
202 
203   PetscFunctionBeginUser;
204   // Print relative error
205   if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
206     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
207   }
208 
209   // Print final time and number of steps
210   PetscCall(TSGetStepNumber(ts, &steps));
211   PetscCall(TSGetConvergedReason(ts, &reason));
212   if (user->app_ctx->test_type == TESTTYPE_NONE) {
213     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
214                           steps, (double)final_time));
215   }
216 
217   // Output numerical values from command line
218   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
219 
220   // Compare reference solution values with current test run for CI
221   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
222     PetscCall(RegressionTest(user->app_ctx, Q));
223   }
224   PetscFunctionReturn(PETSC_SUCCESS);
225 }
226 
227 const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
228 const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
229 const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
230 
231 // Gather initial Q values in case of continuation of simulation
232 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
233   PetscViewer viewer;
234 
235   PetscFunctionBeginUser;
236   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
237   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
238   PetscCall(PetscViewerDestroy(&viewer));
239   PetscFunctionReturn(PETSC_SUCCESS);
240 }
241 
242 // Free a plain data context that was allocated using PETSc; returning libCEED error codes
243 int FreeContextPetsc(void *data) {
244   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
245   return CEED_ERROR_SUCCESS;
246 }
247 
248 // Return mass qfunction specification for number of components N
249 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
250   PetscFunctionBeginUser;
251   switch (N) {
252     case 1:
253       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
254       break;
255     case 5:
256       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
257       break;
258     case 7:
259       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
260       break;
261     case 9:
262       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
263       break;
264     case 22:
265       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
266       break;
267     default:
268       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
269   }
270 
271   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
272   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
273   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
274   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
275   PetscFunctionReturn(PETSC_SUCCESS);
276 }
277 
278 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
279   PetscFunctionBeginUser;
280   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
281 
282   PetscCall(DMDestroy(&context->dm));
283   PetscCall(KSPDestroy(&context->ksp));
284 
285   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
286 
287   PetscCall(PetscFree(context));
288   PetscFunctionReturn(PETSC_SUCCESS);
289 }
290 
291 /*
292  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
293  *
294  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
295  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
296  *
297  * 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.
298  *
299  * @param[in]  comm           MPI_Comm for the program
300  * @param[in]  path           Path to the file
301  * @param[in]  char_array_len Length of the character array that should contain each line
302  * @param[out] dims           Dimensions of the file, taken from the first line of the file
303  * @param[out] fp File        pointer to the opened file
304  */
305 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
306                                  FILE **fp) {
307   int    ndims;
308   char   line[char_array_len];
309   char **array;
310 
311   PetscFunctionBeginUser;
312   PetscCall(PetscFOpen(comm, path, "r", fp));
313   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
314   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
315   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
316 
317   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
318   PetscCall(PetscStrToArrayDestroy(ndims, array));
319   PetscFunctionReturn(PETSC_SUCCESS);
320 }
321 
322 /*
323  * @brief Get the number of rows for the PHASTA file at path.
324  *
325  * 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.
326  *
327  * @param[in]  comm  MPI_Comm for the program
328  * @param[in]  path  Path to the file
329  * @param[out] nrows Number of rows
330  */
331 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
332   const PetscInt char_array_len = 512;
333   PetscInt       dims[2];
334   FILE          *fp;
335 
336   PetscFunctionBeginUser;
337   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
338   *nrows = dims[0];
339   PetscCall(PetscFClose(comm, fp));
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342 
343 PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
344   PetscInt       dims[2];
345   FILE          *fp;
346   const PetscInt char_array_len = 512;
347   char           line[char_array_len];
348 
349   PetscFunctionBeginUser;
350   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
351 
352   for (PetscInt i = 0; i < dims[0]; i++) {
353     int    ndims;
354     char **row_array;
355 
356     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
357     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
358     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
359                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
360 
361     for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
362     PetscCall(PetscStrToArrayDestroy(ndims, row_array));
363   }
364 
365   PetscCall(PetscFClose(comm, fp));
366   PetscFunctionReturn(PETSC_SUCCESS);
367 }
368 
369 // Print information about the given simulation run
370 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
371   Ceed     ceed = user->ceed;
372   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
373 
374   PetscFunctionBeginUser;
375   // Header and rank
376   char        host_name[PETSC_MAX_PATH_LEN];
377   PetscMPIInt rank, comm_size;
378   PetscCall(PetscGetHostName(host_name, sizeof host_name));
379   PetscCallMPI(MPI_Comm_rank(comm, &rank));
380   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
381   PetscCall(PetscPrintf(comm,
382                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
383                         "  MPI:\n"
384                         "    Host Name                          : %s\n"
385                         "    Total ranks                        : %d\n",
386                         host_name, comm_size));
387 
388   // Problem specific info
389   PetscCall(problem->print_info(user, problem, user->app_ctx));
390 
391   // libCEED
392   const char *used_resource;
393   CeedMemType mem_type_backend;
394   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
395   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
396   PetscCall(PetscPrintf(comm,
397                         "  libCEED:\n"
398                         "    libCEED Backend                    : %s\n"
399                         "    libCEED Backend MemType            : %s\n",
400                         used_resource, CeedMemTypes[mem_type_backend]));
401   // PETSc
402   VecType vec_type;
403   char    box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
404   if (problem->dim == 2) box_faces_str[3] = '\0';
405   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
406   PetscCall(DMGetVecType(user->dm, &vec_type));
407   PetscCall(PetscPrintf(comm,
408                         "  PETSc:\n"
409                         "    Box Faces                          : %s\n"
410                         "    DM VecType                         : %s\n"
411                         "    Time Stepping Scheme               : %s\n",
412                         box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
413   {
414     char           pmat_type_str[PETSC_MAX_PATH_LEN];
415     MatType        amat_type, pmat_type;
416     Mat            Amat, Pmat;
417     TSIJacobianFn *ijacob_function;
418 
419     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
420     PetscCall(MatGetType(Amat, &amat_type));
421     PetscCall(MatGetType(Pmat, &pmat_type));
422 
423     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
424     if (!strcmp(pmat_type, MATCEED)) {
425       MatType pmat_coo_type;
426       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
427 
428       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
429       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
430       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
431     }
432     if (ijacob_function) {
433       PetscCall(PetscPrintf(comm,
434                             "    IJacobian A MatType                : %s\n"
435                             "    IJacobian P MatType                : %s\n",
436                             amat_type, pmat_type_str));
437     }
438   }
439   if (user->app_ctx->cont_steps) {
440     PetscCall(PetscPrintf(comm,
441                           "  Continue:\n"
442                           "    Filename:                          : %s\n"
443                           "    Step:                              : %" PetscInt_FMT "\n"
444                           "    Time:                              : %g\n",
445                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
446   }
447   // Mesh
448   const PetscInt num_comp_q = 5;
449   PetscInt       glob_dofs, owned_dofs, local_dofs;
450   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
451   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
452   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
453   PetscCall(PetscPrintf(comm,
454                         "  Mesh:\n"
455                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
456                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
457                         "    Global DoFs                        : %" PetscInt_FMT "\n"
458                         "    DoFs per node                      : %" PetscInt_FMT "\n"
459                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
460                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
461   // -- Get Partition Statistics
462   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
463   {
464     PetscInt *gather_buffer = NULL;
465     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
466     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
467     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
468 
469     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
470     if (!rank) {
471       PetscCall(PetscSortInt(comm_size, gather_buffer));
472       part_owned_dofs[0]             = gather_buffer[0];              // min
473       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
474       part_owned_dofs[2]             = gather_buffer[median_index];   // median
475       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
476       PetscCall(PetscPrintf(
477           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
478           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));
479     }
480 
481     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
482     if (!rank) {
483       PetscCall(PetscSortInt(comm_size, gather_buffer));
484       part_local_dofs[0]             = gather_buffer[0];              // min
485       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
486       part_local_dofs[2]             = gather_buffer[median_index];   // median
487       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
488       PetscCall(PetscPrintf(
489           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
490           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));
491     }
492 
493     if (comm_size != 1) {
494       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
495       {
496         PetscSF            sf;
497         PetscInt           nrranks, niranks;
498         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
499         const PetscMPIInt *rranks, *iranks;
500         PetscCall(DMGetSectionSF(user->dm, &sf));
501         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
502         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
503         for (PetscInt i = 0; i < nrranks; i++) {
504           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
505           num_remote_roots_total += roffset[i + 1] - roffset[i];
506           num_ghost_interface_ranks++;
507         }
508         for (PetscInt i = 0; i < niranks; i++) {
509           if (iranks[i] == rank) continue;
510           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
511           num_owned_interface_ranks++;
512         }
513       }
514       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
515       if (!rank) {
516         PetscCall(PetscSortInt(comm_size, gather_buffer));
517         part_boundary_dofs[0]           = gather_buffer[0];              // min
518         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
519         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
520         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
521         PetscCall(PetscPrintf(
522             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
523             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
524             part_shared_dof_ratio));
525       }
526 
527       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
528       if (!rank) {
529         PetscCall(PetscSortInt(comm_size, gather_buffer));
530         part_neighbors[0]              = gather_buffer[0];              // min
531         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
532         part_neighbors[2]              = gather_buffer[median_index];   // median
533         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
534         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
535                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
536       }
537 
538       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
539       if (!rank) {
540         PetscCall(PetscSortInt(comm_size, gather_buffer));
541         part_boundary_dofs[0]           = gather_buffer[0];              // min
542         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
543         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
544         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
545         PetscCall(PetscPrintf(
546             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
547             num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
548             part_shared_dof_ratio));
549       }
550 
551       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
552       if (!rank) {
553         PetscCall(PetscSortInt(comm_size, gather_buffer));
554         part_neighbors[0]              = gather_buffer[0];              // min
555         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
556         part_neighbors[2]              = gather_buffer[median_index];   // median
557         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
558         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
559                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
560       }
561     }
562 
563     if (!rank) PetscCall(PetscFree(gather_buffer));
564   }
565   PetscFunctionReturn(PETSC_SUCCESS);
566 }
567