1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Miscellaneous utility functions 6a515125bSLeila Ghaffari 7e419654dSJeremy L Thompson #include <ceed.h> 8e419654dSJeremy L Thompson #include <petscdm.h> 9926a6279SJames Wright #include <petscsf.h> 10e419654dSJeremy L Thompson #include <petscts.h> 11e419654dSJeremy L Thompson 12149fb536SJames Wright #include <navierstokes.h> 13a515125bSLeila Ghaffari 14e3663b90SJames Wright PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time) { 150c373b74SJames Wright Ceed ceed = honee->ceed; 16b2948607SJames Wright CeedVector mult_vec; 17b2948607SJames Wright PetscMemType m_mem_type; 18b2948607SJames Wright Vec Multiplicity, Multiplicity_loc; 19b2948607SJames Wright 20a515125bSLeila Ghaffari PetscFunctionBeginUser; 21e3663b90SJames Wright if (honee->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(honee->op_ics_ctx->op, honee->phys->ics_time_label, &time)); 22e3663b90SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, honee->op_ics_ctx)); 23a515125bSLeila Ghaffari 24e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &mult_vec, NULL)); 25a515125bSLeila Ghaffari 26a515125bSLeila Ghaffari // -- Get multiplicity 27b2948607SJames Wright PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 28a7dac1d5SJames Wright PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 29e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(honee->elem_restr_q, mult_vec)); 30a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 31a515125bSLeila Ghaffari 32b2948607SJames Wright PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 33b2948607SJames Wright PetscCall(VecZeroEntries(Multiplicity)); 34b2948607SJames Wright PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 35a515125bSLeila Ghaffari 36a515125bSLeila Ghaffari // -- Fix multiplicity 37b2948607SJames Wright PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 38b2948607SJames Wright PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 39a515125bSLeila Ghaffari 40b2948607SJames Wright PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 41b2948607SJames Wright PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 42b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 43d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 44a515125bSLeila Ghaffari } 45a515125bSLeila Ghaffari 46c56e8d5bSJames Wright // Record boundary values from initial condition 47c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 48c56e8d5bSJames Wright PetscFunctionBeginUser; 49313f2f1eSJames Wright { // Capture initial condition values in Qbc 50313f2f1eSJames Wright Vec Qbc; 51313f2f1eSJames Wright 52c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 53c56e8d5bSJames Wright PetscCall(VecCopy(Q_loc, Qbc)); 54c56e8d5bSJames Wright PetscCall(VecZeroEntries(Q_loc)); 55c56e8d5bSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 56c56e8d5bSJames Wright PetscCall(VecAXPY(Qbc, -1., Q_loc)); 57c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 58313f2f1eSJames Wright } 59c56e8d5bSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 60c56e8d5bSJames Wright 61313f2f1eSJames Wright { // Set boundary mask to zero out essential BCs 62313f2f1eSJames Wright Vec boundary_mask, ones; 63313f2f1eSJames Wright 64c56e8d5bSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 65313f2f1eSJames Wright PetscCall(DMGetGlobalVector(dm, &ones)); 66c56e8d5bSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 67313f2f1eSJames Wright PetscCall(VecSet(ones, 1.0)); 68313f2f1eSJames Wright PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask)); 69c56e8d5bSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 70313f2f1eSJames Wright PetscCall(DMRestoreGlobalVector(dm, &ones)); 71313f2f1eSJames Wright } 72c56e8d5bSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 73c56e8d5bSJames Wright } 74c56e8d5bSJames Wright 75c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 762b916ea7SJeremy L Thompson Vec grad_FVM) { 779d437337SJames Wright Vec Qbc, boundary_mask; 78a515125bSLeila Ghaffari 7906f41313SJames Wright PetscFunctionBeginUser; 802eb7bf1fSJames Wright // Mask (zero) Strong BC entries 819d437337SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 829d437337SJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 839d437337SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 849d437337SJames Wright 852b916ea7SJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 862b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 872b916ea7SJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 88d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 89a515125bSLeila Ghaffari } 90a515125bSLeila Ghaffari 91a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 92c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 93a515125bSLeila Ghaffari Vec Qref; 94a515125bSLeila Ghaffari PetscViewer viewer; 95a515125bSLeila Ghaffari PetscReal error, Qrefnorm; 96e7754af5SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 97a515125bSLeila Ghaffari 9806f41313SJames Wright PetscFunctionBeginUser; 99a515125bSLeila Ghaffari // Read reference file 1002b916ea7SJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 1014c07ec22SJames Wright PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given"); 102e7754af5SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 103360b37c9SJames Wright PetscCall(HoneeLoadBinaryVec(viewer, Qref, NULL, NULL)); 104a515125bSLeila Ghaffari 105a515125bSLeila Ghaffari // Compute error with respect to reference solution 1062b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1072b916ea7SJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1082b916ea7SJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1092b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 110a515125bSLeila Ghaffari 111a515125bSLeila Ghaffari // Check error 112a515125bSLeila Ghaffari if (error > app_ctx->test_tol) { 1132b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 114a515125bSLeila Ghaffari } 115a515125bSLeila Ghaffari 116a515125bSLeila Ghaffari // Cleanup 1172b916ea7SJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1182b916ea7SJeremy L Thompson PetscCall(VecDestroy(&Qref)); 119d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 120a515125bSLeila Ghaffari } 121a515125bSLeila Ghaffari 122a515125bSLeila Ghaffari // Get error for problems with exact solutions 123e3663b90SJames Wright PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time) { 124a515125bSLeila Ghaffari PetscInt loc_nodes; 125a515125bSLeila Ghaffari Vec Q_exact, Q_exact_loc; 126a515125bSLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 127a515125bSLeila Ghaffari 12806f41313SJames Wright PetscFunctionBeginUser; 129a515125bSLeila Ghaffari // Get exact solution at final time 130b2948607SJames Wright PetscCall(DMGetGlobalVector(dm, &Q_exact)); 1312b916ea7SJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1322b916ea7SJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 133e3663b90SJames Wright PetscCall(ICs_FixMultiplicity(dm, honee, Q_exact_loc, Q_exact, final_time)); 134a515125bSLeila Ghaffari 135a515125bSLeila Ghaffari // Get |exact solution - obtained solution| 1362b916ea7SJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1372b916ea7SJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1382b916ea7SJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 139a515125bSLeila Ghaffari 140a515125bSLeila Ghaffari rel_error = norm_error / norm_exact; 1412b916ea7SJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 1422b916ea7SJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 143b2948607SJames Wright PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 144d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 145a515125bSLeila Ghaffari } 146a515125bSLeila Ghaffari 147a515125bSLeila Ghaffari // Post-processing 148e3663b90SJames Wright PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time) { 149a515125bSLeila Ghaffari PetscInt steps; 150f0784ed3SJed Brown TSConvergedReason reason; 151a515125bSLeila Ghaffari 15206f41313SJames Wright PetscFunctionBeginUser; 153a515125bSLeila Ghaffari // Print relative error 1540c373b74SJames Wright if (problem->compute_exact_solution_error && honee->app_ctx->test_type == TESTTYPE_NONE) { 155e3663b90SJames Wright PetscCall(PrintError(dm, honee, Q, final_time)); 156a515125bSLeila Ghaffari } 157a515125bSLeila Ghaffari 158a515125bSLeila Ghaffari // Print final time and number of steps 1592b916ea7SJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 160f0784ed3SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 1610c373b74SJames Wright if (honee->app_ctx->test_type == TESTTYPE_NONE) { 162f0784ed3SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 163f0784ed3SJed Brown steps, (double)final_time)); 164a515125bSLeila Ghaffari } 165a515125bSLeila Ghaffari 166a515125bSLeila Ghaffari // Output numerical values from command line 1672b916ea7SJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 168a515125bSLeila Ghaffari 169a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI 1700c373b74SJames Wright if (honee->app_ctx->test_type == TESTTYPE_SOLVER) { 1710c373b74SJames Wright PetscCall(RegressionTest(honee->app_ctx, Q)); 172a515125bSLeila Ghaffari } 173d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 174a515125bSLeila Ghaffari } 175a515125bSLeila Ghaffari 17615a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 17715a3537eSJed Brown int FreeContextPetsc(void *data) { 1782b916ea7SJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 17915a3537eSJed Brown return CEED_ERROR_SUCCESS; 18015a3537eSJed Brown } 1819f59f36eSJames Wright 1824ea616f4SJames Wright /** 1834ea616f4SJames Wright @brief Destroy `NodalProjectionData` object 1844ea616f4SJames Wright 1854ea616f4SJames Wright @param[in] context `NodalProjectionData` object to destroy 1864ea616f4SJames Wright **/ 1874ea616f4SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData *context) { 1884ea616f4SJames Wright NodalProjectionData context_ = *context; 189*14bd2a07SJames Wright 190457a5831SJames Wright PetscFunctionBeginUser; 1914ea616f4SJames Wright if (context_ == NULL) PetscFunctionReturn(PETSC_SUCCESS); 1924ea616f4SJames Wright PetscCall(DMDestroy(&context_->dm)); 1934ea616f4SJames Wright PetscCall(KSPDestroy(&context_->ksp)); 1944ea616f4SJames Wright PetscCall(OperatorApplyContextDestroy(context_->l2_rhs_ctx)); 1954ea616f4SJames Wright PetscCall(PetscFree(context_)); 1964ea616f4SJames Wright *context = NULL; 197d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 198457a5831SJames Wright } 199457a5831SJames Wright 20059572072SJames Wright /** 20159572072SJames Wright @brief Sets the value of an option if it is not already set 20259572072SJames Wright 20359572072SJames Wright @param[in,out] options `PetscOptions` database, or `NULL` for default global database 20459572072SJames Wright @param[in] name Name of the option (with `-` prepended to it) 20559572072SJames Wright @param[in] value Value to set the option to 20659572072SJames Wright **/ 20759572072SJames Wright PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]) { 20859572072SJames Wright PetscBool has_option; 20959572072SJames Wright 21059572072SJames Wright PetscFunctionBeginUser; 21159572072SJames Wright PetscCall(PetscOptionsHasName(options, NULL, name, &has_option)); 21259572072SJames Wright if (!has_option) PetscCall(PetscOptionsSetValue(options, name, value)); 21359572072SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21459572072SJames Wright } 21559572072SJames Wright 21625125139SJames Wright PetscErrorCode HoneeCalculateDomainSize(Honee honee, PetscScalar *volume) { 21725125139SJames Wright DM dm = honee->dm, dm_coord; 21825125139SJames Wright Ceed ceed = honee->ceed; 21925125139SJames Wright CeedQFunction qf_mass; 22025125139SJames Wright CeedOperator op_mass; 22125125139SJames Wright OperatorApplyContext op_mass_ctx; 22225125139SJames Wright CeedElemRestriction elem_restr_qd; 22325125139SJames Wright CeedVector qdata; 22425125139SJames Wright CeedInt q_data_size, num_comps_x; 22525125139SJames Wright DMLabel domain_label = 0; 22625125139SJames Wright PetscInt label_value = 0; 22725125139SJames Wright Vec u, v; 22825125139SJames Wright 22925125139SJames Wright PetscFunctionBeginUser; 23025125139SJames Wright PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 23125125139SJames Wright 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)); 23225125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comps_x)); 23325125139SJames Wright 23425125139SJames Wright PetscCall(HoneeMassQFunctionCreate(ceed, num_comps_x, q_data_size, &qf_mass)); 23525125139SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 23625125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE)); 23725125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, qdata)); 23825125139SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", honee->elem_restr_x, honee->basis_x, CEED_VECTOR_ACTIVE)); 23925125139SJames Wright 24025125139SJames Wright PetscCall(OperatorApplyContextCreate(NULL, dm_coord, ceed, op_mass, NULL, NULL, NULL, NULL, &op_mass_ctx)); 24125125139SJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_mass, DMReturnVecType(dm), PETSC_COMM_SELF, &u, NULL)); 24225125139SJames Wright PetscCall(DMCreateGlobalVector(dm_coord, &v)); 24325125139SJames Wright PetscCall(VecSet(u, 1.)); 24425125139SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(u, v, op_mass_ctx)); 24525125139SJames Wright PetscCall(VecSum(v, volume)); 24625125139SJames Wright *volume /= num_comps_x; // Correct for number of components != 1 24725125139SJames Wright 24825125139SJames Wright PetscCall(VecDestroy(&u)); 24925125139SJames Wright PetscCall(VecDestroy(&v)); 25025125139SJames Wright PetscCall(OperatorApplyContextDestroy(op_mass_ctx)); 25125125139SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 25225125139SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&qdata)); 25325125139SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 25425125139SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 25525125139SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 25625125139SJames Wright } 25725125139SJames Wright 258926a6279SJames Wright // Print information about the given simulation run 2590c373b74SJames Wright PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts) { 2600c373b74SJames Wright Ceed ceed = honee->ceed; 2613678fae3SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)ts); 2623678fae3SJames Wright 263926a6279SJames Wright PetscFunctionBeginUser; 264926a6279SJames Wright // Header and rank 265926a6279SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 266926a6279SJames Wright PetscMPIInt rank, comm_size; 267926a6279SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 268926a6279SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 269926a6279SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 270926a6279SJames Wright PetscCall(PetscPrintf(comm, 27115747f0fSJames Wright "\n-- HONEE - High-Order Navier-stokes Equation Evaluator --\n" 272926a6279SJames Wright " MPI:\n" 273926a6279SJames Wright " Host Name : %s\n" 274926a6279SJames Wright " Total ranks : %d\n", 275926a6279SJames Wright host_name, comm_size)); 276926a6279SJames Wright 277926a6279SJames Wright // Problem specific info 2780c373b74SJames Wright PetscCall(problem->print_info(honee, problem, honee->app_ctx)); 279926a6279SJames Wright 280926a6279SJames Wright // libCEED 281926a6279SJames Wright const char *used_resource; 282926a6279SJames Wright CeedMemType mem_type_backend; 2830c373b74SJames Wright PetscCallCeed(ceed, CeedGetResource(honee->ceed, &used_resource)); 2840c373b74SJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(honee->ceed, &mem_type_backend)); 285926a6279SJames Wright PetscCall(PetscPrintf(comm, 286926a6279SJames Wright " libCEED:\n" 287926a6279SJames Wright " libCEED Backend : %s\n" 288926a6279SJames Wright " libCEED Backend MemType : %s\n", 289926a6279SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 290926a6279SJames Wright // PETSc 29166d54740SJames Wright { 2923678fae3SJames Wright VecType vec_type; 293926a6279SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 29466d54740SJames Wright PetscInt dim; 29566d54740SJames Wright 2960c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 29766d54740SJames Wright if (dim == 2) box_faces_str[3] = '\0'; 298926a6279SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 2990c373b74SJames Wright PetscCall(DMGetVecType(honee->dm, &vec_type)); 300926a6279SJames Wright PetscCall(PetscPrintf(comm, 301926a6279SJames Wright " PETSc:\n" 302926a6279SJames Wright " Box Faces : %s\n" 303926a6279SJames Wright " DM VecType : %s\n" 304926a6279SJames Wright " Time Stepping Scheme : %s\n", 3053678fae3SJames Wright box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 30666d54740SJames Wright } 3073678fae3SJames Wright { 3083678fae3SJames Wright char pmat_type_str[PETSC_MAX_PATH_LEN]; 3093678fae3SJames Wright MatType amat_type, pmat_type; 3103678fae3SJames Wright Mat Amat, Pmat; 3113678fae3SJames Wright TSIJacobianFn *ijacob_function; 3123678fae3SJames Wright 3133678fae3SJames Wright PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL)); 3143678fae3SJames Wright PetscCall(MatGetType(Amat, &amat_type)); 3153678fae3SJames Wright PetscCall(MatGetType(Pmat, &pmat_type)); 3163678fae3SJames Wright 3173678fae3SJames Wright PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str))); 3183678fae3SJames Wright if (!strcmp(pmat_type, MATCEED)) { 3193678fae3SJames Wright MatType pmat_coo_type; 3203678fae3SJames Wright char pmat_coo_type_str[PETSC_MAX_PATH_LEN]; 3213678fae3SJames Wright 3223678fae3SJames Wright PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type)); 3233678fae3SJames Wright PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type)); 3243678fae3SJames Wright PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str))); 3253678fae3SJames Wright } 3263678fae3SJames Wright if (ijacob_function) { 3273678fae3SJames Wright PetscCall(PetscPrintf(comm, 3283678fae3SJames Wright " IJacobian A MatType : %s\n" 3293678fae3SJames Wright " IJacobian P MatType : %s\n", 3303678fae3SJames Wright amat_type, pmat_type_str)); 3313678fae3SJames Wright } 3323678fae3SJames Wright } 333481d14cbSJames Wright if (honee->app_ctx->use_continue_file) { 334926a6279SJames Wright PetscCall(PetscPrintf(comm, 335926a6279SJames Wright " Continue:\n" 336926a6279SJames Wright " Filename: : %s\n" 337926a6279SJames Wright " Step: : %" PetscInt_FMT "\n" 338926a6279SJames Wright " Time: : %g\n", 3390c373b74SJames Wright honee->app_ctx->cont_file, honee->app_ctx->cont_steps, honee->app_ctx->cont_time)); 340926a6279SJames Wright } 341926a6279SJames Wright // Mesh 342926a6279SJames Wright const PetscInt num_comp_q = 5; 343926a6279SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 3440c373b74SJames Wright const CeedInt num_P = honee->app_ctx->degree + 1, num_Q = num_P + honee->app_ctx->q_extra; 3450c373b74SJames Wright PetscCall(DMGetGlobalVectorInfo(honee->dm, &owned_dofs, &glob_dofs, NULL)); 3460c373b74SJames Wright PetscCall(DMGetLocalVectorInfo(honee->dm, &local_dofs, NULL, NULL)); 347926a6279SJames Wright PetscCall(PetscPrintf(comm, 348926a6279SJames Wright " Mesh:\n" 349926a6279SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 350926a6279SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 351926a6279SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 352926a6279SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 353dfeb939dSJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 354dfeb939dSJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 355926a6279SJames Wright // -- Get Partition Statistics 356926a6279SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 357926a6279SJames Wright { 358926a6279SJames Wright PetscInt *gather_buffer = NULL; 359dfeb939dSJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 360926a6279SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 361926a6279SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 362926a6279SJames Wright 363dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 364926a6279SJames Wright if (!rank) { 365926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 366926a6279SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 367926a6279SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 368926a6279SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 369926a6279SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 370dfeb939dSJames Wright PetscCall(PetscPrintf( 371dfeb939dSJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 372926a6279SJames 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)); 373926a6279SJames Wright } 374926a6279SJames Wright 375dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 376926a6279SJames Wright if (!rank) { 377926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 378926a6279SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 379926a6279SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 380926a6279SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 381926a6279SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 382dfeb939dSJames Wright PetscCall(PetscPrintf( 383dfeb939dSJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 384926a6279SJames 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)); 385926a6279SJames Wright } 386926a6279SJames Wright 38745abf96eSJames Wright if (comm_size != 1) { 388dfeb939dSJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 389926a6279SJames Wright { 390926a6279SJames Wright PetscSF sf; 3913bbbcc04SJames Wright PetscMPIInt nrranks, niranks; 392dfeb939dSJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 393dfeb939dSJames Wright const PetscMPIInt *rranks, *iranks; 3940c373b74SJames Wright PetscCall(DMGetSectionSF(honee->dm, &sf)); 395926a6279SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 396dfeb939dSJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 397926a6279SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 398926a6279SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 399926a6279SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 400dfeb939dSJames Wright num_ghost_interface_ranks++; 401dfeb939dSJames Wright } 402dfeb939dSJames Wright for (PetscInt i = 0; i < niranks; i++) { 403dfeb939dSJames Wright if (iranks[i] == rank) continue; 404dfeb939dSJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 405dfeb939dSJames Wright num_owned_interface_ranks++; 406926a6279SJames Wright } 407926a6279SJames Wright } 408dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 409926a6279SJames Wright if (!rank) { 410926a6279SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 411dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 412dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 413dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 414dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 415dfeb939dSJames Wright PetscCall(PetscPrintf( 41645abf96eSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 41745abf96eSJames 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, 41845abf96eSJames Wright part_shared_dof_ratio)); 419dfeb939dSJames Wright } 420dfeb939dSJames Wright 421dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 422dfeb939dSJames Wright if (!rank) { 423dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 424dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 425dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 426dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 427dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 428dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 429dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 430dfeb939dSJames Wright } 431dfeb939dSJames Wright 432dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 433dfeb939dSJames Wright if (!rank) { 434dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 435dfeb939dSJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 436dfeb939dSJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 437dfeb939dSJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 438dfeb939dSJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 439dfeb939dSJames Wright PetscCall(PetscPrintf( 44045abf96eSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 44145abf96eSJames 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, 44245abf96eSJames Wright part_shared_dof_ratio)); 443dfeb939dSJames Wright } 444dfeb939dSJames Wright 445dfeb939dSJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 446dfeb939dSJames Wright if (!rank) { 447dfeb939dSJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 448dfeb939dSJames Wright part_neighbors[0] = gather_buffer[0]; // min 449dfeb939dSJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 450dfeb939dSJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 451dfeb939dSJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 452dfeb939dSJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 453dfeb939dSJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 454926a6279SJames Wright } 45545abf96eSJames Wright } 456926a6279SJames Wright 457926a6279SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 458926a6279SJames Wright } 459926a6279SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 460926a6279SJames Wright } 461