xref: /libCEED/examples/fluids/src/misc.c (revision 8535f0aae359c08d55230a1e48a051b83ac1b33c)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Miscellaneous utility functions
1077841947SLeila Ghaffari 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdm.h>
132e31f396SJames Wright #include <petscsf.h>
1449aac155SJeremy L Thompson #include <petscts.h>
1549aac155SJeremy L Thompson 
1677841947SLeila Ghaffari #include "../navierstokes.h"
17ef080ff9SJames Wright #include "../qfunctions/mass.h"
1877841947SLeila Ghaffari 
192b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
20a424bcd0SJames Wright   Ceed         ceed = user->ceed;
21f6ce2b0aSJames Wright   CeedVector   mult_vec;
22f6ce2b0aSJames Wright   PetscMemType m_mem_type;
23f6ce2b0aSJames Wright   Vec          Multiplicity, Multiplicity_loc;
24f6ce2b0aSJames Wright 
2577841947SLeila Ghaffari   PetscFunctionBeginUser;
26a424bcd0SJames Wright   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
275263e9c6SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
2877841947SLeila Ghaffari 
29a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
3077841947SLeila Ghaffari 
3177841947SLeila Ghaffari   // -- Get multiplicity
32f6ce2b0aSJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
33d0593705SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
34a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
35d0593705SJames Wright   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
3677841947SLeila Ghaffari 
37f6ce2b0aSJames Wright   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
38f6ce2b0aSJames Wright   PetscCall(VecZeroEntries(Multiplicity));
39f6ce2b0aSJames Wright   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
4077841947SLeila Ghaffari 
4177841947SLeila Ghaffari   // -- Fix multiplicity
42f6ce2b0aSJames Wright   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
43f6ce2b0aSJames Wright   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
4477841947SLeila Ghaffari 
45f6ce2b0aSJames Wright   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
46f6ce2b0aSJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
47a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
48ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4977841947SLeila Ghaffari }
5077841947SLeila Ghaffari 
513c9e7ad1SJames Wright // Record boundary values from initial condition
523c9e7ad1SJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
533c9e7ad1SJames Wright   PetscFunctionBeginUser;
54251425b7SJames Wright   {  // Capture initial condition values in Qbc
55251425b7SJames Wright     Vec Qbc;
56251425b7SJames Wright 
573c9e7ad1SJames Wright     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
583c9e7ad1SJames Wright     PetscCall(VecCopy(Q_loc, Qbc));
593c9e7ad1SJames Wright     PetscCall(VecZeroEntries(Q_loc));
603c9e7ad1SJames Wright     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
613c9e7ad1SJames Wright     PetscCall(VecAXPY(Qbc, -1., Q_loc));
623c9e7ad1SJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
63251425b7SJames Wright   }
643c9e7ad1SJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
653c9e7ad1SJames Wright 
66251425b7SJames Wright   {  // Set boundary mask to zero out essential BCs
67251425b7SJames Wright     Vec boundary_mask, ones;
68251425b7SJames Wright 
693c9e7ad1SJames Wright     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
70251425b7SJames Wright     PetscCall(DMGetGlobalVector(dm, &ones));
713c9e7ad1SJames Wright     PetscCall(VecZeroEntries(boundary_mask));
72251425b7SJames Wright     PetscCall(VecSet(ones, 1.0));
73251425b7SJames Wright     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
743c9e7ad1SJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
75251425b7SJames Wright     PetscCall(DMRestoreGlobalVector(dm, &ones));
76251425b7SJames Wright   }
773c9e7ad1SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
783c9e7ad1SJames Wright }
793c9e7ad1SJames Wright 
803c9e7ad1SJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
812b730f8bSJeremy L Thompson                                                   Vec grad_FVM) {
825571c6fdSJames Wright   Vec Qbc, boundary_mask;
8377841947SLeila Ghaffari 
84f17d818dSJames Wright   PetscFunctionBeginUser;
853722cd23SJames Wright   // Mask (zero) Strong BC entries
865571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
875571c6fdSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
885571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
895571c6fdSJames Wright 
902b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
912b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
922b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
93ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
9477841947SLeila Ghaffari }
9577841947SLeila Ghaffari 
967ceb6423SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
977ceb6423SJed Brown   PetscFunctionBeginUser;
987ceb6423SJed Brown   if (file_type == PETSC_INT32) {
997ceb6423SJed Brown     PetscInt32 val;
1007ceb6423SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
1017ceb6423SJed Brown     *out = val;
1027ceb6423SJed Brown   } else if (file_type == PETSC_INT64) {
1037ceb6423SJed Brown     PetscInt64 val;
1047ceb6423SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
1057ceb6423SJed Brown     *out = val;
1067ceb6423SJed Brown   } else {
1077ceb6423SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
1087ceb6423SJed Brown   }
1097ceb6423SJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
1107ceb6423SJed Brown }
1117ceb6423SJed Brown 
112530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
113530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
1140de6a49fSJames Wright   PetscInt      file_step_number;
1150de6a49fSJames Wright   PetscInt32    token;
116530ad8c4SKenneth E. Jansen   PetscReal     file_time;
1177ceb6423SJed Brown   PetscDataType file_type = PETSC_INT32;
118530ad8c4SKenneth E. Jansen 
119f17d818dSJames Wright   PetscFunctionBeginUser;
1200de6a49fSJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
1210de6a49fSJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
1220de6a49fSJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
1237ceb6423SJed Brown     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
1247ceb6423SJed Brown     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
1257ceb6423SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
126530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
127530ad8c4SKenneth E. Jansen     if (time) *time = file_time;
128530ad8c4SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
129530ad8c4SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
130530ad8c4SKenneth E. Jansen     PetscInt length, N;
1317ceb6423SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
132530ad8c4SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
133530ad8c4SKenneth 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);
134530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
135530ad8c4SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
136530ad8c4SKenneth E. Jansen 
137530ad8c4SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
138ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
139530ad8c4SKenneth E. Jansen }
140530ad8c4SKenneth E. Jansen 
14177841947SLeila Ghaffari // Compare reference solution values with current test run for CI
1423c9e7ad1SJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
143*8535f0aaSJeremy L Thompson   Vec         Q_ref;
14477841947SLeila Ghaffari   PetscViewer viewer;
145*8535f0aaSJeremy L Thompson   PetscReal   error, norm_Q, norm_Q_ref;
146530ad8c4SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
14777841947SLeila Ghaffari 
148f17d818dSJames Wright   PetscFunctionBeginUser;
14977841947SLeila Ghaffari   // Read reference file
150*8535f0aaSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Q_ref));
151013a5551SJames Wright   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
152530ad8c4SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
153*8535f0aaSJeremy L Thompson   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q_ref, NULL, NULL));
15477841947SLeila Ghaffari 
15577841947SLeila Ghaffari   // Compute error with respect to reference solution
156*8535f0aaSJeremy L Thompson   PetscCall(VecNorm(Q_ref, NORM_MAX, &norm_Q));
157*8535f0aaSJeremy L Thompson   PetscCall(VecNorm(Q_ref, NORM_MAX, &norm_Q_ref));
158*8535f0aaSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_ref));
159*8535f0aaSJeremy L Thompson   PetscCall(VecScale(Q, 1. / norm_Q_ref));
1602b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
16177841947SLeila Ghaffari 
16277841947SLeila Ghaffari   // Check error
16377841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
164*8535f0aaSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\nReference solution max norm: %g Computed solution max norm %g\n",
165*8535f0aaSJeremy L Thompson                           (double)error, (double)norm_Q_ref, (double)norm_Q));
16677841947SLeila Ghaffari   }
16777841947SLeila Ghaffari 
16877841947SLeila Ghaffari   // Cleanup
1692b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
170*8535f0aaSJeremy L Thompson   PetscCall(VecDestroy(&Q_ref));
171ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17277841947SLeila Ghaffari }
17377841947SLeila Ghaffari 
17477841947SLeila Ghaffari // Get error for problems with exact solutions
1753c9e7ad1SJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
17677841947SLeila Ghaffari   PetscInt  loc_nodes;
17777841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
17877841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
17977841947SLeila Ghaffari 
180f17d818dSJames Wright   PetscFunctionBeginUser;
18177841947SLeila Ghaffari   // Get exact solution at final time
182f6ce2b0aSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1832b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1842b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1852b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
18677841947SLeila Ghaffari 
18777841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1882b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1892b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1902b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
19177841947SLeila Ghaffari 
19277841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
1932b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1942b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
195f6ce2b0aSJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
196ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19777841947SLeila Ghaffari }
19877841947SLeila Ghaffari 
19977841947SLeila Ghaffari // Post-processing
200731c13d7SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
20177841947SLeila Ghaffari   PetscInt          steps;
202cf7a0454SJed Brown   TSConvergedReason reason;
20377841947SLeila Ghaffari 
204f17d818dSJames Wright   PetscFunctionBeginUser;
20577841947SLeila Ghaffari   // Print relative error
206de327db4SJames Wright   if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
2073c9e7ad1SJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
20877841947SLeila Ghaffari   }
20977841947SLeila Ghaffari 
21077841947SLeila Ghaffari   // Print final time and number of steps
2112b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
212cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
2138fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
214cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
215cf7a0454SJed Brown                           steps, (double)final_time));
21677841947SLeila Ghaffari   }
21777841947SLeila Ghaffari 
21877841947SLeila Ghaffari   // Output numerical values from command line
2192b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
22077841947SLeila Ghaffari 
22177841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
2228fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2233c9e7ad1SJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
22477841947SLeila Ghaffari   }
225ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
22677841947SLeila Ghaffari }
22777841947SLeila Ghaffari 
2280de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
2290de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
2300de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2314de8550aSJed Brown 
23277841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
23377841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
23477841947SLeila Ghaffari   PetscViewer viewer;
2352b730f8bSJeremy L Thompson 
236f17d818dSJames Wright   PetscFunctionBeginUser;
2372b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
238530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2392b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
240ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
24177841947SLeila Ghaffari }
24277841947SLeila Ghaffari 
243841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
244841e4c73SJed Brown int FreeContextPetsc(void *data) {
2452b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
246841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
247841e4c73SJed Brown }
248ef080ff9SJames Wright 
249ef080ff9SJames Wright // Return mass qfunction specification for number of components N
250ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
251ef080ff9SJames Wright   PetscFunctionBeginUser;
252ef080ff9SJames Wright   switch (N) {
253ef080ff9SJames Wright     case 1:
254a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
255ef080ff9SJames Wright       break;
256ef080ff9SJames Wright     case 5:
257a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
258ef080ff9SJames Wright       break;
259052409adSJames Wright     case 7:
260a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
261052409adSJames Wright       break;
262ef080ff9SJames Wright     case 9:
263a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
264ef080ff9SJames Wright       break;
265ef080ff9SJames Wright     case 22:
266a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
267ef080ff9SJames Wright       break;
268ef080ff9SJames Wright     default:
26983ae4962SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
270ef080ff9SJames Wright   }
271ef080ff9SJames Wright 
272a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
273a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
274a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2757f2a9303SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
276ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
277ef080ff9SJames Wright }
2780df12fefSJames Wright 
279999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
280999ff5c7SJames Wright   PetscFunctionBeginUser;
281ee4ca9cbSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
282999ff5c7SJames Wright 
283999ff5c7SJames Wright   PetscCall(DMDestroy(&context->dm));
284999ff5c7SJames Wright   PetscCall(KSPDestroy(&context->ksp));
285999ff5c7SJames Wright 
286999ff5c7SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
287999ff5c7SJames Wright 
288999ff5c7SJames Wright   PetscCall(PetscFree(context));
289ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
290999ff5c7SJames Wright }
291999ff5c7SJames Wright 
2928b892a05SJames Wright /*
2938b892a05SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
2948b892a05SJames Wright  *
2958b892a05SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
2968b892a05SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
2978b892a05SJames Wright  *
2988b892a05SJames 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.
2998b892a05SJames Wright  *
3008b892a05SJames Wright  * @param[in]  comm           MPI_Comm for the program
3018b892a05SJames Wright  * @param[in]  path           Path to the file
3028b892a05SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
3038b892a05SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
3048b892a05SJames Wright  * @param[out] fp File        pointer to the opened file
3058b892a05SJames Wright  */
306cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
3078b892a05SJames Wright                                  FILE **fp) {
308457e73b2SJames Wright   int    ndims;
3098b892a05SJames Wright   char   line[char_array_len];
3108b892a05SJames Wright   char **array;
3118b892a05SJames Wright 
3128b892a05SJames Wright   PetscFunctionBeginUser;
3138b892a05SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
3148b892a05SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
3158b892a05SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
316457e73b2SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
3178b892a05SJames Wright 
3188b892a05SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
3198b892a05SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
320ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3218b892a05SJames Wright }
3228b892a05SJames Wright 
3238b892a05SJames Wright /*
3248b892a05SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
3258b892a05SJames Wright  *
3268b892a05SJames 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.
3278b892a05SJames Wright  *
3288b892a05SJames Wright  * @param[in]  comm  MPI_Comm for the program
3298b892a05SJames Wright  * @param[in]  path  Path to the file
3308b892a05SJames Wright  * @param[out] nrows Number of rows
3318b892a05SJames Wright  */
332cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
3338b892a05SJames Wright   const PetscInt char_array_len = 512;
3348b892a05SJames Wright   PetscInt       dims[2];
3358b892a05SJames Wright   FILE          *fp;
3368b892a05SJames Wright 
3378b892a05SJames Wright   PetscFunctionBeginUser;
338cbef7084SJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3398b892a05SJames Wright   *nrows = dims[0];
3408b892a05SJames Wright   PetscCall(PetscFClose(comm, fp));
341ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3428b892a05SJames Wright }
3434cc9442bSJames Wright 
344cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
345457e73b2SJames Wright   PetscInt       dims[2];
3464cc9442bSJames Wright   FILE          *fp;
3474cc9442bSJames Wright   const PetscInt char_array_len = 512;
3484cc9442bSJames Wright   char           line[char_array_len];
3494cc9442bSJames Wright 
350f17d818dSJames Wright   PetscFunctionBeginUser;
351cbef7084SJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3524cc9442bSJames Wright 
3534cc9442bSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
35438690fecSJames Wright     int    ndims;
35538690fecSJames Wright     char **row_array;
35638690fecSJames Wright 
3574cc9442bSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
3584cc9442bSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3590e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
360457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
3614cc9442bSJames Wright 
36238690fecSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
36338690fecSJames Wright     PetscCall(PetscStrToArrayDestroy(ndims, row_array));
3644cc9442bSJames Wright   }
3654cc9442bSJames Wright 
3664cc9442bSJames Wright   PetscCall(PetscFClose(comm, fp));
367ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3684cc9442bSJames Wright }
36975d1798cSJames Wright 
3702e31f396SJames Wright // Print information about the given simulation run
3717bc7b61fSJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
372a424bcd0SJames Wright   Ceed     ceed = user->ceed;
3737bc7b61fSJames Wright   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
3747bc7b61fSJames Wright 
3752e31f396SJames Wright   PetscFunctionBeginUser;
3762e31f396SJames Wright   // Header and rank
3772e31f396SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
3782e31f396SJames Wright   PetscMPIInt rank, comm_size;
3792e31f396SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
3802e31f396SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3812e31f396SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
3822e31f396SJames Wright   PetscCall(PetscPrintf(comm,
3832e31f396SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
3842e31f396SJames Wright                         "  MPI:\n"
3852e31f396SJames Wright                         "    Host Name                          : %s\n"
3862e31f396SJames Wright                         "    Total ranks                        : %d\n",
3872e31f396SJames Wright                         host_name, comm_size));
3882e31f396SJames Wright 
3892e31f396SJames Wright   // Problem specific info
390367c849eSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
3912e31f396SJames Wright 
3922e31f396SJames Wright   // libCEED
3932e31f396SJames Wright   const char *used_resource;
3942e31f396SJames Wright   CeedMemType mem_type_backend;
395a424bcd0SJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
396a424bcd0SJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
3972e31f396SJames Wright   PetscCall(PetscPrintf(comm,
3982e31f396SJames Wright                         "  libCEED:\n"
3992e31f396SJames Wright                         "    libCEED Backend                    : %s\n"
4002e31f396SJames Wright                         "    libCEED Backend MemType            : %s\n",
4012e31f396SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
4022e31f396SJames Wright   // PETSc
4037bc7b61fSJames Wright   VecType vec_type;
4042e31f396SJames Wright   char    box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
4052e31f396SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
4062e31f396SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
4072e31f396SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
4082e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4092e31f396SJames Wright                         "  PETSc:\n"
4102e31f396SJames Wright                         "    Box Faces                          : %s\n"
4112e31f396SJames Wright                         "    DM VecType                         : %s\n"
4122e31f396SJames Wright                         "    Time Stepping Scheme               : %s\n",
4137bc7b61fSJames Wright                         box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
4147bc7b61fSJames Wright   {
4157bc7b61fSJames Wright     char           pmat_type_str[PETSC_MAX_PATH_LEN];
4167bc7b61fSJames Wright     MatType        amat_type, pmat_type;
4177bc7b61fSJames Wright     Mat            Amat, Pmat;
4187bc7b61fSJames Wright     TSIJacobianFn *ijacob_function;
4197bc7b61fSJames Wright 
4207bc7b61fSJames Wright     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
4217bc7b61fSJames Wright     PetscCall(MatGetType(Amat, &amat_type));
4227bc7b61fSJames Wright     PetscCall(MatGetType(Pmat, &pmat_type));
4237bc7b61fSJames Wright 
4247bc7b61fSJames Wright     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
4257bc7b61fSJames Wright     if (!strcmp(pmat_type, MATCEED)) {
4267bc7b61fSJames Wright       MatType pmat_coo_type;
4277bc7b61fSJames Wright       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
4287bc7b61fSJames Wright 
4297bc7b61fSJames Wright       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
4307bc7b61fSJames Wright       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
4317bc7b61fSJames Wright       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
4327bc7b61fSJames Wright     }
4337bc7b61fSJames Wright     if (ijacob_function) {
4347bc7b61fSJames Wright       PetscCall(PetscPrintf(comm,
4357bc7b61fSJames Wright                             "    IJacobian A MatType                : %s\n"
4367bc7b61fSJames Wright                             "    IJacobian P MatType                : %s\n",
4377bc7b61fSJames Wright                             amat_type, pmat_type_str));
4387bc7b61fSJames Wright     }
4397bc7b61fSJames Wright   }
4402e31f396SJames Wright   if (user->app_ctx->cont_steps) {
4412e31f396SJames Wright     PetscCall(PetscPrintf(comm,
4422e31f396SJames Wright                           "  Continue:\n"
4432e31f396SJames Wright                           "    Filename:                          : %s\n"
4442e31f396SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
4452e31f396SJames Wright                           "    Time:                              : %g\n",
4462e31f396SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
4472e31f396SJames Wright   }
4482e31f396SJames Wright   // Mesh
4492e31f396SJames Wright   const PetscInt num_comp_q = 5;
4502e31f396SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
4512e31f396SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
4522e31f396SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
4532e31f396SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
4542e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4552e31f396SJames Wright                         "  Mesh:\n"
4562e31f396SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
4572e31f396SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
4582e31f396SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
4592e31f396SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
460ce11f295SJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
461ce11f295SJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
4622e31f396SJames Wright   // -- Get Partition Statistics
4632e31f396SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
4642e31f396SJames Wright   {
4652e31f396SJames Wright     PetscInt *gather_buffer = NULL;
466ce11f295SJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
4672e31f396SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
4682e31f396SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
4692e31f396SJames Wright 
470ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4712e31f396SJames Wright     if (!rank) {
4722e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
4732e31f396SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
4742e31f396SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
4752e31f396SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
4762e31f396SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
477ce11f295SJames Wright       PetscCall(PetscPrintf(
478ce11f295SJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
4792e31f396SJames 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));
4802e31f396SJames Wright     }
4812e31f396SJames Wright 
482ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4832e31f396SJames Wright     if (!rank) {
4842e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
4852e31f396SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
4862e31f396SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
4872e31f396SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
4882e31f396SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
489ce11f295SJames Wright       PetscCall(PetscPrintf(
490ce11f295SJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
4912e31f396SJames 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));
4922e31f396SJames Wright     }
4932e31f396SJames Wright 
49465ba01baSJames Wright     if (comm_size != 1) {
495ce11f295SJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
4962e31f396SJames Wright       {
4972e31f396SJames Wright         PetscSF            sf;
4988a310472SJeremy L Thompson         PetscMPIInt        nrranks, niranks;
499ce11f295SJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
500ce11f295SJames Wright         const PetscMPIInt *rranks, *iranks;
5018a310472SJeremy L Thompson 
5022e31f396SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
5032e31f396SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
504ce11f295SJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
5052e31f396SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
5062e31f396SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
5072e31f396SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
508ce11f295SJames Wright           num_ghost_interface_ranks++;
509ce11f295SJames Wright         }
510ce11f295SJames Wright         for (PetscInt i = 0; i < niranks; i++) {
511ce11f295SJames Wright           if (iranks[i] == rank) continue;
512ce11f295SJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
513ce11f295SJames Wright           num_owned_interface_ranks++;
5142e31f396SJames Wright         }
5152e31f396SJames Wright       }
516ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
5172e31f396SJames Wright       if (!rank) {
5182e31f396SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
519ce11f295SJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
520ce11f295SJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
521ce11f295SJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
522ce11f295SJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
523ce11f295SJames Wright         PetscCall(PetscPrintf(
52465ba01baSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
52565ba01baSJames 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,
52665ba01baSJames Wright             part_shared_dof_ratio));
527ce11f295SJames Wright       }
528ce11f295SJames Wright 
529ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
530ce11f295SJames Wright       if (!rank) {
531ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
532ce11f295SJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
533ce11f295SJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
534ce11f295SJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
535ce11f295SJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
536ce11f295SJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
537ce11f295SJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
538ce11f295SJames Wright       }
539ce11f295SJames Wright 
540ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
541ce11f295SJames Wright       if (!rank) {
542ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
543ce11f295SJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
544ce11f295SJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
545ce11f295SJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
546ce11f295SJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
547ce11f295SJames Wright         PetscCall(PetscPrintf(
54865ba01baSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
54965ba01baSJames 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,
55065ba01baSJames Wright             part_shared_dof_ratio));
551ce11f295SJames Wright       }
552ce11f295SJames Wright 
553ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
554ce11f295SJames Wright       if (!rank) {
555ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
556ce11f295SJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
557ce11f295SJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
558ce11f295SJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
559ce11f295SJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
560ce11f295SJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
561ce11f295SJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
5622e31f396SJames Wright       }
56365ba01baSJames Wright     }
5642e31f396SJames Wright 
5652e31f396SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
5662e31f396SJames Wright   }
5672e31f396SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5682e31f396SJames Wright }
569