13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Miscellaneous utility functions 1077841947SLeila Ghaffari 1149aac155SJeremy L Thompson #include <ceed.h> 1249aac155SJeremy L Thompson #include <petscdm.h> 132e31f396SJames Wright #include <petscsf.h> 1449aac155SJeremy L Thompson #include <petscts.h> 1549aac155SJeremy L Thompson 1677841947SLeila Ghaffari #include "../navierstokes.h" 17ef080ff9SJames Wright #include "../qfunctions/mass.h" 1877841947SLeila Ghaffari 192b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) { 20a424bcd0SJames Wright Ceed ceed = user->ceed; 2177841947SLeila Ghaffari PetscFunctionBeginUser; 2277841947SLeila Ghaffari 2377841947SLeila Ghaffari // --------------------------------------------------------------------------- 24a0add3c9SJed Brown // Update time for evaluation 2577841947SLeila Ghaffari // --------------------------------------------------------------------------- 26a424bcd0SJames Wright if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 2777841947SLeila Ghaffari 2877841947SLeila Ghaffari // --------------------------------------------------------------------------- 2977841947SLeila Ghaffari // ICs 3077841947SLeila Ghaffari // --------------------------------------------------------------------------- 3177841947SLeila Ghaffari // -- CEED Restriction 3277841947SLeila Ghaffari CeedVector q0_ceed; 33a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL)); 3477841947SLeila Ghaffari 3577841947SLeila Ghaffari // -- Place PETSc vector in CEED vector 365263e9c6SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 3777841947SLeila Ghaffari 3877841947SLeila Ghaffari // --------------------------------------------------------------------------- 3977841947SLeila Ghaffari // Fix multiplicity for output of ICs 4077841947SLeila Ghaffari // --------------------------------------------------------------------------- 4177841947SLeila Ghaffari // -- CEED Restriction 4277841947SLeila Ghaffari CeedVector mult_vec; 43a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 4477841947SLeila Ghaffari 4577841947SLeila Ghaffari // -- Place PETSc vector in CEED vector 4677841947SLeila Ghaffari PetscMemType m_mem_type; 4777841947SLeila Ghaffari Vec multiplicity_loc; 482b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &multiplicity_loc)); 49c798d55aSJames Wright PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec)); 5077841947SLeila Ghaffari 5177841947SLeila Ghaffari // -- Get multiplicity 52a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 5377841947SLeila Ghaffari 5477841947SLeila Ghaffari // -- Restore vectors 55c798d55aSJames Wright PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc)); 5677841947SLeila Ghaffari 5777841947SLeila Ghaffari // -- Local-to-Global 5877841947SLeila Ghaffari Vec multiplicity; 592b730f8bSJeremy L Thompson PetscCall(DMGetGlobalVector(dm, &multiplicity)); 602b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(multiplicity)); 612b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity)); 6277841947SLeila Ghaffari 6377841947SLeila Ghaffari // -- Fix multiplicity 642b730f8bSJeremy L Thompson PetscCall(VecPointwiseDivide(Q, Q, multiplicity)); 652b730f8bSJeremy L Thompson PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc)); 6677841947SLeila Ghaffari 6777841947SLeila Ghaffari // -- Restore vectors 682b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc)); 692b730f8bSJeremy L Thompson PetscCall(DMRestoreGlobalVector(dm, &multiplicity)); 7077841947SLeila Ghaffari 7177841947SLeila Ghaffari // Cleanup 72a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 73a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q0_ceed)); 7477841947SLeila Ghaffari 75ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 7677841947SLeila Ghaffari } 7777841947SLeila Ghaffari 782b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 792b730f8bSJeremy L Thompson Vec grad_FVM) { 805571c6fdSJames Wright Vec Qbc, boundary_mask; 8177841947SLeila Ghaffari PetscFunctionBegin; 8277841947SLeila Ghaffari 833722cd23SJames Wright // Mask (zero) Strong BC entries 845571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 855571c6fdSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 865571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 875571c6fdSJames Wright 882b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 892b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 902b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 9177841947SLeila Ghaffari 92ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 9377841947SLeila Ghaffari } 9477841947SLeila Ghaffari 95530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 96530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 970de6a49fSJames Wright PetscInt file_step_number; 980de6a49fSJames Wright PetscInt32 token; 99530ad8c4SKenneth E. Jansen PetscReal file_time; 100530ad8c4SKenneth E. Jansen PetscFunctionBeginUser; 101530ad8c4SKenneth E. Jansen 102530ad8c4SKenneth E. Jansen // Attempt 1030de6a49fSJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 1040de6a49fSJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 1050de6a49fSJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 106530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT)); 107530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 108530ad8c4SKenneth E. Jansen if (time) *time = file_time; 109530ad8c4SKenneth E. Jansen if (step_number) *step_number = file_step_number; 110530ad8c4SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 111530ad8c4SKenneth E. Jansen PetscInt length, N; 112530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT)); 113530ad8c4SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 114530ad8c4SKenneth E. Jansen PetscCheck(length == N, comm, PETSC_ERR_ARG_INCOMP, "File Vec has length %" PetscInt_FMT " but DM has global Vec size %" PetscInt_FMT, length, N); 115530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 116530ad8c4SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 117530ad8c4SKenneth E. Jansen 118530ad8c4SKenneth E. Jansen // Load Q from existent solution 119530ad8c4SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 120530ad8c4SKenneth E. Jansen 121ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 122530ad8c4SKenneth E. Jansen } 123530ad8c4SKenneth E. Jansen 12477841947SLeila Ghaffari // Compare reference solution values with current test run for CI 12577841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) { 12677841947SLeila Ghaffari Vec Qref; 12777841947SLeila Ghaffari PetscViewer viewer; 12877841947SLeila Ghaffari PetscReal error, Qrefnorm; 129530ad8c4SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 13077841947SLeila Ghaffari PetscFunctionBegin; 13177841947SLeila Ghaffari 13277841947SLeila Ghaffari // Read reference file 1332b730f8bSJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 134530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 135530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 13677841947SLeila Ghaffari 13777841947SLeila Ghaffari // Compute error with respect to reference solution 1382b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1392b730f8bSJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1402b730f8bSJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1412b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 14277841947SLeila Ghaffari 14377841947SLeila Ghaffari // Check error 14477841947SLeila Ghaffari if (error > app_ctx->test_tol) { 1452b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 14677841947SLeila Ghaffari } 14777841947SLeila Ghaffari 14877841947SLeila Ghaffari // Cleanup 1492b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1502b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Qref)); 15177841947SLeila Ghaffari 152ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 15377841947SLeila Ghaffari } 15477841947SLeila Ghaffari 15577841947SLeila Ghaffari // Get error for problems with exact solutions 1562b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 15777841947SLeila Ghaffari PetscInt loc_nodes; 15877841947SLeila Ghaffari Vec Q_exact, Q_exact_loc; 15977841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 16077841947SLeila Ghaffari PetscFunctionBegin; 16177841947SLeila Ghaffari 16277841947SLeila Ghaffari // Get exact solution at final time 1632b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &Q_exact)); 1642b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1652b730f8bSJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1662b730f8bSJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 16777841947SLeila Ghaffari 16877841947SLeila Ghaffari // Get |exact solution - obtained solution| 1692b730f8bSJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1702b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1712b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 17277841947SLeila Ghaffari 17377841947SLeila Ghaffari // Compute relative error 17477841947SLeila Ghaffari rel_error = norm_error / norm_exact; 17577841947SLeila Ghaffari 17677841947SLeila Ghaffari // Output relative error 1772b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 17877841947SLeila Ghaffari // Cleanup 1792b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 1802b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Q_exact)); 18177841947SLeila Ghaffari 182ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 18377841947SLeila Ghaffari } 18477841947SLeila Ghaffari 18577841947SLeila Ghaffari // Post-processing 1862b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 18777841947SLeila Ghaffari PetscInt steps; 188cf7a0454SJed Brown TSConvergedReason reason; 18977841947SLeila Ghaffari PetscFunctionBegin; 19077841947SLeila Ghaffari 19177841947SLeila Ghaffari // Print relative error 1928fb33541SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 1932b730f8bSJeremy L Thompson PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time)); 19477841947SLeila Ghaffari } 19577841947SLeila Ghaffari 19677841947SLeila Ghaffari // Print final time and number of steps 1972b730f8bSJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 198cf7a0454SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 1998fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 200cf7a0454SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 201cf7a0454SJed Brown steps, (double)final_time)); 20277841947SLeila Ghaffari } 20377841947SLeila Ghaffari 20477841947SLeila Ghaffari // Output numerical values from command line 2052b730f8bSJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 20677841947SLeila Ghaffari 20777841947SLeila Ghaffari // Compare reference solution values with current test run for CI 2088fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 2092b730f8bSJeremy L Thompson PetscCall(RegressionTests_NS(user->app_ctx, Q)); 21077841947SLeila Ghaffari } 211ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21277841947SLeila Ghaffari } 21377841947SLeila Ghaffari 2140de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 2150de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 2160de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2174de8550aSJed Brown 21877841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation 21977841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 22077841947SLeila Ghaffari PetscViewer viewer; 2212b730f8bSJeremy L Thompson 22277841947SLeila Ghaffari PetscFunctionBegin; 22377841947SLeila Ghaffari 2242b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 225530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2262b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 22777841947SLeila Ghaffari 228ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 22977841947SLeila Ghaffari } 23077841947SLeila Ghaffari 23177841947SLeila Ghaffari // Record boundary values from initial condition 23277841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) { 2335571c6fdSJames Wright Vec Qbc, boundary_mask; 23477841947SLeila Ghaffari PetscFunctionBegin; 23577841947SLeila Ghaffari 2362b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 2372b730f8bSJeremy L Thompson PetscCall(VecCopy(Q_loc, Qbc)); 2382b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(Q_loc)); 2392b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 2402b730f8bSJeremy L Thompson PetscCall(VecAXPY(Qbc, -1., Q_loc)); 2412b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 2422b730f8bSJeremy L Thompson PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS)); 24377841947SLeila Ghaffari 2445571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2455571c6fdSJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 2465571c6fdSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 2475571c6fdSJames Wright PetscCall(VecSet(Q, 1.0)); 2485571c6fdSJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 2495571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 2505571c6fdSJames Wright 251ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 25277841947SLeila Ghaffari } 253841e4c73SJed Brown 254841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 255841e4c73SJed Brown int FreeContextPetsc(void *data) { 2562b730f8bSJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 257841e4c73SJed Brown return CEED_ERROR_SUCCESS; 258841e4c73SJed Brown } 259ef080ff9SJames Wright 260ef080ff9SJames Wright // Return mass qfunction specification for number of components N 261ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 262ef080ff9SJames Wright PetscFunctionBeginUser; 263ef080ff9SJames Wright 264ef080ff9SJames Wright switch (N) { 265ef080ff9SJames Wright case 1: 266a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 267ef080ff9SJames Wright break; 268ef080ff9SJames Wright case 5: 269a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 270ef080ff9SJames Wright break; 271052409adSJames Wright case 7: 272a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 273052409adSJames Wright break; 274ef080ff9SJames Wright case 9: 275a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 276ef080ff9SJames Wright break; 277ef080ff9SJames Wright case 22: 278a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 279ef080ff9SJames Wright break; 280ef080ff9SJames Wright default: 28183ae4962SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 282ef080ff9SJames Wright } 283ef080ff9SJames Wright 284a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 285a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 286a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 287ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 288ef080ff9SJames Wright } 2890df12fefSJames Wright 290999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 291999ff5c7SJames Wright PetscFunctionBeginUser; 292ee4ca9cbSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 293999ff5c7SJames Wright 294999ff5c7SJames Wright PetscCall(DMDestroy(&context->dm)); 295999ff5c7SJames Wright PetscCall(KSPDestroy(&context->ksp)); 296999ff5c7SJames Wright 297999ff5c7SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 298999ff5c7SJames Wright 299999ff5c7SJames Wright PetscCall(PetscFree(context)); 300999ff5c7SJames Wright 301ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 302999ff5c7SJames Wright } 303999ff5c7SJames Wright 3048b892a05SJames Wright /* 3058b892a05SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 3068b892a05SJames Wright * 3078b892a05SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 3088b892a05SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 3098b892a05SJames Wright * 3108b892a05SJames Wright * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 3118b892a05SJames Wright * 3128b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3138b892a05SJames Wright * @param[in] path Path to the file 3148b892a05SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 3158b892a05SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 3168b892a05SJames Wright * @param[out] fp File pointer to the opened file 3178b892a05SJames Wright */ 318*cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 3198b892a05SJames Wright FILE **fp) { 320457e73b2SJames Wright int ndims; 3218b892a05SJames Wright char line[char_array_len]; 3228b892a05SJames Wright char **array; 3238b892a05SJames Wright 3248b892a05SJames Wright PetscFunctionBeginUser; 3258b892a05SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 3268b892a05SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 3278b892a05SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 328457e73b2SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 3298b892a05SJames Wright 3308b892a05SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 3318b892a05SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 3328b892a05SJames Wright 333ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3348b892a05SJames Wright } 3358b892a05SJames Wright 3368b892a05SJames Wright /* 3378b892a05SJames Wright * @brief Get the number of rows for the PHASTA file at path. 3388b892a05SJames Wright * 3398b892a05SJames Wright * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 3408b892a05SJames Wright * 3418b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3428b892a05SJames Wright * @param[in] path Path to the file 3438b892a05SJames Wright * @param[out] nrows Number of rows 3448b892a05SJames Wright */ 345*cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 3468b892a05SJames Wright const PetscInt char_array_len = 512; 3478b892a05SJames Wright PetscInt dims[2]; 3488b892a05SJames Wright FILE *fp; 3498b892a05SJames Wright 3508b892a05SJames Wright PetscFunctionBeginUser; 351*cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3528b892a05SJames Wright *nrows = dims[0]; 3538b892a05SJames Wright PetscCall(PetscFClose(comm, fp)); 3548b892a05SJames Wright 355ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3568b892a05SJames Wright } 3574cc9442bSJames Wright 358*cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 359457e73b2SJames Wright PetscInt dims[2]; 360457e73b2SJames Wright int ndims; 3614cc9442bSJames Wright FILE *fp; 3624cc9442bSJames Wright const PetscInt char_array_len = 512; 3634cc9442bSJames Wright char line[char_array_len]; 3644cc9442bSJames Wright char **row_array; 3654cc9442bSJames Wright PetscFunctionBeginUser; 3664cc9442bSJames Wright 367*cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3684cc9442bSJames Wright 3694cc9442bSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 3704cc9442bSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 3714cc9442bSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3720e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 373457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 3744cc9442bSJames Wright 3754cc9442bSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 3764cc9442bSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 3774cc9442bSJames Wright } 3784cc9442bSJames Wright } 3794cc9442bSJames Wright 3804cc9442bSJames Wright PetscCall(PetscFClose(comm, fp)); 3814cc9442bSJames Wright 382ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3834cc9442bSJames Wright } 38475d1798cSJames Wright 38575d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 38675d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemble; 38775d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 38875d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 38975d1798cSJames Wright static PetscClassId libCEED_classid; 39075d1798cSJames Wright 39175d1798cSJames Wright PetscErrorCode RegisterLogEvents() { 39275d1798cSJames Wright PetscFunctionBeginUser; 39375d1798cSJames Wright PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid)); 39475d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply)); 39575d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble)); 39675d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal)); 39775d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal)); 398ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 39975d1798cSJames Wright } 400457e73b2SJames Wright 401457e73b2SJames Wright /** 402457e73b2SJames Wright @brief Translate array of CeedInt to PetscInt. 403457e73b2SJames Wright If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`. 404457e73b2SJames Wright Caller is responsible for freeing `array_petsc` with `free()`. 405457e73b2SJames Wright 406457e73b2SJames Wright @param[in] num_entries Number of array entries 407457e73b2SJames Wright @param[in,out] array_ceed Array of CeedInts 408457e73b2SJames Wright @param[out] array_petsc Array of PetscInts 409457e73b2SJames Wright **/ 410457e73b2SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) { 411457e73b2SJames Wright CeedInt int_c = 0; 412457e73b2SJames Wright PetscInt int_p = 0; 413457e73b2SJames Wright 414457e73b2SJames Wright PetscFunctionBeginUser; 415457e73b2SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 416457e73b2SJames Wright *array_petsc = (PetscInt *)*array_ceed; 417457e73b2SJames Wright } else { 418457e73b2SJames Wright *array_petsc = malloc(num_entries * sizeof(PetscInt)); 419457e73b2SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i]; 420457e73b2SJames Wright free(*array_ceed); 421457e73b2SJames Wright } 422457e73b2SJames Wright *array_ceed = NULL; 423457e73b2SJames Wright 424457e73b2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 425457e73b2SJames Wright } 426457e73b2SJames Wright 427457e73b2SJames Wright /** 428457e73b2SJames Wright @brief Translate array of PetscInt to CeedInt. 429457e73b2SJames Wright If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`. 430457e73b2SJames Wright Caller is responsible for freeing `array_ceed` with `PetscFree()`. 431457e73b2SJames Wright 432457e73b2SJames Wright @param[in] num_entries Number of array entries 433457e73b2SJames Wright @param[in,out] array_petsc Array of PetscInts 434457e73b2SJames Wright @param[out] array_ceed Array of CeedInts 435457e73b2SJames Wright **/ 436457e73b2SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) { 437457e73b2SJames Wright CeedInt int_c = 0; 438457e73b2SJames Wright PetscInt int_p = 0; 439457e73b2SJames Wright 440457e73b2SJames Wright PetscFunctionBeginUser; 441457e73b2SJames Wright if (sizeof(int_c) == sizeof(int_p)) { 442457e73b2SJames Wright *array_ceed = (CeedInt *)*array_petsc; 443457e73b2SJames Wright } else { 444457e73b2SJames Wright PetscCall(PetscMalloc1(num_entries, array_ceed)); 445457e73b2SJames Wright for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i]; 446457e73b2SJames Wright PetscCall(PetscFree(*array_petsc)); 447457e73b2SJames Wright } 448457e73b2SJames Wright *array_petsc = NULL; 449457e73b2SJames Wright 450457e73b2SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 451457e73b2SJames Wright } 4522e31f396SJames Wright 4532e31f396SJames Wright // Print information about the given simulation run 4542e31f396SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) { 455a424bcd0SJames Wright Ceed ceed = user->ceed; 4562e31f396SJames Wright PetscFunctionBeginUser; 4572e31f396SJames Wright // Header and rank 4582e31f396SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 4592e31f396SJames Wright PetscMPIInt rank, comm_size; 4602e31f396SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 4612e31f396SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 4622e31f396SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 4632e31f396SJames Wright PetscCall(PetscPrintf(comm, 4642e31f396SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 4652e31f396SJames Wright " MPI:\n" 4662e31f396SJames Wright " Host Name : %s\n" 4672e31f396SJames Wright " Total ranks : %d\n", 4682e31f396SJames Wright host_name, comm_size)); 4692e31f396SJames Wright 4702e31f396SJames Wright // Problem specific info 471367c849eSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx)); 4722e31f396SJames Wright 4732e31f396SJames Wright // libCEED 4742e31f396SJames Wright const char *used_resource; 4752e31f396SJames Wright CeedMemType mem_type_backend; 476a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 477a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 4782e31f396SJames Wright PetscCall(PetscPrintf(comm, 4792e31f396SJames Wright " libCEED:\n" 4802e31f396SJames Wright " libCEED Backend : %s\n" 4812e31f396SJames Wright " libCEED Backend MemType : %s\n", 4822e31f396SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 4832e31f396SJames Wright // PETSc 4842e31f396SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 4852e31f396SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0'; 4862e31f396SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 4872e31f396SJames Wright MatType mat_type; 4882e31f396SJames Wright VecType vec_type; 4892e31f396SJames Wright PetscCall(DMGetMatType(user->dm, &mat_type)); 4902e31f396SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type)); 4912e31f396SJames Wright PetscCall(PetscPrintf(comm, 4922e31f396SJames Wright " PETSc:\n" 4932e31f396SJames Wright " Box Faces : %s\n" 4942e31f396SJames Wright " DM MatType : %s\n" 4952e31f396SJames Wright " DM VecType : %s\n" 4962e31f396SJames Wright " Time Stepping Scheme : %s\n", 4972e31f396SJames Wright box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 4982e31f396SJames Wright if (user->app_ctx->cont_steps) { 4992e31f396SJames Wright PetscCall(PetscPrintf(comm, 5002e31f396SJames Wright " Continue:\n" 5012e31f396SJames Wright " Filename: : %s\n" 5022e31f396SJames Wright " Step: : %" PetscInt_FMT "\n" 5032e31f396SJames Wright " Time: : %g\n", 5042e31f396SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 5052e31f396SJames Wright } 5062e31f396SJames Wright // Mesh 5072e31f396SJames Wright const PetscInt num_comp_q = 5; 5082e31f396SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 5092e31f396SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 5102e31f396SJames Wright // -- Get global size 5112e31f396SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 5122e31f396SJames Wright // -- Get local size 5132e31f396SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 5142e31f396SJames Wright PetscCall(PetscPrintf(comm, 5152e31f396SJames Wright " Mesh:\n" 5162e31f396SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 5172e31f396SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 5182e31f396SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 5192e31f396SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 520ce11f295SJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 521ce11f295SJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 5222e31f396SJames Wright // -- Get Partition Statistics 5232e31f396SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 5242e31f396SJames Wright { 5252e31f396SJames Wright PetscInt *gather_buffer = NULL; 526ce11f295SJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 5272e31f396SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 5282e31f396SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 5292e31f396SJames Wright 530ce11f295SJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 5312e31f396SJames Wright if (!rank) { 5322e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 5332e31f396SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 5342e31f396SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 5352e31f396SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 5362e31f396SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 537ce11f295SJames Wright PetscCall(PetscPrintf( 538ce11f295SJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 5392e31f396SJames 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)); 5402e31f396SJames Wright } 5412e31f396SJames Wright 542ce11f295SJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 5432e31f396SJames Wright if (!rank) { 5442e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 5452e31f396SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 5462e31f396SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 5472e31f396SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 5482e31f396SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 549ce11f295SJames Wright PetscCall(PetscPrintf( 550ce11f295SJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 5512e31f396SJames 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)); 5522e31f396SJames Wright } 5532e31f396SJames Wright 55465ba01baSJames Wright if (comm_size != 1) { 555ce11f295SJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 5562e31f396SJames Wright { 5572e31f396SJames Wright PetscSF sf; 558ce11f295SJames Wright PetscInt nrranks, niranks; 559ce11f295SJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 560ce11f295SJames Wright const PetscMPIInt *rranks, *iranks; 5612e31f396SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf)); 5622e31f396SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 563ce11f295SJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 5642e31f396SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 5652e31f396SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 5662e31f396SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 567ce11f295SJames Wright num_ghost_interface_ranks++; 568ce11f295SJames Wright } 569ce11f295SJames Wright for (PetscInt i = 0; i < niranks; i++) { 570ce11f295SJames Wright if (iranks[i] == rank) continue; 571ce11f295SJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 572ce11f295SJames Wright num_owned_interface_ranks++; 5732e31f396SJames Wright } 5742e31f396SJames Wright } 575ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 5762e31f396SJames Wright if (!rank) { 5772e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 578ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 579ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 580ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 581ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 582ce11f295SJames Wright PetscCall(PetscPrintf( 58365ba01baSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 58465ba01baSJames 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, 58565ba01baSJames Wright part_shared_dof_ratio)); 586ce11f295SJames Wright } 587ce11f295SJames Wright 588ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 589ce11f295SJames Wright if (!rank) { 590ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 591ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 592ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 593ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 594ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 595ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 596ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 597ce11f295SJames Wright } 598ce11f295SJames Wright 599ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 600ce11f295SJames Wright if (!rank) { 601ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 602ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 603ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 604ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 605ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 606ce11f295SJames Wright PetscCall(PetscPrintf( 60765ba01baSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 60865ba01baSJames 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, 60965ba01baSJames Wright part_shared_dof_ratio)); 610ce11f295SJames Wright } 611ce11f295SJames Wright 612ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 613ce11f295SJames Wright if (!rank) { 614ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 615ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 616ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 617ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 618ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 619ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 620ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 6212e31f396SJames Wright } 62265ba01baSJames Wright } 6232e31f396SJames Wright 6242e31f396SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 6252e31f396SJames Wright } 6262e31f396SJames Wright 6272e31f396SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6282e31f396SJames Wright } 629