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; 21f6ce2b0aSJames Wright CeedVector mult_vec; 22f6ce2b0aSJames Wright PetscMemType m_mem_type; 23f6ce2b0aSJames Wright Vec Multiplicity, Multiplicity_loc; 24f6ce2b0aSJames Wright 2577841947SLeila Ghaffari PetscFunctionBeginUser; 26a424bcd0SJames Wright if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time)); 275263e9c6SJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx)); 2877841947SLeila Ghaffari 29a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL)); 3077841947SLeila Ghaffari 3177841947SLeila Ghaffari // -- Get multiplicity 32f6ce2b0aSJames Wright PetscCall(DMGetLocalVector(dm, &Multiplicity_loc)); 33d0593705SJames Wright PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec)); 34a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec)); 35d0593705SJames Wright PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc)); 3677841947SLeila Ghaffari 37f6ce2b0aSJames Wright PetscCall(DMGetGlobalVector(dm, &Multiplicity)); 38f6ce2b0aSJames Wright PetscCall(VecZeroEntries(Multiplicity)); 39f6ce2b0aSJames Wright PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity)); 4077841947SLeila Ghaffari 4177841947SLeila Ghaffari // -- Fix multiplicity 42f6ce2b0aSJames Wright PetscCall(VecPointwiseDivide(Q, Q, Multiplicity)); 43f6ce2b0aSJames Wright PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc)); 4477841947SLeila Ghaffari 45f6ce2b0aSJames Wright PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc)); 46f6ce2b0aSJames Wright PetscCall(DMRestoreGlobalVector(dm, &Multiplicity)); 47a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec)); 48ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4977841947SLeila Ghaffari } 5077841947SLeila Ghaffari 513c9e7ad1SJames Wright // Record boundary values from initial condition 523c9e7ad1SJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) { 533c9e7ad1SJames Wright Vec Qbc, boundary_mask; 543c9e7ad1SJames Wright 553c9e7ad1SJames Wright PetscFunctionBeginUser; 563c9e7ad1SJames Wright PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 573c9e7ad1SJames Wright PetscCall(VecCopy(Q_loc, Qbc)); 583c9e7ad1SJames Wright PetscCall(VecZeroEntries(Q_loc)); 593c9e7ad1SJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc)); 603c9e7ad1SJames Wright PetscCall(VecAXPY(Qbc, -1., Q_loc)); 613c9e7ad1SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 623c9e7ad1SJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs)); 633c9e7ad1SJames Wright 643c9e7ad1SJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 653c9e7ad1SJames Wright PetscCall(DMGetGlobalVector(dm, &Q)); 663c9e7ad1SJames Wright PetscCall(VecZeroEntries(boundary_mask)); 673c9e7ad1SJames Wright PetscCall(VecSet(Q, 1.0)); 683c9e7ad1SJames Wright PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask)); 693c9e7ad1SJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 703c9e7ad1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 713c9e7ad1SJames Wright } 723c9e7ad1SJames Wright 733c9e7ad1SJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 742b730f8bSJeremy L Thompson Vec grad_FVM) { 755571c6fdSJames Wright Vec Qbc, boundary_mask; 7677841947SLeila Ghaffari 77f17d818dSJames Wright PetscFunctionBeginUser; 783722cd23SJames Wright // Mask (zero) Strong BC entries 795571c6fdSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 805571c6fdSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 815571c6fdSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 825571c6fdSJames Wright 832b730f8bSJeremy L Thompson PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc)); 842b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q_loc, 1., Qbc)); 852b730f8bSJeremy L Thompson PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc)); 86ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 8777841947SLeila Ghaffari } 8877841947SLeila Ghaffari 89*7ceb6423SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) { 90*7ceb6423SJed Brown PetscFunctionBeginUser; 91*7ceb6423SJed Brown if (file_type == PETSC_INT32) { 92*7ceb6423SJed Brown PetscInt32 val; 93*7ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32)); 94*7ceb6423SJed Brown *out = val; 95*7ceb6423SJed Brown } else if (file_type == PETSC_INT64) { 96*7ceb6423SJed Brown PetscInt64 val; 97*7ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64)); 98*7ceb6423SJed Brown *out = val; 99*7ceb6423SJed Brown } else { 100*7ceb6423SJed Brown PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT)); 101*7ceb6423SJed Brown } 102*7ceb6423SJed Brown PetscFunctionReturn(PETSC_SUCCESS); 103*7ceb6423SJed Brown } 104*7ceb6423SJed Brown 105530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number 106530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) { 1070de6a49fSJames Wright PetscInt file_step_number; 1080de6a49fSJames Wright PetscInt32 token; 109530ad8c4SKenneth E. Jansen PetscReal file_time; 110*7ceb6423SJed Brown PetscDataType file_type = PETSC_INT32; 111530ad8c4SKenneth E. Jansen 112f17d818dSJames Wright PetscFunctionBeginUser; 1130de6a49fSJames Wright PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32)); 1140de6a49fSJames Wright if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 || 1150de6a49fSJames Wright token == FLUIDS_FILE_TOKEN) { // New style format; we're reading a file with step number and time in the header 116*7ceb6423SJed Brown if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32; 117*7ceb6423SJed Brown else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64; 118*7ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type)); 119530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL)); 120530ad8c4SKenneth E. Jansen if (time) *time = file_time; 121530ad8c4SKenneth E. Jansen if (step_number) *step_number = file_step_number; 122530ad8c4SKenneth E. Jansen } else if (token == VEC_FILE_CLASSID) { // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ] 123530ad8c4SKenneth E. Jansen PetscInt length, N; 124*7ceb6423SJed Brown PetscCall(BinaryReadIntoInt(viewer, &length, file_type)); 125530ad8c4SKenneth E. Jansen PetscCall(VecGetSize(Q, &N)); 126530ad8c4SKenneth 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); 127530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE)); 128530ad8c4SKenneth E. Jansen } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file"); 129530ad8c4SKenneth E. Jansen 130530ad8c4SKenneth E. Jansen PetscCall(VecLoad(Q, viewer)); 131ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 132530ad8c4SKenneth E. Jansen } 133530ad8c4SKenneth E. Jansen 13477841947SLeila Ghaffari // Compare reference solution values with current test run for CI 1353c9e7ad1SJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) { 13677841947SLeila Ghaffari Vec Qref; 13777841947SLeila Ghaffari PetscViewer viewer; 13877841947SLeila Ghaffari PetscReal error, Qrefnorm; 139530ad8c4SKenneth E. Jansen MPI_Comm comm = PetscObjectComm((PetscObject)Q); 14077841947SLeila Ghaffari 141f17d818dSJames Wright PetscFunctionBeginUser; 14277841947SLeila Ghaffari // Read reference file 1432b730f8bSJeremy L Thompson PetscCall(VecDuplicate(Q, &Qref)); 144530ad8c4SKenneth E. Jansen PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer)); 145530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL)); 14677841947SLeila Ghaffari 14777841947SLeila Ghaffari // Compute error with respect to reference solution 1482b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Qref)); 1492b730f8bSJeremy L Thompson PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm)); 1502b730f8bSJeremy L Thompson PetscCall(VecScale(Q, 1. / Qrefnorm)); 1512b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_MAX, &error)); 15277841947SLeila Ghaffari 15377841947SLeila Ghaffari // Check error 15477841947SLeila Ghaffari if (error > app_ctx->test_tol) { 1552b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error)); 15677841947SLeila Ghaffari } 15777841947SLeila Ghaffari 15877841947SLeila Ghaffari // Cleanup 1592b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 1602b730f8bSJeremy L Thompson PetscCall(VecDestroy(&Qref)); 161ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 16277841947SLeila Ghaffari } 16377841947SLeila Ghaffari 16477841947SLeila Ghaffari // Get error for problems with exact solutions 1653c9e7ad1SJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) { 16677841947SLeila Ghaffari PetscInt loc_nodes; 16777841947SLeila Ghaffari Vec Q_exact, Q_exact_loc; 16877841947SLeila Ghaffari PetscReal rel_error, norm_error, norm_exact; 16977841947SLeila Ghaffari 170f17d818dSJames Wright PetscFunctionBeginUser; 17177841947SLeila Ghaffari // Get exact solution at final time 172f6ce2b0aSJames Wright PetscCall(DMGetGlobalVector(dm, &Q_exact)); 1732b730f8bSJeremy L Thompson PetscCall(DMGetLocalVector(dm, &Q_exact_loc)); 1742b730f8bSJeremy L Thompson PetscCall(VecGetSize(Q_exact_loc, &loc_nodes)); 1752b730f8bSJeremy L Thompson PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time)); 17677841947SLeila Ghaffari 17777841947SLeila Ghaffari // Get |exact solution - obtained solution| 1782b730f8bSJeremy L Thompson PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact)); 1792b730f8bSJeremy L Thompson PetscCall(VecAXPY(Q, -1.0, Q_exact)); 1802b730f8bSJeremy L Thompson PetscCall(VecNorm(Q, NORM_1, &norm_error)); 18177841947SLeila Ghaffari 18277841947SLeila Ghaffari rel_error = norm_error / norm_exact; 1832b730f8bSJeremy L Thompson PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error)); 1842b730f8bSJeremy L Thompson PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc)); 185f6ce2b0aSJames Wright PetscCall(DMRestoreGlobalVector(dm, &Q_exact)); 186ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 18777841947SLeila Ghaffari } 18877841947SLeila Ghaffari 18977841947SLeila Ghaffari // Post-processing 1903c9e7ad1SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) { 19177841947SLeila Ghaffari PetscInt steps; 192cf7a0454SJed Brown TSConvergedReason reason; 19377841947SLeila Ghaffari 194f17d818dSJames Wright PetscFunctionBeginUser; 19577841947SLeila Ghaffari // Print relative error 1968fb33541SJames Wright if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) { 1973c9e7ad1SJames Wright PetscCall(PrintError(ceed_data, dm, user, Q, final_time)); 19877841947SLeila Ghaffari } 19977841947SLeila Ghaffari 20077841947SLeila Ghaffari // Print final time and number of steps 2012b730f8bSJeremy L Thompson PetscCall(TSGetStepNumber(ts, &steps)); 202cf7a0454SJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 2038fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 204cf7a0454SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason], 205cf7a0454SJed Brown steps, (double)final_time)); 20677841947SLeila Ghaffari } 20777841947SLeila Ghaffari 20877841947SLeila Ghaffari // Output numerical values from command line 2092b730f8bSJeremy L Thompson PetscCall(VecViewFromOptions(Q, NULL, "-vec_view")); 21077841947SLeila Ghaffari 21177841947SLeila Ghaffari // Compare reference solution values with current test run for CI 2128fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_SOLVER) { 2133c9e7ad1SJames Wright PetscCall(RegressionTest(user->app_ctx, Q)); 21477841947SLeila Ghaffari } 215ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 21677841947SLeila Ghaffari } 21777841947SLeila Ghaffari 2180de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN = 0xceedf00; // for backwards compatibility 2190de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32; 2200de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64; 2214de8550aSJed Brown 22277841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation 22377841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) { 22477841947SLeila Ghaffari PetscViewer viewer; 2252b730f8bSJeremy L Thompson 226f17d818dSJames Wright PetscFunctionBeginUser; 2272b730f8bSJeremy L Thompson PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer)); 228530ad8c4SKenneth E. Jansen PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps)); 2292b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&viewer)); 230ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 23177841947SLeila Ghaffari } 23277841947SLeila Ghaffari 233841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes 234841e4c73SJed Brown int FreeContextPetsc(void *data) { 2352b730f8bSJeremy L Thompson if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed"); 236841e4c73SJed Brown return CEED_ERROR_SUCCESS; 237841e4c73SJed Brown } 238ef080ff9SJames Wright 239ef080ff9SJames Wright // Return mass qfunction specification for number of components N 240ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) { 241ef080ff9SJames Wright PetscFunctionBeginUser; 242ef080ff9SJames Wright switch (N) { 243ef080ff9SJames Wright case 1: 244a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf)); 245ef080ff9SJames Wright break; 246ef080ff9SJames Wright case 5: 247a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf)); 248ef080ff9SJames Wright break; 249052409adSJames Wright case 7: 250a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf)); 251052409adSJames Wright break; 252ef080ff9SJames Wright case 9: 253a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf)); 254ef080ff9SJames Wright break; 255ef080ff9SJames Wright case 22: 256a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf)); 257ef080ff9SJames Wright break; 258ef080ff9SJames Wright default: 25983ae4962SJames Wright SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N); 260ef080ff9SJames Wright } 261ef080ff9SJames Wright 262a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP)); 263a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE)); 264a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP)); 2657f2a9303SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N)); 266ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 267ef080ff9SJames Wright } 2680df12fefSJames Wright 269999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) { 270999ff5c7SJames Wright PetscFunctionBeginUser; 271ee4ca9cbSJames Wright if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS); 272999ff5c7SJames Wright 273999ff5c7SJames Wright PetscCall(DMDestroy(&context->dm)); 274999ff5c7SJames Wright PetscCall(KSPDestroy(&context->ksp)); 275999ff5c7SJames Wright 276999ff5c7SJames Wright PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx)); 277999ff5c7SJames Wright 278999ff5c7SJames Wright PetscCall(PetscFree(context)); 279ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 280999ff5c7SJames Wright } 281999ff5c7SJames Wright 2828b892a05SJames Wright /* 2838b892a05SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 2848b892a05SJames Wright * 2858b892a05SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 2868b892a05SJames Wright * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 2878b892a05SJames Wright * 2888b892a05SJames 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. 2898b892a05SJames Wright * 2908b892a05SJames Wright * @param[in] comm MPI_Comm for the program 2918b892a05SJames Wright * @param[in] path Path to the file 2928b892a05SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 2938b892a05SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 2948b892a05SJames Wright * @param[out] fp File pointer to the opened file 2958b892a05SJames Wright */ 296cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 2978b892a05SJames Wright FILE **fp) { 298457e73b2SJames Wright int ndims; 2998b892a05SJames Wright char line[char_array_len]; 3008b892a05SJames Wright char **array; 3018b892a05SJames Wright 3028b892a05SJames Wright PetscFunctionBeginUser; 3038b892a05SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 3048b892a05SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 3058b892a05SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 306457e73b2SJames Wright PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path); 3078b892a05SJames Wright 3088b892a05SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 3098b892a05SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 310ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3118b892a05SJames Wright } 3128b892a05SJames Wright 3138b892a05SJames Wright /* 3148b892a05SJames Wright * @brief Get the number of rows for the PHASTA file at path. 3158b892a05SJames Wright * 3168b892a05SJames 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. 3178b892a05SJames Wright * 3188b892a05SJames Wright * @param[in] comm MPI_Comm for the program 3198b892a05SJames Wright * @param[in] path Path to the file 3208b892a05SJames Wright * @param[out] nrows Number of rows 3218b892a05SJames Wright */ 322cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 3238b892a05SJames Wright const PetscInt char_array_len = 512; 3248b892a05SJames Wright PetscInt dims[2]; 3258b892a05SJames Wright FILE *fp; 3268b892a05SJames Wright 3278b892a05SJames Wright PetscFunctionBeginUser; 328cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3298b892a05SJames Wright *nrows = dims[0]; 3308b892a05SJames Wright PetscCall(PetscFClose(comm, fp)); 331ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3328b892a05SJames Wright } 3334cc9442bSJames Wright 334cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) { 335457e73b2SJames Wright PetscInt dims[2]; 336457e73b2SJames Wright int ndims; 3374cc9442bSJames Wright FILE *fp; 3384cc9442bSJames Wright const PetscInt char_array_len = 512; 3394cc9442bSJames Wright char line[char_array_len]; 3404cc9442bSJames Wright char **row_array; 3414cc9442bSJames Wright 342f17d818dSJames Wright PetscFunctionBeginUser; 343cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 3444cc9442bSJames Wright 3454cc9442bSJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 3464cc9442bSJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 3474cc9442bSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array)); 3480e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 349457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 3504cc9442bSJames Wright 3514cc9442bSJames Wright for (PetscInt j = 0; j < dims[1]; j++) { 3524cc9442bSJames Wright array[i * dims[1] + j] = (PetscReal)atof(row_array[j]); 3534cc9442bSJames Wright } 3544cc9442bSJames Wright } 3554cc9442bSJames Wright 3564cc9442bSJames Wright PetscCall(PetscFClose(comm, fp)); 357ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3584cc9442bSJames Wright } 35975d1798cSJames Wright 36075d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorApply; 36175d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemble; 36275d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal; 36375d1798cSJames Wright PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal; 364844b3d4aSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Init; 365844b3d4aSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Meta; 366844b3d4aSRiccardo Balin PetscLogEvent FLUIDS_SmartRedis_Train; 367844b3d4aSRiccardo Balin PetscLogEvent FLUIDS_TrainDataCompute; 3683737f832SJames Wright PetscLogEvent FLUIDS_DifferentialFilter; 3693737f832SJames Wright PetscLogEvent FLUIDS_VelocityGradientProjection; 3703737f832SJames Wright static PetscClassId libCEED_classid, onlineTrain_classid, misc_classid; 37175d1798cSJames Wright 37275d1798cSJames Wright PetscErrorCode RegisterLogEvents() { 37375d1798cSJames Wright PetscFunctionBeginUser; 37475d1798cSJames Wright PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid)); 37575d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply)); 37675d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble)); 37775d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal)); 37875d1798cSJames Wright PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal)); 379844b3d4aSRiccardo Balin 380844b3d4aSRiccardo Balin PetscCall(PetscClassIdRegister("onlineTrain", &onlineTrain_classid)); 381844b3d4aSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Init", onlineTrain_classid, &FLUIDS_SmartRedis_Init)); 382844b3d4aSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Meta", onlineTrain_classid, &FLUIDS_SmartRedis_Meta)); 383844b3d4aSRiccardo Balin PetscCall(PetscLogEventRegister("SmartRedis_Train", onlineTrain_classid, &FLUIDS_SmartRedis_Train)); 384844b3d4aSRiccardo Balin PetscCall(PetscLogEventRegister("TrainDataCompute", onlineTrain_classid, &FLUIDS_TrainDataCompute)); 3853737f832SJames Wright 3863737f832SJames Wright PetscCall(PetscClassIdRegister("Miscellaneous", &misc_classid)); 3873737f832SJames Wright PetscCall(PetscLogEventRegister("DiffFilter", misc_classid, &FLUIDS_DifferentialFilter)); 3883737f832SJames Wright PetscCall(PetscLogEventRegister("VeloGradProj", misc_classid, &FLUIDS_VelocityGradientProjection)); 389ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 39075d1798cSJames Wright } 391457e73b2SJames Wright 3922e31f396SJames Wright // Print information about the given simulation run 3932e31f396SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) { 394a424bcd0SJames Wright Ceed ceed = user->ceed; 3952e31f396SJames Wright PetscFunctionBeginUser; 3962e31f396SJames Wright // Header and rank 3972e31f396SJames Wright char host_name[PETSC_MAX_PATH_LEN]; 3982e31f396SJames Wright PetscMPIInt rank, comm_size; 3992e31f396SJames Wright PetscCall(PetscGetHostName(host_name, sizeof host_name)); 4002e31f396SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 4012e31f396SJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size)); 4022e31f396SJames Wright PetscCall(PetscPrintf(comm, 4032e31f396SJames Wright "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 4042e31f396SJames Wright " MPI:\n" 4052e31f396SJames Wright " Host Name : %s\n" 4062e31f396SJames Wright " Total ranks : %d\n", 4072e31f396SJames Wright host_name, comm_size)); 4082e31f396SJames Wright 4092e31f396SJames Wright // Problem specific info 410367c849eSJames Wright PetscCall(problem->print_info(user, problem, user->app_ctx)); 4112e31f396SJames Wright 4122e31f396SJames Wright // libCEED 4132e31f396SJames Wright const char *used_resource; 4142e31f396SJames Wright CeedMemType mem_type_backend; 415a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource)); 416a424bcd0SJames Wright PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend)); 4172e31f396SJames Wright PetscCall(PetscPrintf(comm, 4182e31f396SJames Wright " libCEED:\n" 4192e31f396SJames Wright " libCEED Backend : %s\n" 4202e31f396SJames Wright " libCEED Backend MemType : %s\n", 4212e31f396SJames Wright used_resource, CeedMemTypes[mem_type_backend])); 4222e31f396SJames Wright // PETSc 4232e31f396SJames Wright char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 4242e31f396SJames Wright if (problem->dim == 2) box_faces_str[3] = '\0'; 4252e31f396SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 426a04f48a3SJames Wright MatType amat_type = user->app_ctx->amat_type, pmat_type; 4272e31f396SJames Wright VecType vec_type; 428a04f48a3SJames Wright PetscCall(DMGetMatType(user->dm, &pmat_type)); 429a04f48a3SJames Wright if (!amat_type) amat_type = pmat_type; 4302e31f396SJames Wright PetscCall(DMGetVecType(user->dm, &vec_type)); 4312e31f396SJames Wright PetscCall(PetscPrintf(comm, 4322e31f396SJames Wright " PETSc:\n" 4332e31f396SJames Wright " Box Faces : %s\n" 434a04f48a3SJames Wright " A MatType : %s\n" 435a04f48a3SJames Wright " P MatType : %s\n" 4362e31f396SJames Wright " DM VecType : %s\n" 4372e31f396SJames Wright " Time Stepping Scheme : %s\n", 438a04f48a3SJames Wright box_faces_str, amat_type, pmat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 4392e31f396SJames Wright if (user->app_ctx->cont_steps) { 4402e31f396SJames Wright PetscCall(PetscPrintf(comm, 4412e31f396SJames Wright " Continue:\n" 4422e31f396SJames Wright " Filename: : %s\n" 4432e31f396SJames Wright " Step: : %" PetscInt_FMT "\n" 4442e31f396SJames Wright " Time: : %g\n", 4452e31f396SJames Wright user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time)); 4462e31f396SJames Wright } 4472e31f396SJames Wright // Mesh 4482e31f396SJames Wright const PetscInt num_comp_q = 5; 4492e31f396SJames Wright PetscInt glob_dofs, owned_dofs, local_dofs; 4502e31f396SJames Wright const CeedInt num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra; 4512e31f396SJames Wright PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL)); 4522e31f396SJames Wright PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL)); 4532e31f396SJames Wright PetscCall(PetscPrintf(comm, 4542e31f396SJames Wright " Mesh:\n" 4552e31f396SJames Wright " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 4562e31f396SJames Wright " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 4572e31f396SJames Wright " Global DoFs : %" PetscInt_FMT "\n" 4582e31f396SJames Wright " DoFs per node : %" PetscInt_FMT "\n" 459ce11f295SJames Wright " Global %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT "\n", 460ce11f295SJames Wright num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q)); 4612e31f396SJames Wright // -- Get Partition Statistics 4622e31f396SJames Wright PetscCall(PetscPrintf(comm, " Partition: (min,max,median,max/median)\n")); 4632e31f396SJames Wright { 4642e31f396SJames Wright PetscInt *gather_buffer = NULL; 465ce11f295SJames Wright PetscInt part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3]; 4662e31f396SJames Wright PetscInt median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1; 4672e31f396SJames Wright if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer)); 4682e31f396SJames Wright 469ce11f295SJames Wright PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 4702e31f396SJames Wright if (!rank) { 4712e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 4722e31f396SJames Wright part_owned_dofs[0] = gather_buffer[0]; // min 4732e31f396SJames Wright part_owned_dofs[1] = gather_buffer[comm_size - 1]; // max 4742e31f396SJames Wright part_owned_dofs[2] = gather_buffer[median_index]; // median 4752e31f396SJames Wright PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2]; 476ce11f295SJames Wright PetscCall(PetscPrintf( 477ce11f295SJames Wright comm, " Global Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 4782e31f396SJames 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)); 4792e31f396SJames Wright } 4802e31f396SJames Wright 481ce11f295SJames Wright PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 4822e31f396SJames Wright if (!rank) { 4832e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 4842e31f396SJames Wright part_local_dofs[0] = gather_buffer[0]; // min 4852e31f396SJames Wright part_local_dofs[1] = gather_buffer[comm_size - 1]; // max 4862e31f396SJames Wright part_local_dofs[2] = gather_buffer[median_index]; // median 4872e31f396SJames Wright PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2]; 488ce11f295SJames Wright PetscCall(PetscPrintf( 489ce11f295SJames Wright comm, " Local Vector %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q, 4902e31f396SJames 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)); 4912e31f396SJames Wright } 4922e31f396SJames Wright 49365ba01baSJames Wright if (comm_size != 1) { 494ce11f295SJames Wright PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0; 4952e31f396SJames Wright { 4962e31f396SJames Wright PetscSF sf; 497ce11f295SJames Wright PetscInt nrranks, niranks; 498ce11f295SJames Wright const PetscInt *roffset, *rmine, *rremote, *ioffset, *irootloc; 499ce11f295SJames Wright const PetscMPIInt *rranks, *iranks; 5002e31f396SJames Wright PetscCall(DMGetSectionSF(user->dm, &sf)); 5012e31f396SJames Wright PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote)); 502ce11f295SJames Wright PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc)); 5032e31f396SJames Wright for (PetscInt i = 0; i < nrranks; i++) { 5042e31f396SJames Wright if (rranks[i] == rank) continue; // Ignore same-part global->local transfers 5052e31f396SJames Wright num_remote_roots_total += roffset[i + 1] - roffset[i]; 506ce11f295SJames Wright num_ghost_interface_ranks++; 507ce11f295SJames Wright } 508ce11f295SJames Wright for (PetscInt i = 0; i < niranks; i++) { 509ce11f295SJames Wright if (iranks[i] == rank) continue; 510ce11f295SJames Wright num_remote_leaves_total += ioffset[i + 1] - ioffset[i]; 511ce11f295SJames Wright num_owned_interface_ranks++; 5122e31f396SJames Wright } 5132e31f396SJames Wright } 514ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 5152e31f396SJames Wright if (!rank) { 5162e31f396SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 517ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 518ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 519ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 520ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 521ce11f295SJames Wright PetscCall(PetscPrintf( 52265ba01baSJames Wright comm, " Ghost Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 52365ba01baSJames 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, 52465ba01baSJames Wright part_shared_dof_ratio)); 525ce11f295SJames Wright } 526ce11f295SJames Wright 527ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 528ce11f295SJames Wright if (!rank) { 529ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 530ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 531ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 532ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 533ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 534ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Ghost Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 535ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 536ce11f295SJames Wright } 537ce11f295SJames Wright 538ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 539ce11f295SJames Wright if (!rank) { 540ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 541ce11f295SJames Wright part_boundary_dofs[0] = gather_buffer[0]; // min 542ce11f295SJames Wright part_boundary_dofs[1] = gather_buffer[comm_size - 1]; // max 543ce11f295SJames Wright part_boundary_dofs[2] = gather_buffer[median_index]; // median 544ce11f295SJames Wright PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2]; 545ce11f295SJames Wright PetscCall(PetscPrintf( 54665ba01baSJames Wright comm, " Owned Interface %" PetscInt_FMT "-DoF nodes : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 54765ba01baSJames 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, 54865ba01baSJames Wright part_shared_dof_ratio)); 549ce11f295SJames Wright } 550ce11f295SJames Wright 551ce11f295SJames Wright PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm)); 552ce11f295SJames Wright if (!rank) { 553ce11f295SJames Wright PetscCall(PetscSortInt(comm_size, gather_buffer)); 554ce11f295SJames Wright part_neighbors[0] = gather_buffer[0]; // min 555ce11f295SJames Wright part_neighbors[1] = gather_buffer[comm_size - 1]; // max 556ce11f295SJames Wright part_neighbors[2] = gather_buffer[median_index]; // median 557ce11f295SJames Wright PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2]; 558ce11f295SJames Wright PetscCall(PetscPrintf(comm, " Owned Interface Ranks : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", 559ce11f295SJames Wright part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio)); 5602e31f396SJames Wright } 56165ba01baSJames Wright } 5622e31f396SJames Wright 5632e31f396SJames Wright if (!rank) PetscCall(PetscFree(gather_buffer)); 5642e31f396SJames Wright } 5652e31f396SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5662e31f396SJames Wright } 567