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 186 PetscCall(DMDestroy(&context->dm)); 187 PetscCall(KSPDestroy(&context->ksp)); 188 189 PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 190 191 PetscCall(PetscFree(context)); 192 PetscFunctionReturn(PETSC_SUCCESS); 193 } 194 195 /** 196 @brief Sets the value of an option if it is not already set 197 198 @param[in,out] options `PetscOptions` database, or `NULL` for default global database 199 @param[in] name Name of the option (with `-` prepended to it) 200 @param[in] value Value to set the option to 201 **/ 202 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]) { 203 PetscBool has_option; 204 205 PetscFunctionBeginUser; 206 PetscCall(PetscOptionsHasName(options, NULL, name, &has_option)); 207 if (!has_option) PetscCall(PetscOptionsSetValue(options, name, value)); 208 PetscFunctionReturn(PETSC_SUCCESS); 209 } 210 211 // Print information about the given simulation run 212 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts) { 213 Ceed ceed = honee->ceed; 214 MPI_Comm comm = PetscObjectComm((PetscObject)ts); 215 216 PetscFunctionBeginUser; 217 // Header and rank 218 char host_name[PETSC_MAX_PATH_LEN]; 219 PetscMPIInt rank, comm_size; 220 PetscCall(PetscGetHostName(host_name, sizeof host_name)); 221 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 222 PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 223 PetscCall(PetscPrintf(comm, 224 "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 225 " MPI:\n" 226 " Host Name : %s\n" 227 " Total ranks : %d\n", 228 host_name, comm_size)); 229 230 // Problem specific info 231 PetscCall(problem->print_info(honee, problem, honee->app_ctx)); 232 233 // libCEED 234 const char *used_resource; 235 CeedMemType mem_type_backend; 236 PetscCallCeed(ceed, CeedGetResource(honee->ceed, &used_resource)); 237 PetscCallCeed(ceed, CeedGetPreferredMemType(honee->ceed, &mem_type_backend)); 238 PetscCall(PetscPrintf(comm, 239 " libCEED:\n" 240 " libCEED Backend : %s\n" 241 " libCEED Backend MemType : %s\n", 242 used_resource, CeedMemTypes[mem_type_backend])); 243 // PETSc 244 { 245 VecType vec_type; 246 char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 247 PetscInt dim; 248 249 PetscCall(DMGetDimension(honee->dm, &dim)); 250 if (dim == 2) box_faces_str[3] = '\0'; 251 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 252 PetscCall(DMGetVecType(honee->dm, &vec_type)); 253 PetscCall(PetscPrintf(comm, 254 " PETSc:\n" 255 " Box Faces : %s\n" 256 " DM VecType : %s\n" 257 " Time Stepping Scheme : %s\n", 258 box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 259 } 260 { 261 char pmat_type_str[PETSC_MAX_PATH_LEN]; 262 MatType amat_type, pmat_type; 263 Mat Amat, Pmat; 264 TSIJacobianFn *ijacob_function; 265 266 PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 267 PetscCall(MatGetType(Amat, &amat_type)); 268 PetscCall(MatGetType(Pmat, &pmat_type)); 269 270 PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 271 if (!strcmp(pmat_type, MATCEED)) { 272 MatType pmat_coo_type; 273 char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 274 275 PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 276 PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 277 PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 278 } 279 if (ijacob_function) { 280 PetscCall(PetscPrintf(comm, 281 " IJacobian A MatType : %s\n" 282 " IJacobian P MatType : %s\n", 283 amat_type, pmat_type_str)); 284 } 285 } 286 if (honee->app_ctx->use_continue_file) { 287 PetscCall(PetscPrintf(comm, 288 " Continue:\n" 289 " Filename: : %s\n" 290 " Step: : %" PetscInt_FMT "\n" 291 " Time: : %g\n", 292 honee->app_ctx->cont_file, honee->app_ctx->cont_steps, honee->app_ctx->cont_time)); 293 } 294 // Mesh 295 const PetscInt num_comp_q = 5; 296 PetscInt glob_dofs, owned_dofs, local_dofs; 297 const CeedInt num_P = honee->app_ctx->degree + 1, num_Q = num_P + honee->app_ctx->q_extra; 298 PetscCall(DMGetGlobalVectorInfo(honee->dm, &owned_dofs, &glob_dofs, NULL)); 299 PetscCall(DMGetLocalVectorInfo(honee->dm, &local_dofs, NULL, NULL)); 300 PetscCall(PetscPrintf(comm, 301 " Mesh:\n" 302 " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 303 " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 304 " Global DoFs : %" PetscInt_FMT "\n" 305 " DoFs per node : %" PetscInt_FMT "\n" 306 " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 307 num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 308 // -- Get Partition Statistics 309 PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 310 { 311 PetscInt *gather_buffer = NULL; 312 PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 313 PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 314 if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 315 316 PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 317 if (!rank) { 318 PetscCall(PetscSortInt(comm_size, gather_buffer)); 319 part_owned_dofs[0] = gather_buffer[0]; // min 320 part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 321 part_owned_dofs[2] = gather_buffer[median_index]; // median 322 PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 323 PetscCall(PetscPrintf( 324 comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 325 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)); 326 } 327 328 PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 329 if (!rank) { 330 PetscCall(PetscSortInt(comm_size, gather_buffer)); 331 part_local_dofs[0] = gather_buffer[0]; // min 332 part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 333 part_local_dofs[2] = gather_buffer[median_index]; // median 334 PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 335 PetscCall(PetscPrintf( 336 comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 337 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)); 338 } 339 340 if (comm_size != 1) { 341 PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 342 { 343 PetscSF sf; 344 PetscMPIInt nrranks, niranks; 345 const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 346 const PetscMPIInt *rranks, *iranks; 347 PetscCall(DMGetSectionSF(honee->dm, &sf)); 348 PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 349 PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 350 for (PetscInt i = 0; i < nrranks; i++) { 351 if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 352 num_remote_roots_total += roffset[i + 1] - roffset[i]; 353 num_ghost_interface_ranks++; 354 } 355 for (PetscInt i = 0; i < niranks; i++) { 356 if (iranks[i] == rank) continue; 357 num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 358 num_owned_interface_ranks++; 359 } 360 } 361 PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 362 if (!rank) { 363 PetscCall(PetscSortInt(comm_size, gather_buffer)); 364 part_boundary_dofs[0] = gather_buffer[0]; // min 365 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 366 part_boundary_dofs[2] = gather_buffer[median_index]; // median 367 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 368 PetscCall(PetscPrintf( 369 comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 370 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 371 part_shared_dof_ratio)); 372 } 373 374 PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 375 if (!rank) { 376 PetscCall(PetscSortInt(comm_size, gather_buffer)); 377 part_neighbors[0] = gather_buffer[0]; // min 378 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 379 part_neighbors[2] = gather_buffer[median_index]; // median 380 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 381 PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 382 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 383 } 384 385 PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 386 if (!rank) { 387 PetscCall(PetscSortInt(comm_size, gather_buffer)); 388 part_boundary_dofs[0] = gather_buffer[0]; // min 389 part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 390 part_boundary_dofs[2] = gather_buffer[median_index]; // median 391 PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 392 PetscCall(PetscPrintf( 393 comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 394 num_comp_q, part_boundary_dofs[0] / num_comp_q, part_boundary_dofs[1] / num_comp_q, part_boundary_dofs[2] / num_comp_q, 395 part_shared_dof_ratio)); 396 } 397 398 PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 399 if (!rank) { 400 PetscCall(PetscSortInt(comm_size, gather_buffer)); 401 part_neighbors[0] = gather_buffer[0]; // min 402 part_neighbors[1] = gather_buffer[comm_size - 1]; // max 403 part_neighbors[2] = gather_buffer[median_index]; // median 404 PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 405 PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 406 part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 407 } 408 } 409 410 if (!rank) PetscCall(PetscFree(gather_buffer)); 411 } 412 PetscFunctionReturn(PETSC_SUCCESS); 413 } 414