xref: /honee/src/misc.c (revision f7325489fbf7ccd07a1d5c11e14143befdb1198d)
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   // ---------------------------------------------------------------------------
24a5f46a7bSJeremy L Thompson   if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics, 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
34a515125bSLeila Ghaffari   PetscMemType q0_mem_type;
35fd969b44SJames Wright   PetscCall(VecP2C(Q_loc, &q0_mem_type, q0_ceed));
36a515125bSLeila Ghaffari 
37a515125bSLeila Ghaffari   // -- Apply CEED Operator
382b916ea7SJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
39a515125bSLeila Ghaffari 
40a515125bSLeila Ghaffari   // -- Restore vectors
41fd969b44SJames Wright   PetscCall(VecC2P(q0_ceed, q0_mem_type, Q_loc));
42a515125bSLeila Ghaffari 
43a515125bSLeila Ghaffari   // -- Local-to-Global
442b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q));
452b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
46a515125bSLeila Ghaffari 
47a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
48a515125bSLeila Ghaffari   // Fix multiplicity for output of ICs
49a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
50a515125bSLeila Ghaffari   // -- CEED Restriction
51a515125bSLeila Ghaffari   CeedVector mult_vec;
52a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
53a515125bSLeila Ghaffari 
54a515125bSLeila Ghaffari   // -- Place PETSc vector in CEED vector
55a515125bSLeila Ghaffari   PetscMemType m_mem_type;
56a515125bSLeila Ghaffari   Vec          multiplicity_loc;
572b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
58fd969b44SJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari   // -- Get multiplicity
61a515125bSLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
62a515125bSLeila Ghaffari 
63a515125bSLeila Ghaffari   // -- Restore vectors
64fd969b44SJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
65a515125bSLeila Ghaffari 
66a515125bSLeila Ghaffari   // -- Local-to-Global
67a515125bSLeila Ghaffari   Vec multiplicity;
682b916ea7SJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
692b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
702b916ea7SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
71a515125bSLeila Ghaffari 
72a515125bSLeila Ghaffari   // -- Fix multiplicity
732b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
742b916ea7SJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
75a515125bSLeila Ghaffari 
76a515125bSLeila Ghaffari   // -- Restore vectors
772b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
782b916ea7SJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
79a515125bSLeila Ghaffari 
80a515125bSLeila Ghaffari   // Cleanup
81a515125bSLeila Ghaffari   CeedVectorDestroy(&mult_vec);
82a515125bSLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
83a515125bSLeila Ghaffari 
84a515125bSLeila Ghaffari   PetscFunctionReturn(0);
85a515125bSLeila Ghaffari }
86a515125bSLeila Ghaffari 
872b916ea7SJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
882b916ea7SJeremy L Thompson                                              Vec grad_FVM) {
899d437337SJames Wright   Vec Qbc, boundary_mask;
90a515125bSLeila Ghaffari   PetscFunctionBegin;
91a515125bSLeila Ghaffari 
929d437337SJames Wright   // Mask (zero) Dirichlet entries
939d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
949d437337SJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
959d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
969d437337SJames Wright 
972b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
982b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
992b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
100a515125bSLeila Ghaffari 
101a515125bSLeila Ghaffari   PetscFunctionReturn(0);
102a515125bSLeila Ghaffari }
103a515125bSLeila Ghaffari 
104e7754af5SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
105e7754af5SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
106e7754af5SKenneth E. Jansen   PetscInt  token, file_step_number;
107e7754af5SKenneth E. Jansen   PetscReal file_time;
108e7754af5SKenneth E. Jansen   PetscFunctionBeginUser;
109e7754af5SKenneth E. Jansen 
110e7754af5SKenneth E. Jansen   // Attempt
111e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
112e7754af5SKenneth E. Jansen   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
113e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
114e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
115e7754af5SKenneth E. Jansen     if (time) *time = file_time;
116e7754af5SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
117e7754af5SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
118e7754af5SKenneth E. Jansen     PetscInt length, N;
119e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
120e7754af5SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
121e7754af5SKenneth 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);
122e7754af5SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
123e7754af5SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
124e7754af5SKenneth E. Jansen 
125e7754af5SKenneth E. Jansen   // Load Q from existent solution
126e7754af5SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
127e7754af5SKenneth E. Jansen 
128e7754af5SKenneth E. Jansen   PetscFunctionReturn(0);
129e7754af5SKenneth E. Jansen }
130e7754af5SKenneth E. Jansen 
131a515125bSLeila Ghaffari // Compare reference solution values with current test run for CI
132a515125bSLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
133a515125bSLeila Ghaffari   Vec         Qref;
134a515125bSLeila Ghaffari   PetscViewer viewer;
135a515125bSLeila Ghaffari   PetscReal   error, Qrefnorm;
136e7754af5SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
137a515125bSLeila Ghaffari   PetscFunctionBegin;
138a515125bSLeila Ghaffari 
139a515125bSLeila Ghaffari   // Read reference file
1402b916ea7SJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
141e7754af5SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
142e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
143a515125bSLeila Ghaffari 
144a515125bSLeila Ghaffari   // Compute error with respect to reference solution
1452b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1462b916ea7SJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1472b916ea7SJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1482b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
149a515125bSLeila Ghaffari 
150a515125bSLeila Ghaffari   // Check error
151a515125bSLeila Ghaffari   if (error > app_ctx->test_tol) {
1522b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
153a515125bSLeila Ghaffari   }
154a515125bSLeila Ghaffari 
155a515125bSLeila Ghaffari   // Cleanup
1562b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1572b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Qref));
158a515125bSLeila Ghaffari 
159a515125bSLeila Ghaffari   PetscFunctionReturn(0);
160a515125bSLeila Ghaffari }
161a515125bSLeila Ghaffari 
162a515125bSLeila Ghaffari // Get error for problems with exact solutions
1632b916ea7SJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
164a515125bSLeila Ghaffari   PetscInt  loc_nodes;
165a515125bSLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
166a515125bSLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
167a515125bSLeila Ghaffari   PetscFunctionBegin;
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari   // Get exact solution at final time
1702b916ea7SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1712b916ea7SJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1722b916ea7SJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1732b916ea7SJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
174a515125bSLeila Ghaffari 
175a515125bSLeila Ghaffari   // Get |exact solution - obtained solution|
1762b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1772b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1782b916ea7SJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
179a515125bSLeila Ghaffari 
180a515125bSLeila Ghaffari   // Compute relative error
181a515125bSLeila Ghaffari   rel_error = norm_error / norm_exact;
182a515125bSLeila Ghaffari 
183a515125bSLeila Ghaffari   // Output relative error
1842b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
185a515125bSLeila Ghaffari   // Cleanup
1862b916ea7SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1872b916ea7SJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
188a515125bSLeila Ghaffari 
189a515125bSLeila Ghaffari   PetscFunctionReturn(0);
190a515125bSLeila Ghaffari }
191a515125bSLeila Ghaffari 
192a515125bSLeila Ghaffari // Post-processing
1932b916ea7SJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
194a515125bSLeila Ghaffari   PetscInt          steps;
195f0784ed3SJed Brown   TSConvergedReason reason;
196a515125bSLeila Ghaffari   PetscFunctionBegin;
197a515125bSLeila Ghaffari 
198a515125bSLeila Ghaffari   // Print relative error
1990e1e9333SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
2002b916ea7SJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
201a515125bSLeila Ghaffari   }
202a515125bSLeila Ghaffari 
203a515125bSLeila Ghaffari   // Print final time and number of steps
2042b916ea7SJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
205f0784ed3SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
2060e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
207f0784ed3SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
208f0784ed3SJed Brown                           steps, (double)final_time));
209a515125bSLeila Ghaffari   }
210a515125bSLeila Ghaffari 
211a515125bSLeila Ghaffari   // Output numerical values from command line
2122b916ea7SJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
213a515125bSLeila Ghaffari 
214a515125bSLeila Ghaffari   // Compare reference solution values with current test run for CI
2150e1e9333SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2162b916ea7SJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
217a515125bSLeila Ghaffari   }
218a515125bSLeila Ghaffari 
219a515125bSLeila Ghaffari   PetscFunctionReturn(0);
220a515125bSLeila Ghaffari }
221a515125bSLeila Ghaffari 
2229293eaa1SJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
2239293eaa1SJed Brown 
224a515125bSLeila Ghaffari // Gather initial Q values in case of continuation of simulation
225a515125bSLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
226a515125bSLeila Ghaffari   PetscViewer viewer;
2272b916ea7SJeremy L Thompson 
228a515125bSLeila Ghaffari   PetscFunctionBegin;
229a515125bSLeila Ghaffari 
2302b916ea7SJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
231e7754af5SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2322b916ea7SJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
233a515125bSLeila Ghaffari 
234a515125bSLeila Ghaffari   PetscFunctionReturn(0);
235a515125bSLeila Ghaffari }
236a515125bSLeila Ghaffari 
237a515125bSLeila Ghaffari // Record boundary values from initial condition
238a515125bSLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2399d437337SJames Wright   Vec Qbc, boundary_mask;
240a515125bSLeila Ghaffari   PetscFunctionBegin;
241a515125bSLeila Ghaffari 
2422b916ea7SJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2432b916ea7SJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2442b916ea7SJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2452b916ea7SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2462b916ea7SJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2472b916ea7SJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2482b916ea7SJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
249a515125bSLeila Ghaffari 
2509d437337SJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2519d437337SJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2529d437337SJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2539d437337SJames Wright   PetscCall(VecSet(Q, 1.0));
2549d437337SJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2559d437337SJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2569d437337SJames Wright 
257a515125bSLeila Ghaffari   PetscFunctionReturn(0);
258a515125bSLeila Ghaffari }
25915a3537eSJed Brown 
26015a3537eSJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
26115a3537eSJed Brown int FreeContextPetsc(void *data) {
2622b916ea7SJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
26315a3537eSJed Brown   return CEED_ERROR_SUCCESS;
26415a3537eSJed Brown }
2659f59f36eSJames Wright 
2669f59f36eSJames Wright // Return mass qfunction specification for number of components N
2679f59f36eSJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
2689f59f36eSJames Wright   CeedQFunctionUser qfunction_ptr;
2699f59f36eSJames Wright   const char       *qfunction_loc;
2709f59f36eSJames Wright   PetscFunctionBeginUser;
2719f59f36eSJames Wright 
2729f59f36eSJames Wright   switch (N) {
2739f59f36eSJames Wright     case 1:
2749f59f36eSJames Wright       qfunction_ptr = Mass_1;
2759f59f36eSJames Wright       qfunction_loc = Mass_1_loc;
2769f59f36eSJames Wright       break;
2779f59f36eSJames Wright     case 5:
2789f59f36eSJames Wright       qfunction_ptr = Mass_5;
2799f59f36eSJames Wright       qfunction_loc = Mass_5_loc;
2809f59f36eSJames Wright       break;
2819f59f36eSJames Wright     case 9:
2829f59f36eSJames Wright       qfunction_ptr = Mass_9;
2839f59f36eSJames Wright       qfunction_loc = Mass_9_loc;
2849f59f36eSJames Wright       break;
2859f59f36eSJames Wright     case 22:
2869f59f36eSJames Wright       qfunction_ptr = Mass_22;
2879f59f36eSJames Wright       qfunction_loc = Mass_22_loc;
2889f59f36eSJames Wright       break;
2899f59f36eSJames Wright     default:
2909f59f36eSJames Wright       SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N);
2919f59f36eSJames Wright   }
2929f59f36eSJames Wright   CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf);
2939f59f36eSJames Wright 
2949f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
2959f59f36eSJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
2969f59f36eSJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
2979f59f36eSJames Wright   PetscFunctionReturn(0);
2989f59f36eSJames Wright }
299e5e81594SJames Wright 
300e5e81594SJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space
301e5e81594SJames Wright  *
302e5e81594SJames Wright  * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum.
303e5e81594SJames Wright  *
304e5e81594SJames Wright  * @param[in]  source_vec    Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc
305e5e81594SJames Wright  * @param[out] target_vec    Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc
306e5e81594SJames Wright  * @param[in]  rhs_matop_ctx MatopApplyContext for performing the RHS evaluation
307e5e81594SJames Wright  * @param[in]  ksp           KSP for solving the consistent projection problem
308e5e81594SJames Wright  */
309*f7325489SJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, OperatorApplyContext rhs_matop_ctx, KSP ksp) {
310e5e81594SJames Wright   PetscFunctionBeginUser;
311e5e81594SJames Wright 
312*f7325489SJames Wright   PetscCall(ApplyCeedOperatorGlobalToGlobal(source_vec, target_vec, rhs_matop_ctx));
313e5e81594SJames Wright   PetscCall(KSPSolve(ksp, target_vec, target_vec));
314e5e81594SJames Wright 
315e5e81594SJames Wright   PetscFunctionReturn(0);
316e5e81594SJames Wright }
317