xref: /honee/problems/advection.c (revision cf8f12c89a11a3279f97d0071ccf240f989bc25c)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Utility functions for setting up ADVECTION
6 
7 #include "../qfunctions/advection.h"
8 
9 #include <ceed.h>
10 #include <petscdm.h>
11 
12 #include <navierstokes.h>
13 
14 const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL};
15 const char *const AdvDifICTypes[]   = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL};
16 const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL};
17 const char *const AdvDifBubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_",
18                                                    NULL};
19 const char *const StabilizationTauTypes[]       = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL};
20 
21 static PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) {
22   MPI_Comm         comm = honee->comm;
23   Ceed             ceed = honee->ceed;
24   SetupContextAdv  setup_ctx;
25   AdvectionContext advection_ctx;
26   PetscInt         dim;
27 
28   PetscFunctionBeginUser;
29   PetscCall(DMGetDimension(honee->dm, &dim));
30   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx));
31   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
32   PetscCall(PetscPrintf(comm,
33                         "  Problem:\n"
34                         "    Problem Name                       : %s\n"
35                         "    Stabilization                      : %s\n"
36                         "    Stabilization Tau                  : %s\n"
37                         "    Wind Type                          : %s\n",
38                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization],
39                         StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type]));
40 
41   if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) {
42     CeedScalar *wind = setup_ctx->wind;
43     switch (dim) {
44       case 2:
45         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", wind[0], wind[1]));
46         break;
47       case 3:
48         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", wind[0], wind[1], wind[2]));
49         break;
50     }
51   }
52 
53   PetscCall(PetscPrintf(comm, "    Initial Condition Type             : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type]));
54   switch (setup_ctx->initial_condition_type) {
55     case ADVDIF_IC_SKEW:
56     case ADVDIF_IC_COSINE_HILL:
57     case ADVDIF_IC_BOUNDARY_LAYER:
58       break;
59     case ADVDIF_IC_BUBBLE_SPHERE:
60     case ADVDIF_IC_BUBBLE_CYLINDER:
61       PetscCall(PetscPrintf(comm, "    Bubble Continuity                  : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type]));
62       break;
63     case ADVDIF_IC_WAVE:
64       PetscCall(PetscPrintf(comm, "    Wave Type                          : %s\n", AdvDifWaveTypes[setup_ctx->wave_type]));
65       break;
66   }
67 
68   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx));
69   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
70   PetscFunctionReturn(PETSC_SUCCESS);
71 }
72 
73 // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
74 //
75 // Only used for SUPG stabilization
76 PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) {
77   Ceed                 ceed = honee->ceed;
78   CeedInt              num_comp_q, q_data_size;
79   CeedQFunction        qf_mass = NULL;
80   CeedElemRestriction  elem_restr_q, elem_restr_qd;
81   CeedBasis            basis_q;
82   CeedVector           q_data;
83   CeedQFunctionContext qfctx = NULL;
84   PetscInt             dim;
85 
86   PetscFunctionBeginUser;
87   PetscCall(DMGetDimension(honee->dm, &dim));
88   {  // Get restriction and basis from the RHS function
89     CeedOperator     *sub_ops;
90     CeedOperatorField op_field;
91     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
92 
93     PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
94     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
95     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
96     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
97     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
98 
99     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
100   }
101 
102   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
103   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
104 
105   switch (dim) {
106     case 2:
107       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
108       break;
109     case 3:
110       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
111       break;
112   }
113 
114   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
115   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
116   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
117   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
118   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
119   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
120   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
121 
122   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
123   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized"));
124   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
125   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
126   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
127   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
128   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
129 
130   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
131   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
132   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
133   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
134   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
135   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
136   PetscFunctionReturn(PETSC_SUCCESS);
137 }
138 
139 /**
140   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
141 
142   @param[in]  honee          `Honee` context
143   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
144   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
145 **/
146 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
147   Ceed                 ceed       = honee->ceed;
148   NodalProjectionData  projection = diff_flux_proj->projection;
149   PetscInt             dim, num_comp_q;
150   CeedQFunctionContext advection_qfctx = NULL;
151 
152   PetscFunctionBeginUser;
153   // -- Get Pre-requisite things
154   PetscCall(DMGetDimension(projection->dm, &dim));
155   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
156 
157   {  // Get advection-diffusion QF context
158     CeedOperator *sub_ops;
159     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
160 
161     if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
162     else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
163     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
164   }
165   PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs));
166   {  // Add the volume integral CeedOperator
167     CeedQFunction       qf_rhs_volume = NULL;
168     CeedOperator        op_rhs_volume;
169     CeedVector          q_data;
170     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL, elem_restr_q;
171     CeedBasis           basis_diff_flux = NULL, basis_q;
172     CeedInt             q_data_size;
173 
174     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
175     PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
176     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
177     PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord,
178                        &elem_restr_qd, &q_data, &q_data_size));
179     switch (dim) {
180       case 2:
181         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc,
182                                                         &qf_rhs_volume));
183         break;
184       case 3:
185         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc,
186                                                         &qf_rhs_volume));
187         break;
188     }
189     PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
190 
191     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx));
192     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
193     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
194     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
195 
196     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
197     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
198     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
199     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
200 
201     PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume));
202 
203     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
204     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
205     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
206     PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
207     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
208     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
209     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
210     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
211   }
212 
213   {  // Add the boundary integral CeedOperator
214     CeedQFunction qf_rhs_boundary;
215     DMLabel       face_sets_label;
216     PetscInt      num_face_set_values, *face_set_values;
217     CeedInt       q_data_size;
218 
219     // -- Build RHS operator
220     switch (dim) {
221       case 2:
222         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc,
223                                                         &qf_rhs_boundary));
224         break;
225       case 3:
226         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc,
227                                                         &qf_rhs_boundary));
228         break;
229     }
230 
231     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
232     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx));
233     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
234     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
235     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
236 
237     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
238     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
239     for (PetscInt f = 0; f < num_face_set_values; f++) {
240       DMLabel  face_orientation_label;
241       PetscInt num_orientations_values, *orientation_values;
242 
243       {
244         char *face_orientation_label_name;
245 
246         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
247         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
248         PetscCall(PetscFree(face_orientation_label_name));
249       }
250       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
251       for (PetscInt o = 0; o < num_orientations_values; o++) {
252         CeedOperator        op_rhs_boundary;
253         CeedBasis           basis_q, basis_diff_flux_boundary;
254         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
255         CeedVector          q_data;
256         CeedInt             q_data_size;
257         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
258 
259         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
260         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
261         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
262                                                   &elem_restr_diff_flux_boundary));
263         PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
264         PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data,
265                                            &q_data_size));
266 
267         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
268         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
269         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
270         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
271                                                  CEED_VECTOR_ACTIVE));
272 
273         PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary));
274 
275         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
276         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
277         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
278         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
279         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
280         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
281         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
282       }
283       PetscCall(PetscFree(orientation_values));
284     }
285     PetscCall(PetscFree(face_set_values));
286     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
287   }
288 
289   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
290   PetscFunctionReturn(PETSC_SUCCESS);
291 }
292 
293 /**
294   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
295 
296   @param[in]  honee          `Honee` context
297   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
298   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
299 **/
300 static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
301   Ceed                 ceed       = honee->ceed;
302   NodalProjectionData  projection = diff_flux_proj->projection;
303   CeedBasis            basis_diff_flux, basis_q;
304   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd, elem_restr_q;
305   CeedVector           q_data;
306   CeedInt              q_data_size;
307   PetscInt             dim, num_comp_q;
308   PetscInt             height = 0, dm_field = 0;
309   CeedQFunction        qf_rhs          = NULL;
310   CeedQFunctionContext advection_qfctx = NULL;
311 
312   PetscFunctionBeginUser;
313   PetscCall(DMGetDimension(projection->dm, &dim));
314   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
315 
316   {  // Get advection-diffusion QF context
317     CeedOperator *sub_ops;
318     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
319 
320     if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
321     else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
322     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
323   }
324   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
325   PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
326   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux));
327   PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux));
328   PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord,
329                      &elem_restr_qd, &q_data, &q_data_size));
330 
331   switch (dim) {
332     case 2:
333       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs));
334       break;
335     case 3:
336       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs));
337       break;
338   }
339   PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
340 
341   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx));
342   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
343   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
344   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
345 
346   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
347   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
348   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
349   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
350 
351   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
352   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
353   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
354   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
355   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
356   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
357   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
358   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
359   PetscFunctionReturn(PETSC_SUCCESS);
360 }
361 
362 static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
363   HoneeBCStruct honee_bc;
364   DM            dm;
365   PetscInt      dim;
366 
367   PetscFunctionBeginUser;
368   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
369   PetscCall(BCDefinitionGetDM(bc_def, &dm));
370   PetscCall(DMGetDimension(dm, &dim));
371   switch (dim) {
372     case 2:
373       PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf));
374       break;
375     case 3:
376       PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf));
377       break;
378   }
379   PetscFunctionReturn(PETSC_SUCCESS);
380 }
381 
382 static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
383 
384 PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) {
385   AdvDifWindType             wind_type;
386   AdvDifICType               advectionic_type;
387   AdvDifBubbleContinuityType bubble_continuity_type = -1;
388   StabilizationType          stab;
389   StabilizationTauType       stab_tau;
390   SetupContextAdv            setup_context;
391   Honee                      honee = *(Honee *)ctx;
392   MPI_Comm                   comm  = honee->comm;
393   Ceed                       ceed  = honee->ceed;
394   PetscBool                  implicit;
395   AdvectionContext           advection_ctx;
396   CeedQFunctionContext       advection_qfctx, ics_qfctx;
397   PetscInt                   dim;
398 
399   PetscFunctionBeginUser;
400   PetscCall(PetscNew(&setup_context));
401   PetscCall(PetscNew(&advection_ctx));
402   PetscCall(DMGetDimension(dm, &dim));
403 
404   // ------------------------------------------------------
405   //             Create the QFunction context
406   // ------------------------------------------------------
407   CeedScalar     rc               = 1000.;  // m (Radius of bubble)
408   CeedScalar     CtauS            = 0.;     // dimensionless
409   PetscBool      strong_form      = PETSC_FALSE;
410   CeedScalar     E_wind           = 1.e6;  // J
411   CeedScalar     Ctau_a           = PetscPowScalarInt(honee->app_ctx->degree, 2);
412   CeedScalar     Ctau_d           = PetscPowScalarInt(honee->app_ctx->degree, 4);
413   CeedScalar     Ctau_t           = 0.;
414   PetscReal      wind[3]          = {1., 0, 0};  // m/s
415   CeedScalar     diffusion_coeff  = 0.;
416   CeedScalar     wave_frequency   = 2 * M_PI;
417   CeedScalar     wave_phase       = 0;
418   AdvDifWaveType wave_type        = -1;
419   PetscScalar    bl_height_factor = 1.;
420   PetscReal      domain_min[3], domain_max[3], domain_size[3] = {0.};
421   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
422   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
423 
424   // ------------------------------------------------------
425   //              Command line Options
426   // ------------------------------------------------------
427   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
428   // -- Physics
429   PetscBool translation;
430   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION),
431                              (PetscEnum *)&wind_type, &translation));
432   PetscInt  n = dim;
433   PetscBool user_wind;
434   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
435   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
436   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
437   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
438                              NULL));
439   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
440   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
441   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
442                              (PetscEnum *)&stab_tau, NULL));
443   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
444   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL));
445   PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL));
446   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
447   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes,
448                              (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
449   // IC-specific options
450   switch (advectionic_type) {
451     case ADVDIF_IC_WAVE:
452       PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL));
453       PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL));
454       PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL));
455       PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE),
456                                  (PetscEnum *)&wave_type, NULL));
457       PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL));
458       PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL));
459       break;
460     case ADVDIF_IC_BOUNDARY_LAYER:
461       PetscCall(PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor,
462                                    NULL));
463       break;
464     case ADVDIF_IC_BUBBLE_CYLINDER:
465     case ADVDIF_IC_BUBBLE_SPHERE:
466       PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL));
467       PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL));
468       PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
469       bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE;
470       PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes,
471                                  (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL));
472       break;
473     case ADVDIF_IC_SKEW:
474     case ADVDIF_IC_COSINE_HILL:
475       break;
476   }
477 
478   // -- Warnings
479   if (wind_type == ADVDIF_WIND_ROTATION && user_wind) {
480     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
481   }
482   if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) {
483     wind[2] = 0;
484     PetscCall(PetscPrintf(comm,
485                           "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
486   }
487   if (stab == STAB_NONE && CtauS != 0) {
488     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
489   }
490   PetscOptionsEnd();
491 
492   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
493 
494   // ------------------------------------------------------
495   //           Set up the QFunction contexts
496   // ------------------------------------------------------
497   // -- Scale variables to desired units
498   Units units = honee->units;
499   E_wind *= units->Joule;
500   rc = fabs(rc) * units->meter;
501   for (PetscInt i = 0; i < dim; i++) {
502     wind[i] *= (units->meter / units->second);
503     domain_size[i] *= units->meter;
504   }
505 
506   // -- Setup Context
507   setup_context->rc                     = rc;
508   setup_context->lx                     = domain_size[0];
509   setup_context->ly                     = domain_size[1];
510   setup_context->lz                     = dim == 3 ? domain_size[2] : 0.;
511   setup_context->wind[0]                = wind[0];
512   setup_context->wind[1]                = wind[1];
513   setup_context->wind[2]                = dim == 3 ? wind[2] : 0.;
514   setup_context->wind_type              = wind_type;
515   setup_context->initial_condition_type = advectionic_type;
516   setup_context->bubble_continuity_type = bubble_continuity_type;
517   setup_context->time                   = 0;
518   setup_context->wave_frequency         = wave_frequency;
519   setup_context->wave_phase             = wave_phase;
520   setup_context->wave_type              = wave_type;
521   setup_context->bl_height_factor       = bl_height_factor;
522 
523   // -- QFunction Context
524   honee->phys->implicit            = implicit;
525   advection_ctx->CtauS             = CtauS;
526   advection_ctx->E_wind            = E_wind;
527   advection_ctx->implicit          = implicit;
528   advection_ctx->strong_form       = strong_form;
529   advection_ctx->stabilization     = stab;
530   advection_ctx->stabilization_tau = stab_tau;
531   advection_ctx->Ctau_a            = Ctau_a;
532   advection_ctx->Ctau_d            = Ctau_d;
533   advection_ctx->Ctau_t            = Ctau_t;
534   advection_ctx->diffusion_coeff   = diffusion_coeff;
535   advection_ctx->divFdiff_method   = honee->app_ctx->divFdiffproj_method;
536 
537   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx));
538   PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
539   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc));
540 
541   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx));
542   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
543   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc));
544   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
545                                                          "Size of timestep, delta t"));
546 
547   // ------------------------------------------------------
548   //               SET UP ADVECTION
549   // ------------------------------------------------------
550   problem->print_info         = PRINT_ADVECTION;
551   problem->num_comps_jac_data = 0;
552   switch (dim) {
553     case 2:
554       problem->ics                          = (HoneeQFSpec){.qf_func_ptr = ICsAdvection2d, .qf_loc = ICsAdvection2d_loc};
555       problem->apply_vol_rhs                = (HoneeQFSpec){.qf_func_ptr = RHS_Advection2d, .qf_loc = RHS_Advection2d_loc};
556       problem->apply_vol_ifunction          = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection2d, .qf_loc = IFunction_Advection2d_loc};
557       problem->compute_exact_solution_error = PETSC_TRUE;
558       break;
559     case 3:
560       problem->ics                          = (HoneeQFSpec){.qf_func_ptr = ICsAdvection, .qf_loc = ICsAdvection_loc};
561       problem->apply_vol_rhs                = (HoneeQFSpec){.qf_func_ptr = RHS_Advection, .qf_loc = RHS_Advection_loc};
562       problem->apply_vol_ifunction          = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection, .qf_loc = IFunction_Advection_loc};
563       problem->compute_exact_solution_error = PETSC_FALSE;
564       break;
565   }
566   problem->ics.qfctx           = ics_qfctx;
567   problem->apply_vol_rhs.qfctx = advection_qfctx;
568   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx));
569 
570   problem->num_components = 5;
571   PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
572   for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i]));
573 
574   PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj));
575   if (honee->diff_flux_proj) {
576     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
577     NodalProjectionData       projection     = diff_flux_proj->projection;
578     PetscSection              section;
579 
580     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_AdvDif;
581     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif;
582     PetscCall(DMGetLocalSection(projection->dm, &section));
583     switch (honee->diff_flux_proj->method) {
584       case DIV_DIFF_FLUX_PROJ_DIRECT: {
585         PetscCall(PetscSectionSetFieldName(section, 0, ""));
586         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar"));
587       } break;
588       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
589         PetscCall(PetscSectionSetFieldName(section, 0, ""));
590         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX"));
591         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY"));
592         if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ"));
593       } break;
594       case DIV_DIFF_FLUX_PROJ_NONE:
595         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
596                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
597         break;
598     }
599   }
600 
601   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
602     BCDefinition bc_def = problem->bc_defs[b];
603     const char  *name;
604 
605     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
606     if (!strcmp(name, "inflow")) {
607       HoneeBCStruct honee_bc;
608 
609       PetscCall(PetscNew(&honee_bc));
610       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx));
611       honee_bc->honee              = honee;
612       honee_bc->num_comps_jac_data = 0;
613       PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
614 
615       PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
616       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
617     }
618   }
619   PetscFunctionReturn(PETSC_SUCCESS);
620 }
621