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