xref: /libCEED/examples/fluids/src/misc.c (revision d0593705e733b5bdd5e4c173fe0008b11db2ed29)
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));
33*d0593705SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
34a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
35*d0593705SJames 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 
89530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
90530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
910de6a49fSJames Wright   PetscInt   file_step_number;
920de6a49fSJames Wright   PetscInt32 token;
93530ad8c4SKenneth E. Jansen   PetscReal  file_time;
94530ad8c4SKenneth E. Jansen 
95f17d818dSJames Wright   PetscFunctionBeginUser;
960de6a49fSJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
970de6a49fSJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
980de6a49fSJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
99530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
100530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
101530ad8c4SKenneth E. Jansen     if (time) *time = file_time;
102530ad8c4SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
103530ad8c4SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
104530ad8c4SKenneth E. Jansen     PetscInt length, N;
105530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
106530ad8c4SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
107530ad8c4SKenneth 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);
108530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
109530ad8c4SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
110530ad8c4SKenneth E. Jansen 
111530ad8c4SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
112ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
113530ad8c4SKenneth E. Jansen }
114530ad8c4SKenneth E. Jansen 
11577841947SLeila Ghaffari // Compare reference solution values with current test run for CI
1163c9e7ad1SJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
11777841947SLeila Ghaffari   Vec         Qref;
11877841947SLeila Ghaffari   PetscViewer viewer;
11977841947SLeila Ghaffari   PetscReal   error, Qrefnorm;
120530ad8c4SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
12177841947SLeila Ghaffari 
122f17d818dSJames Wright   PetscFunctionBeginUser;
12377841947SLeila Ghaffari   // Read reference file
1242b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
125530ad8c4SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
126530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
12777841947SLeila Ghaffari 
12877841947SLeila Ghaffari   // Compute error with respect to reference solution
1292b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1302b730f8bSJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1312b730f8bSJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1322b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
13377841947SLeila Ghaffari 
13477841947SLeila Ghaffari   // Check error
13577841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
1362b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
13777841947SLeila Ghaffari   }
13877841947SLeila Ghaffari 
13977841947SLeila Ghaffari   // Cleanup
1402b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1412b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Qref));
142ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
14377841947SLeila Ghaffari }
14477841947SLeila Ghaffari 
14577841947SLeila Ghaffari // Get error for problems with exact solutions
1463c9e7ad1SJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
14777841947SLeila Ghaffari   PetscInt  loc_nodes;
14877841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
14977841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
15077841947SLeila Ghaffari 
151f17d818dSJames Wright   PetscFunctionBeginUser;
15277841947SLeila Ghaffari   // Get exact solution at final time
153f6ce2b0aSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1542b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1552b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1562b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
15777841947SLeila Ghaffari 
15877841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1592b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1602b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1612b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
16277841947SLeila Ghaffari 
16377841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
1642b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1652b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
166f6ce2b0aSJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
167ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
16877841947SLeila Ghaffari }
16977841947SLeila Ghaffari 
17077841947SLeila Ghaffari // Post-processing
1713c9e7ad1SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
17277841947SLeila Ghaffari   PetscInt          steps;
173cf7a0454SJed Brown   TSConvergedReason reason;
17477841947SLeila Ghaffari 
175f17d818dSJames Wright   PetscFunctionBeginUser;
17677841947SLeila Ghaffari   // Print relative error
1778fb33541SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1783c9e7ad1SJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
17977841947SLeila Ghaffari   }
18077841947SLeila Ghaffari 
18177841947SLeila Ghaffari   // Print final time and number of steps
1822b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
183cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1848fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
185cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
186cf7a0454SJed Brown                           steps, (double)final_time));
18777841947SLeila Ghaffari   }
18877841947SLeila Ghaffari 
18977841947SLeila Ghaffari   // Output numerical values from command line
1902b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
19177841947SLeila Ghaffari 
19277841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
1938fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
1943c9e7ad1SJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
19577841947SLeila Ghaffari   }
196ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19777841947SLeila Ghaffari }
19877841947SLeila Ghaffari 
1990de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
2000de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
2010de6a49fSJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2024de8550aSJed Brown 
20377841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
20477841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
20577841947SLeila Ghaffari   PetscViewer viewer;
2062b730f8bSJeremy L Thompson 
207f17d818dSJames Wright   PetscFunctionBeginUser;
2082b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
209530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2102b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
211ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21277841947SLeila Ghaffari }
21377841947SLeila Ghaffari 
214841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
215841e4c73SJed Brown int FreeContextPetsc(void *data) {
2162b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
217841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
218841e4c73SJed Brown }
219ef080ff9SJames Wright 
220ef080ff9SJames Wright // Return mass qfunction specification for number of components N
221ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
222ef080ff9SJames Wright   PetscFunctionBeginUser;
223ef080ff9SJames Wright   switch (N) {
224ef080ff9SJames Wright     case 1:
225a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
226ef080ff9SJames Wright       break;
227ef080ff9SJames Wright     case 5:
228a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
229ef080ff9SJames Wright       break;
230052409adSJames Wright     case 7:
231a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
232052409adSJames Wright       break;
233ef080ff9SJames Wright     case 9:
234a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
235ef080ff9SJames Wright       break;
236ef080ff9SJames Wright     case 22:
237a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
238ef080ff9SJames Wright       break;
239ef080ff9SJames Wright     default:
24083ae4962SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
241ef080ff9SJames Wright   }
242ef080ff9SJames Wright 
243a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
244a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
245a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2467f2a9303SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
247ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
248ef080ff9SJames Wright }
2490df12fefSJames Wright 
250999ff5c7SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
251999ff5c7SJames Wright   PetscFunctionBeginUser;
252ee4ca9cbSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
253999ff5c7SJames Wright 
254999ff5c7SJames Wright   PetscCall(DMDestroy(&context->dm));
255999ff5c7SJames Wright   PetscCall(KSPDestroy(&context->ksp));
256999ff5c7SJames Wright 
257999ff5c7SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
258999ff5c7SJames Wright 
259999ff5c7SJames Wright   PetscCall(PetscFree(context));
260ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
261999ff5c7SJames Wright }
262999ff5c7SJames Wright 
2638b892a05SJames Wright /*
2648b892a05SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
2658b892a05SJames Wright  *
2668b892a05SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
2678b892a05SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
2688b892a05SJames Wright  *
2698b892a05SJames 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.
2708b892a05SJames Wright  *
2718b892a05SJames Wright  * @param[in]  comm           MPI_Comm for the program
2728b892a05SJames Wright  * @param[in]  path           Path to the file
2738b892a05SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
2748b892a05SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
2758b892a05SJames Wright  * @param[out] fp File        pointer to the opened file
2768b892a05SJames Wright  */
277cbef7084SJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
2788b892a05SJames Wright                                  FILE **fp) {
279457e73b2SJames Wright   int    ndims;
2808b892a05SJames Wright   char   line[char_array_len];
2818b892a05SJames Wright   char **array;
2828b892a05SJames Wright 
2838b892a05SJames Wright   PetscFunctionBeginUser;
2848b892a05SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
2858b892a05SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
2868b892a05SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
287457e73b2SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
2888b892a05SJames Wright 
2898b892a05SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
2908b892a05SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
291ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2928b892a05SJames Wright }
2938b892a05SJames Wright 
2948b892a05SJames Wright /*
2958b892a05SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
2968b892a05SJames Wright  *
2978b892a05SJames 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.
2988b892a05SJames Wright  *
2998b892a05SJames Wright  * @param[in]  comm  MPI_Comm for the program
3008b892a05SJames Wright  * @param[in]  path  Path to the file
3018b892a05SJames Wright  * @param[out] nrows Number of rows
3028b892a05SJames Wright  */
303cbef7084SJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
3048b892a05SJames Wright   const PetscInt char_array_len = 512;
3058b892a05SJames Wright   PetscInt       dims[2];
3068b892a05SJames Wright   FILE          *fp;
3078b892a05SJames Wright 
3088b892a05SJames Wright   PetscFunctionBeginUser;
309cbef7084SJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3108b892a05SJames Wright   *nrows = dims[0];
3118b892a05SJames Wright   PetscCall(PetscFClose(comm, fp));
312ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3138b892a05SJames Wright }
3144cc9442bSJames Wright 
315cbef7084SJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
316457e73b2SJames Wright   PetscInt       dims[2];
317457e73b2SJames Wright   int            ndims;
3184cc9442bSJames Wright   FILE          *fp;
3194cc9442bSJames Wright   const PetscInt char_array_len = 512;
3204cc9442bSJames Wright   char           line[char_array_len];
3214cc9442bSJames Wright   char         **row_array;
3224cc9442bSJames Wright 
323f17d818dSJames Wright   PetscFunctionBeginUser;
324cbef7084SJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
3254cc9442bSJames Wright 
3264cc9442bSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
3274cc9442bSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
3284cc9442bSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3290e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
330457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
3314cc9442bSJames Wright 
3324cc9442bSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
3334cc9442bSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
3344cc9442bSJames Wright     }
3354cc9442bSJames Wright   }
3364cc9442bSJames Wright 
3374cc9442bSJames Wright   PetscCall(PetscFClose(comm, fp));
338ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3394cc9442bSJames Wright }
34075d1798cSJames Wright 
34175d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
34275d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
34375d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
34475d1798cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
345844b3d4aSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Init;
346844b3d4aSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Meta;
347844b3d4aSRiccardo Balin PetscLogEvent       FLUIDS_SmartRedis_Train;
348844b3d4aSRiccardo Balin PetscLogEvent       FLUIDS_TrainDataCompute;
3493737f832SJames Wright PetscLogEvent       FLUIDS_DifferentialFilter;
3503737f832SJames Wright PetscLogEvent       FLUIDS_VelocityGradientProjection;
3513737f832SJames Wright static PetscClassId libCEED_classid, onlineTrain_classid, misc_classid;
35275d1798cSJames Wright 
35375d1798cSJames Wright PetscErrorCode RegisterLogEvents() {
35475d1798cSJames Wright   PetscFunctionBeginUser;
35575d1798cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
35675d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
35775d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
35875d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
35975d1798cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
360844b3d4aSRiccardo Balin 
361844b3d4aSRiccardo Balin   PetscCall(PetscClassIdRegister("onlineTrain", &onlineTrain_classid));
362844b3d4aSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Init", onlineTrain_classid, &FLUIDS_SmartRedis_Init));
363844b3d4aSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Meta", onlineTrain_classid, &FLUIDS_SmartRedis_Meta));
364844b3d4aSRiccardo Balin   PetscCall(PetscLogEventRegister("SmartRedis_Train", onlineTrain_classid, &FLUIDS_SmartRedis_Train));
365844b3d4aSRiccardo Balin   PetscCall(PetscLogEventRegister("TrainDataCompute", onlineTrain_classid, &FLUIDS_TrainDataCompute));
3663737f832SJames Wright 
3673737f832SJames Wright   PetscCall(PetscClassIdRegister("Miscellaneous", &misc_classid));
3683737f832SJames Wright   PetscCall(PetscLogEventRegister("DiffFilter", misc_classid, &FLUIDS_DifferentialFilter));
3693737f832SJames Wright   PetscCall(PetscLogEventRegister("VeloGradProj", misc_classid, &FLUIDS_VelocityGradientProjection));
370ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
37175d1798cSJames Wright }
372457e73b2SJames Wright 
3732e31f396SJames Wright // Print information about the given simulation run
3742e31f396SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm) {
375a424bcd0SJames Wright   Ceed ceed = user->ceed;
3762e31f396SJames Wright   PetscFunctionBeginUser;
3772e31f396SJames Wright   // Header and rank
3782e31f396SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
3792e31f396SJames Wright   PetscMPIInt rank, comm_size;
3802e31f396SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
3812e31f396SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3822e31f396SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
3832e31f396SJames Wright   PetscCall(PetscPrintf(comm,
3842e31f396SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
3852e31f396SJames Wright                         "  MPI:\n"
3862e31f396SJames Wright                         "    Host Name                          : %s\n"
3872e31f396SJames Wright                         "    Total ranks                        : %d\n",
3882e31f396SJames Wright                         host_name, comm_size));
3892e31f396SJames Wright 
3902e31f396SJames Wright   // Problem specific info
391367c849eSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
3922e31f396SJames Wright 
3932e31f396SJames Wright   // libCEED
3942e31f396SJames Wright   const char *used_resource;
3952e31f396SJames Wright   CeedMemType mem_type_backend;
396a424bcd0SJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
397a424bcd0SJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
3982e31f396SJames Wright   PetscCall(PetscPrintf(comm,
3992e31f396SJames Wright                         "  libCEED:\n"
4002e31f396SJames Wright                         "    libCEED Backend                    : %s\n"
4012e31f396SJames Wright                         "    libCEED Backend MemType            : %s\n",
4022e31f396SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
4032e31f396SJames Wright   // PETSc
4042e31f396SJames Wright   char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
4052e31f396SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
4062e31f396SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
407a04f48a3SJames Wright   MatType amat_type = user->app_ctx->amat_type, pmat_type;
4082e31f396SJames Wright   VecType vec_type;
409a04f48a3SJames Wright   PetscCall(DMGetMatType(user->dm, &pmat_type));
410a04f48a3SJames Wright   if (!amat_type) amat_type = pmat_type;
4112e31f396SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
4122e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4132e31f396SJames Wright                         "  PETSc:\n"
4142e31f396SJames Wright                         "    Box Faces                          : %s\n"
415a04f48a3SJames Wright                         "    A MatType                          : %s\n"
416a04f48a3SJames Wright                         "    P MatType                          : %s\n"
4172e31f396SJames Wright                         "    DM VecType                         : %s\n"
4182e31f396SJames Wright                         "    Time Stepping Scheme               : %s\n",
419a04f48a3SJames Wright                         box_faces_str, amat_type, pmat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
4202e31f396SJames Wright   if (user->app_ctx->cont_steps) {
4212e31f396SJames Wright     PetscCall(PetscPrintf(comm,
4222e31f396SJames Wright                           "  Continue:\n"
4232e31f396SJames Wright                           "    Filename:                          : %s\n"
4242e31f396SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
4252e31f396SJames Wright                           "    Time:                              : %g\n",
4262e31f396SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
4272e31f396SJames Wright   }
4282e31f396SJames Wright   // Mesh
4292e31f396SJames Wright   const PetscInt num_comp_q = 5;
4302e31f396SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
4312e31f396SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
4322e31f396SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
4332e31f396SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
4342e31f396SJames Wright   PetscCall(PetscPrintf(comm,
4352e31f396SJames Wright                         "  Mesh:\n"
4362e31f396SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
4372e31f396SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
4382e31f396SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
4392e31f396SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
440ce11f295SJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
441ce11f295SJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
4422e31f396SJames Wright   // -- Get Partition Statistics
4432e31f396SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
4442e31f396SJames Wright   {
4452e31f396SJames Wright     PetscInt *gather_buffer = NULL;
446ce11f295SJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
4472e31f396SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
4482e31f396SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
4492e31f396SJames Wright 
450ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4512e31f396SJames Wright     if (!rank) {
4522e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
4532e31f396SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
4542e31f396SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
4552e31f396SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
4562e31f396SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
457ce11f295SJames Wright       PetscCall(PetscPrintf(
458ce11f295SJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
4592e31f396SJames 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));
4602e31f396SJames Wright     }
4612e31f396SJames Wright 
462ce11f295SJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4632e31f396SJames Wright     if (!rank) {
4642e31f396SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
4652e31f396SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
4662e31f396SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
4672e31f396SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
4682e31f396SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
469ce11f295SJames Wright       PetscCall(PetscPrintf(
470ce11f295SJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
4712e31f396SJames 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));
4722e31f396SJames Wright     }
4732e31f396SJames Wright 
47465ba01baSJames Wright     if (comm_size != 1) {
475ce11f295SJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
4762e31f396SJames Wright       {
4772e31f396SJames Wright         PetscSF            sf;
478ce11f295SJames Wright         PetscInt           nrranks, niranks;
479ce11f295SJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
480ce11f295SJames Wright         const PetscMPIInt *rranks, *iranks;
4812e31f396SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
4822e31f396SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
483ce11f295SJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
4842e31f396SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
4852e31f396SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
4862e31f396SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
487ce11f295SJames Wright           num_ghost_interface_ranks++;
488ce11f295SJames Wright         }
489ce11f295SJames Wright         for (PetscInt i = 0; i < niranks; i++) {
490ce11f295SJames Wright           if (iranks[i] == rank) continue;
491ce11f295SJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
492ce11f295SJames Wright           num_owned_interface_ranks++;
4932e31f396SJames Wright         }
4942e31f396SJames Wright       }
495ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
4962e31f396SJames Wright       if (!rank) {
4972e31f396SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
498ce11f295SJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
499ce11f295SJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
500ce11f295SJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
501ce11f295SJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
502ce11f295SJames Wright         PetscCall(PetscPrintf(
50365ba01baSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
50465ba01baSJames 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,
50565ba01baSJames Wright             part_shared_dof_ratio));
506ce11f295SJames Wright       }
507ce11f295SJames Wright 
508ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
509ce11f295SJames Wright       if (!rank) {
510ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
511ce11f295SJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
512ce11f295SJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
513ce11f295SJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
514ce11f295SJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
515ce11f295SJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
516ce11f295SJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
517ce11f295SJames Wright       }
518ce11f295SJames Wright 
519ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
520ce11f295SJames Wright       if (!rank) {
521ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
522ce11f295SJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
523ce11f295SJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
524ce11f295SJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
525ce11f295SJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
526ce11f295SJames Wright         PetscCall(PetscPrintf(
52765ba01baSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
52865ba01baSJames 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,
52965ba01baSJames Wright             part_shared_dof_ratio));
530ce11f295SJames Wright       }
531ce11f295SJames Wright 
532ce11f295SJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
533ce11f295SJames Wright       if (!rank) {
534ce11f295SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
535ce11f295SJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
536ce11f295SJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
537ce11f295SJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
538ce11f295SJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
539ce11f295SJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
540ce11f295SJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
5412e31f396SJames Wright       }
54265ba01baSJames Wright     }
5432e31f396SJames Wright 
5442e31f396SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
5452e31f396SJames Wright   }
5462e31f396SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5472e31f396SJames Wright }
548