xref: /honee/src/misc.c (revision b4c37c5c26bd208d7f474cbd9aa0850e82e1a854)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Miscellaneous utility functions
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdm.h>
13926a6279SJames Wright #include <petscsf.h>
14e419654dSJeremy L Thompson #include <petscts.h>
15e419654dSJeremy L Thompson 
16a515125bSLeila Ghaffari #include "../navierstokes.h"
179f59f36eSJames Wright #include "../qfunctions/mass.h"
18a515125bSLeila Ghaffari 
192b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
20*b4c37c5cSJames Wright   Ceed ceed = user->ceed;
21a515125bSLeila Ghaffari   PetscFunctionBeginUser;
22a515125bSLeila Ghaffari 
23a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
24b7f03f12SJed Brown   // Update time for evaluation
25a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
26*b4c37c5cSJames Wright   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
27a515125bSLeila Ghaffari 
28a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
29a515125bSLeila Ghaffari   // ICs
30a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
31a515125bSLeila Ghaffari   // -- CEED Restriction
32a515125bSLeila Ghaffari   CeedVector q0_ceed;
33*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL));
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
368f18bb8bSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
37a515125bSLeila Ghaffari 
38a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
39a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
40a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
41a515125bSLeila Ghaffari   // -- CEED Restriction
42a515125bSLeila Ghaffari   CeedVector mult_vec;
43*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
44a515125bSLeila Ghaffari 
45a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
46a515125bSLeila Ghaffari   PetscMemType m_mem_type;
47a515125bSLeila Ghaffari   Vec          multiplicity_loc;
482b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
49fd969b44SJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
50a515125bSLeila Ghaffari 
51a515125bSLeila Ghaffari   // -- Get multiplicity
52*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
53a515125bSLeila Ghaffari 
54a515125bSLeila Ghaffari   // -- Restore vectors
55fd969b44SJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
56a515125bSLeila Ghaffari 
57a515125bSLeila Ghaffari   // -- Local-to-Global
58a515125bSLeila Ghaffari   Vec multiplicity;
592b916ea7SJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
602b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
612b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
62a515125bSLeila Ghaffari 
63a515125bSLeila Ghaffari   // -- Fix multiplicity
642b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
652b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
66a515125bSLeila Ghaffari 
67a515125bSLeila Ghaffari   // -- Restore vectors
682b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
692b916ea7SJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
70a515125bSLeila Ghaffari 
71a515125bSLeila Ghaffari   // Cleanup
72*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
73*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q0_ceed));
74a515125bSLeila Ghaffari 
75d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
76a515125bSLeila Ghaffari }
77a515125bSLeila Ghaffari 
782b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
792b916ea7SJeremy L Thompson                                              Vec grad_FVM) {
809d437337SJames Wright   Vec Qbc, boundary_mask;
81a515125bSLeila Ghaffari   PetscFunctionBegin;
82a515125bSLeila Ghaffari 
832eb7bf1fSJames Wright   // Mask (zero) Strong BC entries
849d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
859d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
869d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
879d437337SJames Wright 
882b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
892b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
902b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
91a515125bSLeila Ghaffari 
92d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
93a515125bSLeila Ghaffari }
94a515125bSLeila Ghaffari 
95e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
96e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
97e1233009SJames Wright   PetscInt   file_step_number;
98e1233009SJames Wright   PetscInt32 token;
99e7754af5SKenneth E. Jansen   PetscReal  file_time;
100e7754af5SKenneth E. Jansen   PetscFunctionBeginUser;
101e7754af5SKenneth E. Jansen 
102e7754af5SKenneth E. Jansen   // Attempt
103e1233009SJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
104e1233009SJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
105e1233009SJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
106e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
107e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
108e7754af5SKenneth E. Jansen     if (time) *time = file_time;
109e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
110e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
111e7754af5SKenneth E. Jansen     PetscInt length, N;
112e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
113e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
114e7754af5SKenneth E. Jansen     PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N);
115e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
116e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
117e7754af5SKenneth E. Jansen 
118e7754af5SKenneth E. Jansen   // Load Q from existent solution
119e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
120e7754af5SKenneth E. Jansen 
121d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
122e7754af5SKenneth E. Jansen }
123e7754af5SKenneth E. Jansen 
124a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
125a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
126a515125bSLeila Ghaffari   Vec         Qref;
127a515125bSLeila Ghaffari   PetscViewer viewer;
128a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
129e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
130a515125bSLeila Ghaffari   PetscFunctionBegin;
131a515125bSLeila Ghaffari 
132a515125bSLeila Ghaffari   // Read reference file
1332b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
134e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
135e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
136a515125bSLeila Ghaffari 
137a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1382b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1392b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1402b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1412b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
142a515125bSLeila Ghaffari 
143a515125bSLeila Ghaffari   // Check error
144a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1452b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
146a515125bSLeila Ghaffari   }
147a515125bSLeila Ghaffari 
148a515125bSLeila Ghaffari   // Cleanup
1492b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1502b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
151a515125bSLeila Ghaffari 
152d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153a515125bSLeila Ghaffari }
154a515125bSLeila Ghaffari 
155a515125bSLeila Ghaffari // Get error for problems with exact solutions
1562b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
157a515125bSLeila Ghaffari   PetscInt  loc_nodes;
158a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
159a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
160a515125bSLeila Ghaffari   PetscFunctionBegin;
161a515125bSLeila Ghaffari 
162a515125bSLeila Ghaffari   // Get exact solution at final time
1632b916ea7SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1642b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1652b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1662b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
167a515125bSLeila Ghaffari 
168a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1692b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1702b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1712b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
172a515125bSLeila Ghaffari 
173a515125bSLeila Ghaffari   // Compute relative error
174a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
175a515125bSLeila Ghaffari 
176a515125bSLeila Ghaffari   // Output relative error
1772b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
178a515125bSLeila Ghaffari   // Cleanup
1792b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1802b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
181a515125bSLeila Ghaffari 
182d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
183a515125bSLeila Ghaffari }
184a515125bSLeila Ghaffari 
185a515125bSLeila Ghaffari // Post-processing
1862b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
187a515125bSLeila Ghaffari   PetscInt          steps;
188f0784ed3SJed Brown   TSConvergedReason reason;
189a515125bSLeila Ghaffari   PetscFunctionBegin;
190a515125bSLeila Ghaffari 
191a515125bSLeila Ghaffari   // Print relative error
1920e1e9333SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1932b916ea7SJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
194a515125bSLeila Ghaffari   }
195a515125bSLeila Ghaffari 
196a515125bSLeila Ghaffari   // Print final time and number of steps
1972b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
198f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1990e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
200f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
201f0784ed3SJed Brown                           steps, (double)final_time));
202a515125bSLeila Ghaffari   }
203a515125bSLeila Ghaffari 
204a515125bSLeila Ghaffari   // Output numerical values from command line
2052b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
206a515125bSLeila Ghaffari 
207a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
2080e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2092b916ea7SJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
210a515125bSLeila Ghaffari   }
211d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
212a515125bSLeila Ghaffari }
213a515125bSLeila Ghaffari 
214e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
215e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
216e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2179293eaa1SJed Brown 
218a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
219a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
220a515125bSLeila Ghaffari   PetscViewer viewer;
2212b916ea7SJeremy L Thompson 
222a515125bSLeila Ghaffari   PetscFunctionBegin;
223a515125bSLeila Ghaffari 
2242b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
225e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2262b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
227a515125bSLeila Ghaffari 
228d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
229a515125bSLeila Ghaffari }
230a515125bSLeila Ghaffari 
231a515125bSLeila Ghaffari // Record boundary values from initial condition
232a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2339d437337SJames Wright   Vec Qbc, boundary_mask;
234a515125bSLeila Ghaffari   PetscFunctionBegin;
235a515125bSLeila Ghaffari 
2362b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2372b916ea7SJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2382b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2392b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2402b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2412b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2422b916ea7SJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
243a515125bSLeila Ghaffari 
2449d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2459d437337SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2469d437337SJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2479d437337SJames Wright   PetscCall(VecSet(Q, 1.0));
2489d437337SJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2499d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2509d437337SJames Wright 
251d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
252a515125bSLeila Ghaffari }
25315a3537eSJed Brown 
25415a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
25515a3537eSJed Brown int FreeContextPetsc(void *data) {
2562b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
25715a3537eSJed Brown   return CEED_ERROR_SUCCESS;
25815a3537eSJed Brown }
2599f59f36eSJames Wright 
2609f59f36eSJames Wright // Return mass qfunction specification for number of components N
2619f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2629f59f36eSJames Wright   PetscFunctionBeginUser;
2639f59f36eSJames Wright 
2649f59f36eSJames Wright   switch (N) {
2659f59f36eSJames Wright     case 1:
266*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
2679f59f36eSJames Wright       break;
2689f59f36eSJames Wright     case 5:
269*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
2709f59f36eSJames Wright       break;
271c38c977aSJames Wright     case 7:
272*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
273c38c977aSJames Wright       break;
2749f59f36eSJames Wright     case 9:
275*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
2769f59f36eSJames Wright       break;
2779f59f36eSJames Wright     case 22:
278*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
2799f59f36eSJames Wright       break;
2809f59f36eSJames Wright     default:
2816f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2829f59f36eSJames Wright   }
2839f59f36eSJames Wright 
284*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
285*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
286*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
287d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2889f59f36eSJames Wright }
289e5e81594SJames Wright 
290457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
291457a5831SJames Wright   PetscFunctionBeginUser;
292d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
293457a5831SJames Wright 
294457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
295457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
296457a5831SJames Wright 
297457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
298457a5831SJames Wright 
299457a5831SJames Wright   PetscCall(PetscFree(context));
300457a5831SJames Wright 
301d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
302457a5831SJames Wright }
303457a5831SJames Wright 
304676080b4SJames Wright /*
305676080b4SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
306676080b4SJames Wright  *
307676080b4SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
308676080b4SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
309676080b4SJames Wright  *
310676080b4SJames Wright  * 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.
311676080b4SJames Wright  *
312676080b4SJames Wright  * @param[in]  comm           MPI_Comm for the program
313676080b4SJames Wright  * @param[in]  path           Path to the file
314676080b4SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
315676080b4SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
316676080b4SJames Wright  * @param[out] fp File        pointer to the opened file
317676080b4SJames Wright  */
318676080b4SJames Wright PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
319676080b4SJames Wright                                  FILE **fp) {
320defe8520SJames Wright   int    ndims;
321676080b4SJames Wright   char   line[char_array_len];
322676080b4SJames Wright   char **array;
323676080b4SJames Wright 
324676080b4SJames Wright   PetscFunctionBeginUser;
325676080b4SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
326676080b4SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
327676080b4SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
328defe8520SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
329676080b4SJames Wright 
330676080b4SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
331676080b4SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
332676080b4SJames Wright 
333d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
334676080b4SJames Wright }
335676080b4SJames Wright 
336676080b4SJames Wright /*
337676080b4SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
338676080b4SJames Wright  *
339676080b4SJames Wright  * 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.
340676080b4SJames Wright  *
341676080b4SJames Wright  * @param[in]  comm  MPI_Comm for the program
342676080b4SJames Wright  * @param[in]  path  Path to the file
343676080b4SJames Wright  * @param[out] nrows Number of rows
344676080b4SJames Wright  */
345676080b4SJames Wright PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
346676080b4SJames Wright   const PetscInt char_array_len = 512;
347676080b4SJames Wright   PetscInt       dims[2];
348676080b4SJames Wright   FILE          *fp;
349676080b4SJames Wright 
350676080b4SJames Wright   PetscFunctionBeginUser;
351676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
352676080b4SJames Wright   *nrows = dims[0];
353676080b4SJames Wright   PetscCall(PetscFClose(comm, fp));
354676080b4SJames Wright 
355d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
356676080b4SJames Wright }
35762b7942eSJames Wright 
35862b7942eSJames Wright PetscErrorCode PHASTADatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
359defe8520SJames Wright   PetscInt       dims[2];
360defe8520SJames Wright   int            ndims;
36162b7942eSJames Wright   FILE          *fp;
36262b7942eSJames Wright   const PetscInt char_array_len = 512;
36362b7942eSJames Wright   char           line[char_array_len];
36462b7942eSJames Wright   char         **row_array;
36562b7942eSJames Wright   PetscFunctionBeginUser;
36662b7942eSJames Wright 
36762b7942eSJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
36862b7942eSJames Wright 
36962b7942eSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
37062b7942eSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
37162b7942eSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3725d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
373defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
37462b7942eSJames Wright 
37562b7942eSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
37662b7942eSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
37762b7942eSJames Wright     }
37862b7942eSJames Wright   }
37962b7942eSJames Wright 
38062b7942eSJames Wright   PetscCall(PetscFClose(comm, fp));
38162b7942eSJames Wright 
382d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
38362b7942eSJames Wright }
3847eedc94cSJames Wright 
3857eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
3867eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
3877eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
3887eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
3897eedc94cSJames Wright static PetscClassId libCEED_classid;
3907eedc94cSJames Wright 
3917eedc94cSJames Wright PetscErrorCode RegisterLogEvents() {
3927eedc94cSJames Wright   PetscFunctionBeginUser;
3937eedc94cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
3947eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
3957eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
3967eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
3977eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
398d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3997eedc94cSJames Wright }
400defe8520SJames Wright 
401defe8520SJames Wright /**
402defe8520SJames Wright   @brief Translate array of CeedInt to PetscInt.
403defe8520SJames Wright     If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`.
404defe8520SJames Wright     Caller is responsible for freeing `array_petsc` with `free()`.
405defe8520SJames Wright 
406defe8520SJames Wright   @param[in]      num_entries  Number of array entries
407defe8520SJames Wright   @param[in,out]  array_ceed   Array of CeedInts
408defe8520SJames Wright   @param[out]     array_petsc  Array of PetscInts
409defe8520SJames Wright **/
410defe8520SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
411defe8520SJames Wright   CeedInt  int_c = 0;
412defe8520SJames Wright   PetscInt int_p = 0;
413defe8520SJames Wright 
414defe8520SJames Wright   PetscFunctionBeginUser;
415defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
416defe8520SJames Wright     *array_petsc = (PetscInt *)*array_ceed;
417defe8520SJames Wright   } else {
418defe8520SJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
419defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
420defe8520SJames Wright     free(*array_ceed);
421defe8520SJames Wright   }
422defe8520SJames Wright   *array_ceed = NULL;
423defe8520SJames Wright 
424defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
425defe8520SJames Wright }
426defe8520SJames Wright 
427defe8520SJames Wright /**
428defe8520SJames Wright   @brief Translate array of PetscInt to CeedInt.
429defe8520SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
430defe8520SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
431defe8520SJames Wright 
432defe8520SJames Wright   @param[in]      num_entries  Number of array entries
433defe8520SJames Wright   @param[in,out]  array_petsc  Array of PetscInts
434defe8520SJames Wright   @param[out]     array_ceed   Array of CeedInts
435defe8520SJames Wright **/
436defe8520SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
437defe8520SJames Wright   CeedInt  int_c = 0;
438defe8520SJames Wright   PetscInt int_p = 0;
439defe8520SJames Wright 
440defe8520SJames Wright   PetscFunctionBeginUser;
441defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
442defe8520SJames Wright     *array_ceed = (CeedInt *)*array_petsc;
443defe8520SJames Wright   } else {
444defe8520SJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
445defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
446defe8520SJames Wright     PetscCall(PetscFree(*array_petsc));
447defe8520SJames Wright   }
448defe8520SJames Wright   *array_petsc = NULL;
449defe8520SJames Wright 
450defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
451defe8520SJames Wright }
452926a6279SJames Wright 
453926a6279SJames Wright // Print information about the given simulation run
454926a6279SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
455*b4c37c5cSJames Wright   Ceed ceed = user->ceed;
456926a6279SJames Wright   PetscFunctionBeginUser;
457926a6279SJames Wright   // Header and rank
458926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
459926a6279SJames Wright   PetscMPIInt rank, comm_size;
460926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
461926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
462926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
463926a6279SJames Wright   PetscCall(PetscPrintf(comm,
464926a6279SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
465926a6279SJames Wright                         "  MPI:\n"
466926a6279SJames Wright                         "    Host Name                          : %s\n"
467926a6279SJames Wright                         "    Total ranks                        : %d\n",
468926a6279SJames Wright                         host_name, comm_size));
469926a6279SJames Wright 
470926a6279SJames Wright   // Problem specific info
4712d49c0afSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
472926a6279SJames Wright 
473926a6279SJames Wright   // libCEED
474926a6279SJames Wright   const char *used_resource;
475926a6279SJames Wright   CeedMemType mem_type_backend;
476*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
477*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
478926a6279SJames Wright   PetscCall(PetscPrintf(comm,
479926a6279SJames Wright                         "  libCEED:\n"
480926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
481926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
482926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
483926a6279SJames Wright   // PETSc
484926a6279SJames Wright   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
485926a6279SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
486926a6279SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
487926a6279SJames Wright   MatType mat_type;
488926a6279SJames Wright   VecType vec_type;
489926a6279SJames Wright   PetscCall(DMGetMatType(user->dm, &mat_type));
490926a6279SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
491926a6279SJames Wright   PetscCall(PetscPrintf(comm,
492926a6279SJames Wright                         "  PETSc:\n"
493926a6279SJames Wright                         "    Box Faces                          : %s\n"
494926a6279SJames Wright                         "    DM MatType                         : %s\n"
495926a6279SJames Wright                         "    DM VecType                         : %s\n"
496926a6279SJames Wright                         "    Time Stepping Scheme               : %s\n",
497926a6279SJames Wright                         box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
498926a6279SJames Wright   if (user->app_ctx->cont_steps) {
499926a6279SJames Wright     PetscCall(PetscPrintf(comm,
500926a6279SJames Wright                           "  Continue:\n"
501926a6279SJames Wright                           "    Filename:                          : %s\n"
502926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
503926a6279SJames Wright                           "    Time:                              : %g\n",
504926a6279SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
505926a6279SJames Wright   }
506926a6279SJames Wright   // Mesh
507926a6279SJames Wright   const PetscInt num_comp_q = 5;
508926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
509926a6279SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
510926a6279SJames Wright   // -- Get global size
511926a6279SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
512926a6279SJames Wright   // -- Get local size
513926a6279SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
514926a6279SJames Wright   PetscCall(PetscPrintf(comm,
515926a6279SJames Wright                         "  Mesh:\n"
516926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
517926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
518926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
519926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
520dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
521dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
522926a6279SJames Wright   // -- Get Partition Statistics
523926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
524926a6279SJames Wright   {
525926a6279SJames Wright     PetscInt *gather_buffer = NULL;
526dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
527926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
528926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
529926a6279SJames Wright 
530dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
531926a6279SJames Wright     if (!rank) {
532926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
533926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
534926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
535926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
536926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
537dfeb939dSJames Wright       PetscCall(PetscPrintf(
538dfeb939dSJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
539926a6279SJames 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));
540926a6279SJames Wright     }
541926a6279SJames Wright 
542dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
543926a6279SJames Wright     if (!rank) {
544926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
545926a6279SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
546926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
547926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
548926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
549dfeb939dSJames Wright       PetscCall(PetscPrintf(
550dfeb939dSJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
551926a6279SJames 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));
552926a6279SJames Wright     }
553926a6279SJames Wright 
55445abf96eSJames Wright     if (comm_size != 1) {
555dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
556926a6279SJames Wright       {
557926a6279SJames Wright         PetscSF            sf;
558dfeb939dSJames Wright         PetscInt           nrranks, niranks;
559dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
560dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
561926a6279SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
562926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
563dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
564926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
565926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
566926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
567dfeb939dSJames Wright           num_ghost_interface_ranks++;
568dfeb939dSJames Wright         }
569dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
570dfeb939dSJames Wright           if (iranks[i] == rank) continue;
571dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
572dfeb939dSJames Wright           num_owned_interface_ranks++;
573926a6279SJames Wright         }
574926a6279SJames Wright       }
575dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
576926a6279SJames Wright       if (!rank) {
577926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
578dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
579dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
580dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
581dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
582dfeb939dSJames Wright         PetscCall(PetscPrintf(
58345abf96eSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
58445abf96eSJames 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,
58545abf96eSJames Wright             part_shared_dof_ratio));
586dfeb939dSJames Wright       }
587dfeb939dSJames Wright 
588dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
589dfeb939dSJames Wright       if (!rank) {
590dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
591dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
592dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
593dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
594dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
595dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
596dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
597dfeb939dSJames Wright       }
598dfeb939dSJames Wright 
599dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
600dfeb939dSJames Wright       if (!rank) {
601dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
602dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
603dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
604dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
605dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
606dfeb939dSJames Wright         PetscCall(PetscPrintf(
60745abf96eSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
60845abf96eSJames 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,
60945abf96eSJames Wright             part_shared_dof_ratio));
610dfeb939dSJames Wright       }
611dfeb939dSJames Wright 
612dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
613dfeb939dSJames Wright       if (!rank) {
614dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
615dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
616dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
617dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
618dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
619dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
620dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
621926a6279SJames Wright       }
62245abf96eSJames Wright     }
623926a6279SJames Wright 
624926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
625926a6279SJames Wright   }
626926a6279SJames Wright 
627926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
628926a6279SJames Wright }
629