1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
ICs_FixMultiplicity(DM dm,CeedData ceed_data,User user,Vec Q_loc,Vec Q,CeedScalar time)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
SetBCsFromICs(DM dm,Vec Q,Vec Q_loc)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
DMPlexInsertBoundaryValues_FromICs(DM dm,PetscBool insert_essential,Vec Q_loc,PetscReal time,Vec face_geom_FVM,Vec cell_geom_FVM,Vec grad_FVM)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
BinaryReadIntoInt(PetscViewer viewer,PetscInt * out,PetscDataType file_type)967ceb6423SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
977ceb6423SJed Brown PetscFunctionBeginUser;
989a2e771eSJames Wright *out = -13; // appease the overzealous GCC compiler warning Gods
997ceb6423SJed Brown if (file_type == PETSC_INT32) {
1007ceb6423SJed Brown PetscInt32 val;
1017ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
1027ceb6423SJed Brown *out = val;
1037ceb6423SJed Brown } else if (file_type == PETSC_INT64) {
1047ceb6423SJed Brown PetscInt64 val;
1057ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
1067ceb6423SJed Brown *out = val;
1077ceb6423SJed Brown } else {
1087ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
1097ceb6423SJed Brown }
1107ceb6423SJed Brown PetscFunctionReturn(PETSC_SUCCESS);
1117ceb6423SJed Brown }
1127ceb6423SJed Brown
113530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
LoadFluidsBinaryVec(MPI_Comm comm,PetscViewer viewer,Vec Q,PetscReal * time,PetscInt * step_number)114530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
1150de6a49fSJames Wright PetscInt file_step_number;
1160de6a49fSJames Wright PetscInt32 token;
117530ad8c4SKenneth E. Jansen PetscReal file_time;
1187ceb6423SJed Brown PetscDataType file_type = PETSC_INT32;
119530ad8c4SKenneth E. Jansen
120f17d818dSJames Wright PetscFunctionBeginUser;
1210de6a49fSJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
1220de6a49fSJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
1230de6a49fSJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header
1247ceb6423SJed Brown if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
1257ceb6423SJed Brown else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
1267ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
127530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
128530ad8c4SKenneth E. Jansen if (time) *time = file_time;
129530ad8c4SKenneth E. Jansen if (step_number) *step_number = file_step_number;
130530ad8c4SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
131530ad8c4SKenneth E. Jansen PetscInt length, N;
1327ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
133530ad8c4SKenneth E. Jansen PetscCall(VecGetSize(Q, &N));
134530ad8c4SKenneth 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);
135530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
136530ad8c4SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
137530ad8c4SKenneth E. Jansen
138530ad8c4SKenneth E. Jansen PetscCall(VecLoad(Q, viewer));
139ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
140530ad8c4SKenneth E. Jansen }
141530ad8c4SKenneth E. Jansen
14277841947SLeila Ghaffari // Compare reference solution values with current test run for CI
RegressionTest(AppCtx app_ctx,Vec Q)1433c9e7ad1SJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
1448535f0aaSJeremy L Thompson Vec Q_ref;
14577841947SLeila Ghaffari PetscViewer viewer;
1468535f0aaSJeremy L Thompson PetscReal error, norm_Q, norm_Q_ref;
147530ad8c4SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q);
14877841947SLeila Ghaffari
149f17d818dSJames Wright PetscFunctionBeginUser;
15077841947SLeila Ghaffari // Read reference file
1518535f0aaSJeremy L Thompson PetscCall(VecDuplicate(Q, &Q_ref));
152013a5551SJames Wright PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
153530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
1548535f0aaSJeremy L Thompson PetscCall(LoadFluidsBinaryVec(comm, viewer, Q_ref, NULL, NULL));
15577841947SLeila Ghaffari
15677841947SLeila Ghaffari // Compute error with respect to reference solution
1578535f0aaSJeremy L Thompson PetscCall(VecNorm(Q_ref, NORM_MAX, &norm_Q));
1588535f0aaSJeremy L Thompson PetscCall(VecNorm(Q_ref, NORM_MAX, &norm_Q_ref));
1598535f0aaSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_ref));
1608535f0aaSJeremy L Thompson PetscCall(VecScale(Q, 1. / norm_Q_ref));
1612b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error));
16277841947SLeila Ghaffari
16377841947SLeila Ghaffari // Check error
16477841947SLeila Ghaffari if (error > app_ctx->test_tol) {
1658535f0aaSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\nReference solution max norm: %g Computed solution max norm %g\n",
1668535f0aaSJeremy L Thompson (double)error, (double)norm_Q_ref, (double)norm_Q));
16777841947SLeila Ghaffari }
16877841947SLeila Ghaffari
16977841947SLeila Ghaffari // Cleanup
1702b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer));
1718535f0aaSJeremy L Thompson PetscCall(VecDestroy(&Q_ref));
172ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
17377841947SLeila Ghaffari }
17477841947SLeila Ghaffari
17577841947SLeila Ghaffari // Get error for problems with exact solutions
PrintError(CeedData ceed_data,DM dm,User user,Vec Q,PetscScalar final_time)1763c9e7ad1SJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
17777841947SLeila Ghaffari PetscInt loc_nodes;
17877841947SLeila Ghaffari Vec Q_exact, Q_exact_loc;
17977841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact;
18077841947SLeila Ghaffari
181f17d818dSJames Wright PetscFunctionBeginUser;
18277841947SLeila Ghaffari // Get exact solution at final time
183f6ce2b0aSJames Wright PetscCall(DMGetGlobalVector(dm, &Q_exact));
1842b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1852b730f8bSJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1862b730f8bSJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
18777841947SLeila Ghaffari
18877841947SLeila Ghaffari // Get |exact solution - obtained solution|
1892b730f8bSJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1902b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact));
1912b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error));
19277841947SLeila Ghaffari
19377841947SLeila Ghaffari rel_error = norm_error / norm_exact;
1942b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1952b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
196f6ce2b0aSJames Wright PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
197ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
19877841947SLeila Ghaffari }
19977841947SLeila Ghaffari
20077841947SLeila Ghaffari // Post-processing
PostProcess(TS ts,CeedData ceed_data,DM dm,ProblemData problem,User user,Vec Q,PetscScalar final_time)201731c13d7SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
20277841947SLeila Ghaffari PetscInt steps;
203cf7a0454SJed Brown TSConvergedReason reason;
20477841947SLeila Ghaffari
205f17d818dSJames Wright PetscFunctionBeginUser;
20677841947SLeila Ghaffari // Print relative error
207de327db4SJames Wright if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
2083c9e7ad1SJames Wright PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
20977841947SLeila Ghaffari }
21077841947SLeila Ghaffari
21177841947SLeila Ghaffari // Print final time and number of steps
2122b730f8bSJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps));
213cf7a0454SJed Brown PetscCall(TSGetConvergedReason(ts, &reason));
2148fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) {
215cf7a0454SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
216cf7a0454SJed Brown steps, (double)final_time));
21777841947SLeila Ghaffari }
21877841947SLeila Ghaffari
21977841947SLeila Ghaffari // Output numerical values from command line
2202b730f8bSJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
22177841947SLeila Ghaffari
22277841947SLeila Ghaffari // Compare reference solution values with current test run for CI
2238fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2243c9e7ad1SJames Wright PetscCall(RegressionTest(user->app_ctx, Q));
22577841947SLeila Ghaffari }
226ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
22777841947SLeila Ghaffari }
22877841947SLeila Ghaffari
2290de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility
2300de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
2310de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2324de8550aSJed Brown
23377841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
SetupICsFromBinary(MPI_Comm comm,AppCtx app_ctx,Vec Q)23477841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
23577841947SLeila Ghaffari PetscViewer viewer;
2362b730f8bSJeremy L Thompson
237f17d818dSJames Wright PetscFunctionBeginUser;
2382b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
239530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2402b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer));
241ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
24277841947SLeila Ghaffari }
24377841947SLeila Ghaffari
244841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
FreeContextPetsc(void * data)245841e4c73SJed Brown int FreeContextPetsc(void *data) {
2462b730f8bSJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
247841e4c73SJed Brown return CEED_ERROR_SUCCESS;
248841e4c73SJed Brown }
249ef080ff9SJames Wright
250ef080ff9SJames Wright // Return mass qfunction specification for number of components N
CreateMassQFunction(Ceed ceed,CeedInt N,CeedInt q_data_size,CeedQFunction * qf)251ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
252ef080ff9SJames Wright PetscFunctionBeginUser;
253ef080ff9SJames Wright switch (N) {
254ef080ff9SJames Wright case 1:
255a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
256ef080ff9SJames Wright break;
257ef080ff9SJames Wright case 5:
258a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
259ef080ff9SJames Wright break;
260052409adSJames Wright case 7:
261a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
262052409adSJames Wright break;
263ef080ff9SJames Wright case 9:
264a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
265ef080ff9SJames Wright break;
266ef080ff9SJames Wright case 22:
267a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
268ef080ff9SJames Wright break;
269ef080ff9SJames Wright default:
27083ae4962SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
271ef080ff9SJames Wright }
272ef080ff9SJames Wright
273a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
274a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
275a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2767f2a9303SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
277ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
278ef080ff9SJames Wright }
2790df12fefSJames Wright
NodalProjectionDataDestroy(NodalProjectionData context)280999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
281999ff5c7SJames Wright PetscFunctionBeginUser;
282ee4ca9cbSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
283999ff5c7SJames Wright
284999ff5c7SJames Wright PetscCall(DMDestroy(&context->dm));
285999ff5c7SJames Wright PetscCall(KSPDestroy(&context->ksp));
286999ff5c7SJames Wright
287999ff5c7SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
288999ff5c7SJames Wright
289999ff5c7SJames Wright PetscCall(PetscFree(context));
290ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
291999ff5c7SJames Wright }
292999ff5c7SJames Wright
2938b892a05SJames Wright /*
2948b892a05SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
2958b892a05SJames Wright *
2968b892a05SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
2978b892a05SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
2988b892a05SJames Wright *
2998b892a05SJames 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.
3008b892a05SJames Wright *
3018b892a05SJames Wright * @param[in] comm MPI_Comm for the program
3028b892a05SJames Wright * @param[in] path Path to the file
3038b892a05SJames Wright * @param[in] char_array_len Length of the character array that should contain each line
3048b892a05SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file
3058b892a05SJames Wright * @param[out] fp File pointer to the opened file
3068b892a05SJames Wright */
PhastaDatFileOpen(const MPI_Comm comm,const char path[PETSC_MAX_PATH_LEN],const PetscInt char_array_len,PetscInt dims[2],FILE ** fp)307cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
3088b892a05SJames Wright FILE **fp) {
309457e73b2SJames Wright int ndims;
3108b892a05SJames Wright char line[char_array_len];
3118b892a05SJames Wright char **array;
3128b892a05SJames Wright
3138b892a05SJames Wright PetscFunctionBeginUser;
3148b892a05SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp));
3158b892a05SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
3168b892a05SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
317457e73b2SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
3188b892a05SJames Wright
3198b892a05SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
3208b892a05SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array));
321ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3228b892a05SJames Wright }
3238b892a05SJames Wright
3248b892a05SJames Wright /*
3258b892a05SJames Wright * @brief Get the number of rows for the PHASTA file at path.
3268b892a05SJames Wright *
3278b892a05SJames 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.
3288b892a05SJames Wright *
3298b892a05SJames Wright * @param[in] comm MPI_Comm for the program
3308b892a05SJames Wright * @param[in] path Path to the file
3318b892a05SJames Wright * @param[out] nrows Number of rows
3328b892a05SJames Wright */
PhastaDatFileGetNRows(const MPI_Comm comm,const char path[PETSC_MAX_PATH_LEN],PetscInt * nrows)333cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
3348b892a05SJames Wright const PetscInt char_array_len = 512;
3358b892a05SJames Wright PetscInt dims[2];
3368b892a05SJames Wright FILE *fp;
3378b892a05SJames Wright
3388b892a05SJames Wright PetscFunctionBeginUser;
339cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3408b892a05SJames Wright *nrows = dims[0];
3418b892a05SJames Wright PetscCall(PetscFClose(comm, fp));
342ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3438b892a05SJames Wright }
3444cc9442bSJames Wright
PhastaDatFileReadToArrayReal(MPI_Comm comm,const char path[PETSC_MAX_PATH_LEN],PetscReal array[])345cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
346457e73b2SJames Wright PetscInt dims[2];
3474cc9442bSJames Wright FILE *fp;
3484cc9442bSJames Wright const PetscInt char_array_len = 512;
3494cc9442bSJames Wright char line[char_array_len];
3504cc9442bSJames Wright
351f17d818dSJames Wright PetscFunctionBeginUser;
352cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3534cc9442bSJames Wright
3544cc9442bSJames Wright for (PetscInt i = 0; i < dims[0]; i++) {
35538690fecSJames Wright int ndims;
35638690fecSJames Wright char **row_array;
35738690fecSJames Wright
3584cc9442bSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
3594cc9442bSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3600e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
361457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
3624cc9442bSJames Wright
36338690fecSJames Wright for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
36438690fecSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, row_array));
3654cc9442bSJames Wright }
3664cc9442bSJames Wright
3674cc9442bSJames Wright PetscCall(PetscFClose(comm, fp));
368ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
3694cc9442bSJames Wright }
37075d1798cSJames Wright
3712e31f396SJames Wright // Print information about the given simulation run
PrintRunInfo(User user,Physics phys_ctx,ProblemData problem,TS ts)3727bc7b61fSJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
373a424bcd0SJames Wright Ceed ceed = user->ceed;
3747bc7b61fSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts);
3757bc7b61fSJames Wright
3762e31f396SJames Wright PetscFunctionBeginUser;
3772e31f396SJames Wright // Header and rank
3782e31f396SJames Wright char host_name[PETSC_MAX_PATH_LEN];
3792e31f396SJames Wright PetscMPIInt rank, comm_size;
3802e31f396SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name));
3812e31f396SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank));
3822e31f396SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size));
3832e31f396SJames Wright PetscCall(PetscPrintf(comm,
3842e31f396SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
3852e31f396SJames Wright " MPI:\n"
3862e31f396SJames Wright " Host Name : %s\n"
3872e31f396SJames Wright " Total ranks : %d\n",
3882e31f396SJames Wright host_name, comm_size));
3892e31f396SJames Wright
3902e31f396SJames Wright // Problem specific info
391367c849eSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx));
3922e31f396SJames Wright
3932e31f396SJames Wright // libCEED
3942e31f396SJames Wright const char *used_resource;
3952e31f396SJames Wright CeedMemType mem_type_backend;
396a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
397a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
3982e31f396SJames Wright PetscCall(PetscPrintf(comm,
3992e31f396SJames Wright " libCEED:\n"
4002e31f396SJames Wright " libCEED Backend : %s\n"
4012e31f396SJames Wright " libCEED Backend MemType : %s\n",
4022e31f396SJames Wright used_resource, CeedMemTypes[mem_type_backend]));
4032e31f396SJames Wright // PETSc
4047bc7b61fSJames Wright VecType vec_type;
4052e31f396SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
4062e31f396SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0';
4072e31f396SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
4082e31f396SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type));
4092e31f396SJames Wright PetscCall(PetscPrintf(comm,
4102e31f396SJames Wright " PETSc:\n"
4112e31f396SJames Wright " Box Faces : %s\n"
4122e31f396SJames Wright " DM VecType : %s\n"
4132e31f396SJames Wright " Time Stepping Scheme : %s\n",
4147bc7b61fSJames Wright box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
4157bc7b61fSJames Wright {
4167bc7b61fSJames Wright char pmat_type_str[PETSC_MAX_PATH_LEN];
4177bc7b61fSJames Wright MatType amat_type, pmat_type;
4187bc7b61fSJames Wright Mat Amat, Pmat;
4197bc7b61fSJames Wright TSIJacobianFn *ijacob_function;
4207bc7b61fSJames Wright
4217bc7b61fSJames Wright PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
4227bc7b61fSJames Wright PetscCall(MatGetType(Amat, &amat_type));
4237bc7b61fSJames Wright PetscCall(MatGetType(Pmat, &pmat_type));
4247bc7b61fSJames Wright
4257bc7b61fSJames Wright PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
4267bc7b61fSJames Wright if (!strcmp(pmat_type, MATCEED)) {
4277bc7b61fSJames Wright MatType pmat_coo_type;
4287bc7b61fSJames Wright char pmat_coo_type_str[PETSC_MAX_PATH_LEN];
4297bc7b61fSJames Wright
4307bc7b61fSJames Wright PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
4317bc7b61fSJames Wright PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
4327bc7b61fSJames Wright PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
4337bc7b61fSJames Wright }
4347bc7b61fSJames Wright if (ijacob_function) {
4357bc7b61fSJames Wright PetscCall(PetscPrintf(comm,
4367bc7b61fSJames Wright " IJacobian A MatType : %s\n"
4377bc7b61fSJames Wright " IJacobian P MatType : %s\n",
4387bc7b61fSJames Wright amat_type, pmat_type_str));
4397bc7b61fSJames Wright }
4407bc7b61fSJames Wright }
4412e31f396SJames Wright if (user->app_ctx->cont_steps) {
4422e31f396SJames Wright PetscCall(PetscPrintf(comm,
4432e31f396SJames Wright " Continue:\n"
4442e31f396SJames Wright " Filename: : %s\n"
4452e31f396SJames Wright " Step: : %" PetscInt_FMT "\n"
4462e31f396SJames Wright " Time: : %g\n",
4472e31f396SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
4482e31f396SJames Wright }
4492e31f396SJames Wright // Mesh
4502e31f396SJames Wright const PetscInt num_comp_q = 5;
4512e31f396SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs;
4522e31f396SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
4532e31f396SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
4542e31f396SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
4552e31f396SJames Wright PetscCall(PetscPrintf(comm,
4562e31f396SJames Wright " Mesh:\n"
4572e31f396SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n"
4582e31f396SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
4592e31f396SJames Wright " Global DoFs : %" PetscInt_FMT "\n"
4602e31f396SJames Wright " DoFs per node : %" PetscInt_FMT "\n"
461ce11f295SJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n",
462ce11f295SJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
4632e31f396SJames Wright // -- Get Partition Statistics
4642e31f396SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n"));
4652e31f396SJames Wright {
4662e31f396SJames Wright PetscInt *gather_buffer = NULL;
467ce11f295SJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
4682e31f396SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
4692e31f396SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
4702e31f396SJames Wright
471ce11f295SJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4722e31f396SJames Wright if (!rank) {
4732e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
4742e31f396SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min
4752e31f396SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max
4762e31f396SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median
4772e31f396SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
4781a8516d0SJames Wright PetscCall(PetscPrintf(comm,
4791a8516d0SJames Wright " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
4801a8516d0SJames Wright num_comp_q, part_owned_dofs[0] / num_comp_q, part_owned_dofs[1] / num_comp_q, part_owned_dofs[2] / num_comp_q,
4811a8516d0SJames Wright part_owned_dof_ratio));
4822e31f396SJames Wright }
4832e31f396SJames Wright
484ce11f295SJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4852e31f396SJames Wright if (!rank) {
4862e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
4872e31f396SJames Wright part_local_dofs[0] = gather_buffer[0]; // min
4882e31f396SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max
4892e31f396SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median
4902e31f396SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
4911a8516d0SJames Wright PetscCall(PetscPrintf(comm,
4921a8516d0SJames Wright " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
4931a8516d0SJames Wright num_comp_q, part_local_dofs[0] / num_comp_q, part_local_dofs[1] / num_comp_q, part_local_dofs[2] / num_comp_q,
4941a8516d0SJames Wright part_local_dof_ratio));
4952e31f396SJames Wright }
4962e31f396SJames Wright
49765ba01baSJames Wright if (comm_size != 1) {
498ce11f295SJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
4992e31f396SJames Wright {
5002e31f396SJames Wright PetscSF sf;
5018a310472SJeremy L Thompson PetscMPIInt nrranks, niranks;
502ce11f295SJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc;
503ce11f295SJames Wright const PetscMPIInt *rranks, *iranks;
5048a310472SJeremy L Thompson
5052e31f396SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf));
5062e31f396SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
507ce11f295SJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
5082e31f396SJames Wright for (PetscInt i = 0; i < nrranks; i++) {
5092e31f396SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers
5102e31f396SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i];
511ce11f295SJames Wright num_ghost_interface_ranks++;
512ce11f295SJames Wright }
513ce11f295SJames Wright for (PetscInt i = 0; i < niranks; i++) {
514ce11f295SJames Wright if (iranks[i] == rank) continue;
515ce11f295SJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
516ce11f295SJames Wright num_owned_interface_ranks++;
5172e31f396SJames Wright }
5182e31f396SJames Wright }
519ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
5202e31f396SJames Wright if (!rank) {
5212e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
522ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min
523ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max
524ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median
525ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
5261a8516d0SJames Wright PetscCall(PetscPrintf(comm,
5271a8516d0SJames Wright " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT
5281a8516d0SJames Wright ", %f\n",
52965ba01baSJames 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,
53065ba01baSJames Wright part_shared_dof_ratio));
531ce11f295SJames Wright }
532ce11f295SJames Wright
533ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
534ce11f295SJames Wright if (!rank) {
535ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
536ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min
537ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max
538ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median
539ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
540ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
541ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
542ce11f295SJames Wright }
543ce11f295SJames Wright
544ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
545ce11f295SJames Wright if (!rank) {
546ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
547ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min
548ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max
549ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median
550ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
5511a8516d0SJames Wright PetscCall(PetscPrintf(comm,
5521a8516d0SJames Wright " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT
5531a8516d0SJames Wright ", %f\n",
55465ba01baSJames 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,
55565ba01baSJames Wright part_shared_dof_ratio));
556ce11f295SJames Wright }
557ce11f295SJames Wright
558ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
559ce11f295SJames Wright if (!rank) {
560ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer));
561ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min
562ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max
563ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median
564ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
565ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
566ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
5672e31f396SJames Wright }
56865ba01baSJames Wright }
5692e31f396SJames Wright
5702e31f396SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer));
5712e31f396SJames Wright }
5722e31f396SJames Wright PetscFunctionReturn(PETSC_SUCCESS);
5732e31f396SJames Wright }
574