xref: /libCEED/examples/fluids/src/misc.c (revision ce11f295d56c538e2836a7b5476020ae6c83aed9)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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) {
2077841947SLeila Ghaffari   PetscFunctionBeginUser;
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari   // ---------------------------------------------------------------------------
23a0add3c9SJed Brown   // Update time for evaluation
2477841947SLeila Ghaffari   // ---------------------------------------------------------------------------
255263e9c6SJames Wright   if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time);
2677841947SLeila Ghaffari 
2777841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2877841947SLeila Ghaffari   // ICs
2977841947SLeila Ghaffari   // ---------------------------------------------------------------------------
3077841947SLeila Ghaffari   // -- CEED Restriction
3177841947SLeila Ghaffari   CeedVector q0_ceed;
3277841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
355263e9c6SJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
3677841947SLeila Ghaffari 
3777841947SLeila Ghaffari   // ---------------------------------------------------------------------------
3877841947SLeila Ghaffari   // Fix multiplicity for output of ICs
3977841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4077841947SLeila Ghaffari   // -- CEED Restriction
4177841947SLeila Ghaffari   CeedVector mult_vec;
4277841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
4377841947SLeila Ghaffari 
4477841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
4577841947SLeila Ghaffari   PetscMemType m_mem_type;
4677841947SLeila Ghaffari   Vec          multiplicity_loc;
472b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
48c798d55aSJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
4977841947SLeila Ghaffari 
5077841947SLeila Ghaffari   // -- Get multiplicity
5177841947SLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
5277841947SLeila Ghaffari 
5377841947SLeila Ghaffari   // -- Restore vectors
54c798d55aSJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
5577841947SLeila Ghaffari 
5677841947SLeila Ghaffari   // -- Local-to-Global
5777841947SLeila Ghaffari   Vec multiplicity;
582b730f8bSJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
592b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
602b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
6177841947SLeila Ghaffari 
6277841947SLeila Ghaffari   // -- Fix multiplicity
632b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
642b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
6577841947SLeila Ghaffari 
6677841947SLeila Ghaffari   // -- Restore vectors
672b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
682b730f8bSJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
6977841947SLeila Ghaffari 
7077841947SLeila Ghaffari   // Cleanup
7177841947SLeila Ghaffari   CeedVectorDestroy(&mult_vec);
7277841947SLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
7377841947SLeila Ghaffari 
74ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
7577841947SLeila Ghaffari }
7677841947SLeila Ghaffari 
772b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
782b730f8bSJeremy L Thompson                                              Vec grad_FVM) {
795571c6fdSJames Wright   Vec Qbc, boundary_mask;
8077841947SLeila Ghaffari   PetscFunctionBegin;
8177841947SLeila Ghaffari 
823722cd23SJames Wright   // Mask (zero) Strong BC entries
835571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
845571c6fdSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
855571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
865571c6fdSJames Wright 
872b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
882b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
892b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
9077841947SLeila Ghaffari 
91ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
9277841947SLeila Ghaffari }
9377841947SLeila Ghaffari 
94530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
95530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
960de6a49fSJames Wright   PetscInt   file_step_number;
970de6a49fSJames Wright   PetscInt32 token;
98530ad8c4SKenneth E. Jansen   PetscReal  file_time;
99530ad8c4SKenneth E. Jansen   PetscFunctionBeginUser;
100530ad8c4SKenneth E. Jansen 
101530ad8c4SKenneth E. Jansen   // Attempt
1020de6a49fSJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
1030de6a49fSJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
1040de6a49fSJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
105530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
106530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
107530ad8c4SKenneth E. Jansen     if (time) *time = file_time;
108530ad8c4SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
109530ad8c4SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
110530ad8c4SKenneth E. Jansen     PetscInt length, N;
111530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
112530ad8c4SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
113530ad8c4SKenneth 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);
114530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
115530ad8c4SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
116530ad8c4SKenneth E. Jansen 
117530ad8c4SKenneth E. Jansen   // Load Q from existent solution
118530ad8c4SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
119530ad8c4SKenneth E. Jansen 
120ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
121530ad8c4SKenneth E. Jansen }
122530ad8c4SKenneth E. Jansen 
12377841947SLeila Ghaffari // Compare reference solution values with current test run for CI
12477841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
12577841947SLeila Ghaffari   Vec         Qref;
12677841947SLeila Ghaffari   PetscViewer viewer;
12777841947SLeila Ghaffari   PetscReal   error, Qrefnorm;
128530ad8c4SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
12977841947SLeila Ghaffari   PetscFunctionBegin;
13077841947SLeila Ghaffari 
13177841947SLeila Ghaffari   // Read reference file
1322b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
133530ad8c4SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
134530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
13577841947SLeila Ghaffari 
13677841947SLeila Ghaffari   // Compute error with respect to reference solution
1372b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1382b730f8bSJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1392b730f8bSJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1402b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
14177841947SLeila Ghaffari 
14277841947SLeila Ghaffari   // Check error
14377841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
1442b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
14577841947SLeila Ghaffari   }
14677841947SLeila Ghaffari 
14777841947SLeila Ghaffari   // Cleanup
1482b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1492b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Qref));
15077841947SLeila Ghaffari 
151ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
15277841947SLeila Ghaffari }
15377841947SLeila Ghaffari 
15477841947SLeila Ghaffari // Get error for problems with exact solutions
1552b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
15677841947SLeila Ghaffari   PetscInt  loc_nodes;
15777841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
15877841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
15977841947SLeila Ghaffari   PetscFunctionBegin;
16077841947SLeila Ghaffari 
16177841947SLeila Ghaffari   // Get exact solution at final time
1622b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1632b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1642b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1652b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
16677841947SLeila Ghaffari 
16777841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1682b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1692b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1702b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
17177841947SLeila Ghaffari 
17277841947SLeila Ghaffari   // Compute relative error
17377841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
17477841947SLeila Ghaffari 
17577841947SLeila Ghaffari   // Output relative error
1762b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
17777841947SLeila Ghaffari   // Cleanup
1782b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1792b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
18077841947SLeila Ghaffari 
181ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
18277841947SLeila Ghaffari }
18377841947SLeila Ghaffari 
18477841947SLeila Ghaffari // Post-processing
1852b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
18677841947SLeila Ghaffari   PetscInt          steps;
187cf7a0454SJed Brown   TSConvergedReason reason;
18877841947SLeila Ghaffari   PetscFunctionBegin;
18977841947SLeila Ghaffari 
19077841947SLeila Ghaffari   // Print relative error
1918fb33541SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1922b730f8bSJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
19377841947SLeila Ghaffari   }
19477841947SLeila Ghaffari 
19577841947SLeila Ghaffari   // Print final time and number of steps
1962b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
197cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1988fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
199cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
200cf7a0454SJed Brown                           steps, (double)final_time));
20177841947SLeila Ghaffari   }
20277841947SLeila Ghaffari 
20377841947SLeila Ghaffari   // Output numerical values from command line
2042b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
20577841947SLeila Ghaffari 
20677841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
2078fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2082b730f8bSJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
20977841947SLeila Ghaffari   }
210ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21177841947SLeila Ghaffari }
21277841947SLeila Ghaffari 
2130de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
2140de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
2150de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2164de8550aSJed Brown 
21777841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
21877841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
21977841947SLeila Ghaffari   PetscViewer viewer;
2202b730f8bSJeremy L Thompson 
22177841947SLeila Ghaffari   PetscFunctionBegin;
22277841947SLeila Ghaffari 
2232b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
224530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2252b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
22677841947SLeila Ghaffari 
227ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
22877841947SLeila Ghaffari }
22977841947SLeila Ghaffari 
23077841947SLeila Ghaffari // Record boundary values from initial condition
23177841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2325571c6fdSJames Wright   Vec Qbc, boundary_mask;
23377841947SLeila Ghaffari   PetscFunctionBegin;
23477841947SLeila Ghaffari 
2352b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2362b730f8bSJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2372b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2382b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2392b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2402b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2412b730f8bSJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
24277841947SLeila Ghaffari 
2435571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2445571c6fdSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2455571c6fdSJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2465571c6fdSJames Wright   PetscCall(VecSet(Q, 1.0));
2475571c6fdSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2485571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2495571c6fdSJames Wright 
250ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25177841947SLeila Ghaffari }
252841e4c73SJed Brown 
253841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
254841e4c73SJed Brown int FreeContextPetsc(void *data) {
2552b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
256841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
257841e4c73SJed Brown }
258ef080ff9SJames Wright 
259ef080ff9SJames Wright // Return mass qfunction specification for number of components N
260ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
261ef080ff9SJames Wright   PetscFunctionBeginUser;
262ef080ff9SJames Wright 
263ef080ff9SJames Wright   switch (N) {
264ef080ff9SJames Wright     case 1:
26583ae4962SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf);
266ef080ff9SJames Wright       break;
267ef080ff9SJames Wright     case 5:
26883ae4962SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf);
269ef080ff9SJames Wright       break;
270052409adSJames Wright     case 7:
27183ae4962SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf);
272052409adSJames Wright       break;
273ef080ff9SJames Wright     case 9:
27483ae4962SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf);
275ef080ff9SJames Wright       break;
276ef080ff9SJames Wright     case 22:
27783ae4962SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf);
278ef080ff9SJames Wright       break;
279ef080ff9SJames Wright     default:
28083ae4962SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
281ef080ff9SJames Wright   }
282ef080ff9SJames Wright 
283ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
284ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
285ef080ff9SJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
286ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
287ef080ff9SJames Wright }
2880df12fefSJames Wright 
2890df12fefSJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space
2900df12fefSJames Wright  *
2910df12fefSJames Wright  * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum.
2920df12fefSJames Wright  *
2930df12fefSJames Wright  * @param[in]  source_vec    Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc
2940df12fefSJames Wright  * @param[out] target_vec    Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc
2950df12fefSJames Wright  * @param[in]  rhs_matop_ctx MatopApplyContext for performing the RHS evaluation
2960df12fefSJames Wright  * @param[in]  ksp           KSP for solving the consistent projection problem
2970df12fefSJames Wright  */
2984021610dSJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp) {
2990df12fefSJames Wright   PetscFunctionBeginUser;
3000df12fefSJames Wright 
3014021610dSJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(source_vec, target_vec, rhs_matop_ctx));
3020df12fefSJames Wright   PetscCall(KSPSolve(ksp, target_vec, target_vec));
3030df12fefSJames Wright 
304ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3050df12fefSJames Wright }
3068b892a05SJames Wright 
307999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
308999ff5c7SJames Wright   PetscFunctionBeginUser;
309ee4ca9cbSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
310999ff5c7SJames Wright 
311999ff5c7SJames Wright   PetscCall(DMDestroy(&context->dm));
312999ff5c7SJames Wright   PetscCall(KSPDestroy(&context->ksp));
313999ff5c7SJames Wright 
314999ff5c7SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
315999ff5c7SJames Wright 
316999ff5c7SJames Wright   PetscCall(PetscFree(context));
317999ff5c7SJames Wright 
318ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
319999ff5c7SJames Wright }
320999ff5c7SJames Wright 
3218b892a05SJames Wright /*
3228b892a05SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
3238b892a05SJames Wright  *
3248b892a05SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
3258b892a05SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
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[in]  char_array_len Length of the character array that should contain each line
3328b892a05SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
3338b892a05SJames Wright  * @param[out] fp File        pointer to the opened file
3348b892a05SJames Wright  */
3358b892a05SJames Wright PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
3368b892a05SJames Wright                                  FILE **fp) {
337457e73b2SJames Wright   int    ndims;
3388b892a05SJames Wright   char   line[char_array_len];
3398b892a05SJames Wright   char **array;
3408b892a05SJames Wright 
3418b892a05SJames Wright   PetscFunctionBeginUser;
3428b892a05SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
3438b892a05SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
3448b892a05SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
345457e73b2SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
3468b892a05SJames Wright 
3478b892a05SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
3488b892a05SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
3498b892a05SJames Wright 
350ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3518b892a05SJames Wright }
3528b892a05SJames Wright 
3538b892a05SJames Wright /*
3548b892a05SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
3558b892a05SJames Wright  *
3568b892a05SJames 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.
3578b892a05SJames Wright  *
3588b892a05SJames Wright  * @param[in]  comm  MPI_Comm for the program
3598b892a05SJames Wright  * @param[in]  path  Path to the file
3608b892a05SJames Wright  * @param[out] nrows Number of rows
3618b892a05SJames Wright  */
3628b892a05SJames Wright PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
3638b892a05SJames Wright   const PetscInt char_array_len = 512;
3648b892a05SJames Wright   PetscInt       dims[2];
3658b892a05SJames Wright   FILE          *fp;
3668b892a05SJames Wright 
3678b892a05SJames Wright   PetscFunctionBeginUser;
3688b892a05SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
3698b892a05SJames Wright   *nrows = dims[0];
3708b892a05SJames Wright   PetscCall(PetscFClose(comm, fp));
3718b892a05SJames Wright 
372ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3738b892a05SJames Wright }
3744cc9442bSJames Wright 
3754cc9442bSJames Wright PetscErrorCode PHASTADatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
376457e73b2SJames Wright   PetscInt       dims[2];
377457e73b2SJames Wright   int            ndims;
3784cc9442bSJames Wright   FILE          *fp;
3794cc9442bSJames Wright   const PetscInt char_array_len = 512;
3804cc9442bSJames Wright   char           line[char_array_len];
3814cc9442bSJames Wright   char         **row_array;
3824cc9442bSJames Wright   PetscFunctionBeginUser;
3834cc9442bSJames Wright 
3844cc9442bSJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
3854cc9442bSJames Wright 
3864cc9442bSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
3874cc9442bSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
3884cc9442bSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3890e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
390457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
3914cc9442bSJames Wright 
3924cc9442bSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
3934cc9442bSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
3944cc9442bSJames Wright     }
3954cc9442bSJames Wright   }
3964cc9442bSJames Wright 
3974cc9442bSJames Wright   PetscCall(PetscFClose(comm, fp));
3984cc9442bSJames Wright 
399ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4004cc9442bSJames Wright }
40175d1798cSJames Wright 
40275d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
40375d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
40475d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
40575d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
40675d1798cSJames Wright static PetscClassId libCEED_classid;
40775d1798cSJames Wright 
40875d1798cSJames Wright PetscErrorCode RegisterLogEvents() {
40975d1798cSJames Wright   PetscFunctionBeginUser;
41075d1798cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
41175d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
41275d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
41375d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
41475d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
415ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
41675d1798cSJames Wright }
417457e73b2SJames Wright 
418457e73b2SJames Wright /**
419457e73b2SJames Wright   @brief Translate array of CeedInt to PetscInt.
420457e73b2SJames Wright     If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`.
421457e73b2SJames Wright     Caller is responsible for freeing `array_petsc` with `free()`.
422457e73b2SJames Wright 
423457e73b2SJames Wright   @param[in]      num_entries  Number of array entries
424457e73b2SJames Wright   @param[in,out]  array_ceed   Array of CeedInts
425457e73b2SJames Wright   @param[out]     array_petsc  Array of PetscInts
426457e73b2SJames Wright **/
427457e73b2SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
428457e73b2SJames Wright   CeedInt  int_c = 0;
429457e73b2SJames Wright   PetscInt int_p = 0;
430457e73b2SJames Wright 
431457e73b2SJames Wright   PetscFunctionBeginUser;
432457e73b2SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
433457e73b2SJames Wright     *array_petsc = (PetscInt *)*array_ceed;
434457e73b2SJames Wright   } else {
435457e73b2SJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
436457e73b2SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
437457e73b2SJames Wright     free(*array_ceed);
438457e73b2SJames Wright   }
439457e73b2SJames Wright   *array_ceed = NULL;
440457e73b2SJames Wright 
441457e73b2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
442457e73b2SJames Wright }
443457e73b2SJames Wright 
444457e73b2SJames Wright /**
445457e73b2SJames Wright   @brief Translate array of PetscInt to CeedInt.
446457e73b2SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
447457e73b2SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
448457e73b2SJames Wright 
449457e73b2SJames Wright   @param[in]      num_entries  Number of array entries
450457e73b2SJames Wright   @param[in,out]  array_petsc  Array of PetscInts
451457e73b2SJames Wright   @param[out]     array_ceed   Array of CeedInts
452457e73b2SJames Wright **/
453457e73b2SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
454457e73b2SJames Wright   CeedInt  int_c = 0;
455457e73b2SJames Wright   PetscInt int_p = 0;
456457e73b2SJames Wright 
457457e73b2SJames Wright   PetscFunctionBeginUser;
458457e73b2SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
459457e73b2SJames Wright     *array_ceed = (CeedInt *)*array_petsc;
460457e73b2SJames Wright   } else {
461457e73b2SJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
462457e73b2SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
463457e73b2SJames Wright     PetscCall(PetscFree(*array_petsc));
464457e73b2SJames Wright   }
465457e73b2SJames Wright   *array_petsc = NULL;
466457e73b2SJames Wright 
467457e73b2SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
468457e73b2SJames Wright }
4692e31f396SJames Wright 
4702e31f396SJames Wright // Print information about the given simulation run
4712e31f396SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
4722e31f396SJames Wright   PetscFunctionBeginUser;
4732e31f396SJames Wright   // Header and rank
4742e31f396SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
4752e31f396SJames Wright   PetscMPIInt rank, comm_size;
4762e31f396SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
4772e31f396SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4782e31f396SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
4792e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4802e31f396SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
4812e31f396SJames Wright                         "  MPI:\n"
4822e31f396SJames Wright                         "    Host Name                          : %s\n"
4832e31f396SJames Wright                         "    Total ranks                        : %d\n",
4842e31f396SJames Wright                         host_name, comm_size));
4852e31f396SJames Wright 
4862e31f396SJames Wright   // Problem specific info
4872e31f396SJames Wright   PetscCall(problem->print_info(problem, user->app_ctx));
4882e31f396SJames Wright 
4892e31f396SJames Wright   // libCEED
4902e31f396SJames Wright   const char *used_resource;
4912e31f396SJames Wright   CeedMemType mem_type_backend;
4922e31f396SJames Wright   CeedGetResource(user->ceed, &used_resource);
4932e31f396SJames Wright   CeedGetPreferredMemType(user->ceed, &mem_type_backend);
4942e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4952e31f396SJames Wright                         "  libCEED:\n"
4962e31f396SJames Wright                         "    libCEED Backend                    : %s\n"
4972e31f396SJames Wright                         "    libCEED Backend MemType            : %s\n",
4982e31f396SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
4992e31f396SJames Wright   // PETSc
5002e31f396SJames Wright   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
5012e31f396SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
5022e31f396SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
5032e31f396SJames Wright   MatType mat_type;
5042e31f396SJames Wright   VecType vec_type;
5052e31f396SJames Wright   PetscCall(DMGetMatType(user->dm, &mat_type));
5062e31f396SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
5072e31f396SJames Wright   PetscCall(PetscPrintf(comm,
5082e31f396SJames Wright                         "  PETSc:\n"
5092e31f396SJames Wright                         "    Box Faces                          : %s\n"
5102e31f396SJames Wright                         "    DM MatType                         : %s\n"
5112e31f396SJames Wright                         "    DM VecType                         : %s\n"
5122e31f396SJames Wright                         "    Time Stepping Scheme               : %s\n",
5132e31f396SJames Wright                         box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
5142e31f396SJames Wright   if (user->app_ctx->cont_steps) {
5152e31f396SJames Wright     PetscCall(PetscPrintf(comm,
5162e31f396SJames Wright                           "  Continue:\n"
5172e31f396SJames Wright                           "    Filename:                          : %s\n"
5182e31f396SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
5192e31f396SJames Wright                           "    Time:                              : %g\n",
5202e31f396SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
5212e31f396SJames Wright   }
5222e31f396SJames Wright   // Mesh
5232e31f396SJames Wright   const PetscInt num_comp_q = 5;
5242e31f396SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
5252e31f396SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
5262e31f396SJames Wright   // -- Get global size
5272e31f396SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
5282e31f396SJames Wright   // -- Get local size
5292e31f396SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
5302e31f396SJames Wright   PetscCall(PetscPrintf(comm,
5312e31f396SJames Wright                         "  Mesh:\n"
5322e31f396SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
5332e31f396SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
5342e31f396SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
5352e31f396SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
536*ce11f295SJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
537*ce11f295SJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
5382e31f396SJames Wright   // -- Get Partition Statistics
5392e31f396SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
5402e31f396SJames Wright   {
5412e31f396SJames Wright     PetscInt *gather_buffer = NULL;
542*ce11f295SJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
5432e31f396SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
5442e31f396SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
5452e31f396SJames Wright 
546*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
5472e31f396SJames Wright     if (!rank) {
5482e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
5492e31f396SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
5502e31f396SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
5512e31f396SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
5522e31f396SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
553*ce11f295SJames Wright       PetscCall(PetscPrintf(
554*ce11f295SJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
5552e31f396SJames 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));
5562e31f396SJames Wright     }
5572e31f396SJames Wright 
558*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
5592e31f396SJames Wright     if (!rank) {
5602e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
5612e31f396SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
5622e31f396SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
5632e31f396SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
5642e31f396SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
565*ce11f295SJames Wright       PetscCall(PetscPrintf(
566*ce11f295SJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
5672e31f396SJames 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));
5682e31f396SJames Wright     }
5692e31f396SJames Wright 
570*ce11f295SJames Wright     PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
5712e31f396SJames Wright     {
5722e31f396SJames Wright       PetscSF            sf;
573*ce11f295SJames Wright       PetscInt           nrranks, niranks;
574*ce11f295SJames Wright       const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
575*ce11f295SJames Wright       const PetscMPIInt *rranks, *iranks;
5762e31f396SJames Wright       PetscCall(DMGetSectionSF(user->dm, &sf));
5772e31f396SJames Wright       PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
578*ce11f295SJames Wright       PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
5792e31f396SJames Wright       for (PetscInt i = 0; i < nrranks; i++) {
5802e31f396SJames Wright         if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
5812e31f396SJames Wright         num_remote_roots_total += roffset[i + 1] - roffset[i];
582*ce11f295SJames Wright         num_ghost_interface_ranks++;
583*ce11f295SJames Wright       }
584*ce11f295SJames Wright       for (PetscInt i = 0; i < niranks; i++) {
585*ce11f295SJames Wright         if (iranks[i] == rank) continue;
586*ce11f295SJames Wright         num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
587*ce11f295SJames Wright         num_owned_interface_ranks++;
5882e31f396SJames Wright       }
5892e31f396SJames Wright     }
590*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
5912e31f396SJames Wright     if (!rank) {
5922e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
593*ce11f295SJames Wright       part_boundary_dofs[0]           = gather_buffer[0];              // min
594*ce11f295SJames Wright       part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
595*ce11f295SJames Wright       part_boundary_dofs[2]           = gather_buffer[median_index];   // median
596*ce11f295SJames Wright       PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
597*ce11f295SJames Wright       PetscCall(PetscPrintf(
598*ce11f295SJames Wright           comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
599*ce11f295SJames Wright           part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, part_shared_dof_ratio));
600*ce11f295SJames Wright     }
601*ce11f295SJames Wright 
602*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
603*ce11f295SJames Wright     if (!rank) {
604*ce11f295SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
605*ce11f295SJames Wright       part_neighbors[0]              = gather_buffer[0];              // min
606*ce11f295SJames Wright       part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
607*ce11f295SJames Wright       part_neighbors[2]              = gather_buffer[median_index];   // median
608*ce11f295SJames Wright       PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
609*ce11f295SJames Wright       PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
610*ce11f295SJames Wright                             part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
611*ce11f295SJames Wright     }
612*ce11f295SJames Wright 
613*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
614*ce11f295SJames Wright     if (!rank) {
615*ce11f295SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
616*ce11f295SJames Wright       part_boundary_dofs[0]           = gather_buffer[0];              // min
617*ce11f295SJames Wright       part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
618*ce11f295SJames Wright       part_boundary_dofs[2]           = gather_buffer[median_index];   // median
619*ce11f295SJames Wright       PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
620*ce11f295SJames Wright       PetscCall(PetscPrintf(
621*ce11f295SJames Wright           comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
622*ce11f295SJames Wright           part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, part_shared_dof_ratio));
623*ce11f295SJames Wright     }
624*ce11f295SJames Wright 
625*ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
626*ce11f295SJames Wright     if (!rank) {
627*ce11f295SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
628*ce11f295SJames Wright       part_neighbors[0]              = gather_buffer[0];              // min
629*ce11f295SJames Wright       part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
630*ce11f295SJames Wright       part_neighbors[2]              = gather_buffer[median_index];   // median
631*ce11f295SJames Wright       PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
632*ce11f295SJames Wright       PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
633*ce11f295SJames Wright                             part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
6342e31f396SJames Wright     }
6352e31f396SJames Wright 
6362e31f396SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
6372e31f396SJames Wright   }
6382e31f396SJames Wright 
6392e31f396SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6402e31f396SJames Wright }
641