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