xref: /honee/src/misc.c (revision 8c85b83518484fc9eaa5bccfe38a7cce91b34691)
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>
139f59f36eSJames Wright #include "../qfunctions/mass.h"
14a515125bSLeila Ghaffari 
152b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
16b4c37c5cSJames Wright   Ceed         ceed = user->ceed;
17b2948607SJames Wright   CeedVector   mult_vec;
18b2948607SJames Wright   PetscMemType m_mem_type;
19b2948607SJames Wright   Vec          Multiplicity, Multiplicity_loc;
20b2948607SJames Wright 
21a515125bSLeila Ghaffari   PetscFunctionBeginUser;
22b4c37c5cSJames Wright   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
238f18bb8bSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
24a515125bSLeila Ghaffari 
25b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
26a515125bSLeila Ghaffari 
27a515125bSLeila Ghaffari   // -- Get multiplicity
28b2948607SJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
29a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
30b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->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));
44d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
45a515125bSLeila Ghaffari }
46a515125bSLeila Ghaffari 
47c56e8d5bSJames Wright // Record boundary values from initial condition
48c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
49c56e8d5bSJames Wright   PetscFunctionBeginUser;
50313f2f1eSJames Wright   {  // Capture initial condition values in Qbc
51313f2f1eSJames Wright     Vec Qbc;
52313f2f1eSJames Wright 
53c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
54c56e8d5bSJames Wright     PetscCall(VecCopy(Q_loc, Qbc));
55c56e8d5bSJames Wright     PetscCall(VecZeroEntries(Q_loc));
56c56e8d5bSJames Wright     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
57c56e8d5bSJames Wright     PetscCall(VecAXPY(Qbc, -1., Q_loc));
58c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
59313f2f1eSJames Wright   }
60c56e8d5bSJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
61c56e8d5bSJames Wright 
62313f2f1eSJames Wright   {  // Set boundary mask to zero out essential BCs
63313f2f1eSJames Wright     Vec boundary_mask, ones;
64313f2f1eSJames Wright 
65c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
66313f2f1eSJames Wright     PetscCall(DMGetGlobalVector(dm, &ones));
67c56e8d5bSJames Wright     PetscCall(VecZeroEntries(boundary_mask));
68313f2f1eSJames Wright     PetscCall(VecSet(ones, 1.0));
69313f2f1eSJames Wright     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
70c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
71313f2f1eSJames Wright     PetscCall(DMRestoreGlobalVector(dm, &ones));
72313f2f1eSJames Wright   }
73c56e8d5bSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
74c56e8d5bSJames Wright }
75c56e8d5bSJames Wright 
76c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
772b916ea7SJeremy L Thompson                                                   Vec grad_FVM) {
789d437337SJames Wright   Vec Qbc, boundary_mask;
79a515125bSLeila Ghaffari 
8006f41313SJames Wright   PetscFunctionBeginUser;
812eb7bf1fSJames Wright   // Mask (zero) Strong BC entries
829d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
839d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
849d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
859d437337SJames Wright 
862b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
872b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
882b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
89d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
90a515125bSLeila Ghaffari }
91a515125bSLeila Ghaffari 
92a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
93c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
94a515125bSLeila Ghaffari   Vec         Qref;
95a515125bSLeila Ghaffari   PetscViewer viewer;
96a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
97e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
98a515125bSLeila Ghaffari 
9906f41313SJames Wright   PetscFunctionBeginUser;
100a515125bSLeila Ghaffari   // Read reference file
1012b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
1024c07ec22SJames Wright   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
103e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
104360b37c9SJames Wright   PetscCall(HoneeLoadBinaryVec(viewer, Qref, NULL, NULL));
105a515125bSLeila Ghaffari 
106a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1072b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1082b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1092b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1102b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   // Check error
113a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1142b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
115a515125bSLeila Ghaffari   }
116a515125bSLeila Ghaffari 
117a515125bSLeila Ghaffari   // Cleanup
1182b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1192b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
120d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
121a515125bSLeila Ghaffari }
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari // Get error for problems with exact solutions
124c56e8d5bSJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
125a515125bSLeila Ghaffari   PetscInt  loc_nodes;
126a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
127a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
128a515125bSLeila Ghaffari 
12906f41313SJames Wright   PetscFunctionBeginUser;
130a515125bSLeila Ghaffari   // Get exact solution at final time
131b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1322b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1332b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1342b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
135a515125bSLeila Ghaffari 
136a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1372b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1382b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1392b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
140a515125bSLeila Ghaffari 
141a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
1422b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1432b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
144b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
145d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
146a515125bSLeila Ghaffari }
147a515125bSLeila Ghaffari 
148a515125bSLeila Ghaffari // Post-processing
149991aef52SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
150a515125bSLeila Ghaffari   PetscInt          steps;
151f0784ed3SJed Brown   TSConvergedReason reason;
152a515125bSLeila Ghaffari 
15306f41313SJames Wright   PetscFunctionBeginUser;
154a515125bSLeila Ghaffari   // Print relative error
15558ce1233SJames Wright   if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
156c56e8d5bSJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
157a515125bSLeila Ghaffari   }
158a515125bSLeila Ghaffari 
159a515125bSLeila Ghaffari   // Print final time and number of steps
1602b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
161f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1620e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
163f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
164f0784ed3SJed Brown                           steps, (double)final_time));
165a515125bSLeila Ghaffari   }
166a515125bSLeila Ghaffari 
167a515125bSLeila Ghaffari   // Output numerical values from command line
1682b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
169a515125bSLeila Ghaffari 
170a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
1710e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
172c56e8d5bSJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
173a515125bSLeila Ghaffari   }
174d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
175a515125bSLeila Ghaffari }
176a515125bSLeila Ghaffari 
17715a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
17815a3537eSJed Brown int FreeContextPetsc(void *data) {
1792b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
18015a3537eSJed Brown   return CEED_ERROR_SUCCESS;
18115a3537eSJed Brown }
1829f59f36eSJames Wright 
1834cb09fddSJames Wright // @brief Return mass qfunction specification for number of components N
1849f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
1859f59f36eSJames Wright   PetscFunctionBeginUser;
1869f59f36eSJames Wright   switch (N) {
1879f59f36eSJames Wright     case 1:
188b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
1899f59f36eSJames Wright       break;
190*8c85b835SJames Wright     case 4:
191*8c85b835SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 4, Mass_4, Mass_4_loc, qf));
192*8c85b835SJames Wright       break;
1939f59f36eSJames Wright     case 5:
194b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
1959f59f36eSJames Wright       break;
196c38c977aSJames Wright     case 7:
197b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
198c38c977aSJames Wright       break;
1999f59f36eSJames Wright     case 9:
200b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
2019f59f36eSJames Wright       break;
202*8c85b835SJames Wright     case 12:
203*8c85b835SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_12, Mass_12_loc, qf));
204*8c85b835SJames Wright       break;
2059f59f36eSJames Wright     case 22:
206b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
2079f59f36eSJames Wright       break;
2089f59f36eSJames Wright     default:
2096f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2109f59f36eSJames Wright   }
2119f59f36eSJames Wright 
212b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
213b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
214b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2153170c09fSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
216d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2179f59f36eSJames Wright }
218e5e81594SJames Wright 
219457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
220457a5831SJames Wright   PetscFunctionBeginUser;
221d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
222457a5831SJames Wright 
223457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
224457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
225457a5831SJames Wright 
226457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
227457a5831SJames Wright 
228457a5831SJames Wright   PetscCall(PetscFree(context));
229d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
230457a5831SJames Wright }
231457a5831SJames Wright 
23259572072SJames Wright /**
23359572072SJames Wright    @brief Sets the value of an option if it is not already set
23459572072SJames Wright 
23559572072SJames Wright    @param[in,out] options `PetscOptions` database, or `NULL` for default global database
23659572072SJames Wright    @param[in]     name    Name of the option (with `-` prepended to it)
23759572072SJames Wright    @param[in]     value   Value to set the option to
23859572072SJames Wright **/
23959572072SJames Wright PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]) {
24059572072SJames Wright   PetscBool has_option;
24159572072SJames Wright 
24259572072SJames Wright   PetscFunctionBeginUser;
24359572072SJames Wright   PetscCall(PetscOptionsHasName(options, NULL, name, &has_option));
24459572072SJames Wright   if (!has_option) PetscCall(PetscOptionsSetValue(options, name, value));
24559572072SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
24659572072SJames Wright }
24759572072SJames Wright 
248926a6279SJames Wright // Print information about the given simulation run
2493678fae3SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
250b4c37c5cSJames Wright   Ceed     ceed = user->ceed;
2513678fae3SJames Wright   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
2523678fae3SJames Wright 
253926a6279SJames Wright   PetscFunctionBeginUser;
254926a6279SJames Wright   // Header and rank
255926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
256926a6279SJames Wright   PetscMPIInt rank, comm_size;
257926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
258926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
259926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
260926a6279SJames Wright   PetscCall(PetscPrintf(comm,
261926a6279SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
262926a6279SJames Wright                         "  MPI:\n"
263926a6279SJames Wright                         "    Host Name                          : %s\n"
264926a6279SJames Wright                         "    Total ranks                        : %d\n",
265926a6279SJames Wright                         host_name, comm_size));
266926a6279SJames Wright 
267926a6279SJames Wright   // Problem specific info
2682d49c0afSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
269926a6279SJames Wright 
270926a6279SJames Wright   // libCEED
271926a6279SJames Wright   const char *used_resource;
272926a6279SJames Wright   CeedMemType mem_type_backend;
273b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
274b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
275926a6279SJames Wright   PetscCall(PetscPrintf(comm,
276926a6279SJames Wright                         "  libCEED:\n"
277926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
278926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
279926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
280926a6279SJames Wright   // PETSc
28166d54740SJames Wright   {
2823678fae3SJames Wright     VecType  vec_type;
283926a6279SJames Wright     char     box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
28466d54740SJames Wright     PetscInt dim;
28566d54740SJames Wright 
28666d54740SJames Wright     PetscCall(DMGetDimension(user->dm, &dim));
28766d54740SJames Wright     if (dim == 2) box_faces_str[3] = '\0';
288926a6279SJames Wright     PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
289926a6279SJames Wright     PetscCall(DMGetVecType(user->dm, &vec_type));
290926a6279SJames Wright     PetscCall(PetscPrintf(comm,
291926a6279SJames Wright                           "  PETSc:\n"
292926a6279SJames Wright                           "    Box Faces                          : %s\n"
293926a6279SJames Wright                           "    DM VecType                         : %s\n"
294926a6279SJames Wright                           "    Time Stepping Scheme               : %s\n",
2953678fae3SJames Wright                           box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
29666d54740SJames Wright   }
2973678fae3SJames Wright   {
2983678fae3SJames Wright     char           pmat_type_str[PETSC_MAX_PATH_LEN];
2993678fae3SJames Wright     MatType        amat_type, pmat_type;
3003678fae3SJames Wright     Mat            Amat, Pmat;
3013678fae3SJames Wright     TSIJacobianFn *ijacob_function;
3023678fae3SJames Wright 
3033678fae3SJames Wright     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
3043678fae3SJames Wright     PetscCall(MatGetType(Amat, &amat_type));
3053678fae3SJames Wright     PetscCall(MatGetType(Pmat, &pmat_type));
3063678fae3SJames Wright 
3073678fae3SJames Wright     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
3083678fae3SJames Wright     if (!strcmp(pmat_type, MATCEED)) {
3093678fae3SJames Wright       MatType pmat_coo_type;
3103678fae3SJames Wright       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
3113678fae3SJames Wright 
3123678fae3SJames Wright       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
3133678fae3SJames Wright       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
3143678fae3SJames Wright       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
3153678fae3SJames Wright     }
3163678fae3SJames Wright     if (ijacob_function) {
3173678fae3SJames Wright       PetscCall(PetscPrintf(comm,
3183678fae3SJames Wright                             "    IJacobian A MatType                : %s\n"
3193678fae3SJames Wright                             "    IJacobian P MatType                : %s\n",
3203678fae3SJames Wright                             amat_type, pmat_type_str));
3213678fae3SJames Wright     }
3223678fae3SJames Wright   }
323926a6279SJames Wright   if (user->app_ctx->cont_steps) {
324926a6279SJames Wright     PetscCall(PetscPrintf(comm,
325926a6279SJames Wright                           "  Continue:\n"
326926a6279SJames Wright                           "    Filename:                          : %s\n"
327926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
328926a6279SJames Wright                           "    Time:                              : %g\n",
329926a6279SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
330926a6279SJames Wright   }
331926a6279SJames Wright   // Mesh
332926a6279SJames Wright   const PetscInt num_comp_q = 5;
333926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
334926a6279SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
335926a6279SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
336926a6279SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
337926a6279SJames Wright   PetscCall(PetscPrintf(comm,
338926a6279SJames Wright                         "  Mesh:\n"
339926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
340926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
341926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
342926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
343dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
344dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
345926a6279SJames Wright   // -- Get Partition Statistics
346926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
347926a6279SJames Wright   {
348926a6279SJames Wright     PetscInt *gather_buffer = NULL;
349dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
350926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
351926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
352926a6279SJames Wright 
353dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
354926a6279SJames Wright     if (!rank) {
355926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
356926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
357926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
358926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
359926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
360dfeb939dSJames Wright       PetscCall(PetscPrintf(
361dfeb939dSJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
362926a6279SJames Wright           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));
363926a6279SJames Wright     }
364926a6279SJames Wright 
365dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
366926a6279SJames Wright     if (!rank) {
367926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
368926a6279SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
369926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
370926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
371926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
372dfeb939dSJames Wright       PetscCall(PetscPrintf(
373dfeb939dSJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
374926a6279SJames Wright           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));
375926a6279SJames Wright     }
376926a6279SJames Wright 
37745abf96eSJames Wright     if (comm_size != 1) {
378dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
379926a6279SJames Wright       {
380926a6279SJames Wright         PetscSF            sf;
381dfeb939dSJames Wright         PetscInt           nrranks, niranks;
382dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
383dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
384926a6279SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
385926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
386dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
387926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
388926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
389926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
390dfeb939dSJames Wright           num_ghost_interface_ranks++;
391dfeb939dSJames Wright         }
392dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
393dfeb939dSJames Wright           if (iranks[i] == rank) continue;
394dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
395dfeb939dSJames Wright           num_owned_interface_ranks++;
396926a6279SJames Wright         }
397926a6279SJames Wright       }
398dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
399926a6279SJames Wright       if (!rank) {
400926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
401dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
402dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
403dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
404dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
405dfeb939dSJames Wright         PetscCall(PetscPrintf(
40645abf96eSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
40745abf96eSJames 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,
40845abf96eSJames Wright             part_shared_dof_ratio));
409dfeb939dSJames Wright       }
410dfeb939dSJames Wright 
411dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
412dfeb939dSJames Wright       if (!rank) {
413dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
414dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
415dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
416dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
417dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
418dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
419dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
420dfeb939dSJames Wright       }
421dfeb939dSJames Wright 
422dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
423dfeb939dSJames Wright       if (!rank) {
424dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
425dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
426dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
427dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
428dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
429dfeb939dSJames Wright         PetscCall(PetscPrintf(
43045abf96eSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
43145abf96eSJames 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,
43245abf96eSJames Wright             part_shared_dof_ratio));
433dfeb939dSJames Wright       }
434dfeb939dSJames Wright 
435dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
436dfeb939dSJames Wright       if (!rank) {
437dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
438dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
439dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
440dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
441dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
442dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
443dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
444926a6279SJames Wright       }
44545abf96eSJames Wright     }
446926a6279SJames Wright 
447926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
448926a6279SJames Wright   }
449926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
450926a6279SJames Wright }
451