xref: /honee/src/diff_flux_projection.c (revision d9e6962164416f0bf8fc41cae34da27fa35469f8)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 /// @file
4 /// Functions for setting up and projecting the divergence of the diffusive flux
5 
6 #include "../qfunctions/diff_flux_projection.h"
7 
8 #include <petscdmplex.h>
9 
10 #include <navierstokes.h>
11 
12 const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL};
13 
14 /**
15   @brief Create `DivDiffFluxProjectionData` for solution DM in `honee`
16 
17   @param[in]  honee               `Honee` context
18   @param[in]  divFdiffproj_method Method used to perform the divergence of diffusive flux method (can be `DIV_DIFF_FLUX_PROJ_NONE`)
19   @param[in]  num_diff_flux_comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion)
20   @param[out] diff_flux_proj      The `DivDiffFluxProjectionData` object created, or set to `NULL` if `divFdiffproj_method = DIV_DIFF_FLUX_PROJ_NONE`
21 **/
22 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, DivDiffFluxProjectionMethod divFdiffproj_method, PetscInt num_diff_flux_comps,
23                                            DivDiffFluxProjectionData *diff_flux_proj) {
24   PetscInt                  label_value = 0, height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx->q_extra;
25   DMLabel                   domain_label = NULL;
26   DivDiffFluxProjectionData diff_flux_proj_;
27   NodalProjectionData       projection;
28 
29   PetscFunctionBeginUser;
30   if (divFdiffproj_method == DIV_DIFF_FLUX_PROJ_NONE) {
31     *diff_flux_proj = NULL;
32     PetscFunctionReturn(PETSC_SUCCESS);
33   }
34   PetscCall(PetscNew(&diff_flux_proj_));
35   PetscCall(PetscNew(&diff_flux_proj_->projection));
36   projection                           = diff_flux_proj_->projection;
37   diff_flux_proj_->method              = divFdiffproj_method;
38   diff_flux_proj_->num_diff_flux_comps = num_diff_flux_comps;
39 
40   PetscCall(DMClone(honee->dm, &projection->dm));
41   PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE));
42   PetscCall(DMGetDimension(projection->dm, &dim));
43   switch (diff_flux_proj_->method) {
44     case DIV_DIFF_FLUX_PROJ_DIRECT: {
45       projection->num_comp = diff_flux_proj_->num_diff_flux_comps;
46       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj"));
47       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
48 
49       PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, domain_label, label_value, height, dm_field,
50                                                 &diff_flux_proj_->elem_restr_div_diff_flux));
51       PetscCallCeed(honee->ceed,
52                     CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL));
53       PetscCall(CreateBasisFromPlex(honee->ceed, projection->dm, domain_label, label_value, height, dm_field, &diff_flux_proj_->basis_div_diff_flux));
54       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP;
55 
56       {  // Create face labels on projection->dm for boundary integrals
57         DMLabel  face_sets_label;
58         PetscInt num_face_set_values, *face_set_values;
59 
60         PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label));
61         PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_set_values));
62         for (PetscInt f = 0; f < num_face_set_values; f++) {
63           DMLabel face_orientation_label;
64           char   *face_orientation_label_name;
65 
66           PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name));
67           PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label));
68           PetscCall(DMAddLabel(projection->dm, face_orientation_label));
69           PetscCall(PetscFree(face_orientation_label_name));
70         }
71         PetscCall(PetscFree(face_set_values));
72       }
73     } break;
74     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
75       projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim;
76       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
77       PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm));
78 
79       PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, domain_label, label_value, height,
80                                                      diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux));
81       PetscCallCeed(honee->ceed,
82                     CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_diff_flux_ceed, NULL));
83       diff_flux_proj_->basis_div_diff_flux     = CEED_BASIS_NONE;
84       diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE;
85     } break;
86     case DIV_DIFF_FLUX_PROJ_NONE:
87       SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
88               DivDiffFluxProjectionMethods[divFdiffproj_method]);
89       break;
90   }
91   *diff_flux_proj = diff_flux_proj_;
92   PetscFunctionReturn(PETSC_SUCCESS);
93 };
94 
95 /**
96   @brief Return the objects required for the Divergence of Diffusive flux to be read by a `CeedOperator`
97 
98   @param[in]  diff_flux_proj Projection object
99   @param[out] elem_restr     Element restriction for the divergence of diffusive flux, or `NULL`
100   @param[out] basis          Basis for the divergence of diffusive flux, or `NULL`
101   @param[out] vector         Vector for the divergence of diffusive flux, or `NULL`
102   @param[out] eval_mode      Eval mode for the divergence of diffusive flux, or `NULL`
103 **/
104 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
105                                                          CeedVector *vector, CeedEvalMode *eval_mode) {
106   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
107 
108   PetscFunctionBeginUser;
109   if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_div_diff_flux, elem_restr));
110   if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis));
111   if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector));
112   if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux;
113   PetscFunctionReturn(PETSC_SUCCESS);
114 }
115 
116 /**
117   @brief Setup direct projection of divergence of diffusive flux
118 
119   @param[in]     honee          `Honee` context
120   @param[in,out] diff_flux_proj Flux projection object to setup
121 **/
122 static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
123   Ceed                ceed       = honee->ceed;
124   NodalProjectionData projection = diff_flux_proj->projection;
125   MPI_Comm            comm       = PetscObjectComm((PetscObject)projection->dm);
126 
127   PetscFunctionBeginUser;
128   {  // Create Projection RHS OperatorApplyContext
129     CeedOperator op_rhs;
130 
131     PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE,
132                "Must define CreateRHSOperator_Direct to use indirect div_diff_flux projection");
133     PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs));
134     PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
135     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
136     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc,
137                                          &projection->l2_rhs_ctx));
138     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
139   }
140 
141   {  // -- Build Mass operator
142     CeedQFunction       qf_mass;
143     CeedOperator        op_mass;
144     CeedBasis           basis_div_diff_flux             = NULL;
145     CeedElemRestriction elem_restr_div_diff_flux_volume = NULL, elem_restr_qd;
146     CeedVector          q_data;
147     CeedInt             q_data_size;
148     PetscInt            label_value  = 0;
149     DMLabel             domain_label = NULL;
150 
151     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux_volume, &basis_div_diff_flux, NULL, NULL));
152     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
153                        &q_data_size));
154 
155     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
156     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
157     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
158     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
159     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_div_diff_flux_volume, basis_div_diff_flux, CEED_VECTOR_ACTIVE));
160 
161     {  // -- Setup KSP for L^2 projection
162       Mat mat_mass;
163 
164       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
165 
166       PetscCall(KSPCreate(comm, &projection->ksp));
167       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
168       {  // lumped by default
169         PC pc;
170         PetscCall(KSPGetPC(projection->ksp, &pc));
171         PetscCall(PCSetType(pc, PCJACOBI));
172         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
173         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
174       }
175       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
176       PetscCall(MatDestroy(&mat_mass));
177     }
178 
179     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux_volume));
180     PetscCallCeed(ceed, CeedBasisDestroy(&basis_div_diff_flux));
181     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
182     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
183     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
184     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
185   }
186   PetscFunctionReturn(PETSC_SUCCESS);
187 }
188 
189 /**
190   @brief Setup indirect projection of divergence of diffusive flux
191 
192   @param[in]     honee          `Honee` context
193   @param[in,out] diff_flux_proj Flux projection object to setup
194 **/
195 static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
196   Ceed                ceed       = honee->ceed;
197   NodalProjectionData projection = diff_flux_proj->projection;
198   CeedBasis           basis_diff_flux;
199   CeedElemRestriction elem_restr_diff_flux, elem_restr_qd;
200   CeedVector          q_data;
201   CeedInt             q_data_size;
202   MPI_Comm            comm = PetscObjectComm((PetscObject)projection->dm);
203 
204   PetscFunctionBeginUser;
205   {
206     PetscInt label_value = 0, height = 0, dm_field = 0;
207     DMLabel  domain_label = NULL;
208 
209     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
210     PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
211     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
212                        &q_data_size));
213   }
214 
215   {
216     CeedOperator op_rhs;
217 
218     PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE,
219                "Must define CreateRHSOperator_Indirect to use indirect div_diff_flux projection");
220     PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs));
221     PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
222     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
223   }
224 
225   {  // -- Build Mass operator
226     CeedQFunction qf_mass;
227     CeedOperator  op_mass;
228 
229     PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass));
230     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
231     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
232     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
233     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
234 
235     {  // -- Setup KSP for L^2 projection
236       Mat mat_mass;
237 
238       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
239 
240       PetscCall(KSPCreate(comm, &projection->ksp));
241       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
242       {  // lumped by default
243         PC pc;
244         PetscCall(KSPGetPC(projection->ksp, &pc));
245         PetscCall(PCSetType(pc, PCJACOBI));
246         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
247         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
248       }
249       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
250       PetscCall(MatDestroy(&mat_mass));
251     }
252     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
253     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
254   }
255 
256   {  // Create OperatorApplyContext to calculate divergence at quadrature points
257     CeedQFunction       qf_calc_divergence = NULL;
258     CeedOperator        op_calc_divergence;
259     CeedElemRestriction elem_restr_div_diff_flux = NULL;
260     PetscInt            dim;
261 
262     PetscCall(DMGetDimension(projection->dm, &dim));
263     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_div_diff_flux, NULL, NULL, NULL));
264 
265     switch (dim) {
266       case 2:
267         switch (diff_flux_proj->num_diff_flux_comps) {
268           case 1:
269             PetscCallCeed(ceed,
270                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux2D_1, ComputeDivDiffusiveFlux2D_1_loc, &qf_calc_divergence));
271             break;
272         }
273         break;
274       case 3:
275         switch (diff_flux_proj->num_diff_flux_comps) {
276           case 1:
277             PetscCallCeed(ceed,
278                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_1, ComputeDivDiffusiveFlux3D_1_loc, &qf_calc_divergence));
279             break;
280           case 4:
281             PetscCallCeed(ceed,
282                           CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
283             break;
284         }
285         break;
286     }
287     PetscCheck(qf_calc_divergence, comm, PETSC_ERR_SUP,
288                "QFunction for calculating divergence of diffusive flux does not exist for"
289                " %" PetscInt_FMT " dimensional grid and %" PetscInt_FMT
290                " number of components.\nA new qfunction can be easily added; see source code for pattern.",
291                dim, diff_flux_proj->num_diff_flux_comps);
292 
293     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
294     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
295     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_diff_flux_comps, CEED_EVAL_NONE));
296 
297     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
298     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
299     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
300     PetscCallCeed(
301         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
302 
303     PetscCall(
304         OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL, &diff_flux_proj->calc_div_diff_flux));
305 
306     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_div_diff_flux));
307     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
308     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
309   }
310   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
311   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
312   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
313   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
314   PetscFunctionReturn(PETSC_SUCCESS);
315 }
316 
317 /**
318   @brief Setup projection of divergence of diffusive flux
319 
320   @param[in]     honee          `Honee` context
321   @param[in,out] diff_flux_proj Flux projection object to setup
322 **/
323 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj) {
324   PetscFunctionBeginUser;
325   switch (honee->app_ctx->divFdiffproj_method) {
326     case DIV_DIFF_FLUX_PROJ_DIRECT:
327       PetscCall(DivDiffFluxProjectionSetup_Direct(honee, diff_flux_proj));
328       break;
329     case DIV_DIFF_FLUX_PROJ_INDIRECT:
330       PetscCall(DivDiffFluxProjectionSetup_Indirect(honee, diff_flux_proj));
331       break;
332     case DIV_DIFF_FLUX_PROJ_NONE:
333       SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
334               DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
335       break;
336   }
337   PetscFunctionReturn(PETSC_SUCCESS);
338 }
339 
340 /**
341   @brief Project the divergence of diffusive flux
342 
343   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
344 
345   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
346   @param[in]  Q_loc          Localized solution vector
347 **/
348 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
349   NodalProjectionData projection = diff_flux_proj->projection;
350 
351   PetscFunctionBeginUser;
352   PetscCall(PetscLogEventBegin(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
353   switch (diff_flux_proj->method) {
354     case DIV_DIFF_FLUX_PROJ_DIRECT: {
355       Vec DivDiffFlux, RHS;
356 
357       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
358       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
359       if (diff_flux_proj->ceed_vec_has_array) {
360         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
361         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
362       }
363       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
364       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
365 
366       {
367         // Run PCApply manually if using ksp_type preonly -pc_type jacobi
368         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
369         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
370         PC        pc;
371         PetscBool ispreonly, isjacobi;
372         PetscCall(KSPGetPC(projection->ksp, &pc));
373         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
374         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
375         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DivDiffFlux));
376         else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux));
377       }
378       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
379 
380       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
381       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
382       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
383 
384       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
385       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
386       break;
387     }
388     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
389       Vec DiffFlux, RHS;
390 
391       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
392       PetscCall(DMGetGlobalVector(projection->dm, &RHS));
393       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx));
394       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
395 
396       {
397         // Run PCApply manually if using -ksp_type preonly -pc_type jacobi
398         // This is to avoid an AllReduce call in KSPSolve_Preonly, which causes significant slowdowns for lumped mass matrix solves.
399         // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix
400         PC        pc;
401         PetscBool ispreonly, isjacobi;
402         PetscCall(KSPGetPC(projection->ksp, &pc));
403         PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly));
404         PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &isjacobi));
405         if (ispreonly && isjacobi) PetscCall(PCApply(pc, RHS, DiffFlux));
406         else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux));
407       }
408       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
409 
410       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
411       PetscCall(DMRestoreGlobalVector(projection->dm, &RHS));
412       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
413     } break;
414     case DIV_DIFF_FLUX_PROJ_NONE:
415       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
416               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
417       break;
418   }
419   PetscCall(PetscLogEventEnd(HONEE_DivDiffFluxProjection, Q_loc, 0, 0, 0));
420   PetscFunctionReturn(PETSC_SUCCESS);
421 }
422 
423 /**
424   @brief Destroy `DivDiffFluxProjectionData` object
425 
426   @param[in,out] diff_flux_proj Object to destroy
427 **/
428 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
429   PetscFunctionBeginUser;
430   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
431   Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed);
432 
433   PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection));
434   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
435   if (diff_flux_proj->ceed_vec_has_array) {
436     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
437     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
438   }
439   PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
440   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux));
441   PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux));
442   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
443   PetscCall(PetscFree(diff_flux_proj));
444   PetscFunctionReturn(PETSC_SUCCESS);
445 }
446