xref: /honee/src/misc.c (revision cf8f12c89a11a3279f97d0071ccf240f989bc25c)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Miscellaneous utility functions
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdm.h>
9926a6279SJames Wright #include <petscsf.h>
10e419654dSJeremy L Thompson #include <petscts.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
13a515125bSLeila Ghaffari 
14e3663b90SJames Wright PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time) {
150c373b74SJames Wright   Ceed                ceed = honee->ceed;
16b2948607SJames Wright   CeedVector          mult_vec;
17b2948607SJames Wright   PetscMemType        m_mem_type;
18b2948607SJames Wright   Vec                 Multiplicity, Multiplicity_loc;
19*cf8f12c8SJames Wright   CeedElemRestriction elem_restr_q;
20b2948607SJames Wright 
21a515125bSLeila Ghaffari   PetscFunctionBeginUser;
22e3663b90SJames Wright   if (honee->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ics_ctx->op, honee->phys->ics_time_label, &time));
23e3663b90SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, honee->op_ics_ctx));
24a515125bSLeila Ghaffari 
25a515125bSLeila Ghaffari   // -- Get multiplicity
26*cf8f12c8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
27*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q, &mult_vec, NULL));
28b2948607SJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
29a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
30*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q, mult_vec));
31a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
32a515125bSLeila Ghaffari 
33b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
34b2948607SJames Wright   PetscCall(VecZeroEntries(Multiplicity));
35b2948607SJames Wright   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari   // -- Fix multiplicity
38b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
39b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
40a515125bSLeila Ghaffari 
41b2948607SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
42b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
43b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
44*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
45d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
46a515125bSLeila Ghaffari }
47a515125bSLeila Ghaffari 
48c56e8d5bSJames Wright // Record boundary values from initial condition
49c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
50c56e8d5bSJames Wright   PetscFunctionBeginUser;
51313f2f1eSJames Wright   {  // Capture initial condition values in Qbc
52313f2f1eSJames Wright     Vec Qbc;
53313f2f1eSJames Wright 
54c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
55c56e8d5bSJames Wright     PetscCall(VecCopy(Q_loc, Qbc));
56c56e8d5bSJames Wright     PetscCall(VecZeroEntries(Q_loc));
57c56e8d5bSJames Wright     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
58c56e8d5bSJames Wright     PetscCall(VecAXPY(Qbc, -1., Q_loc));
59c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
60313f2f1eSJames Wright   }
61c56e8d5bSJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
62c56e8d5bSJames Wright 
63313f2f1eSJames Wright   {  // Set boundary mask to zero out essential BCs
64313f2f1eSJames Wright     Vec boundary_mask, ones;
65313f2f1eSJames Wright 
66c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
67313f2f1eSJames Wright     PetscCall(DMGetGlobalVector(dm, &ones));
68c56e8d5bSJames Wright     PetscCall(VecZeroEntries(boundary_mask));
69313f2f1eSJames Wright     PetscCall(VecSet(ones, 1.0));
70313f2f1eSJames Wright     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
71c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
72313f2f1eSJames Wright     PetscCall(DMRestoreGlobalVector(dm, &ones));
73313f2f1eSJames Wright   }
74c56e8d5bSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
75c56e8d5bSJames Wright }
76c56e8d5bSJames Wright 
77c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
782b916ea7SJeremy L Thompson                                                   Vec grad_FVM) {
799d437337SJames Wright   Vec Qbc, boundary_mask;
80a515125bSLeila Ghaffari 
8106f41313SJames Wright   PetscFunctionBeginUser;
822eb7bf1fSJames Wright   // Mask (zero) Strong BC entries
839d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
849d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
859d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
869d437337SJames Wright 
872b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
882b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
892b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
90d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
91a515125bSLeila Ghaffari }
92a515125bSLeila Ghaffari 
93a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
94c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
95a515125bSLeila Ghaffari   Vec         Qref;
96a515125bSLeila Ghaffari   PetscViewer viewer;
97a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
98e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
99a515125bSLeila Ghaffari 
10006f41313SJames Wright   PetscFunctionBeginUser;
101a515125bSLeila Ghaffari   // Read reference file
1022b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
1034c07ec22SJames Wright   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
104e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
105360b37c9SJames Wright   PetscCall(HoneeLoadBinaryVec(viewer, Qref, NULL, NULL));
106a515125bSLeila Ghaffari 
107a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1082b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1092b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1102b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1112b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
112a515125bSLeila Ghaffari 
113a515125bSLeila Ghaffari   // Check error
114a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1152b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
116a515125bSLeila Ghaffari   }
117a515125bSLeila Ghaffari 
118a515125bSLeila Ghaffari   // Cleanup
1192b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1202b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
121d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
122a515125bSLeila Ghaffari }
123a515125bSLeila Ghaffari 
124a515125bSLeila Ghaffari // Get error for problems with exact solutions
125e3663b90SJames Wright PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time) {
126a515125bSLeila Ghaffari   PetscInt  loc_nodes;
127a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
128a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
129a515125bSLeila Ghaffari 
13006f41313SJames Wright   PetscFunctionBeginUser;
131a515125bSLeila Ghaffari   // Get exact solution at final time
132b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1332b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1342b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
135e3663b90SJames Wright   PetscCall(ICs_FixMultiplicity(dm, honee, Q_exact_loc, Q_exact, final_time));
136a515125bSLeila Ghaffari 
137a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1382b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1392b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1402b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
141a515125bSLeila Ghaffari 
142a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
1432b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1442b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
145b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
146d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
147a515125bSLeila Ghaffari }
148a515125bSLeila Ghaffari 
149a515125bSLeila Ghaffari // Post-processing
150e3663b90SJames Wright PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time) {
151a515125bSLeila Ghaffari   PetscInt          steps;
152f0784ed3SJed Brown   TSConvergedReason reason;
153a515125bSLeila Ghaffari 
15406f41313SJames Wright   PetscFunctionBeginUser;
155a515125bSLeila Ghaffari   // Print relative error
1560c373b74SJames Wright   if (problem->compute_exact_solution_error && honee->app_ctx->test_type == TESTTYPE_NONE) {
157e3663b90SJames Wright     PetscCall(PrintError(dm, honee, Q, final_time));
158a515125bSLeila Ghaffari   }
159a515125bSLeila Ghaffari 
160a515125bSLeila Ghaffari   // Print final time and number of steps
1612b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
162f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1630c373b74SJames Wright   if (honee->app_ctx->test_type == TESTTYPE_NONE) {
164f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
165f0784ed3SJed Brown                           steps, (double)final_time));
166a515125bSLeila Ghaffari   }
167a515125bSLeila Ghaffari 
168a515125bSLeila Ghaffari   // Output numerical values from command line
1692b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
170a515125bSLeila Ghaffari 
171a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
1720c373b74SJames Wright   if (honee->app_ctx->test_type == TESTTYPE_SOLVER) {
1730c373b74SJames Wright     PetscCall(RegressionTest(honee->app_ctx, Q));
174a515125bSLeila Ghaffari   }
175d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
176a515125bSLeila Ghaffari }
177a515125bSLeila Ghaffari 
17815a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
17915a3537eSJed Brown int FreeContextPetsc(void *data) {
1802b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
18115a3537eSJed Brown   return CEED_ERROR_SUCCESS;
18215a3537eSJed Brown }
1839f59f36eSJames Wright 
1844ea616f4SJames Wright /**
1854ea616f4SJames Wright   @brief Destroy `NodalProjectionData` object
1864ea616f4SJames Wright 
1874ea616f4SJames Wright   @param[in] context `NodalProjectionData` object to destroy
1884ea616f4SJames Wright **/
1894ea616f4SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData *context) {
1904ea616f4SJames Wright   NodalProjectionData context_ = *context;
19114bd2a07SJames Wright 
192457a5831SJames Wright   PetscFunctionBeginUser;
1934ea616f4SJames Wright   if (context_ == NULL) PetscFunctionReturn(PETSC_SUCCESS);
1944ea616f4SJames Wright   PetscCall(DMDestroy(&context_->dm));
1954ea616f4SJames Wright   PetscCall(KSPDestroy(&context_->ksp));
1964ea616f4SJames Wright   PetscCall(OperatorApplyContextDestroy(context_->l2_rhs_ctx));
1974ea616f4SJames Wright   PetscCall(PetscFree(context_));
1984ea616f4SJames Wright   *context = NULL;
199d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
200457a5831SJames Wright }
201457a5831SJames Wright 
20259572072SJames Wright /**
20359572072SJames Wright    @brief Sets the value of an option if it is not already set
20459572072SJames Wright 
20559572072SJames Wright    @param[in,out] options `PetscOptions` database, or `NULL` for default global database
20659572072SJames Wright    @param[in]     name    Name of the option (with `-` prepended to it)
20759572072SJames Wright    @param[in]     value   Value to set the option to
20859572072SJames Wright **/
20959572072SJames Wright PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]) {
21059572072SJames Wright   PetscBool has_option;
21159572072SJames Wright 
21259572072SJames Wright   PetscFunctionBeginUser;
21359572072SJames Wright   PetscCall(PetscOptionsHasName(options, NULL, name, &has_option));
21459572072SJames Wright   if (!has_option) PetscCall(PetscOptionsSetValue(options, name, value));
21559572072SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21659572072SJames Wright }
21759572072SJames Wright 
21825125139SJames Wright PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume) {
21925125139SJames Wright   DM                   dm   = honee->dm, dm_coord;
22025125139SJames Wright   Ceed                 ceed = honee->ceed;
22125125139SJames Wright   CeedQFunction        qf_mass;
22225125139SJames Wright   CeedOperator         op_mass;
22325125139SJames Wright   OperatorApplyContext op_mass_ctx;
22425125139SJames Wright   CeedElemRestriction  elem_restr_qd;
22525125139SJames Wright   CeedVector           qdata;
22625125139SJames Wright   CeedInt              q_data_size, num_comps_x;
22725125139SJames Wright   Vec                  u, v;
22825125139SJames Wright 
22925125139SJames Wright   PetscFunctionBeginUser;
23025125139SJames Wright   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
231e3db12f8SJames Wright   PetscCall(QDataGet(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &qdata,
232e3db12f8SJames Wright                      &q_data_size));
23325125139SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comps_x));
23425125139SJames Wright 
23525125139SJames Wright   PetscCall(HoneeMassQFunctionCreate(ceed, num_comps_x, q_data_size, &qf_mass));
23625125139SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
23725125139SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE));
23825125139SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, qdata));
23925125139SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE));
24025125139SJames Wright 
24125125139SJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm_coord, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx));
24225125139SJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_mass, DMReturnVecType(dm), PETSC_COMM_SELF, &u, NULL));
24325125139SJames Wright   PetscCall(DMCreateGlobalVector(dm_coord, &v));
24425125139SJames Wright   PetscCall(VecSet(u, 1.));
24525125139SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(u, v, op_mass_ctx));
24625125139SJames Wright   PetscCall(VecSum(v, volume));
24725125139SJames Wright   *volume /= num_comps_x;  // Correct for number of components != 1
24825125139SJames Wright 
24925125139SJames Wright   PetscCall(VecDestroy(&u));
25025125139SJames Wright   PetscCall(VecDestroy(&v));
25125125139SJames Wright   PetscCall(OperatorApplyContextDestroy(op_mass_ctx));
25225125139SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
25325125139SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&qdata));
25425125139SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
25525125139SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
25625125139SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25725125139SJames Wright }
25825125139SJames Wright 
259926a6279SJames Wright // Print information about the given simulation run
2600c373b74SJames Wright PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts) {
2610c373b74SJames Wright   Ceed     ceed = honee->ceed;
2623678fae3SJames Wright   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
2633678fae3SJames Wright 
264926a6279SJames Wright   PetscFunctionBeginUser;
265926a6279SJames Wright   // Header and rank
266926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
267926a6279SJames Wright   PetscMPIInt rank, comm_size;
268926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
269926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
270926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
271926a6279SJames Wright   PetscCall(PetscPrintf(comm,
27215747f0fSJames Wright                         "\n-- HONEE - High-Order Navier-stokes Equation Evaluator --\n"
273926a6279SJames Wright                         "  MPI:\n"
274926a6279SJames Wright                         "    Host Name                          : %s\n"
275926a6279SJames Wright                         "    Total ranks                        : %d\n",
276926a6279SJames Wright                         host_name, comm_size));
277926a6279SJames Wright 
278926a6279SJames Wright   // Problem specific info
2790c373b74SJames Wright   PetscCall(problem->print_info(honee, problem, honee->app_ctx));
280926a6279SJames Wright 
281926a6279SJames Wright   // libCEED
282926a6279SJames Wright   const char *used_resource;
283926a6279SJames Wright   CeedMemType mem_type_backend;
2840c373b74SJames Wright   PetscCallCeed(ceed, CeedGetResource(honee->ceed, &used_resource));
2850c373b74SJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(honee->ceed, &mem_type_backend));
286926a6279SJames Wright   PetscCall(PetscPrintf(comm,
287926a6279SJames Wright                         "  libCEED:\n"
288926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
289926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
290926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
291926a6279SJames Wright   // PETSc
29266d54740SJames Wright   {
2933678fae3SJames Wright     VecType  vec_type;
294926a6279SJames Wright     char     box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
29566d54740SJames Wright     PetscInt dim;
29666d54740SJames Wright 
2970c373b74SJames Wright     PetscCall(DMGetDimension(honee->dm, &dim));
29866d54740SJames Wright     if (dim == 2) box_faces_str[3] = '\0';
299926a6279SJames Wright     PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
3000c373b74SJames Wright     PetscCall(DMGetVecType(honee->dm, &vec_type));
301926a6279SJames Wright     PetscCall(PetscPrintf(comm,
302926a6279SJames Wright                           "  PETSc:\n"
303926a6279SJames Wright                           "    Box Faces                          : %s\n"
304926a6279SJames Wright                           "    DM VecType                         : %s\n"
305926a6279SJames Wright                           "    Time Stepping Scheme               : %s\n",
3063678fae3SJames Wright                           box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
30766d54740SJames Wright   }
3083678fae3SJames Wright   {
3093678fae3SJames Wright     char           pmat_type_str[PETSC_MAX_PATH_LEN];
3103678fae3SJames Wright     MatType        amat_type, pmat_type;
3113678fae3SJames Wright     Mat            Amat, Pmat;
3123678fae3SJames Wright     TSIJacobianFn *ijacob_function;
3133678fae3SJames Wright 
3143678fae3SJames Wright     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
3153678fae3SJames Wright     PetscCall(MatGetType(Amat, &amat_type));
3163678fae3SJames Wright     PetscCall(MatGetType(Pmat, &pmat_type));
3173678fae3SJames Wright 
3183678fae3SJames Wright     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
3193678fae3SJames Wright     if (!strcmp(pmat_type, MATCEED)) {
3203678fae3SJames Wright       MatType pmat_coo_type;
3213678fae3SJames Wright       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
3223678fae3SJames Wright 
3233678fae3SJames Wright       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
3243678fae3SJames Wright       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
3253678fae3SJames Wright       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
3263678fae3SJames Wright     }
3273678fae3SJames Wright     if (ijacob_function) {
3283678fae3SJames Wright       PetscCall(PetscPrintf(comm,
3293678fae3SJames Wright                             "    IJacobian A MatType                : %s\n"
3303678fae3SJames Wright                             "    IJacobian P MatType                : %s\n",
3313678fae3SJames Wright                             amat_type, pmat_type_str));
3323678fae3SJames Wright     }
3333678fae3SJames Wright   }
334481d14cbSJames Wright   if (honee->app_ctx->use_continue_file) {
335926a6279SJames Wright     PetscCall(PetscPrintf(comm,
336926a6279SJames Wright                           "  Continue:\n"
337926a6279SJames Wright                           "    Filename:                          : %s\n"
338926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
339926a6279SJames Wright                           "    Time:                              : %g\n",
3400c373b74SJames Wright                           honee->app_ctx->cont_file, honee->app_ctx->cont_steps, honee->app_ctx->cont_time));
341926a6279SJames Wright   }
342926a6279SJames Wright   // Mesh
343926a6279SJames Wright   const PetscInt num_comp_q = 5;
344926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
3450c373b74SJames Wright   const CeedInt  num_P = honee->app_ctx->degree + 1, num_Q = num_P + honee->app_ctx->q_extra;
3460c373b74SJames Wright   PetscCall(DMGetGlobalVectorInfo(honee->dm, &owned_dofs, &glob_dofs, NULL));
3470c373b74SJames Wright   PetscCall(DMGetLocalVectorInfo(honee->dm, &local_dofs, NULL, NULL));
348926a6279SJames Wright   PetscCall(PetscPrintf(comm,
349926a6279SJames Wright                         "  Mesh:\n"
350926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
351926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
352926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
353926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
354dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
355dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
356926a6279SJames Wright   // -- Get Partition Statistics
357926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
358926a6279SJames Wright   {
359926a6279SJames Wright     PetscInt *gather_buffer = NULL;
360dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
361926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
362926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
363926a6279SJames Wright 
364dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
365926a6279SJames Wright     if (!rank) {
366926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
367926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
368926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
369926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
370926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
3719eadbee4SJames Wright       PetscCall(PetscPrintf(comm,
3729eadbee4SJames Wright                             "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
3739eadbee4SJames Wright                             num_comp_q, part_owned_dofs[0] / num_comp_q, part_owned_dofs[1] / num_comp_q, part_owned_dofs[2] / num_comp_q,
3749eadbee4SJames Wright                             part_owned_dof_ratio));
375926a6279SJames Wright     }
376926a6279SJames Wright 
377dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
378926a6279SJames Wright     if (!rank) {
379926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
380926a6279SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
381926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
382926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
383926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
3849eadbee4SJames Wright       PetscCall(PetscPrintf(comm,
3859eadbee4SJames Wright                             "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
3869eadbee4SJames Wright                             num_comp_q, part_local_dofs[0] / num_comp_q, part_local_dofs[1] / num_comp_q, part_local_dofs[2] / num_comp_q,
3879eadbee4SJames Wright                             part_local_dof_ratio));
388926a6279SJames Wright     }
389926a6279SJames Wright 
39045abf96eSJames Wright     if (comm_size != 1) {
391dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
392926a6279SJames Wright       {
393926a6279SJames Wright         PetscSF            sf;
3943bbbcc04SJames Wright         PetscMPIInt        nrranks, niranks;
395dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
396dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
3970c373b74SJames Wright         PetscCall(DMGetSectionSF(honee->dm, &sf));
398926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
399dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
400926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
401926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
402926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
403dfeb939dSJames Wright           num_ghost_interface_ranks++;
404dfeb939dSJames Wright         }
405dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
406dfeb939dSJames Wright           if (iranks[i] == rank) continue;
407dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
408dfeb939dSJames Wright           num_owned_interface_ranks++;
409926a6279SJames Wright         }
410926a6279SJames Wright       }
411dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
412926a6279SJames Wright       if (!rank) {
413926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
414dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
415dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
416dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
417dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
4189eadbee4SJames Wright         PetscCall(PetscPrintf(comm,
4199eadbee4SJames Wright                               "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT
4209eadbee4SJames Wright                               ", %f\n",
42145abf96eSJames Wright                               num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
42245abf96eSJames Wright                               part_shared_dof_ratio));
423dfeb939dSJames Wright       }
424dfeb939dSJames Wright 
425dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
426dfeb939dSJames Wright       if (!rank) {
427dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
428dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
429dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
430dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
431dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
432dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
433dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
434dfeb939dSJames Wright       }
435dfeb939dSJames Wright 
436dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
437dfeb939dSJames Wright       if (!rank) {
438dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
439dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
440dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
441dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
442dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
4439eadbee4SJames Wright         PetscCall(PetscPrintf(comm,
4449eadbee4SJames Wright                               "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT
4459eadbee4SJames Wright                               ", %f\n",
44645abf96eSJames Wright                               num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q,
44745abf96eSJames Wright                               part_shared_dof_ratio));
448dfeb939dSJames Wright       }
449dfeb939dSJames Wright 
450dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
451dfeb939dSJames Wright       if (!rank) {
452dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
453dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
454dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
455dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
456dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
457dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
458dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
459926a6279SJames Wright       }
46045abf96eSJames Wright     }
461926a6279SJames Wright 
462926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
463926a6279SJames Wright   }
464926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
465926a6279SJames Wright }
466