1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Miscellaneous utility functions 10 11 #include <ceed.h> 12 #include <petscdm.h> 13 #include <petscsf.h> 14 #include <petscts.h> 15 16 #include "../navierstokes.h" 17 #include "../qfunctions/mass.h" 18 19 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 20 Ceed ceed = user->ceed; 21 CeedVector mult_vec; 22 PetscMemType m_mem_type; 23 Vec Multiplicity, Multiplicity_loc; 24 25 PetscFunctionBeginUser; 26 if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 27 PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 28 29 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 30 31 // -- Get multiplicity 32 PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 33 PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 34 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 35 PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 36 37 PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 38 PetscCall(VecZeroEntries(Multiplicity)); 39 PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 40 41 // -- Fix multiplicity 42 PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 43 PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 44 45 PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 46 PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 47 PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 48 PetscFunctionReturn(PETSC_SUCCESS); 49 } 50 51 // Record boundary values from initial condition 52 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 53 PetscFunctionBeginUser; 54 { // Capture initial condition values in Qbc 55 Vec Qbc; 56 57 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 58 PetscCall(VecCopy(Q_loc, Qbc)); 59 PetscCall(VecZeroEntries(Q_loc)); 60 PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 61 PetscCall(VecAXPY(Qbc, -1., Q_loc)); 62 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 63 } 64 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 65 66 { // Set boundary mask to zero out essential BCs 67 Vec boundary_mask, ones; 68 69 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 70 PetscCall(DMGetGlobalVector(dm, &ones)); 71 PetscCall(VecZeroEntries(boundary_mask)); 72 PetscCall(VecSet(ones, 1.0)); 73 PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask)); 74 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 75 PetscCall(DMRestoreGlobalVector(dm, &ones)); 76 } 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79 80 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 81 Vec grad_FVM) { 82 Vec Qbc, boundary_mask; 83 84 PetscFunctionBeginUser; 85 // Mask (zero) Strong BC entries 86 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 87 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 88 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 89 90 PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 91 PetscCall(VecAXPY(Q_loc, 1., Qbc)); 92 PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 93 PetscFunctionReturn(PETSC_SUCCESS); 94 } 95 96 static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) { 97 PetscFunctionBeginUser; 98 if (file_type == PETSC_INT32) { 99 PetscInt32 val; 100 PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32)); 101 *out = val; 102 } else if (file_type == PETSC_INT64) { 103 PetscInt64 val; 104 PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64)); 105 *out = val; 106 } else { 107 PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT)); 108 } 109 PetscFunctionReturn(PETSC_SUCCESS); 110 } 111 112 // @brief Load vector from binary file, possibly with embedded solution time and step number 113 PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 114 PetscInt file_step_number; 115 PetscInt32 token; 116 PetscReal file_time; 117 PetscDataType file_type = PETSC_INT32; 118 119 PetscFunctionBeginUser; 120 PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 121 if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 122 token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 123 if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32; 124 else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64; 125 PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type)); 126 PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 127 if (time) *time = file_time; 128 if (step_number) *step_number = file_step_number; 129 } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 130 PetscInt length, N; 131 PetscCall(BinaryReadIntoInt(viewer, &length, file_type)); 132 PetscCall(VecGetSize(Q, &N)); 133 PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N); 134 PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 135 } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 136 137 PetscCall(VecLoad(Q, viewer)); 138 PetscFunctionReturn(PETSC_SUCCESS); 139 } 140 141 // Compare reference solution values with current test run for CI 142 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 143 Vec Qref; 144 PetscViewer viewer; 145 PetscReal error, Qrefnorm; 146 MPI_Comm comm = PetscObjectComm((PetscObject)Q); 147 148 PetscFunctionBeginUser; 149 // Read reference file 150 PetscCall(VecDuplicate(Q, &Qref)); 151 PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given"); 152 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 153 PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 154 155 // Compute error with respect to reference solution 156 PetscCall(VecAXPY(Q, -1.0, Qref)); 157 PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 158 PetscCall(VecScale(Q, 1. / Qrefnorm)); 159 PetscCall(VecNorm(Q, NORM_MAX, &error)); 160 161 // Check error 162 if (error > app_ctx->test_tol) { 163 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 164 } 165 166 // Cleanup 167 PetscCall(PetscViewerDestroy(&viewer)); 168 PetscCall(VecDestroy(&Qref)); 169 PetscFunctionReturn(PETSC_SUCCESS); 170 } 171 172 // Get error for problems with exact solutions 173 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 174 PetscInt loc_nodes; 175 Vec Q_exact, Q_exact_loc; 176 PetscReal rel_error, norm_error, norm_exact; 177 178 PetscFunctionBeginUser; 179 // Get exact solution at final time 180 PetscCall(DMGetGlobalVector(dm, &Q_exact)); 181 PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 182 PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 183 PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 184 185 // Get |exact solution - obtained solution| 186 PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 187 PetscCall(VecAXPY(Q, -1.0, Q_exact)); 188 PetscCall(VecNorm(Q, NORM_1, &norm_error)); 189 190 rel_error = norm_error / norm_exact; 191 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 192 PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 193 PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 // Post-processing 198 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) { 199 PetscInt steps; 200 TSConvergedReason reason; 201 202 PetscFunctionBeginUser; 203 // Print relative error 204 if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) { 205 PetscCall(PrintError(ceed_data, dm, user, Q, final_time)); 206 } 207 208 // Print final time and number of steps 209 PetscCall(TSGetStepNumber(ts, &steps)); 210 PetscCall(TSGetConvergedReason(ts, &reason)); 211 if (user->app_ctx->test_type == TESTTYPE_NONE) { 212 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 213 steps, (double)final_time)); 214 } 215 216 // Output numerical values from command line 217 PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 218 219 // Compare reference solution values with current test run for CI 220 if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 221 PetscCall(RegressionTest(user->app_ctx, Q)); 222 } 223 PetscFunctionReturn(PETSC_SUCCESS); 224 } 225 226 const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 227 const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 228 const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 229 230 // Gather initial Q values in case of continuation of simulation 231 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 232 PetscViewer viewer; 233 234 PetscFunctionBeginUser; 235 PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 236 PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 237 PetscCall(PetscViewerDestroy(&viewer)); 238 PetscFunctionReturn(PETSC_SUCCESS); 239 } 240 241 // Free a plain data context that was allocated using PETSc; returning libCEED error codes 242 int FreeContextPetsc(void *data) { 243 if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 244 return CEED_ERROR_SUCCESS; 245 } 246 247 // Return mass qfunction specification for number of components N 248 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 249 PetscFunctionBeginUser; 250 switch (N) { 251 case 1: 252 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 253 break; 254 case 5: 255 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 256 break; 257 case 7: 258 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 259 break; 260 case 9: 261 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 262 break; 263 case 22: 264 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 265 break; 266 default: 267 SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 268 } 269 270 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 271 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 272 PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 273 PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N)); 274 PetscFunctionReturn(PETSC_SUCCESS); 275 } 276 277 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 278 PetscFunctionBeginUser; 279 if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 280 281 PetscCall(DMDestroy(&context->dm)); 282 PetscCall(KSPDestroy(&context->ksp)); 283 284 PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 285 286 PetscCall(PetscFree(context)); 287 PetscFunctionReturn(PETSC_SUCCESS); 288 } 289 290 /* 291 * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 292 * 293 * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 294 * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 295 * 296 * 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. 297 * 298 * @param[in] comm MPI_Comm for the program 299 * @param[in] path Path to the file 300 * @param[in] char_array_len Length of the character array that should contain each line 301 * @param[out] dims Dimensions of the file, taken from the first line of the file 302 * @param[out] fp File pointer to the opened file 303 */ 304 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 305 FILE **fp) { 306 int ndims; 307 char line[char_array_len]; 308 char **array; 309 310 PetscFunctionBeginUser; 311 PetscCall(PetscFOpen(comm, path, "r", fp)); 312 PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 313 PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 314 PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 315 316 for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 317 PetscCall(PetscStrToArrayDestroy(ndims, array)); 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 /* 322 * @brief Get the number of rows for the PHASTA file at path. 323 * 324 * 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. 325 * 326 * @param[in] comm MPI_Comm for the program 327 * @param[in] path Path to the file 328 * @param[out] nrows Number of rows 329 */ 330 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 331 const PetscInt char_array_len = 512; 332 PetscInt dims[2]; 333 FILE *fp; 334 335 PetscFunctionBeginUser; 336 PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 337 *nrows = dims[0]; 338 PetscCall(PetscFClose(comm, fp)); 339 PetscFunctionReturn(PETSC_SUCCESS); 340 } 341 342 PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 343 PetscInt dims[2]; 344 int ndims; 345 FILE *fp; 346 const PetscInt char_array_len = 512; 347 char line[char_array_len]; 348 char **row_array; 349 350 PetscFunctionBeginUser; 351 PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 352 353 for (PetscInt i = 0; i < dims[0]; i++) { 354 PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 355 PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 356 PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 357 "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 358 359 for (PetscInt j = 0; j < dims[1]; j++) { 360 array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 361 } 362 } 363 364 PetscCall(PetscFClose(comm, fp)); 365 PetscFunctionReturn(PETSC_SUCCESS); 366 } 367 368 // Print information about the given simulation run 369 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) { 370 Ceed ceed = user->ceed; 371 MPI_Comm comm = PetscObjectComm((PetscObject)ts); 372 373 PetscFunctionBeginUser; 374 // Header and rank 375 char host_name[PETSC_MAX_PATH_LEN]; 376 PetscMPIInt rank, comm_size; 377 PetscCall(PetscGetHostName(host_name, sizeof host_name)); 378 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 379 PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 380 PetscCall(PetscPrintf(comm, 381 "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 382 " MPI:\n" 383 " Host Name : %s\n" 384 " Total ranks : %d\n", 385 host_name, comm_size)); 386 387 // Problem specific info 388 PetscCall(problem->print_info(user, problem, user->app_ctx)); 389 390 // libCEED 391 const char *used_resource; 392 CeedMemType mem_type_backend; 393 PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 394 PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 395 PetscCall(PetscPrintf(comm, 396 " libCEED:\n" 397 " libCEED Backend : %s\n" 398 " libCEED Backend MemType : %s\n", 399 used_resource, CeedMemTypes[mem_type_backend])); 400 // PETSc 401 VecType vec_type; 402 char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 403 if (problem->dim == 2) box_faces_str[3] = '\0'; 404 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 405 PetscCall(DMGetVecType(user->dm, &vec_type)); 406 PetscCall(PetscPrintf(comm, 407 " PETSc:\n" 408 " Box Faces : %s\n" 409 " DM VecType : %s\n" 410 " Time Stepping Scheme : %s\n", 411 box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 412 { 413 char pmat_type_str[PETSC_MAX_PATH_LEN]; 414 MatType amat_type, pmat_type; 415 Mat Amat, Pmat; 416 TSIJacobianFn *ijacob_function; 417 418 PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 419 PetscCall(MatGetType(Amat, &amat_type)); 420 PetscCall(MatGetType(Pmat, &pmat_type)); 421 422 PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 423 if (!strcmp(pmat_type, MATCEED)) { 424 MatType pmat_coo_type; 425 char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 426 427 PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 428 PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 429 PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 430 } 431 if (ijacob_function) { 432 PetscCall(PetscPrintf(comm, 433 " IJacobian A MatType : %s\n" 434 " IJacobian P MatType : %s\n", 435 amat_type, pmat_type_str)); 436 } 437 } 438 if (user->app_ctx->cont_steps) { 439 PetscCall(PetscPrintf(comm, 440 " Continue:\n" 441 " Filename: : %s\n" 442 " Step: : %" PetscInt_FMT "\n" 443 " Time: : %g\n", 444 user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 445 } 446 // Mesh 447 const PetscInt num_comp_q = 5; 448 PetscInt glob_dofs, owned_dofs, local_dofs; 449 const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 450 PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 451 PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 452 PetscCall(PetscPrintf(comm, 453 " Mesh:\n" 454 " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 455 " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 456 " Global DoFs : %" PetscInt_FMT "\n" 457 " DoFs per node : %" PetscInt_FMT "\n" 458 " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 459 num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 460 // -- Get Partition Statistics 461 PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 462 { 463 PetscInt *gather_buffer = NULL; 464 PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 465 PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 466 if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 467 468 PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 469 if (!rank) { 470 PetscCall(PetscSortInt(comm_size, gather_buffer)); 471 part_owned_dofs[0] = gather_buffer[0]; // min 472 part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 473 part_owned_dofs[2] = gather_buffer[median_index]; // median 474 PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 475 PetscCall(PetscPrintf( 476 comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 477 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)); 478 } 479 480 PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 481 if (!rank) { 482 PetscCall(PetscSortInt(comm_size, gather_buffer)); 483 part_local_dofs[0] = gather_buffer[0]; // min 484 part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 485 part_local_dofs[2] = gather_buffer[median_index]; // median 486 PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 487 PetscCall(PetscPrintf( 488 comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 489 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)); 490 } 491 492 if (comm_size != 1) { 493 PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 494 { 495 PetscSF sf; 496 PetscInt nrranks, niranks; 497 const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 498 const PetscMPIInt *rranks, *iranks; 499 PetscCall(DMGetSectionSF(user->dm, &sf)); 500 PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 501 PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 502 for (PetscInt i = 0; i < nrranks; i++) { 503 if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 504 num_remote_roots_total += roffset[i + 1] - roffset[i]; 505 num_ghost_interface_ranks++; 506 } 507 for (PetscInt i = 0; i < niranks; i++) { 508 if (iranks[i] == rank) continue; 509 num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 510 num_owned_interface_ranks++; 511 } 512 } 513 PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 514 if (!rank) { 515 PetscCall(PetscSortInt(comm_size, gather_buffer)); 516 part_boundary_dofs[0] = gather_buffer[0]; // min 517 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 518 part_boundary_dofs[2] = gather_buffer[median_index]; // median 519 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 520 PetscCall(PetscPrintf( 521 comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 522 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 523 part_shared_dof_ratio)); 524 } 525 526 PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 527 if (!rank) { 528 PetscCall(PetscSortInt(comm_size, gather_buffer)); 529 part_neighbors[0] = gather_buffer[0]; // min 530 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 531 part_neighbors[2] = gather_buffer[median_index]; // median 532 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 533 PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 534 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 535 } 536 537 PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 538 if (!rank) { 539 PetscCall(PetscSortInt(comm_size, gather_buffer)); 540 part_boundary_dofs[0] = gather_buffer[0]; // min 541 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 542 part_boundary_dofs[2] = gather_buffer[median_index]; // median 543 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 544 PetscCall(PetscPrintf( 545 comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 546 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 547 part_shared_dof_ratio)); 548 } 549 550 PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 551 if (!rank) { 552 PetscCall(PetscSortInt(comm_size, gather_buffer)); 553 part_neighbors[0] = gather_buffer[0]; // min 554 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 555 part_neighbors[2] = gather_buffer[median_index]; // median 556 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 557 PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 558 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 559 } 560 } 561 562 if (!rank) PetscCall(PetscFree(gather_buffer)); 563 } 564 PetscFunctionReturn(PETSC_SUCCESS); 565 } 566