xref: /libCEED/examples/fluids/src/misc.c (revision 171d97d0be1a6a7597730b8cc2bec09fb6a98aac)
1 // Copyright (c) 2017-2022, 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   Vec Qbc, boundary_mask;
54 
55   PetscFunctionBeginUser;
56   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
57   PetscCall(VecCopy(Q_loc, Qbc));
58   PetscCall(VecZeroEntries(Q_loc));
59   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
60   PetscCall(VecAXPY(Qbc, -1., Q_loc));
61   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
62   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
63 
64   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
65   PetscCall(DMGetGlobalVector(dm, &Q));
66   PetscCall(VecZeroEntries(boundary_mask));
67   PetscCall(VecSet(Q, 1.0));
68   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
69   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
70   PetscFunctionReturn(PETSC_SUCCESS);
71 }
72 
73 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
74                                                   Vec grad_FVM) {
75   Vec Qbc, boundary_mask;
76 
77   PetscFunctionBeginUser;
78   // Mask (zero) Strong BC entries
79   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
80   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
81   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
82 
83   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
84   PetscCall(VecAXPY(Q_loc, 1., Qbc));
85   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
86   PetscFunctionReturn(PETSC_SUCCESS);
87 }
88 
89 static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
90   PetscFunctionBeginUser;
91   if (file_type == PETSC_INT32) {
92     PetscInt32 val;
93     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
94     *out = val;
95   } else if (file_type == PETSC_INT64) {
96     PetscInt64 val;
97     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
98     *out = val;
99   } else {
100     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
101   }
102   PetscFunctionReturn(PETSC_SUCCESS);
103 }
104 
105 // @brief Load vector from binary file, possibly with embedded solution time and step number
106 PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
107   PetscInt      file_step_number;
108   PetscInt32    token;
109   PetscReal     file_time;
110   PetscDataType file_type = PETSC_INT32;
111 
112   PetscFunctionBeginUser;
113   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
114   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
115       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
116     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
117     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
118     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
119     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
120     if (time) *time = file_time;
121     if (step_number) *step_number = file_step_number;
122   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
123     PetscInt length, N;
124     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
125     PetscCall(VecGetSize(Q, &N));
126     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
127     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
128   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
129 
130   PetscCall(VecLoad(Q, viewer));
131   PetscFunctionReturn(PETSC_SUCCESS);
132 }
133 
134 // Compare reference solution values with current test run for CI
135 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
136   Vec         Qref;
137   PetscViewer viewer;
138   PetscReal   error, Qrefnorm;
139   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
140 
141   PetscFunctionBeginUser;
142   // Read reference file
143   PetscCall(VecDuplicate(Q, &Qref));
144   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
145   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
146 
147   // Compute error with respect to reference solution
148   PetscCall(VecAXPY(Q, -1.0, Qref));
149   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
150   PetscCall(VecScale(Q, 1. / Qrefnorm));
151   PetscCall(VecNorm(Q, NORM_MAX, &error));
152 
153   // Check error
154   if (error > app_ctx->test_tol) {
155     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
156   }
157 
158   // Cleanup
159   PetscCall(PetscViewerDestroy(&viewer));
160   PetscCall(VecDestroy(&Qref));
161   PetscFunctionReturn(PETSC_SUCCESS);
162 }
163 
164 // Get error for problems with exact solutions
165 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
166   PetscInt  loc_nodes;
167   Vec       Q_exact, Q_exact_loc;
168   PetscReal rel_error, norm_error, norm_exact;
169 
170   PetscFunctionBeginUser;
171   // Get exact solution at final time
172   PetscCall(DMGetGlobalVector(dm, &Q_exact));
173   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
174   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
175   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
176 
177   // Get |exact solution - obtained solution|
178   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
179   PetscCall(VecAXPY(Q, -1.0, Q_exact));
180   PetscCall(VecNorm(Q, NORM_1, &norm_error));
181 
182   rel_error = norm_error / norm_exact;
183   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
184   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
185   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
186   PetscFunctionReturn(PETSC_SUCCESS);
187 }
188 
189 // Post-processing
190 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
191   PetscInt          steps;
192   TSConvergedReason reason;
193 
194   PetscFunctionBeginUser;
195   // Print relative error
196   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
197     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
198   }
199 
200   // Print final time and number of steps
201   PetscCall(TSGetStepNumber(ts, &steps));
202   PetscCall(TSGetConvergedReason(ts, &reason));
203   if (user->app_ctx->test_type == TESTTYPE_NONE) {
204     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
205                           steps, (double)final_time));
206   }
207 
208   // Output numerical values from command line
209   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
210 
211   // Compare reference solution values with current test run for CI
212   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
213     PetscCall(RegressionTest(user->app_ctx, Q));
214   }
215   PetscFunctionReturn(PETSC_SUCCESS);
216 }
217 
218 const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
219 const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
220 const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
221 
222 // Gather initial Q values in case of continuation of simulation
223 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
224   PetscViewer viewer;
225 
226   PetscFunctionBeginUser;
227   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
228   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
229   PetscCall(PetscViewerDestroy(&viewer));
230   PetscFunctionReturn(PETSC_SUCCESS);
231 }
232 
233 // Free a plain data context that was allocated using PETSc; returning libCEED error codes
234 int FreeContextPetsc(void *data) {
235   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
236   return CEED_ERROR_SUCCESS;
237 }
238 
239 // Return mass qfunction specification for number of components N
240 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
241   PetscFunctionBeginUser;
242   switch (N) {
243     case 1:
244       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
245       break;
246     case 5:
247       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
248       break;
249     case 7:
250       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
251       break;
252     case 9:
253       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
254       break;
255     case 22:
256       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
257       break;
258     default:
259       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
260   }
261 
262   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
263   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
264   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
265   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
266   PetscFunctionReturn(PETSC_SUCCESS);
267 }
268 
269 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
270   PetscFunctionBeginUser;
271   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
272 
273   PetscCall(DMDestroy(&context->dm));
274   PetscCall(KSPDestroy(&context->ksp));
275 
276   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
277 
278   PetscCall(PetscFree(context));
279   PetscFunctionReturn(PETSC_SUCCESS);
280 }
281 
282 /*
283  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
284  *
285  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
286  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
287  *
288  * 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.
289  *
290  * @param[in]  comm           MPI_Comm for the program
291  * @param[in]  path           Path to the file
292  * @param[in]  char_array_len Length of the character array that should contain each line
293  * @param[out] dims           Dimensions of the file, taken from the first line of the file
294  * @param[out] fp File        pointer to the opened file
295  */
296 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
297                                  FILE **fp) {
298   int    ndims;
299   char   line[char_array_len];
300   char **array;
301 
302   PetscFunctionBeginUser;
303   PetscCall(PetscFOpen(comm, path, "r", fp));
304   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
305   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
306   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
307 
308   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
309   PetscCall(PetscStrToArrayDestroy(ndims, array));
310   PetscFunctionReturn(PETSC_SUCCESS);
311 }
312 
313 /*
314  * @brief Get the number of rows for the PHASTA file at path.
315  *
316  * 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.
317  *
318  * @param[in]  comm  MPI_Comm for the program
319  * @param[in]  path  Path to the file
320  * @param[out] nrows Number of rows
321  */
322 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
323   const PetscInt char_array_len = 512;
324   PetscInt       dims[2];
325   FILE          *fp;
326 
327   PetscFunctionBeginUser;
328   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
329   *nrows = dims[0];
330   PetscCall(PetscFClose(comm, fp));
331   PetscFunctionReturn(PETSC_SUCCESS);
332 }
333 
334 PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
335   PetscInt       dims[2];
336   int            ndims;
337   FILE          *fp;
338   const PetscInt char_array_len = 512;
339   char           line[char_array_len];
340   char         **row_array;
341 
342   PetscFunctionBeginUser;
343   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
344 
345   for (PetscInt i = 0; i < dims[0]; i++) {
346     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
347     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
348     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
349                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
350 
351     for (PetscInt j = 0; j < dims[1]; j++) {
352       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
353     }
354   }
355 
356   PetscCall(PetscFClose(comm, fp));
357   PetscFunctionReturn(PETSC_SUCCESS);
358 }
359 
360 PetscLogEvent       FLUIDS_CeedOperatorApply;
361 PetscLogEvent       FLUIDS_CeedOperatorAssemble;
362 PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
363 PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
364 PetscLogEvent       FLUIDS_SmartRedis_Init;
365 PetscLogEvent       FLUIDS_SmartRedis_Meta;
366 PetscLogEvent       FLUIDS_SmartRedis_Train;
367 PetscLogEvent       FLUIDS_TrainDataCompute;
368 PetscLogEvent       FLUIDS_DifferentialFilter;
369 PetscLogEvent       FLUIDS_VelocityGradientProjection;
370 static PetscClassId libCEED_classid, onlineTrain_classid, misc_classid;
371 
372 PetscErrorCode RegisterLogEvents() {
373   PetscFunctionBeginUser;
374   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
375   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
376   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
377   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
378   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
379 
380   PetscCall(PetscClassIdRegister("onlineTrain", &onlineTrain_classid));
381   PetscCall(PetscLogEventRegister("SmartRedis_Init", onlineTrain_classid, &FLUIDS_SmartRedis_Init));
382   PetscCall(PetscLogEventRegister("SmartRedis_Meta", onlineTrain_classid, &FLUIDS_SmartRedis_Meta));
383   PetscCall(PetscLogEventRegister("SmartRedis_Train", onlineTrain_classid, &FLUIDS_SmartRedis_Train));
384   PetscCall(PetscLogEventRegister("TrainDataCompute", onlineTrain_classid, &FLUIDS_TrainDataCompute));
385 
386   PetscCall(PetscClassIdRegister("Miscellaneous", &misc_classid));
387   PetscCall(PetscLogEventRegister("DiffFilter", misc_classid, &FLUIDS_DifferentialFilter));
388   PetscCall(PetscLogEventRegister("VeloGradProj", misc_classid, &FLUIDS_VelocityGradientProjection));
389   PetscFunctionReturn(PETSC_SUCCESS);
390 }
391 
392 // Print information about the given simulation run
393 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
394   Ceed ceed = user->ceed;
395   PetscFunctionBeginUser;
396   // Header and rank
397   char        host_name[PETSC_MAX_PATH_LEN];
398   PetscMPIInt rank, comm_size;
399   PetscCall(PetscGetHostName(host_name, sizeof host_name));
400   PetscCallMPI(MPI_Comm_rank(comm, &rank));
401   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
402   PetscCall(PetscPrintf(comm,
403                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
404                         "  MPI:\n"
405                         "    Host Name                          : %s\n"
406                         "    Total ranks                        : %d\n",
407                         host_name, comm_size));
408 
409   // Problem specific info
410   PetscCall(problem->print_info(user, problem, user->app_ctx));
411 
412   // libCEED
413   const char *used_resource;
414   CeedMemType mem_type_backend;
415   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
416   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
417   PetscCall(PetscPrintf(comm,
418                         "  libCEED:\n"
419                         "    libCEED Backend                    : %s\n"
420                         "    libCEED Backend MemType            : %s\n",
421                         used_resource, CeedMemTypes[mem_type_backend]));
422   // PETSc
423   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
424   if (problem->dim == 2) box_faces_str[3] = '\0';
425   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
426   MatType amat_type = user->app_ctx->amat_type, pmat_type;
427   VecType vec_type;
428   PetscCall(DMGetMatType(user->dm, &pmat_type));
429   if (!amat_type) amat_type = pmat_type;
430   PetscCall(DMGetVecType(user->dm, &vec_type));
431   PetscCall(PetscPrintf(comm,
432                         "  PETSc:\n"
433                         "    Box Faces                          : %s\n"
434                         "    A MatType                          : %s\n"
435                         "    P MatType                          : %s\n"
436                         "    DM VecType                         : %s\n"
437                         "    Time Stepping Scheme               : %s\n",
438                         box_faces_str, amat_type, pmat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
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