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