xref: /honee/src/misc.c (revision defe8520f6a758707c23a8f3c4b0156b231420b2)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Miscellaneous utility functions
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdm.h>
13e419654dSJeremy L Thompson #include <petscts.h>
14e419654dSJeremy L Thompson 
15a515125bSLeila Ghaffari #include "../navierstokes.h"
169f59f36eSJames Wright #include "../qfunctions/mass.h"
17a515125bSLeila Ghaffari 
182b916ea7SJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
19a515125bSLeila Ghaffari   PetscFunctionBeginUser;
20a515125bSLeila Ghaffari 
21a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
22b7f03f12SJed Brown   // Update time for evaluation
23a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
248f18bb8bSJames Wright   if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics_ctx->op, user->phys->ics_time_label, &time);
25a515125bSLeila Ghaffari 
26a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
27a515125bSLeila Ghaffari   // ICs
28a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
29a515125bSLeila Ghaffari   // -- CEED Restriction
30a515125bSLeila Ghaffari   CeedVector q0_ceed;
31a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
32a515125bSLeila Ghaffari 
33a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
348f18bb8bSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(NULL, Q, ceed_data->op_ics_ctx));
35a515125bSLeila Ghaffari 
36a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
37a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
38a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
39a515125bSLeila Ghaffari   // -- CEED Restriction
40a515125bSLeila Ghaffari   CeedVector mult_vec;
41a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
42a515125bSLeila Ghaffari 
43a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
44a515125bSLeila Ghaffari   PetscMemType m_mem_type;
45a515125bSLeila Ghaffari   Vec          multiplicity_loc;
462b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
47fd969b44SJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   // -- Get multiplicity
50a515125bSLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
51a515125bSLeila Ghaffari 
52a515125bSLeila Ghaffari   // -- Restore vectors
53fd969b44SJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
54a515125bSLeila Ghaffari 
55a515125bSLeila Ghaffari   // -- Local-to-Global
56a515125bSLeila Ghaffari   Vec multiplicity;
572b916ea7SJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
582b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
592b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
60a515125bSLeila Ghaffari 
61a515125bSLeila Ghaffari   // -- Fix multiplicity
622b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
632b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
64a515125bSLeila Ghaffari 
65a515125bSLeila Ghaffari   // -- Restore vectors
662b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
672b916ea7SJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
68a515125bSLeila Ghaffari 
69a515125bSLeila Ghaffari   // Cleanup
70a515125bSLeila Ghaffari   CeedVectorDestroy(&mult_vec);
71a515125bSLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
72a515125bSLeila Ghaffari 
73d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
74a515125bSLeila Ghaffari }
75a515125bSLeila Ghaffari 
762b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(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   PetscFunctionBegin;
80a515125bSLeila Ghaffari 
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));
89a515125bSLeila Ghaffari 
90d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
91a515125bSLeila Ghaffari }
92a515125bSLeila Ghaffari 
93e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
94e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
95e7754af5SKenneth E. Jansen   PetscInt  token, file_step_number;
96e7754af5SKenneth E. Jansen   PetscReal file_time;
97e7754af5SKenneth E. Jansen   PetscFunctionBeginUser;
98e7754af5SKenneth E. Jansen 
99e7754af5SKenneth E. Jansen   // Attempt
100e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
101e7754af5SKenneth E. Jansen   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
102e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
103e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
104e7754af5SKenneth E. Jansen     if (time) *time = file_time;
105e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
106e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
107e7754af5SKenneth E. Jansen     PetscInt length, N;
108e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
109e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
110e7754af5SKenneth 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);
111e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
112e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
113e7754af5SKenneth E. Jansen 
114e7754af5SKenneth E. Jansen   // Load Q from existent solution
115e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
116e7754af5SKenneth E. Jansen 
117d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
118e7754af5SKenneth E. Jansen }
119e7754af5SKenneth E. Jansen 
120a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
121a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
122a515125bSLeila Ghaffari   Vec         Qref;
123a515125bSLeila Ghaffari   PetscViewer viewer;
124a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
125e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
126a515125bSLeila Ghaffari   PetscFunctionBegin;
127a515125bSLeila Ghaffari 
128a515125bSLeila Ghaffari   // Read reference file
1292b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
130e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
131e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
132a515125bSLeila Ghaffari 
133a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1342b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1352b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1362b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1372b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
138a515125bSLeila Ghaffari 
139a515125bSLeila Ghaffari   // Check error
140a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1412b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
142a515125bSLeila Ghaffari   }
143a515125bSLeila Ghaffari 
144a515125bSLeila Ghaffari   // Cleanup
1452b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1462b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
147a515125bSLeila Ghaffari 
148d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
149a515125bSLeila Ghaffari }
150a515125bSLeila Ghaffari 
151a515125bSLeila Ghaffari // Get error for problems with exact solutions
1522b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
153a515125bSLeila Ghaffari   PetscInt  loc_nodes;
154a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
155a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
156a515125bSLeila Ghaffari   PetscFunctionBegin;
157a515125bSLeila Ghaffari 
158a515125bSLeila Ghaffari   // Get exact solution at final time
1592b916ea7SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1602b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1612b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1622b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
163a515125bSLeila Ghaffari 
164a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1652b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1662b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1672b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari   // Compute relative error
170a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
171a515125bSLeila Ghaffari 
172a515125bSLeila Ghaffari   // Output relative error
1732b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
174a515125bSLeila Ghaffari   // Cleanup
1752b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1762b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
177a515125bSLeila Ghaffari 
178d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
179a515125bSLeila Ghaffari }
180a515125bSLeila Ghaffari 
181a515125bSLeila Ghaffari // Post-processing
1822b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
183a515125bSLeila Ghaffari   PetscInt          steps;
184f0784ed3SJed Brown   TSConvergedReason reason;
185a515125bSLeila Ghaffari   PetscFunctionBegin;
186a515125bSLeila Ghaffari 
187a515125bSLeila Ghaffari   // Print relative error
1880e1e9333SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1892b916ea7SJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
190a515125bSLeila Ghaffari   }
191a515125bSLeila Ghaffari 
192a515125bSLeila Ghaffari   // Print final time and number of steps
1932b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
194f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
1950e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
196f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
197f0784ed3SJed Brown                           steps, (double)final_time));
198a515125bSLeila Ghaffari   }
199a515125bSLeila Ghaffari 
200a515125bSLeila Ghaffari   // Output numerical values from command line
2012b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
202a515125bSLeila Ghaffari 
203a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
2040e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2052b916ea7SJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
206a515125bSLeila Ghaffari   }
207a515125bSLeila Ghaffari 
208d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
209a515125bSLeila Ghaffari }
210a515125bSLeila Ghaffari 
2119293eaa1SJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
2129293eaa1SJed Brown 
213a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
214a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
215a515125bSLeila Ghaffari   PetscViewer viewer;
2162b916ea7SJeremy L Thompson 
217a515125bSLeila Ghaffari   PetscFunctionBegin;
218a515125bSLeila Ghaffari 
2192b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
220e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2212b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
222a515125bSLeila Ghaffari 
223d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
224a515125bSLeila Ghaffari }
225a515125bSLeila Ghaffari 
226a515125bSLeila Ghaffari // Record boundary values from initial condition
227a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2289d437337SJames Wright   Vec Qbc, boundary_mask;
229a515125bSLeila Ghaffari   PetscFunctionBegin;
230a515125bSLeila Ghaffari 
2312b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2322b916ea7SJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2332b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2342b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2352b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2362b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2372b916ea7SJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
238a515125bSLeila Ghaffari 
2399d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2409d437337SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2419d437337SJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2429d437337SJames Wright   PetscCall(VecSet(Q, 1.0));
2439d437337SJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2449d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2459d437337SJames Wright 
246d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
247a515125bSLeila Ghaffari }
24815a3537eSJed Brown 
24915a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
25015a3537eSJed Brown int FreeContextPetsc(void *data) {
2512b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
25215a3537eSJed Brown   return CEED_ERROR_SUCCESS;
25315a3537eSJed Brown }
2549f59f36eSJames Wright 
2559f59f36eSJames Wright // Return mass qfunction specification for number of components N
2569f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2579f59f36eSJames Wright   PetscFunctionBeginUser;
2589f59f36eSJames Wright 
2599f59f36eSJames Wright   switch (N) {
2609f59f36eSJames Wright     case 1:
2616f539f71SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_1, Mass_1_loc, qf);
2629f59f36eSJames Wright       break;
2639f59f36eSJames Wright     case 5:
2646f539f71SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_5, Mass_5_loc, qf);
2659f59f36eSJames Wright       break;
266c38c977aSJames Wright     case 7:
2676f539f71SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_7, Mass_7_loc, qf);
268c38c977aSJames Wright       break;
2699f59f36eSJames Wright     case 9:
2706f539f71SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_9, Mass_9_loc, qf);
2719f59f36eSJames Wright       break;
2729f59f36eSJames Wright     case 22:
2736f539f71SJames Wright       CeedQFunctionCreateInterior(ceed, 1, Mass_22, Mass_22_loc, qf);
2749f59f36eSJames Wright       break;
2759f59f36eSJames Wright     default:
2766f539f71SJames Wright       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Could not find mass qfunction of size %d", N);
2779f59f36eSJames Wright   }
2789f59f36eSJames Wright 
2799f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
2809f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
2819f59f36eSJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
282d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2839f59f36eSJames Wright }
284e5e81594SJames Wright 
285e5e81594SJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space
286e5e81594SJames Wright  *
287e5e81594SJames Wright  * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum.
288e5e81594SJames Wright  *
289e5e81594SJames Wright  * @param[in]  source_vec    Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc
290e5e81594SJames Wright  * @param[out] target_vec    Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc
291e5e81594SJames Wright  * @param[in]  rhs_matop_ctx MatopApplyContext for performing the RHS evaluation
292e5e81594SJames Wright  * @param[in]  ksp           KSP for solving the consistent projection problem
293e5e81594SJames Wright  */
294f7325489SJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp) {
295e5e81594SJames Wright   PetscFunctionBeginUser;
296e5e81594SJames Wright 
297f7325489SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(source_vec, target_vec, rhs_matop_ctx));
298e5e81594SJames Wright   PetscCall(KSPSolve(ksp, target_vec, target_vec));
299e5e81594SJames Wright 
300d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
301e5e81594SJames Wright }
302676080b4SJames Wright 
303457a5831SJames Wright PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context) {
304457a5831SJames Wright   PetscFunctionBeginUser;
305d949ddfcSJames Wright   if (context == NULL) PetscFunctionReturn(PETSC_SUCCESS);
306457a5831SJames Wright 
307457a5831SJames Wright   PetscCall(DMDestroy(&context->dm));
308457a5831SJames Wright   PetscCall(KSPDestroy(&context->ksp));
309457a5831SJames Wright 
310457a5831SJames Wright   PetscCall(OperatorApplyContextDestroy(context->l2_rhs_ctx));
311457a5831SJames Wright 
312457a5831SJames Wright   PetscCall(PetscFree(context));
313457a5831SJames Wright 
314d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
315457a5831SJames Wright }
316457a5831SJames Wright 
317676080b4SJames Wright /*
318676080b4SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
319676080b4SJames Wright  *
320676080b4SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
321676080b4SJames Wright  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
322676080b4SJames Wright  *
323676080b4SJames 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.
324676080b4SJames Wright  *
325676080b4SJames Wright  * @param[in]  comm           MPI_Comm for the program
326676080b4SJames Wright  * @param[in]  path           Path to the file
327676080b4SJames Wright  * @param[in]  char_array_len Length of the character array that should contain each line
328676080b4SJames Wright  * @param[out] dims           Dimensions of the file, taken from the first line of the file
329676080b4SJames Wright  * @param[out] fp File        pointer to the opened file
330676080b4SJames Wright  */
331676080b4SJames Wright PetscErrorCode PHASTADatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
332676080b4SJames Wright                                  FILE **fp) {
333*defe8520SJames Wright   int    ndims;
334676080b4SJames Wright   char   line[char_array_len];
335676080b4SJames Wright   char **array;
336676080b4SJames Wright 
337676080b4SJames Wright   PetscFunctionBeginUser;
338676080b4SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
339676080b4SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
340676080b4SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
341*defe8520SJames Wright   PetscCheck(ndims == 2, comm, PETSC_ERR_FILE_UNEXPECTED, "Found %d dimensions instead of 2 on the first line of %s", ndims, path);
342676080b4SJames Wright 
343676080b4SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
344676080b4SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
345676080b4SJames Wright 
346d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
347676080b4SJames Wright }
348676080b4SJames Wright 
349676080b4SJames Wright /*
350676080b4SJames Wright  * @brief Get the number of rows for the PHASTA file at path.
351676080b4SJames Wright  *
352676080b4SJames 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.
353676080b4SJames Wright  *
354676080b4SJames Wright  * @param[in]  comm  MPI_Comm for the program
355676080b4SJames Wright  * @param[in]  path  Path to the file
356676080b4SJames Wright  * @param[out] nrows Number of rows
357676080b4SJames Wright  */
358676080b4SJames Wright PetscErrorCode PHASTADatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
359676080b4SJames Wright   const PetscInt char_array_len = 512;
360676080b4SJames Wright   PetscInt       dims[2];
361676080b4SJames Wright   FILE          *fp;
362676080b4SJames Wright 
363676080b4SJames Wright   PetscFunctionBeginUser;
364676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
365676080b4SJames Wright   *nrows = dims[0];
366676080b4SJames Wright   PetscCall(PetscFClose(comm, fp));
367676080b4SJames Wright 
368d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
369676080b4SJames Wright }
37062b7942eSJames Wright 
37162b7942eSJames Wright PetscErrorCode PHASTADatFileReadToArrayReal(MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]) {
372*defe8520SJames Wright   PetscInt       dims[2];
373*defe8520SJames Wright   int            ndims;
37462b7942eSJames Wright   FILE          *fp;
37562b7942eSJames Wright   const PetscInt char_array_len = 512;
37662b7942eSJames Wright   char           line[char_array_len];
37762b7942eSJames Wright   char         **row_array;
37862b7942eSJames Wright   PetscFunctionBeginUser;
37962b7942eSJames Wright 
38062b7942eSJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
38162b7942eSJames Wright 
38262b7942eSJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
38362b7942eSJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
38462b7942eSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &row_array));
3855d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
386*defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain enough columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
38762b7942eSJames Wright 
38862b7942eSJames Wright     for (PetscInt j = 0; j < dims[1]; j++) {
38962b7942eSJames Wright       array[i * dims[1] + j] = (PetscReal)atof(row_array[j]);
39062b7942eSJames Wright     }
39162b7942eSJames Wright   }
39262b7942eSJames Wright 
39362b7942eSJames Wright   PetscCall(PetscFClose(comm, fp));
39462b7942eSJames Wright 
395d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
39662b7942eSJames Wright }
3977eedc94cSJames Wright 
3987eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorApply;
3997eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemble;
4007eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssembleDiagonal;
4017eedc94cSJames Wright PetscLogEvent       FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
4027eedc94cSJames Wright static PetscClassId libCEED_classid;
4037eedc94cSJames Wright 
4047eedc94cSJames Wright PetscErrorCode RegisterLogEvents() {
4057eedc94cSJames Wright   PetscFunctionBeginUser;
4067eedc94cSJames Wright   PetscCall(PetscClassIdRegister("libCEED", &libCEED_classid));
4077eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpApply", libCEED_classid, &FLUIDS_CeedOperatorApply));
4087eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsm", libCEED_classid, &FLUIDS_CeedOperatorAssemble));
4097eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmD", libCEED_classid, &FLUIDS_CeedOperatorAssembleDiagonal));
4107eedc94cSJames Wright   PetscCall(PetscLogEventRegister("CeedOpAsmPBD", libCEED_classid, &FLUIDS_CeedOperatorAssemblePointBlockDiagonal));
411d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4127eedc94cSJames Wright }
413*defe8520SJames Wright 
414*defe8520SJames Wright /**
415*defe8520SJames Wright   @brief Translate array of CeedInt to PetscInt.
416*defe8520SJames Wright     If the types differ, `array_ceed` is freed with `free()` and `array_petsc` is allocated with `malloc()`.
417*defe8520SJames Wright     Caller is responsible for freeing `array_petsc` with `free()`.
418*defe8520SJames Wright 
419*defe8520SJames Wright   @param[in]      num_entries  Number of array entries
420*defe8520SJames Wright   @param[in,out]  array_ceed   Array of CeedInts
421*defe8520SJames Wright   @param[out]     array_petsc  Array of PetscInts
422*defe8520SJames Wright **/
423*defe8520SJames Wright PetscErrorCode IntArrayC2P(PetscInt num_entries, CeedInt **array_ceed, PetscInt **array_petsc) {
424*defe8520SJames Wright   CeedInt  int_c = 0;
425*defe8520SJames Wright   PetscInt int_p = 0;
426*defe8520SJames Wright 
427*defe8520SJames Wright   PetscFunctionBeginUser;
428*defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
429*defe8520SJames Wright     *array_petsc = (PetscInt *)*array_ceed;
430*defe8520SJames Wright   } else {
431*defe8520SJames Wright     *array_petsc = malloc(num_entries * sizeof(PetscInt));
432*defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_petsc)[i] = (*array_ceed)[i];
433*defe8520SJames Wright     free(*array_ceed);
434*defe8520SJames Wright   }
435*defe8520SJames Wright   *array_ceed = NULL;
436*defe8520SJames Wright 
437*defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
438*defe8520SJames Wright }
439*defe8520SJames Wright 
440*defe8520SJames Wright /**
441*defe8520SJames Wright   @brief Translate array of PetscInt to CeedInt.
442*defe8520SJames Wright     If the types differ, `array_petsc` is freed with `PetscFree()` and `array_ceed` is allocated with `PetscMalloc1()`.
443*defe8520SJames Wright     Caller is responsible for freeing `array_ceed` with `PetscFree()`.
444*defe8520SJames Wright 
445*defe8520SJames Wright   @param[in]      num_entries  Number of array entries
446*defe8520SJames Wright   @param[in,out]  array_petsc  Array of PetscInts
447*defe8520SJames Wright   @param[out]     array_ceed   Array of CeedInts
448*defe8520SJames Wright **/
449*defe8520SJames Wright PetscErrorCode IntArrayP2C(PetscInt num_entries, PetscInt **array_petsc, CeedInt **array_ceed) {
450*defe8520SJames Wright   CeedInt  int_c = 0;
451*defe8520SJames Wright   PetscInt int_p = 0;
452*defe8520SJames Wright 
453*defe8520SJames Wright   PetscFunctionBeginUser;
454*defe8520SJames Wright   if (sizeof(int_c) == sizeof(int_p)) {
455*defe8520SJames Wright     *array_ceed = (CeedInt *)*array_petsc;
456*defe8520SJames Wright   } else {
457*defe8520SJames Wright     PetscCall(PetscMalloc1(num_entries, array_ceed));
458*defe8520SJames Wright     for (PetscInt i = 0; i < num_entries; i++) (*array_ceed)[i] = (*array_petsc)[i];
459*defe8520SJames Wright     PetscCall(PetscFree(*array_petsc));
460*defe8520SJames Wright   }
461*defe8520SJames Wright   *array_petsc = NULL;
462*defe8520SJames Wright 
463*defe8520SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
464*defe8520SJames Wright }
465