xref: /libCEED/examples/fluids/src/misc.c (revision 530ad8c4fcac036810a1c230b0a577d30882de6e)
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 
1177841947SLeila Ghaffari #include "../navierstokes.h"
12ef080ff9SJames Wright #include "../qfunctions/mass.h"
1377841947SLeila Ghaffari 
142b730f8bSJeremy L Thompson PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time) {
1577841947SLeila Ghaffari   PetscFunctionBeginUser;
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari   // ---------------------------------------------------------------------------
18a0add3c9SJed Brown   // Update time for evaluation
1977841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2017b0d5c6SJeremy L Thompson   if (user->phys->ics_time_label) CeedOperatorSetContextDouble(ceed_data->op_ics, user->phys->ics_time_label, &time);
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2377841947SLeila Ghaffari   // ICs
2477841947SLeila Ghaffari   // ---------------------------------------------------------------------------
2577841947SLeila Ghaffari   // -- CEED Restriction
2677841947SLeila Ghaffari   CeedVector q0_ceed;
2777841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &q0_ceed, NULL);
2877841947SLeila Ghaffari 
2977841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
3077841947SLeila Ghaffari   PetscMemType q0_mem_type;
31c798d55aSJames Wright   PetscCall(VecP2C(Q_loc, &q0_mem_type, q0_ceed));
3277841947SLeila Ghaffari 
3377841947SLeila Ghaffari   // -- Apply CEED Operator
342b730f8bSJeremy L Thompson   CeedOperatorApply(ceed_data->op_ics, ceed_data->x_coord, q0_ceed, CEED_REQUEST_IMMEDIATE);
3577841947SLeila Ghaffari 
3677841947SLeila Ghaffari   // -- Restore vectors
37c798d55aSJames Wright   PetscCall(VecC2P(q0_ceed, q0_mem_type, Q_loc));
3877841947SLeila Ghaffari 
3977841947SLeila Ghaffari   // -- Local-to-Global
402b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q));
412b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, Q_loc, ADD_VALUES, Q));
4277841947SLeila Ghaffari 
4377841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4477841947SLeila Ghaffari   // Fix multiplicity for output of ICs
4577841947SLeila Ghaffari   // ---------------------------------------------------------------------------
4677841947SLeila Ghaffari   // -- CEED Restriction
4777841947SLeila Ghaffari   CeedVector mult_vec;
4877841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &mult_vec, NULL);
4977841947SLeila Ghaffari 
5077841947SLeila Ghaffari   // -- Place PETSc vector in CEED vector
5177841947SLeila Ghaffari   PetscMemType m_mem_type;
5277841947SLeila Ghaffari   Vec          multiplicity_loc;
532b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &multiplicity_loc));
54c798d55aSJames Wright   PetscCall(VecP2C(multiplicity_loc, &m_mem_type, mult_vec));
5577841947SLeila Ghaffari 
5677841947SLeila Ghaffari   // -- Get multiplicity
5777841947SLeila Ghaffari   CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, mult_vec);
5877841947SLeila Ghaffari 
5977841947SLeila Ghaffari   // -- Restore vectors
60c798d55aSJames Wright   PetscCall(VecC2P(mult_vec, m_mem_type, multiplicity_loc));
6177841947SLeila Ghaffari 
6277841947SLeila Ghaffari   // -- Local-to-Global
6377841947SLeila Ghaffari   Vec multiplicity;
642b730f8bSJeremy L Thompson   PetscCall(DMGetGlobalVector(dm, &multiplicity));
652b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(multiplicity));
662b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, multiplicity_loc, ADD_VALUES, multiplicity));
6777841947SLeila Ghaffari 
6877841947SLeila Ghaffari   // -- Fix multiplicity
692b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q, Q, multiplicity));
702b730f8bSJeremy L Thompson   PetscCall(VecPointwiseDivide(Q_loc, Q_loc, multiplicity_loc));
7177841947SLeila Ghaffari 
7277841947SLeila Ghaffari   // -- Restore vectors
732b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &multiplicity_loc));
742b730f8bSJeremy L Thompson   PetscCall(DMRestoreGlobalVector(dm, &multiplicity));
7577841947SLeila Ghaffari 
7677841947SLeila Ghaffari   // Cleanup
7777841947SLeila Ghaffari   CeedVectorDestroy(&mult_vec);
7877841947SLeila Ghaffari   CeedVectorDestroy(&q0_ceed);
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari   PetscFunctionReturn(0);
8177841947SLeila Ghaffari }
8277841947SLeila Ghaffari 
832b730f8bSJeremy L Thompson PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
842b730f8bSJeremy L Thompson                                              Vec grad_FVM) {
855571c6fdSJames Wright   Vec Qbc, boundary_mask;
8677841947SLeila Ghaffari   PetscFunctionBegin;
8777841947SLeila Ghaffari 
885571c6fdSJames Wright   // Mask (zero) Dirichlet entries
895571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
905571c6fdSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
915571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
925571c6fdSJames Wright 
932b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
942b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q_loc, 1., Qbc));
952b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
9677841947SLeila Ghaffari 
9777841947SLeila Ghaffari   PetscFunctionReturn(0);
9877841947SLeila Ghaffari }
9977841947SLeila Ghaffari 
100*530ad8c4SKenneth E. Jansen // @brief Load vector from binary file, possibly with embedded solution time and step number
101*530ad8c4SKenneth E. Jansen PetscErrorCode LoadFluidsBinaryVec(MPI_Comm comm, PetscViewer viewer, Vec Q, PetscReal *time, PetscInt *step_number) {
102*530ad8c4SKenneth E. Jansen   PetscInt  token, file_step_number;
103*530ad8c4SKenneth E. Jansen   PetscReal file_time;
104*530ad8c4SKenneth E. Jansen   PetscFunctionBeginUser;
105*530ad8c4SKenneth E. Jansen 
106*530ad8c4SKenneth E. Jansen   // Attempt
107*530ad8c4SKenneth E. Jansen   PetscCall(PetscViewerBinaryRead(viewer, &token, 1, NULL, PETSC_INT));
108*530ad8c4SKenneth E. Jansen   if (token == FLUIDS_FILE_TOKEN) {  // New style format; we're reading a file with step number and time in the header
109*530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_step_number, 1, NULL, PETSC_INT));
110*530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &file_time, 1, NULL, PETSC_REAL));
111*530ad8c4SKenneth E. Jansen     if (time) *time = file_time;
112*530ad8c4SKenneth E. Jansen     if (step_number) *step_number = file_step_number;
113*530ad8c4SKenneth E. Jansen   } else if (token == VEC_FILE_CLASSID) {  // Legacy format of just the vector, encoded as [VEC_FILE_CLASSID, length, ]
114*530ad8c4SKenneth E. Jansen     PetscInt length, N;
115*530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinaryRead(viewer, &length, 1, NULL, PETSC_INT));
116*530ad8c4SKenneth E. Jansen     PetscCall(VecGetSize(Q, &N));
117*530ad8c4SKenneth 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);
118*530ad8c4SKenneth E. Jansen     PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
119*530ad8c4SKenneth E. Jansen   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Not a fluids header token or a PETSc Vec in file");
120*530ad8c4SKenneth E. Jansen 
121*530ad8c4SKenneth E. Jansen   // Load Q from existent solution
122*530ad8c4SKenneth E. Jansen   PetscCall(VecLoad(Q, viewer));
123*530ad8c4SKenneth E. Jansen 
124*530ad8c4SKenneth E. Jansen   PetscFunctionReturn(0);
125*530ad8c4SKenneth E. Jansen }
126*530ad8c4SKenneth E. Jansen 
12777841947SLeila Ghaffari // Compare reference solution values with current test run for CI
12877841947SLeila Ghaffari PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q) {
12977841947SLeila Ghaffari   Vec         Qref;
13077841947SLeila Ghaffari   PetscViewer viewer;
13177841947SLeila Ghaffari   PetscReal   error, Qrefnorm;
132*530ad8c4SKenneth E. Jansen   MPI_Comm    comm = PetscObjectComm((PetscObject)Q);
13377841947SLeila Ghaffari   PetscFunctionBegin;
13477841947SLeila Ghaffari 
13577841947SLeila Ghaffari   // Read reference file
1362b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(Q, &Qref));
137*530ad8c4SKenneth E. Jansen   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->test_file_path, FILE_MODE_READ, &viewer));
138*530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Qref, NULL, NULL));
13977841947SLeila Ghaffari 
14077841947SLeila Ghaffari   // Compute error with respect to reference solution
1412b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Qref));
1422b730f8bSJeremy L Thompson   PetscCall(VecNorm(Qref, NORM_MAX, &Qrefnorm));
1432b730f8bSJeremy L Thompson   PetscCall(VecScale(Q, 1. / Qrefnorm));
1442b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_MAX, &error));
14577841947SLeila Ghaffari 
14677841947SLeila Ghaffari   // Check error
14777841947SLeila Ghaffari   if (error > app_ctx->test_tol) {
1482b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test failed with error norm %g\n", (double)error));
14977841947SLeila Ghaffari   }
15077841947SLeila Ghaffari 
15177841947SLeila Ghaffari   // Cleanup
1522b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
1532b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Qref));
15477841947SLeila Ghaffari 
15577841947SLeila Ghaffari   PetscFunctionReturn(0);
15677841947SLeila Ghaffari }
15777841947SLeila Ghaffari 
15877841947SLeila Ghaffari // Get error for problems with exact solutions
1592b730f8bSJeremy L Thompson PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time) {
16077841947SLeila Ghaffari   PetscInt  loc_nodes;
16177841947SLeila Ghaffari   Vec       Q_exact, Q_exact_loc;
16277841947SLeila Ghaffari   PetscReal rel_error, norm_error, norm_exact;
16377841947SLeila Ghaffari   PetscFunctionBegin;
16477841947SLeila Ghaffari 
16577841947SLeila Ghaffari   // Get exact solution at final time
1662b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &Q_exact));
1672b730f8bSJeremy L Thompson   PetscCall(DMGetLocalVector(dm, &Q_exact_loc));
1682b730f8bSJeremy L Thompson   PetscCall(VecGetSize(Q_exact_loc, &loc_nodes));
1692b730f8bSJeremy L Thompson   PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, Q_exact_loc, Q_exact, final_time));
17077841947SLeila Ghaffari 
17177841947SLeila Ghaffari   // Get |exact solution - obtained solution|
1722b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q_exact, NORM_1, &norm_exact));
1732b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Q, -1.0, Q_exact));
1742b730f8bSJeremy L Thompson   PetscCall(VecNorm(Q, NORM_1, &norm_error));
17577841947SLeila Ghaffari 
17677841947SLeila Ghaffari   // Compute relative error
17777841947SLeila Ghaffari   rel_error = norm_error / norm_exact;
17877841947SLeila Ghaffari 
17977841947SLeila Ghaffari   // Output relative error
1802b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Relative Error: %g\n", (double)rel_error));
18177841947SLeila Ghaffari   // Cleanup
1822b730f8bSJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm, &Q_exact_loc));
1832b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&Q_exact));
18477841947SLeila Ghaffari 
18577841947SLeila Ghaffari   PetscFunctionReturn(0);
18677841947SLeila Ghaffari }
18777841947SLeila Ghaffari 
18877841947SLeila Ghaffari // Post-processing
1892b730f8bSJeremy L Thompson PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time) {
19077841947SLeila Ghaffari   PetscInt          steps;
191cf7a0454SJed Brown   TSConvergedReason reason;
19277841947SLeila Ghaffari   PetscFunctionBegin;
19377841947SLeila Ghaffari 
19477841947SLeila Ghaffari   // Print relative error
1958fb33541SJames Wright   if (problem->non_zero_time && user->app_ctx->test_type == TESTTYPE_NONE) {
1962b730f8bSJeremy L Thompson     PetscCall(GetError_NS(ceed_data, dm, user, Q, final_time));
19777841947SLeila Ghaffari   }
19877841947SLeila Ghaffari 
19977841947SLeila Ghaffari   // Print final time and number of steps
2002b730f8bSJeremy L Thompson   PetscCall(TSGetStepNumber(ts, &steps));
201cf7a0454SJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
2028fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_NONE) {
203cf7a0454SJed Brown     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Time integrator %s on time step %" PetscInt_FMT " with final time %g\n", TSConvergedReasons[reason],
204cf7a0454SJed Brown                           steps, (double)final_time));
20577841947SLeila Ghaffari   }
20677841947SLeila Ghaffari 
20777841947SLeila Ghaffari   // Output numerical values from command line
2082b730f8bSJeremy L Thompson   PetscCall(VecViewFromOptions(Q, NULL, "-vec_view"));
20977841947SLeila Ghaffari 
21077841947SLeila Ghaffari   // Compare reference solution values with current test run for CI
2118fb33541SJames Wright   if (user->app_ctx->test_type == TESTTYPE_SOLVER) {
2122b730f8bSJeremy L Thompson     PetscCall(RegressionTests_NS(user->app_ctx, Q));
21377841947SLeila Ghaffari   }
21477841947SLeila Ghaffari 
21577841947SLeila Ghaffari   PetscFunctionReturn(0);
21677841947SLeila Ghaffari }
21777841947SLeila Ghaffari 
2184de8550aSJed Brown const PetscInt FLUIDS_FILE_TOKEN = 0xceedf00;
2194de8550aSJed Brown 
22077841947SLeila Ghaffari // Gather initial Q values in case of continuation of simulation
22177841947SLeila Ghaffari PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q) {
22277841947SLeila Ghaffari   PetscViewer viewer;
2232b730f8bSJeremy L Thompson 
22477841947SLeila Ghaffari   PetscFunctionBegin;
22577841947SLeila Ghaffari 
2262b730f8bSJeremy L Thompson   PetscCall(PetscViewerBinaryOpen(comm, app_ctx->cont_file, FILE_MODE_READ, &viewer));
227*530ad8c4SKenneth E. Jansen   PetscCall(LoadFluidsBinaryVec(comm, viewer, Q, &app_ctx->cont_time, &app_ctx->cont_steps));
2282b730f8bSJeremy L Thompson   PetscCall(PetscViewerDestroy(&viewer));
22977841947SLeila Ghaffari 
23077841947SLeila Ghaffari   PetscFunctionReturn(0);
23177841947SLeila Ghaffari }
23277841947SLeila Ghaffari 
23377841947SLeila Ghaffari // Record boundary values from initial condition
23477841947SLeila Ghaffari PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc) {
2355571c6fdSJames Wright   Vec Qbc, boundary_mask;
23677841947SLeila Ghaffari   PetscFunctionBegin;
23777841947SLeila Ghaffari 
2382b730f8bSJeremy L Thompson   PetscCall(DMGetNamedLocalVector(dm, "Qbc", &Qbc));
2392b730f8bSJeremy L Thompson   PetscCall(VecCopy(Q_loc, Qbc));
2402b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(Q_loc));
2412b730f8bSJeremy L Thompson   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, Q_loc));
2422b730f8bSJeremy L Thompson   PetscCall(VecAXPY(Qbc, -1., Q_loc));
2432b730f8bSJeremy L Thompson   PetscCall(DMRestoreNamedLocalVector(dm, "Qbc", &Qbc));
2442b730f8bSJeremy L Thompson   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS));
24577841947SLeila Ghaffari 
2465571c6fdSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2475571c6fdSJames Wright   PetscCall(DMGetGlobalVector(dm, &Q));
2485571c6fdSJames Wright   PetscCall(VecZeroEntries(boundary_mask));
2495571c6fdSJames Wright   PetscCall(VecSet(Q, 1.0));
2505571c6fdSJames Wright   PetscCall(DMGlobalToLocal(dm, Q, INSERT_VALUES, boundary_mask));
2515571c6fdSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2525571c6fdSJames Wright 
25377841947SLeila Ghaffari   PetscFunctionReturn(0);
25477841947SLeila Ghaffari }
255841e4c73SJed Brown 
256841e4c73SJed Brown // Free a plain data context that was allocated using PETSc; returning libCEED error codes
257841e4c73SJed Brown int FreeContextPetsc(void *data) {
2582b730f8bSJeremy L Thompson   if (PetscFree(data)) return CeedError(NULL, CEED_ERROR_ACCESS, "PetscFree failed");
259841e4c73SJed Brown   return CEED_ERROR_SUCCESS;
260841e4c73SJed Brown }
261ef080ff9SJames Wright 
262ef080ff9SJames Wright // Return mass qfunction specification for number of components N
263ef080ff9SJames Wright PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf) {
264ef080ff9SJames Wright   CeedQFunctionUser qfunction_ptr;
265ef080ff9SJames Wright   const char       *qfunction_loc;
266ef080ff9SJames Wright   PetscFunctionBeginUser;
267ef080ff9SJames Wright 
268ef080ff9SJames Wright   switch (N) {
269ef080ff9SJames Wright     case 1:
270ef080ff9SJames Wright       qfunction_ptr = Mass_1;
271ef080ff9SJames Wright       qfunction_loc = Mass_1_loc;
272ef080ff9SJames Wright       break;
273ef080ff9SJames Wright     case 5:
274ef080ff9SJames Wright       qfunction_ptr = Mass_5;
275ef080ff9SJames Wright       qfunction_loc = Mass_5_loc;
276ef080ff9SJames Wright       break;
277ef080ff9SJames Wright     case 9:
278ef080ff9SJames Wright       qfunction_ptr = Mass_9;
279ef080ff9SJames Wright       qfunction_loc = Mass_9_loc;
280ef080ff9SJames Wright       break;
281ef080ff9SJames Wright     case 22:
282ef080ff9SJames Wright       qfunction_ptr = Mass_22;
283ef080ff9SJames Wright       qfunction_loc = Mass_22_loc;
284ef080ff9SJames Wright       break;
285ef080ff9SJames Wright     default:
286ef080ff9SJames Wright       SETERRQ(PETSC_COMM_WORLD, -1, "Could not find mass qfunction of size %d", N);
287ef080ff9SJames Wright   }
288ef080ff9SJames Wright   CeedQFunctionCreateInterior(ceed, 1, qfunction_ptr, qfunction_loc, qf);
289ef080ff9SJames Wright 
290ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "u", N, CEED_EVAL_INTERP);
291ef080ff9SJames Wright   CeedQFunctionAddInput(*qf, "qdata", q_data_size, CEED_EVAL_NONE);
292ef080ff9SJames Wright   CeedQFunctionAddOutput(*qf, "v", N, CEED_EVAL_INTERP);
293ef080ff9SJames Wright   PetscFunctionReturn(0);
294ef080ff9SJames Wright }
2950df12fefSJames Wright 
2960df12fefSJames Wright /* @brief L^2 Projection of a source FEM function to a target FEM space
2970df12fefSJames Wright  *
2980df12fefSJames Wright  * To solve system using a lumped mass matrix, pass a KSP object with ksp_type=preonly, pc_type=jacobi, pc_jacobi_type=rowsum.
2990df12fefSJames Wright  *
3000df12fefSJames Wright  * @param[in]  source_vec    Global Vec of the source FEM function. NULL indicates using rhs_matop_ctx->X_loc
3010df12fefSJames Wright  * @param[out] target_vec    Global Vec of the target (result) FEM function. NULL indicates using rhs_matop_ctx->Y_loc
3020df12fefSJames Wright  * @param[in]  rhs_matop_ctx MatopApplyContext for performing the RHS evaluation
3030df12fefSJames Wright  * @param[in]  ksp           KSP for solving the consistent projection problem
3040df12fefSJames Wright  */
3050df12fefSJames Wright PetscErrorCode ComputeL2Projection(Vec source_vec, Vec target_vec, MatopApplyContext rhs_matop_ctx, KSP ksp) {
3060df12fefSJames Wright   PetscFunctionBeginUser;
3070df12fefSJames Wright 
3080df12fefSJames Wright   PetscCall(ApplyLocal_Ceed(source_vec, target_vec, rhs_matop_ctx));
3090df12fefSJames Wright   PetscCall(KSPSolve(ksp, target_vec, target_vec));
3100df12fefSJames Wright 
3110df12fefSJames Wright   PetscFunctionReturn(0);
3120df12fefSJames Wright }
313