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