1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Miscellaneous utility functions 6 7 #include <ceed.h> 8 #include <petscdm.h> 9 #include <petscsf.h> 10 #include <petscts.h> 11 12 #include <navierstokes.h> 13 14 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time) { 15 Ceed ceed = honee->ceed; 16 CeedVector mult_vec; 17 PetscMemType m_mem_type; 18 Vec Multiplicity, Multiplicity_loc; 19 20 PetscFunctionBeginUser; 21 if (honee->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ics_ctx->op, honee->phys->ics_time_label, &time)); 22 PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, honee->op_ics_ctx)); 23 24 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &mult_vec, NULL)); 25 26 // -- Get multiplicity 27 PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 28 PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 29 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(honee->elem_restr_q, mult_vec)); 30 PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 31 32 PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 33 PetscCall(VecZeroEntries(Multiplicity)); 34 PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 35 36 // -- Fix multiplicity 37 PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 38 PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 39 40 PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 41 PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 42 PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 43 PetscFunctionReturn(PETSC_SUCCESS); 44 } 45 46 // Record boundary values from initial condition 47 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 48 PetscFunctionBeginUser; 49 { // Capture initial condition values in Qbc 50 Vec Qbc; 51 52 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 53 PetscCall(VecCopy(Q_loc, Qbc)); 54 PetscCall(VecZeroEntries(Q_loc)); 55 PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 56 PetscCall(VecAXPY(Qbc, -1., Q_loc)); 57 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 58 } 59 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 60 61 { // Set boundary mask to zero out essential BCs 62 Vec boundary_mask, ones; 63 64 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 65 PetscCall(DMGetGlobalVector(dm, &ones)); 66 PetscCall(VecZeroEntries(boundary_mask)); 67 PetscCall(VecSet(ones, 1.0)); 68 PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask)); 69 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 70 PetscCall(DMRestoreGlobalVector(dm, &ones)); 71 } 72 PetscFunctionReturn(PETSC_SUCCESS); 73 } 74 75 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 76 Vec grad_FVM) { 77 Vec Qbc, boundary_mask; 78 79 PetscFunctionBeginUser; 80 // Mask (zero) Strong BC entries 81 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 82 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 83 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 84 85 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 86 PetscCall(VecAXPY(Q_loc, 1., Qbc)); 87 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 88 PetscFunctionReturn(PETSC_SUCCESS); 89 } 90 91 // Compare reference solution values with current test run for CI 92 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 93 Vec Qref; 94 PetscViewer viewer; 95 PetscReal error, Qrefnorm; 96 MPI_Comm comm = PetscObjectComm((PetscObject)Q); 97 98 PetscFunctionBeginUser; 99 // Read reference file 100 PetscCall(VecDuplicate(Q, &Qref)); 101 PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given"); 102 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 103 PetscCall(HoneeLoadBinaryVec(viewer, Qref, NULL, NULL)); 104 105 // Compute error with respect to reference solution 106 PetscCall(VecAXPY(Q, -1.0, Qref)); 107 PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 108 PetscCall(VecScale(Q, 1. / Qrefnorm)); 109 PetscCall(VecNorm(Q, NORM_MAX, &error)); 110 111 // Check error 112 if (error > app_ctx->test_tol) { 113 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 114 } 115 116 // Cleanup 117 PetscCall(PetscViewerDestroy(&viewer)); 118 PetscCall(VecDestroy(&Qref)); 119 PetscFunctionReturn(PETSC_SUCCESS); 120 } 121 122 // Get error for problems with exact solutions 123 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time) { 124 PetscInt loc_nodes; 125 Vec Q_exact, Q_exact_loc; 126 PetscReal rel_error, norm_error, norm_exact; 127 128 PetscFunctionBeginUser; 129 // Get exact solution at final time 130 PetscCall(DMGetGlobalVector(dm, &Q_exact)); 131 PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 132 PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 133 PetscCall(ICs_FixMultiplicity(dm, honee, Q_exact_loc, Q_exact, final_time)); 134 135 // Get |exact solution - obtained solution| 136 PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 137 PetscCall(VecAXPY(Q, -1.0, Q_exact)); 138 PetscCall(VecNorm(Q, NORM_1, &norm_error)); 139 140 rel_error = norm_error / norm_exact; 141 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 142 PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 143 PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 144 PetscFunctionReturn(PETSC_SUCCESS); 145 } 146 147 // Post-processing 148 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time) { 149 PetscInt steps; 150 TSConvergedReason reason; 151 152 PetscFunctionBeginUser; 153 // Print relative error 154 if (problem->compute_exact_solution_error && honee->app_ctx->test_type == TESTTYPE_NONE) { 155 PetscCall(PrintError(dm, honee, Q, final_time)); 156 } 157 158 // Print final time and number of steps 159 PetscCall(TSGetStepNumber(ts, &steps)); 160 PetscCall(TSGetConvergedReason(ts, &reason)); 161 if (honee->app_ctx->test_type == TESTTYPE_NONE) { 162 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 163 steps, (double)final_time)); 164 } 165 166 // Output numerical values from command line 167 PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 168 169 // Compare reference solution values with current test run for CI 170 if (honee->app_ctx->test_type == TESTTYPE_SOLVER) { 171 PetscCall(RegressionTest(honee->app_ctx, Q)); 172 } 173 PetscFunctionReturn(PETSC_SUCCESS); 174 } 175 176 // Free a plain data context that was allocated using PETSc; returning libCEED error codes 177 int FreeContextPetsc(void *data) { 178 if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 179 return CEED_ERROR_SUCCESS; 180 } 181 182 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 183 PetscFunctionBeginUser; 184 if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 185 PetscCall(DMDestroy(&context->dm)); 186 PetscCall(KSPDestroy(&context->ksp)); 187 PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 188 PetscCall(PetscFree(context)); 189 PetscFunctionReturn(PETSC_SUCCESS); 190 } 191 192 /** 193 @brief Sets the value of an option if it is not already set 194 195 @param[in,out] options `PetscOptions` database, or `NULL` for default global database 196 @param[in] name Name of the option (with `-` prepended to it) 197 @param[in] value Value to set the option to 198 **/ 199 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]) { 200 PetscBool has_option; 201 202 PetscFunctionBeginUser; 203 PetscCall(PetscOptionsHasName(options, NULL, name, &has_option)); 204 if (!has_option) PetscCall(PetscOptionsSetValue(options, name, value)); 205 PetscFunctionReturn(PETSC_SUCCESS); 206 } 207 208 PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume) { 209 DM dm = honee->dm, dm_coord; 210 Ceed ceed = honee->ceed; 211 CeedQFunction qf_mass; 212 CeedOperator op_mass; 213 OperatorApplyContext op_mass_ctx; 214 CeedElemRestriction elem_restr_qd; 215 CeedVector qdata; 216 CeedInt q_data_size, num_comps_x; 217 DMLabel domain_label = 0; 218 PetscInt label_value = 0; 219 Vec u, v; 220 221 PetscFunctionBeginUser; 222 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 223 PetscCall(QDataGet(ceed, dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &qdata, &q_data_size)); 224 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comps_x)); 225 226 PetscCall(HoneeMassQFunctionCreate(ceed, num_comps_x, q_data_size, &qf_mass)); 227 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 228 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE)); 229 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, qdata)); 230 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE)); 231 232 PetscCall(OperatorApplyContextCreate(NULL, dm_coord, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx)); 233 PetscCall(CeedOperatorCreateLocalVecs(op_mass, DMReturnVecType(dm), PETSC_COMM_SELF, &u, NULL)); 234 PetscCall(DMCreateGlobalVector(dm_coord, &v)); 235 PetscCall(VecSet(u, 1.)); 236 PetscCall(ApplyCeedOperatorLocalToGlobal(u, v, op_mass_ctx)); 237 PetscCall(VecSum(v, volume)); 238 *volume /= num_comps_x; // Correct for number of components != 1 239 240 PetscCall(VecDestroy(&u)); 241 PetscCall(VecDestroy(&v)); 242 PetscCall(OperatorApplyContextDestroy(op_mass_ctx)); 243 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 244 PetscCallCeed(ceed, CeedVectorDestroy(&qdata)); 245 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 246 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 247 PetscFunctionReturn(PETSC_SUCCESS); 248 } 249 250 // Print information about the given simulation run 251 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts) { 252 Ceed ceed = honee->ceed; 253 MPI_Comm comm = PetscObjectComm((PetscObject)ts); 254 255 PetscFunctionBeginUser; 256 // Header and rank 257 char host_name[PETSC_MAX_PATH_LEN]; 258 PetscMPIInt rank, comm_size; 259 PetscCall(PetscGetHostName(host_name, sizeof host_name)); 260 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 261 PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 262 PetscCall(PetscPrintf(comm, 263 "\n-- HONEE - High-Order Navier-stokes Equation Evaluator --\n" 264 " MPI:\n" 265 " Host Name : %s\n" 266 " Total ranks : %d\n", 267 host_name, comm_size)); 268 269 // Problem specific info 270 PetscCall(problem->print_info(honee, problem, honee->app_ctx)); 271 272 // libCEED 273 const char *used_resource; 274 CeedMemType mem_type_backend; 275 PetscCallCeed(ceed, CeedGetResource(honee->ceed, &used_resource)); 276 PetscCallCeed(ceed, CeedGetPreferredMemType(honee->ceed, &mem_type_backend)); 277 PetscCall(PetscPrintf(comm, 278 " libCEED:\n" 279 " libCEED Backend : %s\n" 280 " libCEED Backend MemType : %s\n", 281 used_resource, CeedMemTypes[mem_type_backend])); 282 // PETSc 283 { 284 VecType vec_type; 285 char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 286 PetscInt dim; 287 288 PetscCall(DMGetDimension(honee->dm, &dim)); 289 if (dim == 2) box_faces_str[3] = '\0'; 290 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 291 PetscCall(DMGetVecType(honee->dm, &vec_type)); 292 PetscCall(PetscPrintf(comm, 293 " PETSc:\n" 294 " Box Faces : %s\n" 295 " DM VecType : %s\n" 296 " Time Stepping Scheme : %s\n", 297 box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 298 } 299 { 300 char pmat_type_str[PETSC_MAX_PATH_LEN]; 301 MatType amat_type, pmat_type; 302 Mat Amat, Pmat; 303 TSIJacobianFn *ijacob_function; 304 305 PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 306 PetscCall(MatGetType(Amat, &amat_type)); 307 PetscCall(MatGetType(Pmat, &pmat_type)); 308 309 PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 310 if (!strcmp(pmat_type, MATCEED)) { 311 MatType pmat_coo_type; 312 char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 313 314 PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 315 PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 316 PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 317 } 318 if (ijacob_function) { 319 PetscCall(PetscPrintf(comm, 320 " IJacobian A MatType : %s\n" 321 " IJacobian P MatType : %s\n", 322 amat_type, pmat_type_str)); 323 } 324 } 325 if (honee->app_ctx->use_continue_file) { 326 PetscCall(PetscPrintf(comm, 327 " Continue:\n" 328 " Filename: : %s\n" 329 " Step: : %" PetscInt_FMT "\n" 330 " Time: : %g\n", 331 honee->app_ctx->cont_file, honee->app_ctx->cont_steps, honee->app_ctx->cont_time)); 332 } 333 // Mesh 334 const PetscInt num_comp_q = 5; 335 PetscInt glob_dofs, owned_dofs, local_dofs; 336 const CeedInt num_P = honee->app_ctx->degree + 1, num_Q = num_P + honee->app_ctx->q_extra; 337 PetscCall(DMGetGlobalVectorInfo(honee->dm, &owned_dofs, &glob_dofs, NULL)); 338 PetscCall(DMGetLocalVectorInfo(honee->dm, &local_dofs, NULL, NULL)); 339 PetscCall(PetscPrintf(comm, 340 " Mesh:\n" 341 " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 342 " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 343 " Global DoFs : %" PetscInt_FMT "\n" 344 " DoFs per node : %" PetscInt_FMT "\n" 345 " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 346 num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 347 // -- Get Partition Statistics 348 PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 349 { 350 PetscInt *gather_buffer = NULL; 351 PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 352 PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 353 if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 354 355 PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 356 if (!rank) { 357 PetscCall(PetscSortInt(comm_size, gather_buffer)); 358 part_owned_dofs[0] = gather_buffer[0]; // min 359 part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 360 part_owned_dofs[2] = gather_buffer[median_index]; // median 361 PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 362 PetscCall(PetscPrintf( 363 comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 364 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)); 365 } 366 367 PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 368 if (!rank) { 369 PetscCall(PetscSortInt(comm_size, gather_buffer)); 370 part_local_dofs[0] = gather_buffer[0]; // min 371 part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 372 part_local_dofs[2] = gather_buffer[median_index]; // median 373 PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 374 PetscCall(PetscPrintf( 375 comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 376 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)); 377 } 378 379 if (comm_size != 1) { 380 PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 381 { 382 PetscSF sf; 383 PetscMPIInt nrranks, niranks; 384 const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 385 const PetscMPIInt *rranks, *iranks; 386 PetscCall(DMGetSectionSF(honee->dm, &sf)); 387 PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 388 PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 389 for (PetscInt i = 0; i < nrranks; i++) { 390 if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 391 num_remote_roots_total += roffset[i + 1] - roffset[i]; 392 num_ghost_interface_ranks++; 393 } 394 for (PetscInt i = 0; i < niranks; i++) { 395 if (iranks[i] == rank) continue; 396 num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 397 num_owned_interface_ranks++; 398 } 399 } 400 PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 401 if (!rank) { 402 PetscCall(PetscSortInt(comm_size, gather_buffer)); 403 part_boundary_dofs[0] = gather_buffer[0]; // min 404 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 405 part_boundary_dofs[2] = gather_buffer[median_index]; // median 406 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 407 PetscCall(PetscPrintf( 408 comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 409 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 410 part_shared_dof_ratio)); 411 } 412 413 PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 414 if (!rank) { 415 PetscCall(PetscSortInt(comm_size, gather_buffer)); 416 part_neighbors[0] = gather_buffer[0]; // min 417 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 418 part_neighbors[2] = gather_buffer[median_index]; // median 419 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 420 PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 421 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 422 } 423 424 PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 425 if (!rank) { 426 PetscCall(PetscSortInt(comm_size, gather_buffer)); 427 part_boundary_dofs[0] = gather_buffer[0]; // min 428 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 429 part_boundary_dofs[2] = gather_buffer[median_index]; // median 430 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 431 PetscCall(PetscPrintf( 432 comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 433 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 434 part_shared_dof_ratio)); 435 } 436 437 PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 438 if (!rank) { 439 PetscCall(PetscSortInt(comm_size, gather_buffer)); 440 part_neighbors[0] = gather_buffer[0]; // min 441 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 442 part_neighbors[2] = gather_buffer[median_index]; // median 443 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 444 PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 445 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 446 } 447 } 448 449 if (!rank) PetscCall(PetscFree(gather_buffer)); 450 } 451 PetscFunctionReturn(PETSC_SUCCESS); 452 } 453