1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Miscellaneous utility functions 10a515125bSLeila Ghaffari 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdm.h> 13926a6279SJames Wright #include <petscsf.h> 14e419654dSJeremy L Thompson #include <petscts.h> 15e419654dSJeremy L Thompson 16a515125bSLeila Ghaffari #include "../navierstokes.h" 179f59f36eSJames Wright #include "../qfunctions/mass.h" 18a515125bSLeila Ghaffari 192b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 20b4c37c5cSJames Wright Ceed ceed = user->ceed; 21b2948607SJames Wright CeedVector mult_vec; 22b2948607SJames Wright PetscMemType m_mem_type; 23b2948607SJames Wright Vec Multiplicity, Multiplicity_loc; 24b2948607SJames Wright 25a515125bSLeila Ghaffari PetscFunctionBeginUser; 26b4c37c5cSJames Wright if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 278f18bb8bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 28a515125bSLeila Ghaffari 29b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 30a515125bSLeila Ghaffari 31a515125bSLeila Ghaffari // -- Get multiplicity 32b2948607SJames Wright PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 33a7dac1d5SJames Wright PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 34b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 35a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 36a515125bSLeila Ghaffari 37b2948607SJames Wright PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 38b2948607SJames Wright PetscCall(VecZeroEntries(Multiplicity)); 39b2948607SJames Wright PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 40a515125bSLeila Ghaffari 41a515125bSLeila Ghaffari // -- Fix multiplicity 42b2948607SJames Wright PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 43b2948607SJames Wright PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 44a515125bSLeila Ghaffari 45b2948607SJames Wright PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 46b2948607SJames Wright PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 47b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 48d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 49a515125bSLeila Ghaffari } 50a515125bSLeila Ghaffari 51c56e8d5bSJames Wright // Record boundary values from initial condition 52c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 53c56e8d5bSJames Wright PetscFunctionBeginUser; 54313f2f1eSJames Wright { // Capture initial condition values in Qbc 55313f2f1eSJames Wright Vec Qbc; 56313f2f1eSJames Wright 57c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 58c56e8d5bSJames Wright PetscCall(VecCopy(Q_loc, Qbc)); 59c56e8d5bSJames Wright PetscCall(VecZeroEntries(Q_loc)); 60c56e8d5bSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 61c56e8d5bSJames Wright PetscCall(VecAXPY(Qbc, -1., Q_loc)); 62c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 63313f2f1eSJames Wright } 64c56e8d5bSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 65c56e8d5bSJames Wright 66313f2f1eSJames Wright { // Set boundary mask to zero out essential BCs 67313f2f1eSJames Wright Vec boundary_mask, ones; 68313f2f1eSJames Wright 69c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 70313f2f1eSJames Wright PetscCall(DMGetGlobalVector(dm, &ones)); 71c56e8d5bSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 72313f2f1eSJames Wright PetscCall(VecSet(ones, 1.0)); 73313f2f1eSJames Wright PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask)); 74c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 75313f2f1eSJames Wright PetscCall(DMRestoreGlobalVector(dm, &ones)); 76313f2f1eSJames Wright } 77c56e8d5bSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 78c56e8d5bSJames Wright } 79c56e8d5bSJames Wright 80c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 812b916ea7SJeremy L Thompson Vec grad_FVM) { 829d437337SJames Wright Vec Qbc, boundary_mask; 83a515125bSLeila Ghaffari 8406f41313SJames Wright PetscFunctionBeginUser; 852eb7bf1fSJames Wright // Mask (zero) Strong BC entries 869d437337SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 879d437337SJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 889d437337SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 899d437337SJames Wright 902b916ea7SJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 912b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 922b916ea7SJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 93d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 94a515125bSLeila Ghaffari } 95a515125bSLeila Ghaffari 96990f1db0SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) { 97990f1db0SJed Brown PetscFunctionBeginUser; 98990f1db0SJed Brown if (file_type == PETSC_INT32) { 99990f1db0SJed Brown PetscInt32 val; 100990f1db0SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32)); 101990f1db0SJed Brown *out = val; 102990f1db0SJed Brown } else if (file_type == PETSC_INT64) { 103990f1db0SJed Brown PetscInt64 val; 104990f1db0SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64)); 105990f1db0SJed Brown *out = val; 106990f1db0SJed Brown } else { 107990f1db0SJed Brown PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT)); 108990f1db0SJed Brown } 109990f1db0SJed Brown PetscFunctionReturn(PETSC_SUCCESS); 110990f1db0SJed Brown } 111990f1db0SJed Brown 112e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 113e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 114e1233009SJames Wright PetscInt file_step_number; 115e1233009SJames Wright PetscInt32 token; 116e7754af5SKenneth E. Jansen PetscReal file_time; 117990f1db0SJed Brown PetscDataType file_type = PETSC_INT32; 118e7754af5SKenneth E. Jansen 11906f41313SJames Wright PetscFunctionBeginUser; 120e1233009SJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 121e1233009SJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 122e1233009SJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 123990f1db0SJed Brown if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32; 124990f1db0SJed Brown else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64; 125990f1db0SJed Brown PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type)); 126e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 127e7754af5SKenneth E. Jansen if (time) *time = file_time; 128e7754af5SKenneth E. Jansen if (step_number) *step_number = file_step_number; 129e7754af5SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 130e7754af5SKenneth E. Jansen PetscInt length, N; 131990f1db0SJed Brown PetscCall(BinaryReadIntoInt(viewer, &length, file_type)); 132e7754af5SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 133e7754af5SKenneth 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); 134e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 135e7754af5SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 136e7754af5SKenneth E. Jansen 137e7754af5SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 138d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 139e7754af5SKenneth E. Jansen } 140e7754af5SKenneth E. Jansen 141a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 142c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 143a515125bSLeila Ghaffari Vec Qref; 144a515125bSLeila Ghaffari PetscViewer viewer; 145a515125bSLeila Ghaffari PetscReal error, Qrefnorm; 146e7754af5SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 147a515125bSLeila Ghaffari 14806f41313SJames Wright PetscFunctionBeginUser; 149a515125bSLeila Ghaffari // Read reference file 1502b916ea7SJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 151*4c07ec22SJames Wright PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given"); 152e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 153e7754af5SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 154a515125bSLeila Ghaffari 155a515125bSLeila Ghaffari // Compute error with respect to reference solution 1562b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1572b916ea7SJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1582b916ea7SJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1592b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 160a515125bSLeila Ghaffari 161a515125bSLeila Ghaffari // Check error 162a515125bSLeila Ghaffari if (error > app_ctx->test_tol) { 1632b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 164a515125bSLeila Ghaffari } 165a515125bSLeila Ghaffari 166a515125bSLeila Ghaffari // Cleanup 1672b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1682b916ea7SJeremy L Thompson PetscCall(VecDestroy(&Qref)); 169d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 170a515125bSLeila Ghaffari } 171a515125bSLeila Ghaffari 172a515125bSLeila Ghaffari // Get error for problems with exact solutions 173c56e8d5bSJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 174a515125bSLeila Ghaffari PetscInt loc_nodes; 175a515125bSLeila Ghaffari Vec Q_exact, Q_exact_loc; 176a515125bSLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 177a515125bSLeila Ghaffari 17806f41313SJames Wright PetscFunctionBeginUser; 179a515125bSLeila Ghaffari // Get exact solution at final time 180b2948607SJames Wright PetscCall(DMGetGlobalVector(dm, &Q_exact)); 1812b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1822b916ea7SJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1832b916ea7SJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 184a515125bSLeila Ghaffari 185a515125bSLeila Ghaffari // Get |exact solution - obtained solution| 1862b916ea7SJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1872b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1882b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 189a515125bSLeila Ghaffari 190a515125bSLeila Ghaffari rel_error = norm_error / norm_exact; 1912b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 1922b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 193b2948607SJames Wright PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 194d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 195a515125bSLeila Ghaffari } 196a515125bSLeila Ghaffari 197a515125bSLeila Ghaffari // Post-processing 198991aef52SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) { 199a515125bSLeila Ghaffari PetscInt steps; 200f0784ed3SJed Brown TSConvergedReason reason; 201a515125bSLeila Ghaffari 20206f41313SJames Wright PetscFunctionBeginUser; 203a515125bSLeila Ghaffari // Print relative error 2040e1e9333SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 205c56e8d5bSJames Wright PetscCall(PrintError(ceed_data, dm, user, Q, final_time)); 206a515125bSLeila Ghaffari } 207a515125bSLeila Ghaffari 208a515125bSLeila Ghaffari // Print final time and number of steps 2092b916ea7SJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 210f0784ed3SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 2110e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 212f0784ed3SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 213f0784ed3SJed Brown steps, (double)final_time)); 214a515125bSLeila Ghaffari } 215a515125bSLeila Ghaffari 216a515125bSLeila Ghaffari // Output numerical values from command line 2172b916ea7SJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 218a515125bSLeila Ghaffari 219a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 2200e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 221c56e8d5bSJames Wright PetscCall(RegressionTest(user->app_ctx, Q)); 222a515125bSLeila Ghaffari } 223d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 224a515125bSLeila Ghaffari } 225a515125bSLeila Ghaffari 226e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 227e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 228e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2299293eaa1SJed Brown 230a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation 231a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 232a515125bSLeila Ghaffari PetscViewer viewer; 2332b916ea7SJeremy L Thompson 23406f41313SJames Wright PetscFunctionBeginUser; 2352b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 236e7754af5SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2372b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 238d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 239a515125bSLeila Ghaffari } 240a515125bSLeila Ghaffari 24115a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 24215a3537eSJed Brown int FreeContextPetsc(void *data) { 2432b916ea7SJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 24415a3537eSJed Brown return CEED_ERROR_SUCCESS; 24515a3537eSJed Brown } 2469f59f36eSJames Wright 2479f59f36eSJames Wright // Return mass qfunction specification for number of components N 2489f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 2499f59f36eSJames Wright PetscFunctionBeginUser; 2509f59f36eSJames Wright switch (N) { 2519f59f36eSJames Wright case 1: 252b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 2539f59f36eSJames Wright break; 2549f59f36eSJames Wright case 5: 255b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 2569f59f36eSJames Wright break; 257c38c977aSJames Wright case 7: 258b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 259c38c977aSJames Wright break; 2609f59f36eSJames Wright case 9: 261b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 2629f59f36eSJames Wright break; 2639f59f36eSJames Wright case 22: 264b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 2659f59f36eSJames Wright break; 2669f59f36eSJames Wright default: 2676f539f71SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 2689f59f36eSJames Wright } 2699f59f36eSJames Wright 270b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 271b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 272b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 2733170c09fSJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N)); 274d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2759f59f36eSJames Wright } 276e5e81594SJames Wright 277457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 278457a5831SJames Wright PetscFunctionBeginUser; 279d949ddfcSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 280457a5831SJames Wright 281457a5831SJames Wright PetscCall(DMDestroy(&context->dm)); 282457a5831SJames Wright PetscCall(KSPDestroy(&context->ksp)); 283457a5831SJames Wright 284457a5831SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 285457a5831SJames Wright 286457a5831SJames Wright PetscCall(PetscFree(context)); 287d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 288457a5831SJames Wright } 289457a5831SJames Wright 290676080b4SJames Wright /* 291676080b4SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 292676080b4SJames Wright * 293676080b4SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 294676080b4SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 295676080b4SJames Wright * 296676080b4SJames 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. 297676080b4SJames Wright * 298676080b4SJames Wright * @param[in] comm MPI_Comm for the program 299676080b4SJames Wright * @param[in] path Path to the file 300676080b4SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 301676080b4SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 302676080b4SJames Wright * @param[out] fp File pointer to the opened file 303676080b4SJames Wright */ 30442454adaSJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 305676080b4SJames Wright FILE **fp) { 306defe8520SJames Wright int ndims; 307676080b4SJames Wright char line[char_array_len]; 308676080b4SJames Wright char **array; 309676080b4SJames Wright 310676080b4SJames Wright PetscFunctionBeginUser; 311676080b4SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 312676080b4SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 313676080b4SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 314defe8520SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 315676080b4SJames Wright 316676080b4SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 317676080b4SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 318d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 319676080b4SJames Wright } 320676080b4SJames Wright 321676080b4SJames Wright /* 322676080b4SJames Wright * @brief Get the number of rows for the PHASTA file at path. 323676080b4SJames Wright * 324676080b4SJames 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. 325676080b4SJames Wright * 326676080b4SJames Wright * @param[in] comm MPI_Comm for the program 327676080b4SJames Wright * @param[in] path Path to the file 328676080b4SJames Wright * @param[out] nrows Number of rows 329676080b4SJames Wright */ 33042454adaSJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 331676080b4SJames Wright const PetscInt char_array_len = 512; 332676080b4SJames Wright PetscInt dims[2]; 333676080b4SJames Wright FILE *fp; 334676080b4SJames Wright 335676080b4SJames Wright PetscFunctionBeginUser; 33642454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 337676080b4SJames Wright *nrows = dims[0]; 338676080b4SJames Wright PetscCall(PetscFClose(comm, fp)); 339d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 340676080b4SJames Wright } 34162b7942eSJames Wright 34242454adaSJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 343defe8520SJames Wright PetscInt dims[2]; 344defe8520SJames Wright int ndims; 34562b7942eSJames Wright FILE *fp; 34662b7942eSJames Wright const PetscInt char_array_len = 512; 34762b7942eSJames Wright char line[char_array_len]; 34862b7942eSJames Wright char **row_array; 34962b7942eSJames Wright 35006f41313SJames Wright PetscFunctionBeginUser; 35142454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 35262b7942eSJames Wright 35362b7942eSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 35462b7942eSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 35562b7942eSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3565d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 357defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 35862b7942eSJames Wright 35962b7942eSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 36062b7942eSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 36162b7942eSJames Wright } 36262b7942eSJames Wright } 36362b7942eSJames Wright 36462b7942eSJames Wright PetscCall(PetscFClose(comm, fp)); 365d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 36662b7942eSJames Wright } 3677eedc94cSJames Wright 3687eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 369ad2e713eSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Init; 370ad2e713eSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Meta; 371ad2e713eSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Train; 372ad2e713eSRiccardo Balin PetscLogEvent FLUIDS_TrainDataCompute; 373855536edSJames Wright PetscLogEvent FLUIDS_DifferentialFilter; 374855536edSJames Wright PetscLogEvent FLUIDS_VelocityGradientProjection; 375855536edSJames Wright static PetscClassId libCEED_classid, onlineTrain_classid, misc_classid; 3767eedc94cSJames Wright 3777eedc94cSJames Wright PetscErrorCode RegisterLogEvents() { 3787eedc94cSJames Wright PetscFunctionBeginUser; 3797eedc94cSJames Wright PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid)); 3807eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply)); 381ad2e713eSRiccardo Balin 382ad2e713eSRiccardo Balin PetscCall(PetscClassIdRegister("onlineTrain", &onlineTrain_classid)); 383ad2e713eSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Init", onlineTrain_classid, &FLUIDS_SmartRedis_Init)); 384ad2e713eSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Meta", onlineTrain_classid, &FLUIDS_SmartRedis_Meta)); 385ad2e713eSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Train", onlineTrain_classid, &FLUIDS_SmartRedis_Train)); 386ad2e713eSRiccardo Balin PetscCall(PetscLogEventRegister("TrainDataCompute", onlineTrain_classid, &FLUIDS_TrainDataCompute)); 387855536edSJames Wright 388855536edSJames Wright PetscCall(PetscClassIdRegister("Miscellaneous", &misc_classid)); 389855536edSJames Wright PetscCall(PetscLogEventRegister("DiffFilter", misc_classid, &FLUIDS_DifferentialFilter)); 390855536edSJames Wright PetscCall(PetscLogEventRegister("VeloGradProj", misc_classid, &FLUIDS_VelocityGradientProjection)); 391d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3927eedc94cSJames Wright } 393defe8520SJames Wright 394926a6279SJames Wright // Print information about the given simulation run 3953678fae3SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) { 396b4c37c5cSJames Wright Ceed ceed = user->ceed; 3973678fae3SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 3983678fae3SJames Wright 399926a6279SJames Wright PetscFunctionBeginUser; 400926a6279SJames Wright // Header and rank 401926a6279SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 402926a6279SJames Wright PetscMPIInt rank, comm_size; 403926a6279SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 404926a6279SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 405926a6279SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 406926a6279SJames Wright PetscCall(PetscPrintf(comm, 407926a6279SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 408926a6279SJames Wright " MPI:\n" 409926a6279SJames Wright " Host Name : %s\n" 410926a6279SJames Wright " Total ranks : %d\n", 411926a6279SJames Wright host_name, comm_size)); 412926a6279SJames Wright 413926a6279SJames Wright // Problem specific info 4142d49c0afSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx)); 415926a6279SJames Wright 416926a6279SJames Wright // libCEED 417926a6279SJames Wright const char *used_resource; 418926a6279SJames Wright CeedMemType mem_type_backend; 419b4c37c5cSJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 420b4c37c5cSJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 421926a6279SJames Wright PetscCall(PetscPrintf(comm, 422926a6279SJames Wright " libCEED:\n" 423926a6279SJames Wright " libCEED Backend : %s\n" 424926a6279SJames Wright " libCEED Backend MemType : %s\n", 425926a6279SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 426926a6279SJames Wright // PETSc 4273678fae3SJames Wright VecType vec_type; 428926a6279SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 429926a6279SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0'; 430926a6279SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 431926a6279SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type)); 432926a6279SJames Wright PetscCall(PetscPrintf(comm, 433926a6279SJames Wright " PETSc:\n" 434926a6279SJames Wright " Box Faces : %s\n" 435926a6279SJames Wright " DM VecType : %s\n" 436926a6279SJames Wright " Time Stepping Scheme : %s\n", 4373678fae3SJames Wright box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 4383678fae3SJames Wright { 4393678fae3SJames Wright char pmat_type_str[PETSC_MAX_PATH_LEN]; 4403678fae3SJames Wright MatType amat_type, pmat_type; 4413678fae3SJames Wright Mat Amat, Pmat; 4423678fae3SJames Wright TSIJacobianFn *ijacob_function; 4433678fae3SJames Wright 4443678fae3SJames Wright PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 4453678fae3SJames Wright PetscCall(MatGetType(Amat, &amat_type)); 4463678fae3SJames Wright PetscCall(MatGetType(Pmat, &pmat_type)); 4473678fae3SJames Wright 4483678fae3SJames Wright PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 4493678fae3SJames Wright if (!strcmp(pmat_type, MATCEED)) { 4503678fae3SJames Wright MatType pmat_coo_type; 4513678fae3SJames Wright char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 4523678fae3SJames Wright 4533678fae3SJames Wright PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 4543678fae3SJames Wright PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 4553678fae3SJames Wright PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 4563678fae3SJames Wright } 4573678fae3SJames Wright if (ijacob_function) { 4583678fae3SJames Wright PetscCall(PetscPrintf(comm, 4593678fae3SJames Wright " IJacobian A MatType : %s\n" 4603678fae3SJames Wright " IJacobian P MatType : %s\n", 4613678fae3SJames Wright amat_type, pmat_type_str)); 4623678fae3SJames Wright } 4633678fae3SJames Wright } 464926a6279SJames Wright if (user->app_ctx->cont_steps) { 465926a6279SJames Wright PetscCall(PetscPrintf(comm, 466926a6279SJames Wright " Continue:\n" 467926a6279SJames Wright " Filename: : %s\n" 468926a6279SJames Wright " Step: : %" PetscInt_FMT "\n" 469926a6279SJames Wright " Time: : %g\n", 470926a6279SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 471926a6279SJames Wright } 472926a6279SJames Wright // Mesh 473926a6279SJames Wright const PetscInt num_comp_q = 5; 474926a6279SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 475926a6279SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 476926a6279SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 477926a6279SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 478926a6279SJames Wright PetscCall(PetscPrintf(comm, 479926a6279SJames Wright " Mesh:\n" 480926a6279SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 481926a6279SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 482926a6279SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 483926a6279SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 484dfeb939dSJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 485dfeb939dSJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 486926a6279SJames Wright // -- Get Partition Statistics 487926a6279SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 488926a6279SJames Wright { 489926a6279SJames Wright PetscInt *gather_buffer = NULL; 490dfeb939dSJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 491926a6279SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 492926a6279SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 493926a6279SJames Wright 494dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 495926a6279SJames Wright if (!rank) { 496926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 497926a6279SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 498926a6279SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 499926a6279SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 500926a6279SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 501dfeb939dSJames Wright PetscCall(PetscPrintf( 502dfeb939dSJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 503926a6279SJames 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)); 504926a6279SJames Wright } 505926a6279SJames Wright 506dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 507926a6279SJames Wright if (!rank) { 508926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 509926a6279SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 510926a6279SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 511926a6279SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 512926a6279SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 513dfeb939dSJames Wright PetscCall(PetscPrintf( 514dfeb939dSJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 515926a6279SJames 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)); 516926a6279SJames Wright } 517926a6279SJames Wright 51845abf96eSJames Wright if (comm_size != 1) { 519dfeb939dSJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 520926a6279SJames Wright { 521926a6279SJames Wright PetscSF sf; 522dfeb939dSJames Wright PetscInt nrranks, niranks; 523dfeb939dSJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 524dfeb939dSJames Wright const PetscMPIInt *rranks, *iranks; 525926a6279SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf)); 526926a6279SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 527dfeb939dSJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 528926a6279SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 529926a6279SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 530926a6279SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 531dfeb939dSJames Wright num_ghost_interface_ranks++; 532dfeb939dSJames Wright } 533dfeb939dSJames Wright for (PetscInt i = 0; i < niranks; i++) { 534dfeb939dSJames Wright if (iranks[i] == rank) continue; 535dfeb939dSJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 536dfeb939dSJames Wright num_owned_interface_ranks++; 537926a6279SJames Wright } 538926a6279SJames Wright } 539dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 540926a6279SJames Wright if (!rank) { 541926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 542dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 543dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 544dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 545dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 546dfeb939dSJames Wright PetscCall(PetscPrintf( 54745abf96eSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 54845abf96eSJames 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, 54945abf96eSJames Wright part_shared_dof_ratio)); 550dfeb939dSJames Wright } 551dfeb939dSJames Wright 552dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 553dfeb939dSJames Wright if (!rank) { 554dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 555dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 556dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 557dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 558dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 559dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 560dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 561dfeb939dSJames Wright } 562dfeb939dSJames Wright 563dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 564dfeb939dSJames Wright if (!rank) { 565dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 566dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 567dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 568dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 569dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 570dfeb939dSJames Wright PetscCall(PetscPrintf( 57145abf96eSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 57245abf96eSJames 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, 57345abf96eSJames Wright part_shared_dof_ratio)); 574dfeb939dSJames Wright } 575dfeb939dSJames Wright 576dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 577dfeb939dSJames Wright if (!rank) { 578dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 579dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 580dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 581dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 582dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 583dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 584dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 585926a6279SJames Wright } 58645abf96eSJames Wright } 587926a6279SJames Wright 588926a6279SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 589926a6279SJames Wright } 590926a6279SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 591926a6279SJames Wright } 592