xref: /honee/src/diff_flux_projection.c (revision 8c85b83518484fc9eaa5bccfe38a7cce91b34691)
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 Initialize projection of divergence of diffusive flux
14 
15   Creates underlying `DM` for the projection operation and creates the restriction and basis to use with the CeedOperator
16 
17   @param[in]  user                     `User` context
18   @param[out] elem_restr_div_diff_flux `CeedElemRestriction` of the divergence of diffusive flux vector
19   @param[out] basis_div_diff_flux      `CeedBasis` of the divergence of diffusive flux vector
20   @param[out] eval_mode_diff_flux      `CeedEvalMode` for the divergence of the diffusive flux
21 **/
22 PetscErrorCode DiffFluxProjectionInitialize(User user, CeedElemRestriction *elem_restr_div_diff_flux, CeedBasis *basis_div_diff_flux,
23                                             CeedEvalMode *eval_mode_diff_flux) {
24   PetscSection              section;
25   PetscInt                  label_value = 0, height = 0, dm_field = 0, dim;
26   DMLabel                   domain_label = NULL;
27   DivDiffFluxProjectionData diff_flux_proj;
28   NodalProjectionData       projection;
29 
30   PetscFunctionBeginUser;
31   PetscCall(PetscNew(&user->diff_flux_proj));
32   diff_flux_proj = user->diff_flux_proj;
33   PetscCall(PetscNew(&user->diff_flux_proj->projection));
34   projection                          = user->diff_flux_proj->projection;
35   diff_flux_proj->method              = user->app_ctx->divFdiffproj_method;
36   diff_flux_proj->num_diff_flux_comps = 4;
37 
38   PetscCall(DMClone(user->dm, &projection->dm));
39   PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE));
40   PetscCall(DMGetDimension(projection->dm, &dim));
41   switch (diff_flux_proj->method) {
42     case DIV_DIFF_FLUX_PROJ_DIRECT: {
43       projection->num_comp = diff_flux_proj->num_diff_flux_comps;
44       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj"));
45       PetscCall(
46           DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &projection->num_comp, projection->dm));
47 
48       PetscCall(DMGetLocalSection(projection->dm, &section));
49       PetscCall(PetscSectionSetFieldName(section, 0, ""));
50       PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_MomentumX"));
51       PetscCall(PetscSectionSetComponentName(section, 0, 1, "DivDiffusiveFlux_MomentumY"));
52       PetscCall(PetscSectionSetComponentName(section, 0, 2, "DivDiffusiveFlux_MomentumZ"));
53       PetscCall(PetscSectionSetComponentName(section, 0, 3, "DivDiffusiveFlux_Energy"));
54 
55       PetscCall(DMPlexCeedElemRestrictionCreate(user->ceed, projection->dm, domain_label, label_value, height, dm_field, elem_restr_div_diff_flux));
56       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(*elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
57       PetscCall(CreateBasisFromPlex(user->ceed, projection->dm, domain_label, label_value, height, dm_field, basis_div_diff_flux));
58       *eval_mode_diff_flux = CEED_EVAL_INTERP;
59     } break;
60     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
61       projection->num_comp = diff_flux_proj->num_diff_flux_comps * dim;
62       PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj"));
63       PetscCall(
64           DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &projection->num_comp, projection->dm));
65 
66       PetscCall(DMGetLocalSection(projection->dm, &section));
67       PetscCall(PetscSectionSetFieldName(section, 0, ""));
68       PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_MomentumXX"));
69       PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_MomentumXY"));
70       PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_MomentumXZ"));
71       PetscCall(PetscSectionSetComponentName(section, 0, 3, "DiffusiveFlux_MomentumYX"));
72       PetscCall(PetscSectionSetComponentName(section, 0, 4, "DiffusiveFlux_MomentumYY"));
73       PetscCall(PetscSectionSetComponentName(section, 0, 5, "DiffusiveFlux_MomentumYZ"));
74       PetscCall(PetscSectionSetComponentName(section, 0, 6, "DiffusiveFlux_MomentumZX"));
75       PetscCall(PetscSectionSetComponentName(section, 0, 7, "DiffusiveFlux_MomentumZY"));
76       PetscCall(PetscSectionSetComponentName(section, 0, 8, "DiffusiveFlux_MomentumZZ"));
77       PetscCall(PetscSectionSetComponentName(section, 0, 9, "DiffusiveFlux_EnergyX"));
78       PetscCall(PetscSectionSetComponentName(section, 0, 10, "DiffusiveFlux_EnergyY"));
79       PetscCall(PetscSectionSetComponentName(section, 0, 11, "DiffusiveFlux_EnergyZ"));
80 
81       PetscCall(DMPlexCeedElemRestrictionQDataCreate(user->ceed, projection->dm, domain_label, label_value, height,
82                                                      diff_flux_proj->num_diff_flux_comps, elem_restr_div_diff_flux));
83       PetscCallCeed(user->ceed, CeedElemRestrictionCreateVector(*elem_restr_div_diff_flux, &diff_flux_proj->div_diff_flux_ceed, NULL));
84       *basis_div_diff_flux = CEED_BASIS_NONE;
85       *eval_mode_diff_flux = CEED_EVAL_NONE;
86     } break;
87     case DIV_DIFF_FLUX_PROJ_NONE:
88       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
89               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
90       break;
91   }
92   PetscFunctionReturn(PETSC_SUCCESS);
93 };
94 
95 /**
96   @brief Setup direct projection of divergence of diffusive flux
97 
98   @param[in] ceed      `Ceed` context
99   @param[in] user      `User` context
100   @param[in] ceed_data `CeedData` context
101   @param[in] problem   `ProblemData` context
102 **/
103 static PetscErrorCode DivDiffFluxProjectionSetup_Direct(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
104   DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj;
105   NodalProjectionData       projection     = diff_flux_proj->projection;
106   CeedOperator              op_rhs;
107   CeedBasis                 basis_diff_flux;
108   CeedElemRestriction       elem_restr_diff_flux_volume, elem_restr_qd;
109   CeedVector                q_data;
110   CeedInt                   num_comp_q, q_data_size;
111   PetscInt                  dim, label_value = 0;
112   DMLabel                   domain_label = NULL;
113 
114   PetscFunctionBeginUser;
115   // -- Get Pre-requisite things
116   PetscCall(DMGetDimension(projection->dm, &dim));
117   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
118 
119   {  // Get elem_restr_diff_flux and basis_diff_flux
120     CeedOperator     *sub_ops;
121     CeedOperatorField op_field;
122     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
123 
124     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops));
125     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "div F_diff", &op_field));
126     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_diff_flux_volume));
127     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_diff_flux));
128   }
129   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &elem_restr_qd,
130                      &q_data, &q_data_size));
131 
132   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
133   {  // Add the volume integral CeedOperator
134     CeedQFunction qf_rhs_volume;
135     CeedOperator  op_rhs_volume;
136 
137     switch (user->phys->state_var) {
138       case STATEVAR_PRIMITIVE:
139         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Prim, DivDiffusiveFluxVolumeRHS_Prim_loc, &qf_rhs_volume));
140         break;
141       case STATEVAR_CONSERVATIVE:
142         PetscCallCeed(ceed,
143                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Conserv, DivDiffusiveFluxVolumeRHS_Conserv_loc, &qf_rhs_volume));
144         break;
145       case STATEVAR_ENTROPY:
146         PetscCallCeed(ceed,
147                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_Entropy, DivDiffusiveFluxVolumeRHS_Entropy_loc, &qf_rhs_volume));
148         break;
149     }
150 
151     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, problem->apply_vol_ifunction.qfctx));
152     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "q", num_comp_q, CEED_EVAL_INTERP));
153     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
154     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
155     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
156 
157     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
158     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
159     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
160     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
161     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
162 
163     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_volume));
164 
165     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
166     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
167   }
168 
169   {  // Add the boundary integral CeedOperator
170     CeedQFunction qf_rhs_boundary;
171     DMLabel       face_sets_label;
172     PetscInt      num_face_set_values, *face_set_values;
173     CeedInt       q_data_size;
174 
175     // -- Build RHS operator
176     switch (user->phys->state_var) {
177       case STATEVAR_PRIMITIVE:
178         PetscCallCeed(ceed,
179                       CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Prim, DivDiffusiveFluxBoundaryRHS_Prim_loc, &qf_rhs_boundary));
180         break;
181       case STATEVAR_CONSERVATIVE:
182         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Conserv, DivDiffusiveFluxBoundaryRHS_Conserv_loc,
183                                                         &qf_rhs_boundary));
184         break;
185       case STATEVAR_ENTROPY:
186         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_Entropy, DivDiffusiveFluxBoundaryRHS_Entropy_loc,
187                                                         &qf_rhs_boundary));
188         break;
189     }
190 
191     PetscCall(QDataBoundaryGradientGetNumComponents(user->dm, &q_data_size));
192     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, problem->apply_vol_ifunction.qfctx));
193     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "q", num_comp_q, CEED_EVAL_INTERP));
194     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
195     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
196     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
197 
198     PetscCall(DMGetLabel(user->dm, "Face Sets", &face_sets_label));
199     PetscCall(DMLabelCreateGlobalValueArray(user->dm, face_sets_label, &num_face_set_values, &face_set_values));
200     for (PetscInt f = 0; f < num_face_set_values; f++) {
201       DMLabel  face_orientation_label;
202       PetscInt num_orientations_values, *orientation_values;
203 
204       {
205         char *face_orientation_label_name;
206 
207         PetscCall(DMPlexCreateFaceLabel(user->dm, face_set_values[f], &face_orientation_label_name));
208         PetscCall(DMGetLabel(user->dm, face_orientation_label_name, &face_orientation_label));
209         PetscCall(DMAddLabel(projection->dm, face_orientation_label));
210         PetscCall(PetscFree(face_orientation_label_name));
211       }
212       PetscCall(DMLabelCreateGlobalValueArray(user->dm, face_orientation_label, &num_orientations_values, &orientation_values));
213       for (PetscInt o = 0; o < num_orientations_values; o++) {
214         CeedOperator        op_rhs_boundary;
215         CeedBasis           basis_q, basis_diff_flux_boundary;
216         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
217         CeedVector          q_data;
218         CeedInt             q_data_size;
219         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
220 
221         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, user->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
222         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
223                                                   &elem_restr_diff_flux_boundary));
224         PetscCall(QDataBoundaryGradientGet(ceed, user->dm, face_orientation_label, orientation, ceed_data->x_coord, &elem_restr_qdata, &q_data,
225                                            &q_data_size));
226         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, user->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
227         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
228 
229         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
230         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
231         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
232         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
233         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
234                                                  CEED_VECTOR_ACTIVE));
235 
236         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_boundary));
237 
238         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
239         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
240         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
241         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
242         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
243         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
244         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
245       }
246       PetscCall(PetscFree(orientation_values));
247     }
248     PetscCall(PetscFree(face_set_values));
249     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
250   }
251 
252   PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc));
253   diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
254   PetscCall(
255       OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, diff_flux_proj->DivDiffFlux_loc, &projection->l2_rhs_ctx));
256 
257   {  // -- Build Mass operator
258     CeedQFunction qf_mass;
259     CeedOperator  op_mass;
260 
261     PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass));
262     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
263     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
264     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
265     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
266 
267     {  // -- Setup KSP for L^2 projection
268       Mat      mat_mass;
269       MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm);
270 
271       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
272 
273       PetscCall(KSPCreate(comm, &projection->ksp));
274       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
275       {  // lumped by default
276         PC pc;
277         PetscCall(KSPGetPC(projection->ksp, &pc));
278         PetscCall(PCSetType(pc, PCJACOBI));
279         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
280         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
281       }
282       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
283       PetscCall(MatDestroy(&mat_mass));
284     }
285     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
286     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
287   }
288   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
289   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
290   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
291   PetscFunctionReturn(PETSC_SUCCESS);
292 }
293 
294 /**
295   @brief Setup indirect projection of divergence of diffusive flux
296 
297   @param[in]     ceed      `Ceed` context
298   @param[in,out] user      `User` context
299   @param[in]     ceed_data `CeedData` context
300   @param[in]     problem   `ProblemData` context
301 **/
302 static PetscErrorCode DivDiffFluxProjectionSetup_Indirect(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
303   DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj;
304   NodalProjectionData       projection     = diff_flux_proj->projection;
305   CeedBasis                 basis_diff_flux;
306   CeedElemRestriction       elem_restr_diff_flux, elem_restr_qd;
307   CeedVector                q_data;
308   CeedInt                   num_comp_q, q_data_size;
309   PetscInt                  dim;
310   PetscInt                  label_value = 0, height = 0, dm_field = 0;
311   DMLabel                   domain_label = NULL;
312 
313   PetscFunctionBeginUser;
314   PetscCall(DMGetDimension(projection->dm, &dim));
315   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
316 
317   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
318   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
319   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &elem_restr_qd,
320                      &q_data, &q_data_size));
321 
322   {  // Create RHS CeedOperator for L^2 projection
323     CeedQFunction qf_rhs;
324     CeedOperator  op_rhs;
325 
326     switch (user->phys->state_var) {
327       case STATEVAR_PRIMITIVE:
328         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Prim, DiffusiveFluxRHS_Prim_loc, &qf_rhs));
329         break;
330       case STATEVAR_CONSERVATIVE:
331         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Conserv, DiffusiveFluxRHS_Conserv_loc, &qf_rhs));
332         break;
333       case STATEVAR_ENTROPY:
334         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_Entropy, DiffusiveFluxRHS_Entropy_loc, &qf_rhs));
335         break;
336     }
337 
338     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, problem->apply_vol_ifunction.qfctx));
339     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "q", num_comp_q, CEED_EVAL_INTERP));
340     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
341     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
342     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
343 
344     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, &op_rhs));
345     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
346     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
347     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
348     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
349 
350     PetscCall(OperatorApplyContextCreate(user->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL, &projection->l2_rhs_ctx));
351 
352     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
353     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
354   }
355 
356   {  // -- Build Mass operator
357     CeedQFunction qf_mass;
358     CeedOperator  op_mass;
359 
360     PetscCall(CreateMassQFunction(ceed, projection->num_comp, q_data_size, &qf_mass));
361     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
362     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
363     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
364     PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
365 
366     {  // -- Setup KSP for L^2 projection
367       Mat      mat_mass;
368       MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm);
369 
370       PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass));
371 
372       PetscCall(KSPCreate(comm, &projection->ksp));
373       PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_"));
374       {  // lumped by default
375         PC pc;
376         PetscCall(KSPGetPC(projection->ksp, &pc));
377         PetscCall(PCSetType(pc, PCJACOBI));
378         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
379         PetscCall(KSPSetType(projection->ksp, KSPPREONLY));
380       }
381       PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass));
382       PetscCall(MatDestroy(&mat_mass));
383     }
384     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
385     PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
386   }
387 
388   {  // Create OperatorApplyContext to calculate divergence at quadrature points
389     CeedQFunction       qf_calc_divergence;
390     CeedOperator        op_calc_divergence;
391     CeedElemRestriction elem_restr_div_diff_flux;
392 
393     {  // Get elem_restr_div_diff_flux
394       CeedOperator     *sub_ops;
395       CeedOperatorField op_field;
396       PetscInt          sub_op_index = 0;  // will be 0 for the volume op
397 
398       PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops));
399       PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "div F_diff", &op_field));
400       PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_div_diff_flux));
401     }
402 
403     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeDivDiffusiveFlux3D_4, ComputeDivDiffusiveFlux3D_4_loc, &qf_calc_divergence));
404 
405     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp * dim, CEED_EVAL_GRAD));
406     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "qdata", q_data_size, CEED_EVAL_NONE));
407     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", 4, CEED_EVAL_NONE));
408 
409     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_calc_divergence, NULL, NULL, &op_calc_divergence));
410     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "Grad F_diff", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
411     PetscCallCeed(ceed, CeedOperatorSetField(op_calc_divergence, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
412     PetscCallCeed(
413         ceed, CeedOperatorSetField(op_calc_divergence, "Div F_diff", elem_restr_div_diff_flux, CEED_BASIS_NONE, diff_flux_proj->div_diff_flux_ceed));
414 
415     PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, NULL, NULL, NULL,
416                                          &user->diff_flux_proj->calc_div_diff_flux));
417 
418     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_calc_divergence));
419     PetscCallCeed(ceed, CeedOperatorDestroy(&op_calc_divergence));
420   }
421   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
422   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
423   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
424   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
425   PetscFunctionReturn(PETSC_SUCCESS);
426 }
427 
428 /**
429   @brief Setup projection of divergence of diffusive flux
430 
431   @param[in]     ceed      `Ceed` context
432   @param[in,out] user      `User` context
433   @param[in]     ceed_data `CeedData` context
434   @param[in]     problem   `ProblemData` context
435 **/
436 PetscErrorCode DivDiffFluxProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
437   PetscFunctionBeginUser;
438   switch (user->app_ctx->divFdiffproj_method) {
439     case DIV_DIFF_FLUX_PROJ_DIRECT:
440       PetscCall(DivDiffFluxProjectionSetup_Direct(ceed, user, ceed_data, problem));
441       break;
442     case DIV_DIFF_FLUX_PROJ_INDIRECT:
443       PetscCall(DivDiffFluxProjectionSetup_Indirect(ceed, user, ceed_data, problem));
444       break;
445     case DIV_DIFF_FLUX_PROJ_NONE:
446       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
447               DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
448       break;
449   }
450   PetscFunctionReturn(PETSC_SUCCESS);
451 }
452 
453 /**
454   @brief Project the divergence of diffusive flux
455 
456   This implicitly sets the `CeedVector` input (`div_diff_flux_ceed`) to the divergence of diffusive flux.
457 
458   @param[in]  diff_flux_proj `NodalProjectionData` for the projection
459   @param[in]  Q_loc          Localized solution vector
460 **/
461 PetscErrorCode DiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc) {
462   NodalProjectionData projection = diff_flux_proj->projection;
463 
464   PetscFunctionBeginUser;
465   PetscCall(PetscLogEventBegin(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
466   switch (diff_flux_proj->method) {
467     case DIV_DIFF_FLUX_PROJ_DIRECT: {
468       Vec DivDiffFlux;
469 
470       PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux));
471       if (diff_flux_proj->ceed_vec_has_array) {
472         PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
473         diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
474       }
475       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DivDiffFlux, projection->l2_rhs_ctx));
476       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
477 
478       PetscCall(KSPSolve(projection->ksp, DivDiffFlux, DivDiffFlux));
479       PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view"));
480 
481       PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_loc));
482       PetscCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->div_diff_flux_ceed));
483       diff_flux_proj->ceed_vec_has_array = PETSC_TRUE;
484 
485       PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux));
486       break;
487     }
488     case DIV_DIFF_FLUX_PROJ_INDIRECT: {
489       Vec DiffFlux;
490 
491       PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux));
492       PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, DiffFlux, projection->l2_rhs_ctx));
493       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view"));
494 
495       PetscCall(KSPSolve(projection->ksp, DiffFlux, DiffFlux));
496       PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view"));
497 
498       PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux));
499       PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux));
500     } break;
501     case DIV_DIFF_FLUX_PROJ_NONE:
502       SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
503               DivDiffFluxProjectionMethods[diff_flux_proj->method]);
504       break;
505   }
506   PetscCall(PetscLogEventEnd(FLUIDS_DivDiffFluxProjection, Q_loc, 0, 0, 0));
507   PetscFunctionReturn(PETSC_SUCCESS);
508 }
509 
510 /**
511   @brief Destroy `DivDiffFluxProjectionData` object
512 
513   @param[in,out] diff_flux_proj Object to destroy
514 **/
515 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj) {
516   PetscFunctionBeginUser;
517   if (diff_flux_proj == NULL) PetscFunctionReturn(PETSC_SUCCESS);
518   PetscCall(NodalProjectionDataDestroy(diff_flux_proj->projection));
519   PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux));
520   if (diff_flux_proj->ceed_vec_has_array) {
521     PetscCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, diff_flux_proj->DivDiffFlux_loc));
522     diff_flux_proj->ceed_vec_has_array = PETSC_FALSE;
523   }
524   PetscCall(CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed));
525   PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc));
526   PetscCall(PetscFree(diff_flux_proj));
527   PetscFunctionReturn(PETSC_SUCCESS);
528 }
529