xref: /honee/problems/advection.c (revision 337840fc9f8b699b98947f09e30443f8c0c7f962)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Utility functions for setting up ADVECTION
6a515125bSLeila Ghaffari 
72b916ea7SJeremy L Thompson #include "../qfunctions/advection.h"
82b916ea7SJeremy L Thompson 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson 
12149fb536SJames Wright #include <navierstokes.h>
13a515125bSLeila Ghaffari 
14e5a8cae0SJames Wright const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL};
15e5a8cae0SJames Wright const char *const AdvDifICTypes[]   = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "WAVE", "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_", NULL};
16e5a8cae0SJames Wright const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL};
17e5a8cae0SJames Wright const char *const AdvDifBubbleContinuityTypes[] = {"SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_",
18e5a8cae0SJames Wright                                                    NULL};
19e5a8cae0SJames Wright const char *const StabilizationTauTypes[]       = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL};
20e5a8cae0SJames Wright 
PRINT_ADVECTION(Honee honee,ProblemData problem,AppCtx app_ctx)219b05e62eSJames Wright static PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) {
229b05e62eSJames Wright   MPI_Comm         comm = honee->comm;
239b05e62eSJames Wright   Ceed             ceed = honee->ceed;
249b05e62eSJames Wright   SetupContextAdv  setup_ctx;
259b05e62eSJames Wright   AdvectionContext advection_ctx;
269b05e62eSJames Wright   PetscInt         dim;
279b05e62eSJames Wright 
289b05e62eSJames Wright   PetscFunctionBeginUser;
299b05e62eSJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
309b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx));
319b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
329b05e62eSJames Wright   PetscCall(PetscPrintf(comm,
339b05e62eSJames Wright                         "  Problem:\n"
349b05e62eSJames Wright                         "    Problem Name                       : %s\n"
359b05e62eSJames Wright                         "    Stabilization                      : %s\n"
369b05e62eSJames Wright                         "    Stabilization Tau                  : %s\n"
379b05e62eSJames Wright                         "    Wind Type                          : %s\n",
389b05e62eSJames Wright                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization],
399b05e62eSJames Wright                         StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type]));
409b05e62eSJames Wright 
419b05e62eSJames Wright   if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) {
429b05e62eSJames Wright     CeedScalar *wind = setup_ctx->wind;
439b05e62eSJames Wright     switch (dim) {
449b05e62eSJames Wright       case 2:
459b05e62eSJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", wind[0], wind[1]));
469b05e62eSJames Wright         break;
479b05e62eSJames Wright       case 3:
489b05e62eSJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", wind[0], wind[1], wind[2]));
499b05e62eSJames Wright         break;
509b05e62eSJames Wright     }
519b05e62eSJames Wright   }
529b05e62eSJames Wright 
539b05e62eSJames Wright   PetscCall(PetscPrintf(comm, "    Initial Condition Type             : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type]));
549b05e62eSJames Wright   switch (setup_ctx->initial_condition_type) {
559b05e62eSJames Wright     case ADVDIF_IC_SKEW:
569b05e62eSJames Wright     case ADVDIF_IC_COSINE_HILL:
579b05e62eSJames Wright     case ADVDIF_IC_BOUNDARY_LAYER:
589b05e62eSJames Wright       break;
599b05e62eSJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
609b05e62eSJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER:
619b05e62eSJames Wright       PetscCall(PetscPrintf(comm, "    Bubble Continuity                  : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type]));
629b05e62eSJames Wright       break;
639b05e62eSJames Wright     case ADVDIF_IC_WAVE:
649b05e62eSJames Wright       PetscCall(PetscPrintf(comm, "    Wave Type                          : %s\n", AdvDifWaveTypes[setup_ctx->wave_type]));
659b05e62eSJames Wright       break;
669b05e62eSJames Wright   }
679b05e62eSJames Wright 
689b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx));
699b05e62eSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
709b05e62eSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
719b05e62eSJames Wright }
729b05e62eSJames Wright 
73a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
74a78efa86SJames Wright //
75a78efa86SJames Wright // Only used for SUPG stabilization
CreateKSPMassOperator_AdvectionStabilized(Honee honee,CeedOperator * op_mass)760c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) {
770c373b74SJames Wright   Ceed                 ceed = honee->ceed;
78a78efa86SJames Wright   CeedInt              num_comp_q, q_data_size;
7987d3884fSJeremy L Thompson   CeedQFunction        qf_mass = NULL;
8001e19bfaSJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd;
81a78efa86SJames Wright   CeedBasis            basis_q;
82a78efa86SJames Wright   CeedVector           q_data;
83236a8ba6SJames Wright   CeedQFunctionContext qfctx = NULL;
84a78efa86SJames Wright   PetscInt             dim;
85a78efa86SJames Wright 
86a78efa86SJames Wright   PetscFunctionBeginUser;
870c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
88a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
89a78efa86SJames Wright     CeedOperator     *sub_ops;
9001e19bfaSJames Wright     CeedOperatorField op_field;
91a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
92a78efa86SJames Wright 
93da7f3a07SJames Wright     PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
9401e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
9501e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
9601e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
9701e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
98a78efa86SJames Wright 
99236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
100a78efa86SJames Wright   }
101a78efa86SJames Wright 
102a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
10301e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
104a78efa86SJames Wright 
105a78efa86SJames Wright   switch (dim) {
106a78efa86SJames Wright     case 2:
107a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
108a78efa86SJames Wright       break;
109a78efa86SJames Wright     case 3:
110a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
111a78efa86SJames Wright       break;
112a78efa86SJames Wright   }
113a78efa86SJames Wright 
114236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
115a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
116a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
117a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
118a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
119a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
120a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
121a78efa86SJames Wright 
122a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
12371f2ed29SJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized"));
124a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
1250c373b74SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
12601e19bfaSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
127a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
128a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
129a78efa86SJames Wright 
130fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
13101e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
132fff85bd3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
133fff85bd3SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
134236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
135a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
136a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
137a78efa86SJames Wright }
138a78efa86SJames Wright 
139c2d90829SJames Wright /**
140c2d90829SJames Wright   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
141c2d90829SJames Wright 
1420c373b74SJames Wright   @param[in]  honee          `Honee` context
143c2d90829SJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
144c2d90829SJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
145c2d90829SJames Wright **/
DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee,DivDiffFluxProjectionData diff_flux_proj,CeedOperator * op_rhs)146e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
1470c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
148c2d90829SJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
149cf8f12c8SJames Wright   PetscInt             dim, num_comp_q;
150c2d90829SJames Wright   CeedQFunctionContext advection_qfctx = NULL;
151c2d90829SJames Wright 
152c2d90829SJames Wright   PetscFunctionBeginUser;
153c2d90829SJames Wright   // -- Get Pre-requisite things
154c2d90829SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
155cf8f12c8SJames Wright   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
156c2d90829SJames Wright 
15740b78511SJames Wright   {  // Get advection-diffusion QF context
158c2d90829SJames Wright     CeedOperator *sub_ops;
159c2d90829SJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
160c2d90829SJames Wright 
161da7f3a07SJames Wright     if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
162da7f3a07SJames Wright     else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
163c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
164c2d90829SJames Wright   }
165da7f3a07SJames Wright   PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, op_rhs));
166c2d90829SJames Wright   {  // Add the volume integral CeedOperator
167c2d90829SJames Wright     CeedQFunction       qf_rhs_volume = NULL;
168c2d90829SJames Wright     CeedOperator        op_rhs_volume;
169c2d90829SJames Wright     CeedVector          q_data;
170cf8f12c8SJames Wright     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL, elem_restr_q;
171cf8f12c8SJames Wright     CeedBasis           basis_diff_flux = NULL, basis_q;
172c2d90829SJames Wright     CeedInt             q_data_size;
173c2d90829SJames Wright 
174cf8f12c8SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
175cf8f12c8SJames Wright     PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
176c2d90829SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
1779018c49aSJames Wright     PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
178c2d90829SJames Wright     switch (dim) {
179c2d90829SJames Wright       case 2:
1809eadbee4SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc,
1819eadbee4SJames Wright                                                         &qf_rhs_volume));
182c2d90829SJames Wright         break;
183c2d90829SJames Wright       case 3:
1849eadbee4SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc,
1859eadbee4SJames Wright                                                         &qf_rhs_volume));
186c2d90829SJames Wright         break;
187c2d90829SJames Wright     }
1880c373b74SJames Wright     PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
189c2d90829SJames Wright 
190c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx));
191c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
192c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
193c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
194c2d90829SJames Wright 
195c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
196cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
197c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
198c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
199c2d90829SJames Wright 
200da7f3a07SJames Wright     PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_volume));
201c2d90829SJames Wright 
202c2d90829SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
203c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
204cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
205cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
206c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
207c2d90829SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
208c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
209c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
210c2d90829SJames Wright   }
211c2d90829SJames Wright 
212c2d90829SJames Wright   {  // Add the boundary integral CeedOperator
213c2d90829SJames Wright     CeedQFunction qf_rhs_boundary;
214c2d90829SJames Wright     DMLabel       face_sets_label;
215c2d90829SJames Wright     PetscInt      num_face_set_values, *face_set_values;
216c2d90829SJames Wright     CeedInt       q_data_size;
217c2d90829SJames Wright 
218c2d90829SJames Wright     // -- Build RHS operator
219c2d90829SJames Wright     switch (dim) {
220c2d90829SJames Wright       case 2:
221c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc,
222c2d90829SJames Wright                                                         &qf_rhs_boundary));
223c2d90829SJames Wright         break;
224c2d90829SJames Wright       case 3:
225c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc,
226c2d90829SJames Wright                                                         &qf_rhs_boundary));
227c2d90829SJames Wright         break;
228c2d90829SJames Wright     }
229c2d90829SJames Wright 
2300c373b74SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
231c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx));
232c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
233c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
234c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
235c2d90829SJames Wright 
236c2d90829SJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
237c2d90829SJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
238c2d90829SJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
239c2d90829SJames Wright       DMLabel  face_orientation_label;
240c2d90829SJames Wright       PetscInt num_orientations_values, *orientation_values;
241c2d90829SJames Wright 
242c2d90829SJames Wright       {
243c2d90829SJames Wright         char *face_orientation_label_name;
244c2d90829SJames Wright 
245c2d90829SJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
246c2d90829SJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
247c2d90829SJames Wright         PetscCall(PetscFree(face_orientation_label_name));
248c2d90829SJames Wright       }
249c2d90829SJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
250c2d90829SJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
251c2d90829SJames Wright         CeedOperator        op_rhs_boundary;
252c2d90829SJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
253c2d90829SJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
254c2d90829SJames Wright         CeedVector          q_data;
255c2d90829SJames Wright         CeedInt             q_data_size;
256c2d90829SJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
257c2d90829SJames Wright 
2580c373b74SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
2590c373b74SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
260c2d90829SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
261c2d90829SJames Wright                                                   &elem_restr_diff_flux_boundary));
2626b9fd993SJames Wright         PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
263d510ce3cSJames Wright         PetscCall(QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, &elem_restr_qdata, &q_data, &q_data_size));
264c2d90829SJames Wright 
265c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
266c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
267c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
268c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
269c2d90829SJames Wright                                                  CEED_VECTOR_ACTIVE));
270c2d90829SJames Wright 
271da7f3a07SJames Wright         PetscCallCeed(ceed, CeedOperatorCompositeAddSub(*op_rhs, op_rhs_boundary));
272c2d90829SJames Wright 
273c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
274c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
275c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
276c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
277c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
278c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
279c2d90829SJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
280c2d90829SJames Wright       }
281c2d90829SJames Wright       PetscCall(PetscFree(orientation_values));
282c2d90829SJames Wright     }
283c2d90829SJames Wright     PetscCall(PetscFree(face_set_values));
284c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
285c2d90829SJames Wright   }
286c2d90829SJames Wright 
287c2d90829SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
288c2d90829SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
289c2d90829SJames Wright }
290c2d90829SJames Wright 
29140b78511SJames Wright /**
29240b78511SJames Wright   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
29340b78511SJames Wright 
2940c373b74SJames Wright   @param[in]  honee          `Honee` context
29540b78511SJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
29640b78511SJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
29740b78511SJames Wright **/
DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee,DivDiffFluxProjectionData diff_flux_proj,CeedOperator * op_rhs)298e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
2990c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
30040b78511SJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
301cf8f12c8SJames Wright   CeedBasis            basis_diff_flux, basis_q;
302cf8f12c8SJames Wright   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd, elem_restr_q;
30340b78511SJames Wright   CeedVector           q_data;
304cf8f12c8SJames Wright   CeedInt              q_data_size;
305cf8f12c8SJames Wright   PetscInt             dim, num_comp_q;
306e3db12f8SJames Wright   PetscInt             height = 0, dm_field = 0;
30740b78511SJames Wright   CeedQFunction        qf_rhs          = NULL;
30840b78511SJames Wright   CeedQFunctionContext advection_qfctx = NULL;
30940b78511SJames Wright 
31040b78511SJames Wright   PetscFunctionBeginUser;
31140b78511SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
312cf8f12c8SJames Wright   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
31340b78511SJames Wright 
31440b78511SJames Wright   {  // Get advection-diffusion QF context
31540b78511SJames Wright     CeedOperator *sub_ops;
31640b78511SJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
31740b78511SJames Wright 
318da7f3a07SJames Wright     if (honee->op_ifunction) PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_ifunction, &sub_ops));
319da7f3a07SJames Wright     else PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(honee->op_rhs_ctx->op, &sub_ops));
32040b78511SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
32140b78511SJames Wright   }
322cf8f12c8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
323cf8f12c8SJames Wright   PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
324e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_diff_flux));
325e3db12f8SJames Wright   PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_diff_flux));
3269018c49aSJames Wright   PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &q_data, &q_data_size));
32740b78511SJames Wright 
32840b78511SJames Wright   switch (dim) {
32940b78511SJames Wright     case 2:
33040b78511SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs));
33140b78511SJames Wright       break;
33240b78511SJames Wright     case 3:
33340b78511SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs));
33440b78511SJames Wright       break;
33540b78511SJames Wright   }
3360c373b74SJames Wright   PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
33740b78511SJames Wright 
33840b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx));
33940b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
34040b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
34140b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
34240b78511SJames Wright 
34340b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
344cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
34540b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
34640b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
34740b78511SJames Wright 
34840b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
34940b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
350cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
351cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
35240b78511SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
35340b78511SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
35440b78511SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
355cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
35640b78511SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
35740b78511SJames Wright }
35840b78511SJames Wright 
AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def,CeedQFunction * qf)359d6cac220SJames Wright static PetscErrorCode AdvectionInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
360d6cac220SJames Wright   HoneeBCStruct honee_bc;
361d6cac220SJames Wright   DM            dm;
362d6cac220SJames Wright   PetscInt      dim;
363d6cac220SJames Wright 
364d6cac220SJames Wright   PetscFunctionBeginUser;
365d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
366d6cac220SJames Wright   PetscCall(BCDefinitionGetDM(bc_def, &dm));
367d6cac220SJames Wright   PetscCall(DMGetDimension(dm, &dim));
368d6cac220SJames Wright   switch (dim) {
369d6cac220SJames Wright     case 2:
370d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection2d_InOutFlow, Advection2d_InOutFlow_loc, honee_bc->qfctx, qf));
371d6cac220SJames Wright       break;
372d6cac220SJames Wright     case 3:
373d6cac220SJames Wright       PetscCall(HoneeBCCreateIFunctionQF(bc_def, Advection_InOutFlow, Advection_InOutFlow_loc, honee_bc->qfctx, qf));
374d6cac220SJames Wright       break;
375d6cac220SJames Wright   }
376d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
377d6cac220SJames Wright }
378d6cac220SJames Wright 
379f27c2204SJames Wright static const char *const component_names[] = {"Density", "MomentumX", "MomentumY", "MomentumZ", "TotalEnergy"};
380f27c2204SJames Wright 
NS_ADVECTION(ProblemData problem,DM dm,void * ctx)381d3c60affSJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx) {
3825f636aeaSJames Wright   AdvDifWindType             wind_type;
3835f636aeaSJames Wright   AdvDifICType               advectionic_type;
384108c6689SJames Wright   AdvDifBubbleContinuityType bubble_continuity_type = -1;
385a515125bSLeila Ghaffari   StabilizationType          stab;
38657272ee0SJames Wright   StabilizationTauType       stab_tau;
3873636f6a4SJames Wright   SetupContextAdv            setup_context;
3880c373b74SJames Wright   Honee                      honee = *(Honee *)ctx;
3890c373b74SJames Wright   MPI_Comm                   comm  = honee->comm;
3900c373b74SJames Wright   Ceed                       ceed  = honee->ceed;
391a515125bSLeila Ghaffari   PetscBool                  implicit;
39215a3537eSJed Brown   AdvectionContext           advection_ctx;
393f5dc303cSJames Wright   CeedQFunctionContext       advection_qfctx, ics_qfctx;
3949529d636SJames Wright   PetscInt                   dim;
395a515125bSLeila Ghaffari 
39615a3537eSJed Brown   PetscFunctionBeginUser;
3972d898fa6SJames Wright   PetscCall(PetscNew(&setup_context));
3982d898fa6SJames Wright   PetscCall(PetscNew(&advection_ctx));
3999529d636SJames Wright   PetscCall(DMGetDimension(dm, &dim));
400a515125bSLeila Ghaffari 
401a515125bSLeila Ghaffari   // ------------------------------------------------------
402ea615d4cSJames Wright   //             Create the QFunction context
403a515125bSLeila Ghaffari   // ------------------------------------------------------
404a515125bSLeila Ghaffari   CeedScalar     rc               = 1000.;  // m (Radius of bubble)
405a515125bSLeila Ghaffari   CeedScalar     CtauS            = 0.;     // dimensionless
406059d1c40SJames Wright   PetscBool      strong_form      = PETSC_FALSE;
407a515125bSLeila Ghaffari   CeedScalar     E_wind           = 1.e6;  // J
4080c373b74SJames Wright   CeedScalar     Ctau_a           = PetscPowScalarInt(honee->app_ctx->degree, 2);
4090c373b74SJames Wright   CeedScalar     Ctau_d           = PetscPowScalarInt(honee->app_ctx->degree, 4);
41057272ee0SJames Wright   CeedScalar     Ctau_t           = 0.;
411a515125bSLeila Ghaffari   PetscReal      wind[3]          = {1., 0, 0};  // m/s
412c8d249deSJames Wright   CeedScalar     diffusion_coeff  = 0.;
413a62be6baSJames Wright   CeedScalar     wave_frequency   = 2 * M_PI;
414a62be6baSJames Wright   CeedScalar     wave_phase       = 0;
415108c6689SJames Wright   AdvDifWaveType wave_type        = -1;
416b4fd18dfSJames Wright   PetscScalar    bl_height_factor = 1.;
41787d3884fSJeremy L Thompson   PetscReal      domain_min[3], domain_max[3], domain_size[3] = {0.};
4182b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
41966d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
42005a512bdSLeila Ghaffari 
421a515125bSLeila Ghaffari   // ------------------------------------------------------
422a515125bSLeila Ghaffari   //              Command line Options
423a515125bSLeila Ghaffari   // ------------------------------------------------------
4241485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
425a515125bSLeila Ghaffari   // -- Physics
426a515125bSLeila Ghaffari   PetscBool translation;
4275f636aeaSJames Wright   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION),
4285f636aeaSJames Wright                              (PetscEnum *)&wind_type, &translation));
42966d54740SJames Wright   PetscInt  n = dim;
430a515125bSLeila Ghaffari   PetscBool user_wind;
4312b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
432c8d249deSJames Wright   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
4332b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
434059d1c40SJames Wright   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
435059d1c40SJames Wright                              NULL));
4362b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
4372b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
43857272ee0SJames Wright   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
43957272ee0SJames Wright                              (PetscEnum *)&stab_tau, NULL));
44057272ee0SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
441fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL));
442fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL));
4432b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
44480e9ac5bSJames Wright   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes,
44580e9ac5bSJames Wright                              (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
44657bb6b22SJames Wright   // IC-specific options
44757bb6b22SJames Wright   switch (advectionic_type) {
44857bb6b22SJames Wright     case ADVDIF_IC_WAVE:
44980e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL));
45080e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL));
45180e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL));
45280e9ac5bSJames Wright       PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE),
453108c6689SJames Wright                                  (PetscEnum *)&wave_type, NULL));
45480e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL));
45580e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL));
45657bb6b22SJames Wright       break;
45757bb6b22SJames Wright     case ADVDIF_IC_BOUNDARY_LAYER:
4589eadbee4SJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor,
4599eadbee4SJames Wright                                    NULL));
46057bb6b22SJames Wright       break;
46157bb6b22SJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER:
46257bb6b22SJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
46380e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL));
46480e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL));
46580e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
466108c6689SJames Wright       bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE;
46780e9ac5bSJames Wright       PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes,
468108c6689SJames Wright                                  (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL));
46957bb6b22SJames Wright       break;
47057bb6b22SJames Wright     case ADVDIF_IC_SKEW:
47157bb6b22SJames Wright     case ADVDIF_IC_COSINE_HILL:
47257bb6b22SJames Wright       break;
473108c6689SJames Wright   }
474a62be6baSJames Wright 
475a515125bSLeila Ghaffari   // -- Warnings
4765f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_ROTATION && user_wind) {
4772b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
478a515125bSLeila Ghaffari   }
4795f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) {
480a515125bSLeila Ghaffari     wind[2] = 0;
4819eadbee4SJames Wright     PetscCall(PetscPrintf(comm,
4829eadbee4SJames Wright                           "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
483a515125bSLeila Ghaffari   }
484a515125bSLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
4852b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
486a515125bSLeila Ghaffari   }
4871485969bSJeremy L Thompson   PetscOptionsEnd();
488a515125bSLeila Ghaffari 
489a78efa86SJames Wright   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
490a78efa86SJames Wright 
491a515125bSLeila Ghaffari   // ------------------------------------------------------
492ea615d4cSJames Wright   //           Set up the QFunction contexts
493a515125bSLeila Ghaffari   // ------------------------------------------------------
494a515125bSLeila Ghaffari   // -- Scale variables to desired units
495c9f37605SMohammed Amin   Units units = honee->units;
496c9f37605SMohammed Amin   E_wind *= units->Joule;
497c9f37605SMohammed Amin   rc = fabs(rc) * units->meter;
49866d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) {
499c9f37605SMohammed Amin     wind[i] *= (units->meter / units->second);
50005a512bdSLeila Ghaffari   }
501a515125bSLeila Ghaffari 
502a515125bSLeila Ghaffari   // -- Setup Context
503a515125bSLeila Ghaffari   setup_context->rc                     = rc;
50405a512bdSLeila Ghaffari   setup_context->lx                     = domain_size[0];
50505a512bdSLeila Ghaffari   setup_context->ly                     = domain_size[1];
50666d54740SJames Wright   setup_context->lz                     = dim == 3 ? domain_size[2] : 0.;
507a515125bSLeila Ghaffari   setup_context->wind[0]                = wind[0];
508a515125bSLeila Ghaffari   setup_context->wind[1]                = wind[1];
50966d54740SJames Wright   setup_context->wind[2]                = dim == 3 ? wind[2] : 0.;
510a515125bSLeila Ghaffari   setup_context->wind_type              = wind_type;
511c51f031aSJames Wright   setup_context->initial_condition_type = advectionic_type;
512a515125bSLeila Ghaffari   setup_context->bubble_continuity_type = bubble_continuity_type;
513a515125bSLeila Ghaffari   setup_context->time                   = 0;
514a62be6baSJames Wright   setup_context->wave_frequency         = wave_frequency;
515a62be6baSJames Wright   setup_context->wave_phase             = wave_phase;
516a62be6baSJames Wright   setup_context->wave_type              = wave_type;
517b4fd18dfSJames Wright   setup_context->bl_height_factor       = bl_height_factor;
518a515125bSLeila Ghaffari 
519a515125bSLeila Ghaffari   // -- QFunction Context
5200c373b74SJames Wright   honee->phys->implicit            = implicit;
52115a3537eSJed Brown   advection_ctx->CtauS             = CtauS;
52215a3537eSJed Brown   advection_ctx->E_wind            = E_wind;
52315a3537eSJed Brown   advection_ctx->implicit          = implicit;
52415a3537eSJed Brown   advection_ctx->strong_form       = strong_form;
52515a3537eSJed Brown   advection_ctx->stabilization     = stab;
52657272ee0SJames Wright   advection_ctx->stabilization_tau = stab_tau;
52757272ee0SJames Wright   advection_ctx->Ctau_a            = Ctau_a;
528fbabb365SJames Wright   advection_ctx->Ctau_d            = Ctau_d;
52957272ee0SJames Wright   advection_ctx->Ctau_t            = Ctau_t;
530c8d249deSJames Wright   advection_ctx->diffusion_coeff   = diffusion_coeff;
5310c373b74SJames Wright   advection_ctx->divFdiff_method   = honee->app_ctx->divFdiffproj_method;
532a515125bSLeila Ghaffari 
533f5dc303cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &ics_qfctx));
534f5dc303cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(ics_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
535f5dc303cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(ics_qfctx, CEED_MEM_HOST, FreeContextPetsc));
536a515125bSLeila Ghaffari 
5370c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx));
538e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
539e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc));
540e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
54157272ee0SJames Wright                                                          "Size of timestep, delta t"));
542f5dc303cSJames Wright 
543f5dc303cSJames Wright   // ------------------------------------------------------
544f5dc303cSJames Wright   //               SET UP ADVECTION
545f5dc303cSJames Wright   // ------------------------------------------------------
546f5dc303cSJames Wright   problem->print_info         = PRINT_ADVECTION;
547f5dc303cSJames Wright   problem->num_comps_jac_data = 0;
548f5dc303cSJames Wright   switch (dim) {
549f5dc303cSJames Wright     case 2:
550f5dc303cSJames Wright       problem->ics                          = (HoneeQFSpec){.qf_func_ptr = ICsAdvection2d, .qf_loc = ICsAdvection2d_loc};
551f5dc303cSJames Wright       problem->apply_vol_rhs                = (HoneeQFSpec){.qf_func_ptr = RHS_Advection2d, .qf_loc = RHS_Advection2d_loc};
552f5dc303cSJames Wright       problem->apply_vol_ifunction          = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection2d, .qf_loc = IFunction_Advection2d_loc};
553f5dc303cSJames Wright       problem->compute_exact_solution_error = PETSC_TRUE;
554f5dc303cSJames Wright       break;
555f5dc303cSJames Wright     case 3:
556f5dc303cSJames Wright       problem->ics                          = (HoneeQFSpec){.qf_func_ptr = ICsAdvection, .qf_loc = ICsAdvection_loc};
557f5dc303cSJames Wright       problem->apply_vol_rhs                = (HoneeQFSpec){.qf_func_ptr = RHS_Advection, .qf_loc = RHS_Advection_loc};
558f5dc303cSJames Wright       problem->apply_vol_ifunction          = (HoneeQFSpec){.qf_func_ptr = IFunction_Advection, .qf_loc = IFunction_Advection_loc};
559f5dc303cSJames Wright       problem->compute_exact_solution_error = PETSC_FALSE;
560f5dc303cSJames Wright       break;
561f5dc303cSJames Wright   }
562f5dc303cSJames Wright   problem->ics.qfctx           = ics_qfctx;
563e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = advection_qfctx;
564e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx));
565d6cac220SJames Wright 
566f5dc303cSJames Wright   problem->num_components = 5;
567f5dc303cSJames Wright   PetscCall(PetscMalloc1(problem->num_components, &problem->component_names));
568f5dc303cSJames Wright   for (PetscInt i = 0; i < 5; i++) PetscCall(PetscStrallocpy(component_names[i], &problem->component_names[i]));
569f5dc303cSJames Wright 
570f5dc303cSJames Wright   PetscCall(DivDiffFluxProjectionCreate(honee, honee->app_ctx->divFdiffproj_method, 1, &honee->diff_flux_proj));
571f5dc303cSJames Wright   if (honee->diff_flux_proj) {
572f5dc303cSJames Wright     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
573f5dc303cSJames Wright     NodalProjectionData       projection     = diff_flux_proj->projection;
574ae2e5212SJames Wright     PetscSection              section;
575f5dc303cSJames Wright 
576f5dc303cSJames Wright     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_AdvDif;
577f5dc303cSJames Wright     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif;
578ae2e5212SJames Wright     PetscCall(DMGetLocalSection(projection->dm, &section));
579f5dc303cSJames Wright     switch (honee->diff_flux_proj->method) {
580f5dc303cSJames Wright       case DIV_DIFF_FLUX_PROJ_DIRECT: {
581f5dc303cSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
582f5dc303cSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar"));
583f5dc303cSJames Wright       } break;
584f5dc303cSJames Wright       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
585f5dc303cSJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
586f5dc303cSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX"));
587f5dc303cSJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY"));
588f5dc303cSJames Wright         if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ"));
589f5dc303cSJames Wright       } break;
590f5dc303cSJames Wright       case DIV_DIFF_FLUX_PROJ_NONE:
591f5dc303cSJames Wright         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
592f5dc303cSJames Wright                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
593f5dc303cSJames Wright         break;
594f5dc303cSJames Wright     }
595f5dc303cSJames Wright   }
596f5dc303cSJames Wright 
597d6cac220SJames Wright   for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
598d6cac220SJames Wright     BCDefinition bc_def = problem->bc_defs[b];
599d6cac220SJames Wright     const char  *name;
600d6cac220SJames Wright 
601d6cac220SJames Wright     PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
602d6cac220SJames Wright     if (!strcmp(name, "inflow")) {
603d6cac220SJames Wright       HoneeBCStruct honee_bc;
604d6cac220SJames Wright 
605d6cac220SJames Wright       PetscCall(PetscNew(&honee_bc));
606d6cac220SJames Wright       PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &honee_bc->qfctx));
607d6cac220SJames Wright       honee_bc->honee              = honee;
6081abc2837SJames Wright       honee_bc->num_comps_jac_data = 0;
609*26d401f3SJames Wright       PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc));
610d6cac220SJames Wright 
611d6cac220SJames Wright       PetscCall(BCDefinitionSetIFunction(bc_def, AdvectionInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
612d6cac220SJames Wright       PetscCall(BCDefinitionSetIJacobian(bc_def, NULL, NULL));
613d6cac220SJames Wright     }
614d6cac220SJames Wright   }
615d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
616a515125bSLeila Ghaffari }
617