xref: /honee/src/misc.c (revision ae2b091fac884a554e48acc4b4c187524c2a2818)
1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2*ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Miscellaneous utility functions
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdm.h>
9926a6279SJames Wright #include <petscsf.h>
10e419654dSJeremy L Thompson #include <petscts.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
139f59f36eSJames Wright #include "../qfunctions/mass.h"
14a515125bSLeila Ghaffari 
152b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
16b4c37c5cSJames Wright   Ceed         ceed = user->ceed;
17b2948607SJames Wright   CeedVector   mult_vec;
18b2948607SJames Wright   PetscMemType m_mem_type;
19b2948607SJames Wright   Vec          Multiplicity, Multiplicity_loc;
20b2948607SJames Wright 
21a515125bSLeila Ghaffari   PetscFunctionBeginUser;
22b4c37c5cSJames Wright   if (user->phys->ics_time_label) PetscCallCeed(ceed, CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time));
238f18bb8bSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
24a515125bSLeila Ghaffari 
25b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL));
26a515125bSLeila Ghaffari 
27a515125bSLeila Ghaffari   // -- Get multiplicity
28b2948607SJames Wright   PetscCall(DMGetLocalVector(dm, &Multiplicity_loc));
29a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Multiplicity_loc, &m_mem_type, mult_vec));
30b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec));
31a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(mult_vec, m_mem_type, Multiplicity_loc));
32a515125bSLeila Ghaffari 
33b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Multiplicity));
34b2948607SJames Wright   PetscCall(VecZeroEntries(Multiplicity));
35b2948607SJames Wright   PetscCall(DMLocalToGlobal(dm, Multiplicity_loc, ADD_VALUES, Multiplicity));
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari   // -- Fix multiplicity
38b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q, Q, Multiplicity));
39b2948607SJames Wright   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, Multiplicity_loc));
40a515125bSLeila Ghaffari 
41b2948607SJames Wright   PetscCall(DMRestoreLocalVector(dm, &Multiplicity_loc));
42b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Multiplicity));
43b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&mult_vec));
44d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
45a515125bSLeila Ghaffari }
46a515125bSLeila Ghaffari 
47c56e8d5bSJames Wright // Record boundary values from initial condition
48c56e8d5bSJames Wright PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc) {
49c56e8d5bSJames Wright   PetscFunctionBeginUser;
50313f2f1eSJames Wright   {  // Capture initial condition values in Qbc
51313f2f1eSJames Wright     Vec Qbc;
52313f2f1eSJames Wright 
53c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
54c56e8d5bSJames Wright     PetscCall(VecCopy(Q_loc, Qbc));
55c56e8d5bSJames Wright     PetscCall(VecZeroEntries(Q_loc));
56c56e8d5bSJames Wright     PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
57c56e8d5bSJames Wright     PetscCall(VecAXPY(Qbc, -1., Q_loc));
58c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
59313f2f1eSJames Wright   }
60c56e8d5bSJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_FromICs));
61c56e8d5bSJames Wright 
62313f2f1eSJames Wright   {  // Set boundary mask to zero out essential BCs
63313f2f1eSJames Wright     Vec boundary_mask, ones;
64313f2f1eSJames Wright 
65c56e8d5bSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
66313f2f1eSJames Wright     PetscCall(DMGetGlobalVector(dm, &ones));
67c56e8d5bSJames Wright     PetscCall(VecZeroEntries(boundary_mask));
68313f2f1eSJames Wright     PetscCall(VecSet(ones, 1.0));
69313f2f1eSJames Wright     PetscCall(DMGlobalToLocal(dm, ones, INSERT_VALUES, boundary_mask));
70c56e8d5bSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
71313f2f1eSJames Wright     PetscCall(DMRestoreGlobalVector(dm, &ones));
72313f2f1eSJames Wright   }
73c56e8d5bSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
74c56e8d5bSJames Wright }
75c56e8d5bSJames Wright 
76c56e8d5bSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
772b916ea7SJeremy L Thompson                                                   Vec grad_FVM) {
789d437337SJames Wright   Vec Qbc, boundary_mask;
79a515125bSLeila Ghaffari 
8006f41313SJames Wright   PetscFunctionBeginUser;
812eb7bf1fSJames Wright   // Mask (zero) Strong BC entries
829d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
839d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
849d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
859d437337SJames Wright 
862b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
872b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
882b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
89d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
90a515125bSLeila Ghaffari }
91a515125bSLeila Ghaffari 
92990f1db0SJed Brown static PetscErrorCode BinaryReadIntoInt(PetscViewer viewer, PetscInt *out, PetscDataType file_type) {
93990f1db0SJed Brown   PetscFunctionBeginUser;
9430ff0608SJames Wright   *out = -13;  // appease the overzealous GCC compiler warning Gods
95990f1db0SJed Brown   if (file_type == PETSC_INT32) {
96990f1db0SJed Brown     PetscInt32 val;
97990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT32));
98990f1db0SJed Brown     *out = val;
99990f1db0SJed Brown   } else if (file_type == PETSC_INT64) {
100990f1db0SJed Brown     PetscInt64 val;
101990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, &val, 1, NULL, PETSC_INT64));
102990f1db0SJed Brown     *out = val;
103990f1db0SJed Brown   } else {
104990f1db0SJed Brown     PetscCall(PetscViewerBinaryRead(viewer, out, 1, NULL, PETSC_INT));
105990f1db0SJed Brown   }
106990f1db0SJed Brown   PetscFunctionReturn(PETSC_SUCCESS);
107990f1db0SJed Brown }
108990f1db0SJed Brown 
109e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
110e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
111e1233009SJames Wright   PetscInt      file_step_number;
112e1233009SJames Wright   PetscInt32    token;
113e7754af5SKenneth E. Jansen   PetscReal     file_time;
114990f1db0SJed Brown   PetscDataType file_type = PETSC_INT32;
115e7754af5SKenneth E. Jansen 
11606f41313SJames Wright   PetscFunctionBeginUser;
117e1233009SJames Wright   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT32));
118e1233009SJames Wright   if (token == FLUIDS_FILE_TOKEN_32 || token == FLUIDS_FILE_TOKEN_64 ||
119e1233009SJames Wright       token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
120990f1db0SJed Brown     if (token == FLUIDS_FILE_TOKEN_32) file_type = PETSC_INT32;
121990f1db0SJed Brown     else if (token == FLUIDS_FILE_TOKEN_64) file_type = PETSC_INT64;
122990f1db0SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &file_step_number, file_type));
123e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
124e7754af5SKenneth E. Jansen     if (time) *time = file_time;
125e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
126e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
127e7754af5SKenneth E. Jansen     PetscInt length, N;
128990f1db0SJed Brown     PetscCall(BinaryReadIntoInt(viewer, &length, file_type));
129e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
130e7754af5SKenneth 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);
131e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
132e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
133e7754af5SKenneth E. Jansen 
134e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
135d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
136e7754af5SKenneth E. Jansen }
137e7754af5SKenneth E. Jansen 
138a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
139c56e8d5bSJames Wright PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q) {
140a515125bSLeila Ghaffari   Vec         Qref;
141a515125bSLeila Ghaffari   PetscViewer viewer;
142a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
143e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
144a515125bSLeila Ghaffari 
14506f41313SJames Wright   PetscFunctionBeginUser;
146a515125bSLeila Ghaffari   // Read reference file
1472b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
1484c07ec22SJames Wright   PetscCheck(strcmp(app_ctx->test_file_path, "") != 0, comm, PETSC_ERR_FILE_READ, "File for regression test not given");
149e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
150e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
151a515125bSLeila Ghaffari 
152a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1532b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1542b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1552b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1562b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari   // Check error
159a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1602b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
161a515125bSLeila Ghaffari   }
162a515125bSLeila Ghaffari 
163a515125bSLeila Ghaffari   // Cleanup
1642b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1652b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
166d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
167a515125bSLeila Ghaffari }
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari // Get error for problems with exact solutions
170c56e8d5bSJames Wright PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
171a515125bSLeila Ghaffari   PetscInt  loc_nodes;
172a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
173a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
174a515125bSLeila Ghaffari 
17506f41313SJames Wright   PetscFunctionBeginUser;
176a515125bSLeila Ghaffari   // Get exact solution at final time
177b2948607SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q_exact));
1782b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1792b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1802b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
181a515125bSLeila Ghaffari 
182a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1832b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1842b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1852b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
186a515125bSLeila Ghaffari 
187a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
1882b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
1892b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
190b2948607SJames Wright   PetscCall(DMRestoreGlobalVector(dm, &Q_exact));
191d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
192a515125bSLeila Ghaffari }
193a515125bSLeila Ghaffari 
194a515125bSLeila Ghaffari // Post-processing
195991aef52SJames Wright PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time) {
196a515125bSLeila Ghaffari   PetscInt          steps;
197f0784ed3SJed Brown   TSConvergedReason reason;
198a515125bSLeila Ghaffari 
19906f41313SJames Wright   PetscFunctionBeginUser;
200a515125bSLeila Ghaffari   // Print relative error
20158ce1233SJames Wright   if (problem->compute_exact_solution_error && user->app_ctx->test_type == TESTTYPE_NONE) {
202c56e8d5bSJames Wright     PetscCall(PrintError(ceed_data, dm, user, Q, final_time));
203a515125bSLeila Ghaffari   }
204a515125bSLeila Ghaffari 
205a515125bSLeila Ghaffari   // Print final time and number of steps
2062b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
207f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
2080e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
209f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
210f0784ed3SJed Brown                           steps, (double)final_time));
211a515125bSLeila Ghaffari   }
212a515125bSLeila Ghaffari 
213a515125bSLeila Ghaffari   // Output numerical values from command line
2142b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
215a515125bSLeila Ghaffari 
216a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
2170e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
218c56e8d5bSJames Wright     PetscCall(RegressionTest(user->app_ctx, Q));
219a515125bSLeila Ghaffari   }
220d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
221a515125bSLeila Ghaffari }
222a515125bSLeila Ghaffari 
223e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN    = 0xceedf00;  // for backwards compatibility
224e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_32 = 0xceedf32;
225e1233009SJames Wright const PetscInt32 FLUIDS_FILE_TOKEN_64 = 0xceedf64;
2269293eaa1SJed Brown 
227a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
228a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
229a515125bSLeila Ghaffari   PetscViewer viewer;
2302b916ea7SJeremy L Thompson 
23106f41313SJames Wright   PetscFunctionBeginUser;
2322b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
233e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2342b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
235d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
236a515125bSLeila Ghaffari }
237a515125bSLeila Ghaffari 
23815a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
23915a3537eSJed Brown int FreeContextPetsc(void *data) {
2402b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
24115a3537eSJed Brown   return CEED_ERROR_SUCCESS;
24215a3537eSJed Brown }
2439f59f36eSJames Wright 
2449f59f36eSJames Wright // Return mass qfunction specification for number of components N
2459f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2469f59f36eSJames Wright   PetscFunctionBeginUser;
2479f59f36eSJames Wright   switch (N) {
2489f59f36eSJames Wright     case 1:
249b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf));
2509f59f36eSJames Wright       break;
2519f59f36eSJames Wright     case 5:
252b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf));
2539f59f36eSJames Wright       break;
254c38c977aSJames Wright     case 7:
255b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf));
256c38c977aSJames Wright       break;
2579f59f36eSJames Wright     case 9:
258b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf));
2599f59f36eSJames Wright       break;
2609f59f36eSJames Wright     case 22:
261b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf));
2629f59f36eSJames Wright       break;
2639f59f36eSJames Wright     default:
2646f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2659f59f36eSJames Wright   }
2669f59f36eSJames Wright 
267b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP));
268b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE));
269b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP));
2703170c09fSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf, N));
271d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2729f59f36eSJames Wright }
273e5e81594SJames Wright 
274457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
275457a5831SJames Wright   PetscFunctionBeginUser;
276d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
277457a5831SJames Wright 
278457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
279457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
280457a5831SJames Wright 
281457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
282457a5831SJames Wright 
283457a5831SJames Wright   PetscCall(PetscFree(context));
284d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
285457a5831SJames Wright }
286457a5831SJames Wright 
287676080b4SJames Wright /*
288676080b4SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
289676080b4SJames Wright  *
290676080b4SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
291676080b4SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
292676080b4SJames Wright  *
293676080b4SJames 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.
294676080b4SJames Wright  *
295676080b4SJames Wright  * @param[in]  comm           MPI_Comm for the program
296676080b4SJames Wright  * @param[in]  path           Path to the file
297676080b4SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
298676080b4SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
299676080b4SJames Wright  * @param[out] fp File        pointer to the opened file
300676080b4SJames Wright  */
30142454adaSJames Wright PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
302676080b4SJames Wright                                  FILE **fp) {
303defe8520SJames Wright   int    ndims;
304676080b4SJames Wright   char   line[char_array_len];
305676080b4SJames Wright   char **array;
306676080b4SJames Wright 
307676080b4SJames Wright   PetscFunctionBeginUser;
308676080b4SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
309676080b4SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
310676080b4SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
311defe8520SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
312676080b4SJames Wright 
313676080b4SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
314676080b4SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
315d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
316676080b4SJames Wright }
317676080b4SJames Wright 
318676080b4SJames Wright /*
319676080b4SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
320676080b4SJames Wright  *
321676080b4SJames 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.
322676080b4SJames Wright  *
323676080b4SJames Wright  * @param[in]  comm  MPI_Comm for the program
324676080b4SJames Wright  * @param[in]  path  Path to the file
325676080b4SJames Wright  * @param[out] nrows Number of rows
326676080b4SJames Wright  */
32742454adaSJames Wright PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
328676080b4SJames Wright   const PetscInt char_array_len = 512;
329676080b4SJames Wright   PetscInt       dims[2];
330676080b4SJames Wright   FILE          *fp;
331676080b4SJames Wright 
332676080b4SJames Wright   PetscFunctionBeginUser;
33342454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
334676080b4SJames Wright   *nrows = dims[0];
335676080b4SJames Wright   PetscCall(PetscFClose(comm, fp));
336d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
337676080b4SJames Wright }
33862b7942eSJames Wright 
33942454adaSJames Wright PetscErrorCode PhastaDatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
340defe8520SJames Wright   PetscInt       dims[2];
34162b7942eSJames Wright   FILE          *fp;
34262b7942eSJames Wright   const PetscInt char_array_len = 512;
34362b7942eSJames Wright   char           line[char_array_len];
34462b7942eSJames Wright 
34506f41313SJames Wright   PetscFunctionBeginUser;
34642454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
34762b7942eSJames Wright 
34862b7942eSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
349039caf0dSJames Wright     int    ndims;
350039caf0dSJames Wright     char **row_array;
351039caf0dSJames Wright 
35262b7942eSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
35362b7942eSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3545d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
355defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
35662b7942eSJames Wright 
357039caf0dSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
358039caf0dSJames Wright     PetscCall(PetscStrToArrayDestroy(ndims, row_array));
35962b7942eSJames Wright   }
36062b7942eSJames Wright 
36162b7942eSJames Wright   PetscCall(PetscFClose(comm, fp));
362d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
36362b7942eSJames Wright }
3647eedc94cSJames Wright 
365926a6279SJames Wright // Print information about the given simulation run
3663678fae3SJames Wright PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts) {
367b4c37c5cSJames Wright   Ceed     ceed = user->ceed;
3683678fae3SJames Wright   MPI_Comm comm = PetscObjectComm((PetscObject)ts);
3693678fae3SJames Wright 
370926a6279SJames Wright   PetscFunctionBeginUser;
371926a6279SJames Wright   // Header and rank
372926a6279SJames Wright   char        host_name[PETSC_MAX_PATH_LEN];
373926a6279SJames Wright   PetscMPIInt rank, comm_size;
374926a6279SJames Wright   PetscCall(PetscGetHostName(host_name, sizeof host_name));
375926a6279SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
376926a6279SJames Wright   PetscCallMPI(MPI_Comm_size(comm, &comm_size));
377926a6279SJames Wright   PetscCall(PetscPrintf(comm,
378926a6279SJames Wright                         "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
379926a6279SJames Wright                         "  MPI:\n"
380926a6279SJames Wright                         "    Host Name                          : %s\n"
381926a6279SJames Wright                         "    Total ranks                        : %d\n",
382926a6279SJames Wright                         host_name, comm_size));
383926a6279SJames Wright 
384926a6279SJames Wright   // Problem specific info
3852d49c0afSJames Wright   PetscCall(problem->print_info(user, problem, user->app_ctx));
386926a6279SJames Wright 
387926a6279SJames Wright   // libCEED
388926a6279SJames Wright   const char *used_resource;
389926a6279SJames Wright   CeedMemType mem_type_backend;
390b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetResource(user->ceed, &used_resource));
391b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedGetPreferredMemType(user->ceed, &mem_type_backend));
392926a6279SJames Wright   PetscCall(PetscPrintf(comm,
393926a6279SJames Wright                         "  libCEED:\n"
394926a6279SJames Wright                         "    libCEED Backend                    : %s\n"
395926a6279SJames Wright                         "    libCEED Backend MemType            : %s\n",
396926a6279SJames Wright                         used_resource, CeedMemTypes[mem_type_backend]));
397926a6279SJames Wright   // PETSc
3983678fae3SJames Wright   VecType vec_type;
399926a6279SJames Wright   char    box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3";
400926a6279SJames Wright   if (problem->dim == 2) box_faces_str[3] = '\0';
401926a6279SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL));
402926a6279SJames Wright   PetscCall(DMGetVecType(user->dm, &vec_type));
403926a6279SJames Wright   PetscCall(PetscPrintf(comm,
404926a6279SJames Wright                         "  PETSc:\n"
405926a6279SJames Wright                         "    Box Faces                          : %s\n"
406926a6279SJames Wright                         "    DM VecType                         : %s\n"
407926a6279SJames Wright                         "    Time Stepping Scheme               : %s\n",
4083678fae3SJames Wright                         box_faces_str, vec_type, phys_ctx->implicit ? "implicit" : "explicit"));
4093678fae3SJames Wright   {
4103678fae3SJames Wright     char           pmat_type_str[PETSC_MAX_PATH_LEN];
4113678fae3SJames Wright     MatType        amat_type, pmat_type;
4123678fae3SJames Wright     Mat            Amat, Pmat;
4133678fae3SJames Wright     TSIJacobianFn *ijacob_function;
4143678fae3SJames Wright 
4153678fae3SJames Wright     PetscCall(TSGetIJacobian(ts, &Amat, &Pmat, &ijacob_function, NULL));
4163678fae3SJames Wright     PetscCall(MatGetType(Amat, &amat_type));
4173678fae3SJames Wright     PetscCall(MatGetType(Pmat, &pmat_type));
4183678fae3SJames Wright 
4193678fae3SJames Wright     PetscCall(PetscStrncpy(pmat_type_str, pmat_type, sizeof(pmat_type_str)));
4203678fae3SJames Wright     if (!strcmp(pmat_type, MATCEED)) {
4213678fae3SJames Wright       MatType pmat_coo_type;
4223678fae3SJames Wright       char    pmat_coo_type_str[PETSC_MAX_PATH_LEN];
4233678fae3SJames Wright 
4243678fae3SJames Wright       PetscCall(MatCeedGetCOOMatType(Pmat, &pmat_coo_type));
4253678fae3SJames Wright       PetscCall(PetscSNPrintf(pmat_coo_type_str, sizeof(pmat_coo_type_str), " (COO MatType: %s)", pmat_coo_type));
4263678fae3SJames Wright       PetscCall(PetscStrlcat(pmat_type_str, pmat_coo_type_str, sizeof(pmat_type_str)));
4273678fae3SJames Wright     }
4283678fae3SJames Wright     if (ijacob_function) {
4293678fae3SJames Wright       PetscCall(PetscPrintf(comm,
4303678fae3SJames Wright                             "    IJacobian A MatType                : %s\n"
4313678fae3SJames Wright                             "    IJacobian P MatType                : %s\n",
4323678fae3SJames Wright                             amat_type, pmat_type_str));
4333678fae3SJames Wright     }
4343678fae3SJames Wright   }
435926a6279SJames Wright   if (user->app_ctx->cont_steps) {
436926a6279SJames Wright     PetscCall(PetscPrintf(comm,
437926a6279SJames Wright                           "  Continue:\n"
438926a6279SJames Wright                           "    Filename:                          : %s\n"
439926a6279SJames Wright                           "    Step:                              : %" PetscInt_FMT "\n"
440926a6279SJames Wright                           "    Time:                              : %g\n",
441926a6279SJames Wright                           user->app_ctx->cont_file, user->app_ctx->cont_steps, user->app_ctx->cont_time));
442926a6279SJames Wright   }
443926a6279SJames Wright   // Mesh
444926a6279SJames Wright   const PetscInt num_comp_q = 5;
445926a6279SJames Wright   PetscInt       glob_dofs, owned_dofs, local_dofs;
446926a6279SJames Wright   const CeedInt  num_P = user->app_ctx->degree + 1, num_Q = num_P + user->app_ctx->q_extra;
447926a6279SJames Wright   PetscCall(DMGetGlobalVectorInfo(user->dm, &owned_dofs, &glob_dofs, NULL));
448926a6279SJames Wright   PetscCall(DMGetLocalVectorInfo(user->dm, &local_dofs, NULL, NULL));
449926a6279SJames Wright   PetscCall(PetscPrintf(comm,
450926a6279SJames Wright                         "  Mesh:\n"
451926a6279SJames Wright                         "    Number of 1D Basis Nodes (P)       : %" CeedInt_FMT "\n"
452926a6279SJames Wright                         "    Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n"
453926a6279SJames Wright                         "    Global DoFs                        : %" PetscInt_FMT "\n"
454926a6279SJames Wright                         "    DoFs per node                      : %" PetscInt_FMT "\n"
455dfeb939dSJames Wright                         "    Global %" PetscInt_FMT "-DoF nodes                 : %" PetscInt_FMT "\n",
456dfeb939dSJames Wright                         num_P, num_Q, glob_dofs, num_comp_q, num_comp_q, glob_dofs / num_comp_q));
457926a6279SJames Wright   // -- Get Partition Statistics
458926a6279SJames Wright   PetscCall(PetscPrintf(comm, "  Partition:                             (min,max,median,max/median)\n"));
459926a6279SJames Wright   {
460926a6279SJames Wright     PetscInt *gather_buffer = NULL;
461dfeb939dSJames Wright     PetscInt  part_owned_dofs[3], part_local_dofs[3], part_boundary_dofs[3], part_neighbors[3];
462926a6279SJames Wright     PetscInt  median_index = comm_size % 2 ? comm_size / 2 : comm_size / 2 - 1;
463926a6279SJames Wright     if (!rank) PetscCall(PetscMalloc1(comm_size, &gather_buffer));
464926a6279SJames Wright 
465dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&owned_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
466926a6279SJames Wright     if (!rank) {
467926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
468926a6279SJames Wright       part_owned_dofs[0]             = gather_buffer[0];              // min
469926a6279SJames Wright       part_owned_dofs[1]             = gather_buffer[comm_size - 1];  // max
470926a6279SJames Wright       part_owned_dofs[2]             = gather_buffer[median_index];   // median
471926a6279SJames Wright       PetscReal part_owned_dof_ratio = (PetscReal)part_owned_dofs[1] / (PetscReal)part_owned_dofs[2];
472dfeb939dSJames Wright       PetscCall(PetscPrintf(
473dfeb939dSJames Wright           comm, "    Global Vector %" PetscInt_FMT "-DoF nodes          : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
474926a6279SJames 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));
475926a6279SJames Wright     }
476926a6279SJames Wright 
477dfeb939dSJames Wright     PetscCallMPI(MPI_Gather(&local_dofs, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
478926a6279SJames Wright     if (!rank) {
479926a6279SJames Wright       PetscCall(PetscSortInt(comm_size, gather_buffer));
480926a6279SJames Wright       part_local_dofs[0]             = gather_buffer[0];              // min
481926a6279SJames Wright       part_local_dofs[1]             = gather_buffer[comm_size - 1];  // max
482926a6279SJames Wright       part_local_dofs[2]             = gather_buffer[median_index];   // median
483926a6279SJames Wright       PetscReal part_local_dof_ratio = (PetscReal)part_local_dofs[1] / (PetscReal)part_local_dofs[2];
484dfeb939dSJames Wright       PetscCall(PetscPrintf(
485dfeb939dSJames Wright           comm, "    Local Vector %" PetscInt_FMT "-DoF nodes           : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n", num_comp_q,
486926a6279SJames 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));
487926a6279SJames Wright     }
488926a6279SJames Wright 
48945abf96eSJames Wright     if (comm_size != 1) {
490dfeb939dSJames Wright       PetscInt num_remote_roots_total = 0, num_remote_leaves_total = 0, num_ghost_interface_ranks = 0, num_owned_interface_ranks = 0;
491926a6279SJames Wright       {
492926a6279SJames Wright         PetscSF            sf;
493dfeb939dSJames Wright         PetscInt           nrranks, niranks;
494dfeb939dSJames Wright         const PetscInt    *roffset, *rmine, *rremote, *ioffset, *irootloc;
495dfeb939dSJames Wright         const PetscMPIInt *rranks, *iranks;
496926a6279SJames Wright         PetscCall(DMGetSectionSF(user->dm, &sf));
497926a6279SJames Wright         PetscCall(PetscSFGetRootRanks(sf, &nrranks, &rranks, &roffset, &rmine, &rremote));
498dfeb939dSJames Wright         PetscCall(PetscSFGetLeafRanks(sf, &niranks, &iranks, &ioffset, &irootloc));
499926a6279SJames Wright         for (PetscInt i = 0; i < nrranks; i++) {
500926a6279SJames Wright           if (rranks[i] == rank) continue;  // Ignore same-part global->local transfers
501926a6279SJames Wright           num_remote_roots_total += roffset[i + 1] - roffset[i];
502dfeb939dSJames Wright           num_ghost_interface_ranks++;
503dfeb939dSJames Wright         }
504dfeb939dSJames Wright         for (PetscInt i = 0; i < niranks; i++) {
505dfeb939dSJames Wright           if (iranks[i] == rank) continue;
506dfeb939dSJames Wright           num_remote_leaves_total += ioffset[i + 1] - ioffset[i];
507dfeb939dSJames Wright           num_owned_interface_ranks++;
508926a6279SJames Wright         }
509926a6279SJames Wright       }
510dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_roots_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
511926a6279SJames Wright       if (!rank) {
512926a6279SJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
513dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
514dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
515dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
516dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
517dfeb939dSJames Wright         PetscCall(PetscPrintf(
51845abf96eSJames Wright             comm, "    Ghost Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
51945abf96eSJames 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,
52045abf96eSJames Wright             part_shared_dof_ratio));
521dfeb939dSJames Wright       }
522dfeb939dSJames Wright 
523dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_ghost_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
524dfeb939dSJames Wright       if (!rank) {
525dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
526dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
527dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
528dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
529dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
530dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Ghost Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
531dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
532dfeb939dSJames Wright       }
533dfeb939dSJames Wright 
534dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_remote_leaves_total, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
535dfeb939dSJames Wright       if (!rank) {
536dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
537dfeb939dSJames Wright         part_boundary_dofs[0]           = gather_buffer[0];              // min
538dfeb939dSJames Wright         part_boundary_dofs[1]           = gather_buffer[comm_size - 1];  // max
539dfeb939dSJames Wright         part_boundary_dofs[2]           = gather_buffer[median_index];   // median
540dfeb939dSJames Wright         PetscReal part_shared_dof_ratio = (PetscReal)part_boundary_dofs[1] / (PetscReal)part_boundary_dofs[2];
541dfeb939dSJames Wright         PetscCall(PetscPrintf(
54245abf96eSJames Wright             comm, "    Owned Interface %" PetscInt_FMT "-DoF nodes        : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
54345abf96eSJames 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,
54445abf96eSJames Wright             part_shared_dof_ratio));
545dfeb939dSJames Wright       }
546dfeb939dSJames Wright 
547dfeb939dSJames Wright       PetscCallMPI(MPI_Gather(&num_owned_interface_ranks, 1, MPIU_INT, gather_buffer, 1, MPIU_INT, 0, comm));
548dfeb939dSJames Wright       if (!rank) {
549dfeb939dSJames Wright         PetscCall(PetscSortInt(comm_size, gather_buffer));
550dfeb939dSJames Wright         part_neighbors[0]              = gather_buffer[0];              // min
551dfeb939dSJames Wright         part_neighbors[1]              = gather_buffer[comm_size - 1];  // max
552dfeb939dSJames Wright         part_neighbors[2]              = gather_buffer[median_index];   // median
553dfeb939dSJames Wright         PetscReal part_neighbors_ratio = (PetscReal)part_neighbors[1] / (PetscReal)part_neighbors[2];
554dfeb939dSJames Wright         PetscCall(PetscPrintf(comm, "    Owned Interface Ranks              : %" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ", %f\n",
555dfeb939dSJames Wright                               part_neighbors[0], part_neighbors[1], part_neighbors[2], part_neighbors_ratio));
556926a6279SJames Wright       }
55745abf96eSJames Wright     }
558926a6279SJames Wright 
559926a6279SJames Wright     if (!rank) PetscCall(PetscFree(gather_buffer));
560926a6279SJames Wright   }
561926a6279SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
562926a6279SJames Wright }
563