xref: /honee/problems/advection.c (revision 71f2ed299fae2990cd472d07f295294bae9f840c)
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 
14a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
15a78efa86SJames Wright //
16a78efa86SJames Wright // Only used for SUPG stabilization
170c373b74SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(Honee honee, CeedOperator *op_mass) {
180c373b74SJames Wright   Ceed                 ceed = honee->ceed;
19a78efa86SJames Wright   CeedInt              num_comp_q, q_data_size;
2087d3884fSJeremy L Thompson   CeedQFunction        qf_mass = NULL;
2101e19bfaSJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd;
22a78efa86SJames Wright   CeedBasis            basis_q;
23a78efa86SJames Wright   CeedVector           q_data;
24236a8ba6SJames Wright   CeedQFunctionContext qfctx = NULL;
25a78efa86SJames Wright   PetscInt             dim;
26a78efa86SJames Wright 
27a78efa86SJames Wright   PetscFunctionBeginUser;
280c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
29a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
30a78efa86SJames Wright     CeedOperator     *sub_ops;
3101e19bfaSJames Wright     CeedOperatorField op_field;
32a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
33a78efa86SJames Wright 
340c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
3501e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
3601e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
3701e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
3801e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
39a78efa86SJames Wright 
40236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
41a78efa86SJames Wright   }
42a78efa86SJames Wright 
43a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
4401e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
45a78efa86SJames Wright 
46a78efa86SJames Wright   switch (dim) {
47a78efa86SJames Wright     case 2:
48a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
49a78efa86SJames Wright       break;
50a78efa86SJames Wright     case 3:
51a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
52a78efa86SJames Wright       break;
53a78efa86SJames Wright   }
54a78efa86SJames Wright 
55236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
56a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
57a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
58a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
59a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
60a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
61a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
62a78efa86SJames Wright 
63a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
64*71f2ed29SJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator, Advection-Diffusion Stabilized"));
65a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
660c373b74SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, honee->q_ceed));
6701e19bfaSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
68a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
69a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
70a78efa86SJames Wright 
71fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
7201e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
73fff85bd3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
74fff85bd3SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
75236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
76a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
77a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
78a78efa86SJames Wright }
79a78efa86SJames Wright 
80c2d90829SJames Wright /**
81c2d90829SJames Wright   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
82c2d90829SJames Wright 
830c373b74SJames Wright   @param[in]  honee          `Honee` context
84c2d90829SJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
85c2d90829SJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
86c2d90829SJames Wright **/
87e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
880c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
89c2d90829SJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
90c2d90829SJames Wright   CeedInt              num_comp_q;
91c2d90829SJames Wright   PetscInt             dim, label_value = 0;
92c2d90829SJames Wright   DMLabel              domain_label    = NULL;
93c2d90829SJames Wright   CeedQFunctionContext advection_qfctx = NULL;
94c2d90829SJames Wright 
95c2d90829SJames Wright   PetscFunctionBeginUser;
96c2d90829SJames Wright   // -- Get Pre-requisite things
97c2d90829SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
98e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
99c2d90829SJames Wright 
10040b78511SJames Wright   {  // Get advection-diffusion QF context
101c2d90829SJames Wright     CeedOperator *sub_ops;
102c2d90829SJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
103c2d90829SJames Wright 
1040c373b74SJames Wright     if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops));
1050c373b74SJames Wright     else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
106c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
107c2d90829SJames Wright   }
108c2d90829SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs));
109c2d90829SJames Wright   {  // Add the volume integral CeedOperator
110c2d90829SJames Wright     CeedQFunction       qf_rhs_volume = NULL;
111c2d90829SJames Wright     CeedOperator        op_rhs_volume;
112c2d90829SJames Wright     CeedVector          q_data;
113c2d90829SJames Wright     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL;
114c2d90829SJames Wright     CeedBasis           basis_diff_flux = NULL;
115c2d90829SJames Wright     CeedInt             q_data_size;
116c2d90829SJames Wright 
117c2d90829SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
118e3663b90SJames Wright     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
119e3663b90SJames Wright                        &q_data_size));
120c2d90829SJames Wright     switch (dim) {
121c2d90829SJames Wright       case 2:
122c2d90829SJames Wright         PetscCallCeed(
123c2d90829SJames Wright             ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume));
124c2d90829SJames Wright         break;
125c2d90829SJames Wright       case 3:
126c2d90829SJames Wright         PetscCallCeed(
127c2d90829SJames Wright             ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume));
128c2d90829SJames Wright         break;
129c2d90829SJames Wright     }
1300c373b74SJames Wright     PetscCheck(qf_rhs_volume, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
131c2d90829SJames Wright 
132c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx));
133c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
134c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
135c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
136c2d90829SJames Wright 
137c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
138e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
139c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
140c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
141c2d90829SJames Wright 
142c2d90829SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume));
143c2d90829SJames Wright 
144c2d90829SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
145c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
146c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
147c2d90829SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
148c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
149c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
150c2d90829SJames Wright   }
151c2d90829SJames Wright 
152c2d90829SJames Wright   {  // Add the boundary integral CeedOperator
153c2d90829SJames Wright     CeedQFunction qf_rhs_boundary;
154c2d90829SJames Wright     DMLabel       face_sets_label;
155c2d90829SJames Wright     PetscInt      num_face_set_values, *face_set_values;
156c2d90829SJames Wright     CeedInt       q_data_size;
157c2d90829SJames Wright 
158c2d90829SJames Wright     // -- Build RHS operator
159c2d90829SJames Wright     switch (dim) {
160c2d90829SJames Wright       case 2:
161c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc,
162c2d90829SJames Wright                                                         &qf_rhs_boundary));
163c2d90829SJames Wright         break;
164c2d90829SJames Wright       case 3:
165c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc,
166c2d90829SJames Wright                                                         &qf_rhs_boundary));
167c2d90829SJames Wright         break;
168c2d90829SJames Wright     }
169c2d90829SJames Wright 
1700c373b74SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(honee->dm, &q_data_size));
171c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx));
172c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
173c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
174c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
175c2d90829SJames Wright 
176c2d90829SJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
177c2d90829SJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
178c2d90829SJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
179c2d90829SJames Wright       DMLabel  face_orientation_label;
180c2d90829SJames Wright       PetscInt num_orientations_values, *orientation_values;
181c2d90829SJames Wright 
182c2d90829SJames Wright       {
183c2d90829SJames Wright         char *face_orientation_label_name;
184c2d90829SJames Wright 
185c2d90829SJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
186c2d90829SJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
187c2d90829SJames Wright         PetscCall(PetscFree(face_orientation_label_name));
188c2d90829SJames Wright       }
189c2d90829SJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
190c2d90829SJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
191c2d90829SJames Wright         CeedOperator        op_rhs_boundary;
192c2d90829SJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
193c2d90829SJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
194c2d90829SJames Wright         CeedVector          q_data;
195c2d90829SJames Wright         CeedInt             q_data_size;
196c2d90829SJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
197c2d90829SJames Wright 
1980c373b74SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
1990c373b74SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, honee->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
200c2d90829SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
201c2d90829SJames Wright                                                   &elem_restr_diff_flux_boundary));
202c2d90829SJames Wright         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
203e3663b90SJames Wright         PetscCall(
204e3663b90SJames Wright             QDataBoundaryGradientGet(ceed, honee->dm, face_orientation_label, orientation, honee->x_coord, &elem_restr_qdata, &q_data, &q_data_size));
205c2d90829SJames Wright 
206c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
207c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
208c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
209c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
210c2d90829SJames Wright                                                  CEED_VECTOR_ACTIVE));
211c2d90829SJames Wright 
212c2d90829SJames Wright         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary));
213c2d90829SJames Wright 
214c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
215c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
216c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
217c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
218c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
219c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
220c2d90829SJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
221c2d90829SJames Wright       }
222c2d90829SJames Wright       PetscCall(PetscFree(orientation_values));
223c2d90829SJames Wright     }
224c2d90829SJames Wright     PetscCall(PetscFree(face_set_values));
225c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
226c2d90829SJames Wright   }
227c2d90829SJames Wright 
228c2d90829SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
229c2d90829SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
230c2d90829SJames Wright }
231c2d90829SJames Wright 
23240b78511SJames Wright /**
23340b78511SJames Wright   @brief Create RHS CeedOperator for indirect projection of divergence of diffusive flux
23440b78511SJames Wright 
2350c373b74SJames Wright   @param[in]  honee          `Honee` context
23640b78511SJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
23740b78511SJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
23840b78511SJames Wright **/
239e3663b90SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Indirect_AdvDif(Honee honee, DivDiffFluxProjectionData diff_flux_proj, CeedOperator *op_rhs) {
2400c373b74SJames Wright   Ceed                 ceed       = honee->ceed;
24140b78511SJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
24240b78511SJames Wright   CeedBasis            basis_diff_flux;
24340b78511SJames Wright   CeedElemRestriction  elem_restr_diff_flux, elem_restr_qd;
24440b78511SJames Wright   CeedVector           q_data;
24540b78511SJames Wright   CeedInt              num_comp_q, q_data_size;
24640b78511SJames Wright   PetscInt             dim;
24740b78511SJames Wright   PetscInt             label_value = 0, height = 0, dm_field = 0;
24840b78511SJames Wright   DMLabel              domain_label    = NULL;
24940b78511SJames Wright   CeedQFunction        qf_rhs          = NULL;
25040b78511SJames Wright   CeedQFunctionContext advection_qfctx = NULL;
25140b78511SJames Wright 
25240b78511SJames Wright   PetscFunctionBeginUser;
25340b78511SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
254e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q));
25540b78511SJames Wright 
25640b78511SJames Wright   {  // Get advection-diffusion QF context
25740b78511SJames Wright     CeedOperator *sub_ops;
25840b78511SJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
25940b78511SJames Wright 
2600c373b74SJames Wright     if (honee->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_ifunction, &sub_ops));
2610c373b74SJames Wright     else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
26240b78511SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
26340b78511SJames Wright   }
26440b78511SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, domain_label, label_value, height, dm_field, &elem_restr_diff_flux));
26540b78511SJames Wright   PetscCall(CreateBasisFromPlex(ceed, projection->dm, domain_label, label_value, height, dm_field, &basis_diff_flux));
266e3663b90SJames Wright   PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
267e3663b90SJames Wright                      &q_data_size));
26840b78511SJames Wright 
26940b78511SJames Wright   switch (dim) {
27040b78511SJames Wright     case 2:
27140b78511SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_2D, DiffusiveFluxRHS_AdvDif_2D_loc, &qf_rhs));
27240b78511SJames Wright       break;
27340b78511SJames Wright     case 3:
27440b78511SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DiffusiveFluxRHS_AdvDif_3D, DiffusiveFluxRHS_AdvDif_3D_loc, &qf_rhs));
27540b78511SJames Wright       break;
27640b78511SJames Wright   }
2770c373b74SJames Wright   PetscCheck(qf_rhs, honee->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
27840b78511SJames Wright 
27940b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs, advection_qfctx));
28040b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
28140b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs, "qdata", q_data_size, CEED_EVAL_NONE));
28240b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs, "F_diff RHS", projection->num_comp, CEED_EVAL_INTERP));
28340b78511SJames Wright 
28440b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs, NULL, NULL, op_rhs));
285e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
28640b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
28740b78511SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_rhs, "F_diff RHS", elem_restr_diff_flux, basis_diff_flux, CEED_VECTOR_ACTIVE));
28840b78511SJames Wright 
28940b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs));
29040b78511SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
29140b78511SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
29240b78511SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
29340b78511SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
29440b78511SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
29540b78511SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
29640b78511SJames Wright }
29740b78511SJames Wright 
298991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
2995f636aeaSJames Wright   AdvDifWindType             wind_type;
3005f636aeaSJames Wright   AdvDifICType               advectionic_type;
301108c6689SJames Wright   AdvDifBubbleContinuityType bubble_continuity_type = -1;
302a515125bSLeila Ghaffari   StabilizationType          stab;
30357272ee0SJames Wright   StabilizationTauType       stab_tau;
3043636f6a4SJames Wright   SetupContextAdv            setup_context;
3050c373b74SJames Wright   Honee                      honee = *(Honee *)ctx;
3060c373b74SJames Wright   MPI_Comm                   comm  = honee->comm;
3070c373b74SJames Wright   Ceed                       ceed  = honee->ceed;
308a515125bSLeila Ghaffari   PetscBool                  implicit;
30915a3537eSJed Brown   AdvectionContext           advection_ctx;
310e07531f7SJames Wright   CeedQFunctionContext       advection_qfctx;
3119529d636SJames Wright   PetscInt                   dim;
312a515125bSLeila Ghaffari 
31315a3537eSJed Brown   PetscFunctionBeginUser;
3142b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
3152b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &advection_ctx));
3169529d636SJames Wright   PetscCall(DMGetDimension(dm, &dim));
317a515125bSLeila Ghaffari 
318a515125bSLeila Ghaffari   // ------------------------------------------------------
319a515125bSLeila Ghaffari   //               SET UP ADVECTION
320a515125bSLeila Ghaffari   // ------------------------------------------------------
32128160fc2SJames Wright   problem->print_info        = PRINT_ADVECTION;
32228160fc2SJames Wright   problem->jac_data_size_vol = 0;
32328160fc2SJames Wright   problem->jac_data_size_sur = 0;
3249529d636SJames Wright   switch (dim) {
3259529d636SJames Wright     case 2:
326e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsAdvection2d;
327e07531f7SJames Wright       problem->ics.qf_loc                      = ICsAdvection2d_loc;
328e07531f7SJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHS_Advection2d;
329e07531f7SJames Wright       problem->apply_vol_rhs.qf_loc            = RHS_Advection2d_loc;
330e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d;
331e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Advection2d_loc;
332e07531f7SJames Wright       problem->apply_inflow.qf_func_ptr        = Advection2d_InOutFlow;
333e07531f7SJames Wright       problem->apply_inflow.qf_loc             = Advection2d_InOutFlow_loc;
33458ce1233SJames Wright       problem->compute_exact_solution_error    = PETSC_TRUE;
3359529d636SJames Wright       break;
3369529d636SJames Wright     case 3:
337e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsAdvection;
338e07531f7SJames Wright       problem->ics.qf_loc                      = ICsAdvection_loc;
339e07531f7SJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHS_Advection;
340e07531f7SJames Wright       problem->apply_vol_rhs.qf_loc            = RHS_Advection_loc;
341e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection;
342e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Advection_loc;
343e07531f7SJames Wright       problem->apply_inflow.qf_func_ptr        = Advection_InOutFlow;
344e07531f7SJames Wright       problem->apply_inflow.qf_loc             = Advection_InOutFlow_loc;
34558ce1233SJames Wright       problem->compute_exact_solution_error    = PETSC_FALSE;
3469529d636SJames Wright       break;
3479529d636SJames Wright   }
348a515125bSLeila Ghaffari 
3490c373b74SJames Wright   PetscCall(DivDiffFluxProjectionCreate(honee, 1, &honee->diff_flux_proj));
3500c373b74SJames Wright   if (honee->diff_flux_proj) {
3510c373b74SJames Wright     DivDiffFluxProjectionData diff_flux_proj = honee->diff_flux_proj;
352c2d90829SJames Wright     NodalProjectionData       projection     = diff_flux_proj->projection;
353c2d90829SJames Wright 
354c2d90829SJames Wright     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_AdvDif;
35540b78511SJames Wright     diff_flux_proj->CreateRHSOperator_Indirect = DivDiffFluxProjectionCreateRHS_Indirect_AdvDif;
356c2d90829SJames Wright 
3570c373b74SJames Wright     switch (honee->diff_flux_proj->method) {
358c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_DIRECT: {
359c2d90829SJames Wright         PetscSection section;
360c2d90829SJames Wright 
361c2d90829SJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
362c2d90829SJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
363c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar"));
364c2d90829SJames Wright       } break;
365c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
366c2d90829SJames Wright         PetscSection section;
367c2d90829SJames Wright 
368c2d90829SJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
369c2d90829SJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
370c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX"));
371c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY"));
37240b78511SJames Wright         if (dim >= 3) PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ"));
373c2d90829SJames Wright       } break;
374c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_NONE:
3750c373b74SJames Wright         SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
3760c373b74SJames Wright                 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]);
377c2d90829SJames Wright         break;
378c2d90829SJames Wright     }
379c2d90829SJames Wright   }
380c2d90829SJames Wright 
381a515125bSLeila Ghaffari   // ------------------------------------------------------
382ea615d4cSJames Wright   //             Create the QFunction context
383a515125bSLeila Ghaffari   // ------------------------------------------------------
384a515125bSLeila Ghaffari   CeedScalar     rc               = 1000.;  // m (Radius of bubble)
385a515125bSLeila Ghaffari   CeedScalar     CtauS            = 0.;     // dimensionless
386059d1c40SJames Wright   PetscBool      strong_form      = PETSC_FALSE;
387a515125bSLeila Ghaffari   CeedScalar     E_wind           = 1.e6;  // J
3880c373b74SJames Wright   CeedScalar     Ctau_a           = PetscPowScalarInt(honee->app_ctx->degree, 2);
3890c373b74SJames Wright   CeedScalar     Ctau_d           = PetscPowScalarInt(honee->app_ctx->degree, 4);
39057272ee0SJames Wright   CeedScalar     Ctau_t           = 0.;
391a515125bSLeila Ghaffari   PetscReal      wind[3]          = {1., 0, 0};  // m/s
392c8d249deSJames Wright   CeedScalar     diffusion_coeff  = 0.;
393a62be6baSJames Wright   CeedScalar     wave_frequency   = 2 * M_PI;
394a62be6baSJames Wright   CeedScalar     wave_phase       = 0;
395108c6689SJames Wright   AdvDifWaveType wave_type        = -1;
396b4fd18dfSJames Wright   PetscScalar    bl_height_factor = 1.;
39787d3884fSJeremy L Thompson   PetscReal      domain_min[3], domain_max[3], domain_size[3] = {0.};
3982b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
39966d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
40005a512bdSLeila Ghaffari 
401a515125bSLeila Ghaffari   // ------------------------------------------------------
402a515125bSLeila Ghaffari   //             Create the PETSc context
403a515125bSLeila Ghaffari   // ------------------------------------------------------
404a515125bSLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
405a515125bSLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
406a515125bSLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
407a515125bSLeila Ghaffari   PetscScalar Joule;
408a515125bSLeila Ghaffari 
409a515125bSLeila Ghaffari   // ------------------------------------------------------
410a515125bSLeila Ghaffari   //              Command line Options
411a515125bSLeila Ghaffari   // ------------------------------------------------------
4121485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
413a515125bSLeila Ghaffari   // -- Physics
414a515125bSLeila Ghaffari   PetscBool translation;
4155f636aeaSJames Wright   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION),
4165f636aeaSJames Wright                              (PetscEnum *)&wind_type, &translation));
41766d54740SJames Wright   PetscInt  n = dim;
418a515125bSLeila Ghaffari   PetscBool user_wind;
4192b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
420c8d249deSJames Wright   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
4212b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
422059d1c40SJames Wright   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
423059d1c40SJames Wright                              NULL));
4242b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
4252b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
42657272ee0SJames Wright   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
42757272ee0SJames Wright                              (PetscEnum *)&stab_tau, NULL));
42857272ee0SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
429fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL));
430fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL));
4312b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
43280e9ac5bSJames Wright   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes,
43380e9ac5bSJames Wright                              (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
43457bb6b22SJames Wright   // IC-specific options
43557bb6b22SJames Wright   switch (advectionic_type) {
43657bb6b22SJames Wright     case ADVDIF_IC_WAVE:
43780e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_type", "-advection_ic_wave_type", "HONEE 0.0", NULL));
43880e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_frequency", "-advection_ic_wave_frequency", "HONEE 0.0", NULL));
43980e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-wave_phase", "-advection_ic_wave_phase", "HONEE 0.0", NULL));
44080e9ac5bSJames Wright       PetscCall(PetscOptionsEnum("-advection_ic_wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE),
441108c6689SJames Wright                                  (PetscEnum *)&wave_type, NULL));
44280e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL));
44380e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL));
44457bb6b22SJames Wright       break;
44557bb6b22SJames Wright     case ADVDIF_IC_BOUNDARY_LAYER:
44680e9ac5bSJames Wright       PetscCall(
44780e9ac5bSJames Wright           PetscOptionsScalar("-advection_ic_bl_height_factor", "Height of boundary layer in IC", NULL, bl_height_factor, &bl_height_factor, NULL));
44857bb6b22SJames Wright       break;
44957bb6b22SJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER:
45057bb6b22SJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
45180e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-rc", "-advection_ic_bubble_rc", "HONEE 0.0", NULL));
45280e9ac5bSJames Wright       PetscCall(PetscOptionsDeprecated("-bubble_continuity", "-advection_ic_bubble_continuity", "HONEE 0.0", NULL));
45380e9ac5bSJames Wright       PetscCall(PetscOptionsScalar("-advection_ic_bubble_rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
454108c6689SJames Wright       bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE;
45580e9ac5bSJames Wright       PetscCall(PetscOptionsEnum("-advection_ic_bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes,
456108c6689SJames Wright                                  (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL));
45757bb6b22SJames Wright       break;
45857bb6b22SJames Wright     case ADVDIF_IC_SKEW:
45957bb6b22SJames Wright     case ADVDIF_IC_COSINE_HILL:
46057bb6b22SJames Wright       break;
461108c6689SJames Wright   }
462a62be6baSJames Wright 
463a515125bSLeila Ghaffari   // -- Units
4642b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
465a515125bSLeila Ghaffari   meter = fabs(meter);
4662b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
467a515125bSLeila Ghaffari   kilogram = fabs(kilogram);
4682b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
469a515125bSLeila Ghaffari   second = fabs(second);
470a515125bSLeila Ghaffari 
471a515125bSLeila Ghaffari   // -- Warnings
4725f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_ROTATION && user_wind) {
4732b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
474a515125bSLeila Ghaffari   }
4755f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) {
476a515125bSLeila Ghaffari     wind[2] = 0;
477c51f031aSJames Wright     PetscCall(
478c51f031aSJames Wright         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
479a515125bSLeila Ghaffari   }
480a515125bSLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
4812b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
482a515125bSLeila Ghaffari   }
4831485969bSJeremy L Thompson   PetscOptionsEnd();
484a515125bSLeila Ghaffari 
485a78efa86SJames Wright   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
486a78efa86SJames Wright 
487a515125bSLeila Ghaffari   // ------------------------------------------------------
488a515125bSLeila Ghaffari   //           Set up the PETSc context
489a515125bSLeila Ghaffari   // ------------------------------------------------------
490a515125bSLeila Ghaffari   // -- Define derived units
491a515125bSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
492a515125bSLeila Ghaffari 
4930c373b74SJames Wright   honee->units->meter    = meter;
4940c373b74SJames Wright   honee->units->kilogram = kilogram;
4950c373b74SJames Wright   honee->units->second   = second;
4960c373b74SJames Wright   honee->units->Joule    = Joule;
497a515125bSLeila Ghaffari 
498a515125bSLeila Ghaffari   // ------------------------------------------------------
499ea615d4cSJames Wright   //           Set up the QFunction contexts
500a515125bSLeila Ghaffari   // ------------------------------------------------------
501a515125bSLeila Ghaffari   // -- Scale variables to desired units
502a515125bSLeila Ghaffari   E_wind *= Joule;
503a515125bSLeila Ghaffari   rc = fabs(rc) * meter;
50466d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) {
50505a512bdSLeila Ghaffari     wind[i] *= (meter / second);
50605a512bdSLeila Ghaffari     domain_size[i] *= meter;
50705a512bdSLeila Ghaffari   }
508a515125bSLeila Ghaffari 
509a515125bSLeila Ghaffari   // -- Setup Context
510a515125bSLeila Ghaffari   setup_context->rc                     = rc;
51105a512bdSLeila Ghaffari   setup_context->lx                     = domain_size[0];
51205a512bdSLeila Ghaffari   setup_context->ly                     = domain_size[1];
51366d54740SJames Wright   setup_context->lz                     = dim == 3 ? domain_size[2] : 0.;
514a515125bSLeila Ghaffari   setup_context->wind[0]                = wind[0];
515a515125bSLeila Ghaffari   setup_context->wind[1]                = wind[1];
51666d54740SJames Wright   setup_context->wind[2]                = dim == 3 ? wind[2] : 0.;
517a515125bSLeila Ghaffari   setup_context->wind_type              = wind_type;
518c51f031aSJames Wright   setup_context->initial_condition_type = advectionic_type;
519a515125bSLeila Ghaffari   setup_context->bubble_continuity_type = bubble_continuity_type;
520a515125bSLeila Ghaffari   setup_context->time                   = 0;
521a62be6baSJames Wright   setup_context->wave_frequency         = wave_frequency;
522a62be6baSJames Wright   setup_context->wave_phase             = wave_phase;
523a62be6baSJames Wright   setup_context->wave_type              = wave_type;
524b4fd18dfSJames Wright   setup_context->bl_height_factor       = bl_height_factor;
525a515125bSLeila Ghaffari 
526a515125bSLeila Ghaffari   // -- QFunction Context
5270c373b74SJames Wright   honee->phys->implicit            = implicit;
52815a3537eSJed Brown   advection_ctx->CtauS             = CtauS;
52915a3537eSJed Brown   advection_ctx->E_wind            = E_wind;
53015a3537eSJed Brown   advection_ctx->implicit          = implicit;
53115a3537eSJed Brown   advection_ctx->strong_form       = strong_form;
53215a3537eSJed Brown   advection_ctx->stabilization     = stab;
53357272ee0SJames Wright   advection_ctx->stabilization_tau = stab_tau;
53457272ee0SJames Wright   advection_ctx->Ctau_a            = Ctau_a;
535fbabb365SJames Wright   advection_ctx->Ctau_d            = Ctau_d;
53657272ee0SJames Wright   advection_ctx->Ctau_t            = Ctau_t;
537c8d249deSJames Wright   advection_ctx->diffusion_coeff   = diffusion_coeff;
5380c373b74SJames Wright   advection_ctx->divFdiff_method   = honee->app_ctx->divFdiffproj_method;
539a515125bSLeila Ghaffari 
5400c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &problem->ics.qfctx));
541e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
542e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
543a515125bSLeila Ghaffari 
5440c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &advection_qfctx));
545e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
546e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc));
547e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
54857272ee0SJames Wright                                                          "Size of timestep, delta t"));
549e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = advection_qfctx;
550e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx));
551e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx));
552d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
553a515125bSLeila Ghaffari }
554a515125bSLeila Ghaffari 
5550c373b74SJames Wright PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx) {
5560c373b74SJames Wright   MPI_Comm         comm = honee->comm;
5570c373b74SJames Wright   Ceed             ceed = honee->ceed;
5583636f6a4SJames Wright   SetupContextAdv  setup_ctx;
55915a3537eSJed Brown   AdvectionContext advection_ctx;
56066d54740SJames Wright   PetscInt         dim;
561a515125bSLeila Ghaffari 
56215a3537eSJed Brown   PetscFunctionBeginUser;
5630c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
564e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx));
565e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
5662b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
567a515125bSLeila Ghaffari                         "  Problem:\n"
568a515125bSLeila Ghaffari                         "    Problem Name                       : %s\n"
569a515125bSLeila Ghaffari                         "    Stabilization                      : %s\n"
5705bcc2007SJames Wright                         "    Stabilization Tau                  : %s\n"
571a515125bSLeila Ghaffari                         "    Wind Type                          : %s\n",
5725bcc2007SJames Wright                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization],
5735f636aeaSJames Wright                         StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type]));
574a515125bSLeila Ghaffari 
5755f636aeaSJames Wright   if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) {
576c781aee8SJames Wright     CeedScalar *wind = setup_ctx->wind;
57766d54740SJames Wright     switch (dim) {
5789529d636SJames Wright       case 2:
579c781aee8SJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", wind[0], wind[1]));
5809529d636SJames Wright         break;
5819529d636SJames Wright       case 3:
582c781aee8SJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", wind[0], wind[1], wind[2]));
5839529d636SJames Wright         break;
5849529d636SJames Wright     }
585a515125bSLeila Ghaffari   }
586c781aee8SJames Wright 
5875f636aeaSJames Wright   PetscCall(PetscPrintf(comm, "    Initial Condition Type             : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type]));
588c781aee8SJames Wright   switch (setup_ctx->initial_condition_type) {
5895f636aeaSJames Wright     case ADVDIF_IC_SKEW:
5905f636aeaSJames Wright     case ADVDIF_IC_COSINE_HILL:
5913d1afcc1SJames Wright     case ADVDIF_IC_BOUNDARY_LAYER:
592c781aee8SJames Wright       break;
5935f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
5945f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER:
5955f636aeaSJames Wright       PetscCall(PetscPrintf(comm, "    Bubble Continuity                  : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type]));
596c781aee8SJames Wright       break;
5975f636aeaSJames Wright     case ADVDIF_IC_WAVE:
598c781aee8SJames Wright       PetscCall(PetscPrintf(comm, "    Wave Type                          : %s\n", AdvDifWaveTypes[setup_ctx->wave_type]));
599c781aee8SJames Wright       break;
600c781aee8SJames Wright   }
601c781aee8SJames Wright 
602e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx));
603e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
604d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
605a515125bSLeila Ghaffari }
606