1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3a515125bSLeila Ghaffari // 4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5a515125bSLeila Ghaffari // 6727da7e7SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 7a515125bSLeila Ghaffari 8a515125bSLeila Ghaffari /// @file 9a515125bSLeila Ghaffari /// Miscellaneous utility functions 10a515125bSLeila Ghaffari 11e419654dSJeremy L Thompson #include <ceed.h> 12e419654dSJeremy L Thompson #include <petscdm.h> 13926a6279SJames Wright #include <petscsf.h> 14e419654dSJeremy L Thompson #include <petscts.h> 15e419654dSJeremy L Thompson 16a515125bSLeila Ghaffari #include "../navierstokes.h" 179f59f36eSJames Wright #include "../qfunctions/mass.h" 18a515125bSLeila Ghaffari 192b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 20b4c37c5cSJames Wright Ceed ceed = user->ceed; 21a515125bSLeila Ghaffari PetscFunctionBeginUser; 22a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 23b7f03f12SJed Brown // Update time for evaluation 24a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 25b4c37c5cSJames Wright if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 26a515125bSLeila Ghaffari 27a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 28a515125bSLeila Ghaffari // ICs 29a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 30a515125bSLeila Ghaffari // -- CEED Restriction 31a515125bSLeila Ghaffari CeedVector q0_ceed; 32b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL)); 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari // -- Place PETSc vector in CEED vector 358f18bb8bSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 36a515125bSLeila Ghaffari 37a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 38a515125bSLeila Ghaffari // Fix multiplicity for output of ICs 39a515125bSLeila Ghaffari // --------------------------------------------------------------------------- 40a515125bSLeila Ghaffari // -- CEED Restriction 41a515125bSLeila Ghaffari CeedVector mult_vec; 42b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 43a515125bSLeila Ghaffari 44a515125bSLeila Ghaffari // -- Place PETSc vector in CEED vector 45a515125bSLeila Ghaffari PetscMemType m_mem_type; 46a515125bSLeila Ghaffari Vec multiplicity_loc; 472b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &multiplicity_loc)); 48fd969b44SJames Wright PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec)); 49a515125bSLeila Ghaffari 50a515125bSLeila Ghaffari // -- Get multiplicity 51b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 52a515125bSLeila Ghaffari 53a515125bSLeila Ghaffari // -- Restore vectors 54fd969b44SJames Wright PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc)); 55a515125bSLeila Ghaffari 56a515125bSLeila Ghaffari // -- Local-to-Global 57a515125bSLeila Ghaffari Vec multiplicity; 582b916ea7SJeremy L Thompson PetscCall(DMGetGlobalVector(dm, &multiplicity)); 592b916ea7SJeremy L Thompson PetscCall(VecZeroEntries(multiplicity)); 602b916ea7SJeremy L Thompson PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity)); 61a515125bSLeila Ghaffari 62a515125bSLeila Ghaffari // -- Fix multiplicity 632b916ea7SJeremy L Thompson PetscCall(VecPointwiseDivide(Q, Q, multiplicity)); 642b916ea7SJeremy L Thompson PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc)); 65a515125bSLeila Ghaffari 66a515125bSLeila Ghaffari // -- Restore vectors 672b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc)); 682b916ea7SJeremy L Thompson PetscCall(DMRestoreGlobalVector(dm, &multiplicity)); 69a515125bSLeila Ghaffari 70a515125bSLeila Ghaffari // Cleanup 71b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 72b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q0_ceed)); 73d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 74a515125bSLeila Ghaffari } 75a515125bSLeila Ghaffari 76*c56e8d5bSJames Wright // Record boundary values from initial condition 77*c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 78*c56e8d5bSJames Wright Vec Qbc, boundary_mask; 79*c56e8d5bSJames Wright 80*c56e8d5bSJames Wright PetscFunctionBeginUser; 81*c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 82*c56e8d5bSJames Wright PetscCall(VecCopy(Q_loc, Qbc)); 83*c56e8d5bSJames Wright PetscCall(VecZeroEntries(Q_loc)); 84*c56e8d5bSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 85*c56e8d5bSJames Wright PetscCall(VecAXPY(Qbc, -1., Q_loc)); 86*c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 87*c56e8d5bSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 88*c56e8d5bSJames Wright 89*c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 90*c56e8d5bSJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 91*c56e8d5bSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 92*c56e8d5bSJames Wright PetscCall(VecSet(Q, 1.0)); 93*c56e8d5bSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 94*c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 95*c56e8d5bSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 96*c56e8d5bSJames Wright } 97*c56e8d5bSJames Wright 98*c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 992b916ea7SJeremy L Thompson Vec grad_FVM) { 1009d437337SJames Wright Vec Qbc, boundary_mask; 101a515125bSLeila Ghaffari 10206f41313SJames Wright PetscFunctionBeginUser; 1032eb7bf1fSJames Wright // Mask (zero) Strong BC entries 1049d437337SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 1059d437337SJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 1069d437337SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 1079d437337SJames Wright 1082b916ea7SJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 1092b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 1102b916ea7SJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 111d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 112a515125bSLeila Ghaffari } 113a515125bSLeila Ghaffari 114e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 115e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 116e1233009SJames Wright PetscInt file_step_number; 117e1233009SJames Wright PetscInt32 token; 118e7754af5SKenneth E. Jansen PetscReal file_time; 119e7754af5SKenneth E. Jansen 12006f41313SJames Wright PetscFunctionBeginUser; 121e7754af5SKenneth E. Jansen // Attempt 122e1233009SJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 123e1233009SJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 124e1233009SJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 125e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT)); 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; 131e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT)); 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 // Load Q from existent solution 138e7754af5SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 139d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 140e7754af5SKenneth E. Jansen } 141e7754af5SKenneth E. Jansen 142a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 143*c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 144a515125bSLeila Ghaffari Vec Qref; 145a515125bSLeila Ghaffari PetscViewer viewer; 146a515125bSLeila Ghaffari PetscReal error, Qrefnorm; 147e7754af5SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 148a515125bSLeila Ghaffari 14906f41313SJames Wright PetscFunctionBeginUser; 150a515125bSLeila Ghaffari // Read reference file 1512b916ea7SJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 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 173*c56e8d5bSJames 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 1802b916ea7SJeremy L Thompson PetscCall(DMCreateGlobalVector(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 // Compute relative error 191a515125bSLeila Ghaffari rel_error = norm_error / norm_exact; 192a515125bSLeila Ghaffari 193a515125bSLeila Ghaffari // Output relative error 1942b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 195a515125bSLeila Ghaffari // Cleanup 1962b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 1972b916ea7SJeremy L Thompson PetscCall(VecDestroy(&Q_exact)); 198d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 199a515125bSLeila Ghaffari } 200a515125bSLeila Ghaffari 201a515125bSLeila Ghaffari // Post-processing 202*c56e8d5bSJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 203a515125bSLeila Ghaffari PetscInt steps; 204f0784ed3SJed Brown TSConvergedReason reason; 205a515125bSLeila Ghaffari 20606f41313SJames Wright PetscFunctionBeginUser; 207a515125bSLeila Ghaffari // Print relative error 2080e1e9333SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 209*c56e8d5bSJames Wright PetscCall(PrintError(ceed_data, dm, user, Q, final_time)); 210a515125bSLeila Ghaffari } 211a515125bSLeila Ghaffari 212a515125bSLeila Ghaffari // Print final time and number of steps 2132b916ea7SJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 214f0784ed3SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 2150e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 216f0784ed3SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 217f0784ed3SJed Brown steps, (double)final_time)); 218a515125bSLeila Ghaffari } 219a515125bSLeila Ghaffari 220a515125bSLeila Ghaffari // Output numerical values from command line 2212b916ea7SJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 222a515125bSLeila Ghaffari 223a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 2240e1e9333SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 225*c56e8d5bSJames Wright PetscCall(RegressionTest(user->app_ctx, Q)); 226a515125bSLeila Ghaffari } 227d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 228a515125bSLeila Ghaffari } 229a515125bSLeila Ghaffari 230e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 231e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 232e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2339293eaa1SJed Brown 234a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation 235a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 236a515125bSLeila Ghaffari PetscViewer viewer; 2372b916ea7SJeremy L Thompson 23806f41313SJames Wright PetscFunctionBeginUser; 2392b916ea7SJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 240e7754af5SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2412b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 242d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 243a515125bSLeila Ghaffari } 244a515125bSLeila Ghaffari 24515a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 24615a3537eSJed Brown int FreeContextPetsc(void *data) { 2472b916ea7SJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 24815a3537eSJed Brown return CEED_ERROR_SUCCESS; 24915a3537eSJed Brown } 2509f59f36eSJames Wright 2519f59f36eSJames Wright // Return mass qfunction specification for number of components N 2529f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 2539f59f36eSJames Wright PetscFunctionBeginUser; 2549f59f36eSJames Wright switch (N) { 2559f59f36eSJames Wright case 1: 256b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 2579f59f36eSJames Wright break; 2589f59f36eSJames Wright case 5: 259b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 2609f59f36eSJames Wright break; 261c38c977aSJames Wright case 7: 262b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 263c38c977aSJames Wright break; 2649f59f36eSJames Wright case 9: 265b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 2669f59f36eSJames Wright break; 2679f59f36eSJames Wright case 22: 268b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 2699f59f36eSJames Wright break; 2709f59f36eSJames Wright default: 2716f539f71SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 2729f59f36eSJames Wright } 2739f59f36eSJames Wright 274b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 275b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 276b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 277d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2789f59f36eSJames Wright } 279e5e81594SJames Wright 280457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 281457a5831SJames Wright PetscFunctionBeginUser; 282d949ddfcSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 283457a5831SJames Wright 284457a5831SJames Wright PetscCall(DMDestroy(&context->dm)); 285457a5831SJames Wright PetscCall(KSPDestroy(&context->ksp)); 286457a5831SJames Wright 287457a5831SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 288457a5831SJames Wright 289457a5831SJames Wright PetscCall(PetscFree(context)); 290d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 291457a5831SJames Wright } 292457a5831SJames Wright 293676080b4SJames Wright /* 294676080b4SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 295676080b4SJames Wright * 296676080b4SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 297676080b4SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 298676080b4SJames Wright * 299676080b4SJames 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. 300676080b4SJames Wright * 301676080b4SJames Wright * @param[in] comm MPI_Comm for the program 302676080b4SJames Wright * @param[in] path Path to the file 303676080b4SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 304676080b4SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 305676080b4SJames Wright * @param[out] fp File pointer to the opened file 306676080b4SJames Wright */ 30742454adaSJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 308676080b4SJames Wright FILE **fp) { 309defe8520SJames Wright int ndims; 310676080b4SJames Wright char line[char_array_len]; 311676080b4SJames Wright char **array; 312676080b4SJames Wright 313676080b4SJames Wright PetscFunctionBeginUser; 314676080b4SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 315676080b4SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 316676080b4SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 317defe8520SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 318676080b4SJames Wright 319676080b4SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 320676080b4SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 321d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 322676080b4SJames Wright } 323676080b4SJames Wright 324676080b4SJames Wright /* 325676080b4SJames Wright * @brief Get the number of rows for the PHASTA file at path. 326676080b4SJames Wright * 327676080b4SJames 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. 328676080b4SJames Wright * 329676080b4SJames Wright * @param[in] comm MPI_Comm for the program 330676080b4SJames Wright * @param[in] path Path to the file 331676080b4SJames Wright * @param[out] nrows Number of rows 332676080b4SJames Wright */ 33342454adaSJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 334676080b4SJames Wright const PetscInt char_array_len = 512; 335676080b4SJames Wright PetscInt dims[2]; 336676080b4SJames Wright FILE *fp; 337676080b4SJames Wright 338676080b4SJames Wright PetscFunctionBeginUser; 33942454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 340676080b4SJames Wright *nrows = dims[0]; 341676080b4SJames Wright PetscCall(PetscFClose(comm, fp)); 342d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 343676080b4SJames Wright } 34462b7942eSJames Wright 34542454adaSJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 346defe8520SJames Wright PetscInt dims[2]; 347defe8520SJames Wright int ndims; 34862b7942eSJames Wright FILE *fp; 34962b7942eSJames Wright const PetscInt char_array_len = 512; 35062b7942eSJames Wright char line[char_array_len]; 35162b7942eSJames Wright char **row_array; 35262b7942eSJames Wright 35306f41313SJames Wright PetscFunctionBeginUser; 35442454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 35562b7942eSJames Wright 35662b7942eSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 35762b7942eSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 35862b7942eSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3595d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 360defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 36162b7942eSJames Wright 36262b7942eSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 36362b7942eSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 36462b7942eSJames Wright } 36562b7942eSJames Wright } 36662b7942eSJames Wright 36762b7942eSJames Wright PetscCall(PetscFClose(comm, fp)); 368d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 36962b7942eSJames Wright } 3707eedc94cSJames Wright 3717eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 3727eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemble; 3737eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 3747eedc94cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 3757eedc94cSJames Wright static PetscClassId libCEED_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)); 3817eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble)); 3827eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal)); 3837eedc94cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal)); 384d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3857eedc94cSJames Wright } 386defe8520SJames Wright 387defe8520SJames Wright /** 388defe8520SJames Wright @brief Translate array of CeedInt to PetscInt. 389defe8520SJames Wright If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`. 390defe8520SJames Wright Caller is responsible for freeing `array_petsc` with `free()`. 391defe8520SJames Wright 392defe8520SJames Wright @param[in] num_entries Number of array entries 393defe8520SJames Wright @param[in,out] array_ceed Array of CeedInts 394defe8520SJames Wright @param[out] array_petsc Array of PetscInts 395defe8520SJames Wright **/ 396defe8520SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) { 397defe8520SJames Wright CeedInt int_c = 0; 398defe8520SJames Wright PetscInt int_p = 0; 399defe8520SJames Wright 400defe8520SJames Wright PetscFunctionBeginUser; 401defe8520SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 402defe8520SJames Wright *array_petsc = (PetscInt *)*array_ceed; 403defe8520SJames Wright } else { 404defe8520SJames Wright *array_petsc = malloc(num_entries * sizeof(PetscInt)); 405defe8520SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i]; 406defe8520SJames Wright free(*array_ceed); 407defe8520SJames Wright } 408defe8520SJames Wright *array_ceed = NULL; 409defe8520SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 410defe8520SJames Wright } 411defe8520SJames Wright 412defe8520SJames Wright /** 413defe8520SJames Wright @brief Translate array of PetscInt to CeedInt. 414defe8520SJames Wright If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`. 415defe8520SJames Wright Caller is responsible for freeing `array_ceed` with `PetscFree()`. 416defe8520SJames Wright 417defe8520SJames Wright @param[in] num_entries Number of array entries 418defe8520SJames Wright @param[in,out] array_petsc Array of PetscInts 419defe8520SJames Wright @param[out] array_ceed Array of CeedInts 420defe8520SJames Wright **/ 421defe8520SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) { 422defe8520SJames Wright CeedInt int_c = 0; 423defe8520SJames Wright PetscInt int_p = 0; 424defe8520SJames Wright 425defe8520SJames Wright PetscFunctionBeginUser; 426defe8520SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 427defe8520SJames Wright *array_ceed = (CeedInt *)*array_petsc; 428defe8520SJames Wright } else { 429defe8520SJames Wright PetscCall(PetscMalloc1(num_entries, array_ceed)); 430defe8520SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i]; 431defe8520SJames Wright PetscCall(PetscFree(*array_petsc)); 432defe8520SJames Wright } 433defe8520SJames Wright *array_petsc = NULL; 434defe8520SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 435defe8520SJames Wright } 436926a6279SJames Wright 437926a6279SJames Wright // Print information about the given simulation run 438926a6279SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) { 439b4c37c5cSJames Wright Ceed ceed = user->ceed; 440926a6279SJames Wright PetscFunctionBeginUser; 441926a6279SJames Wright // Header and rank 442926a6279SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 443926a6279SJames Wright PetscMPIInt rank, comm_size; 444926a6279SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 445926a6279SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 446926a6279SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 447926a6279SJames Wright PetscCall(PetscPrintf(comm, 448926a6279SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 449926a6279SJames Wright " MPI:\n" 450926a6279SJames Wright " Host Name : %s\n" 451926a6279SJames Wright " Total ranks : %d\n", 452926a6279SJames Wright host_name, comm_size)); 453926a6279SJames Wright 454926a6279SJames Wright // Problem specific info 4552d49c0afSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx)); 456926a6279SJames Wright 457926a6279SJames Wright // libCEED 458926a6279SJames Wright const char *used_resource; 459926a6279SJames Wright CeedMemType mem_type_backend; 460b4c37c5cSJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 461b4c37c5cSJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 462926a6279SJames Wright PetscCall(PetscPrintf(comm, 463926a6279SJames Wright " libCEED:\n" 464926a6279SJames Wright " libCEED Backend : %s\n" 465926a6279SJames Wright " libCEED Backend MemType : %s\n", 466926a6279SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 467926a6279SJames Wright // PETSc 468926a6279SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 469926a6279SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0'; 470926a6279SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 471926a6279SJames Wright MatType mat_type; 472926a6279SJames Wright VecType vec_type; 473926a6279SJames Wright PetscCall(DMGetMatType(user->dm, &mat_type)); 474926a6279SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type)); 475926a6279SJames Wright PetscCall(PetscPrintf(comm, 476926a6279SJames Wright " PETSc:\n" 477926a6279SJames Wright " Box Faces : %s\n" 478926a6279SJames Wright " DM MatType : %s\n" 479926a6279SJames Wright " DM VecType : %s\n" 480926a6279SJames Wright " Time Stepping Scheme : %s\n", 481926a6279SJames Wright box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 482926a6279SJames Wright if (user->app_ctx->cont_steps) { 483926a6279SJames Wright PetscCall(PetscPrintf(comm, 484926a6279SJames Wright " Continue:\n" 485926a6279SJames Wright " Filename: : %s\n" 486926a6279SJames Wright " Step: : %" PetscInt_FMT "\n" 487926a6279SJames Wright " Time: : %g\n", 488926a6279SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 489926a6279SJames Wright } 490926a6279SJames Wright // Mesh 491926a6279SJames Wright const PetscInt num_comp_q = 5; 492926a6279SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 493926a6279SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 494926a6279SJames Wright // -- Get global size 495926a6279SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 496926a6279SJames Wright // -- Get local size 497926a6279SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 498926a6279SJames Wright PetscCall(PetscPrintf(comm, 499926a6279SJames Wright " Mesh:\n" 500926a6279SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 501926a6279SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 502926a6279SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 503926a6279SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 504dfeb939dSJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 505dfeb939dSJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 506926a6279SJames Wright // -- Get Partition Statistics 507926a6279SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 508926a6279SJames Wright { 509926a6279SJames Wright PetscInt *gather_buffer = NULL; 510dfeb939dSJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 511926a6279SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 512926a6279SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 513926a6279SJames Wright 514dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 515926a6279SJames Wright if (!rank) { 516926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 517926a6279SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 518926a6279SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 519926a6279SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 520926a6279SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 521dfeb939dSJames Wright PetscCall(PetscPrintf( 522dfeb939dSJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 523926a6279SJames 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)); 524926a6279SJames Wright } 525926a6279SJames Wright 526dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 527926a6279SJames Wright if (!rank) { 528926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 529926a6279SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 530926a6279SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 531926a6279SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 532926a6279SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 533dfeb939dSJames Wright PetscCall(PetscPrintf( 534dfeb939dSJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 535926a6279SJames 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)); 536926a6279SJames Wright } 537926a6279SJames Wright 53845abf96eSJames Wright if (comm_size != 1) { 539dfeb939dSJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 540926a6279SJames Wright { 541926a6279SJames Wright PetscSF sf; 542dfeb939dSJames Wright PetscInt nrranks, niranks; 543dfeb939dSJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 544dfeb939dSJames Wright const PetscMPIInt *rranks, *iranks; 545926a6279SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf)); 546926a6279SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 547dfeb939dSJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 548926a6279SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 549926a6279SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 550926a6279SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 551dfeb939dSJames Wright num_ghost_interface_ranks++; 552dfeb939dSJames Wright } 553dfeb939dSJames Wright for (PetscInt i = 0; i < niranks; i++) { 554dfeb939dSJames Wright if (iranks[i] == rank) continue; 555dfeb939dSJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 556dfeb939dSJames Wright num_owned_interface_ranks++; 557926a6279SJames Wright } 558926a6279SJames Wright } 559dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 560926a6279SJames Wright if (!rank) { 561926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 562dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 563dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 564dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 565dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 566dfeb939dSJames Wright PetscCall(PetscPrintf( 56745abf96eSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 56845abf96eSJames 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, 56945abf96eSJames Wright part_shared_dof_ratio)); 570dfeb939dSJames Wright } 571dfeb939dSJames Wright 572dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 573dfeb939dSJames Wright if (!rank) { 574dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 575dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 576dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 577dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 578dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 579dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 580dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 581dfeb939dSJames Wright } 582dfeb939dSJames Wright 583dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 584dfeb939dSJames Wright if (!rank) { 585dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 586dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 587dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 588dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 589dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 590dfeb939dSJames Wright PetscCall(PetscPrintf( 59145abf96eSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 59245abf96eSJames 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, 59345abf96eSJames Wright part_shared_dof_ratio)); 594dfeb939dSJames Wright } 595dfeb939dSJames Wright 596dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 597dfeb939dSJames Wright if (!rank) { 598dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 599dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 600dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 601dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 602dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 603dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 604dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 605926a6279SJames Wright } 60645abf96eSJames Wright } 607926a6279SJames Wright 608926a6279SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 609926a6279SJames Wright } 610926a6279SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 611926a6279SJames Wright } 612