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