xref: /honee/problems/advection.c (revision a78efa869d5e60cbf1c38b3b4d784572ed49ea21)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Utility functions for setting up ADVECTION
10a515125bSLeila Ghaffari 
112b916ea7SJeremy L Thompson #include "../qfunctions/advection.h"
122b916ea7SJeremy L Thompson 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson 
16a515125bSLeila Ghaffari #include "../navierstokes.h"
17a515125bSLeila Ghaffari #include "../qfunctions/setupgeo.h"
189529d636SJames Wright #include "../qfunctions/setupgeo2d.h"
19a515125bSLeila Ghaffari 
20*a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
21*a78efa86SJames Wright //
22*a78efa86SJames Wright // Only used for SUPG stabilization
23*a78efa86SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) {
24*a78efa86SJames Wright   Ceed                 ceed = user->ceed;
25*a78efa86SJames Wright   CeedInt              num_comp_q, q_data_size;
26*a78efa86SJames Wright   CeedQFunction        qf_mass;
27*a78efa86SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
28*a78efa86SJames Wright   CeedBasis            basis_q;
29*a78efa86SJames Wright   CeedVector           q_data;
30*a78efa86SJames Wright   CeedQFunctionContext qf_ctx = NULL;
31*a78efa86SJames Wright   PetscInt             dim;
32*a78efa86SJames Wright 
33*a78efa86SJames Wright   PetscFunctionBeginUser;
34*a78efa86SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
35*a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
36*a78efa86SJames Wright     CeedOperator     *sub_ops;
37*a78efa86SJames Wright     CeedOperatorField field;
38*a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
39*a78efa86SJames Wright 
40*a78efa86SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
41*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
42*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
43*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
44*a78efa86SJames Wright 
45*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
46*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
47*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
48*a78efa86SJames Wright 
49*a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
50*a78efa86SJames Wright   }
51*a78efa86SJames Wright 
52*a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
53*a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
54*a78efa86SJames Wright 
55*a78efa86SJames Wright   switch (dim) {
56*a78efa86SJames Wright     case 2:
57*a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
58*a78efa86SJames Wright       break;
59*a78efa86SJames Wright     case 3:
60*a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
61*a78efa86SJames Wright       break;
62*a78efa86SJames Wright   }
63*a78efa86SJames Wright 
64*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
65*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
66*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
67*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
68*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
69*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
70*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
71*a78efa86SJames Wright 
72*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
73*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
74*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
75*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
76*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
77*a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
78*a78efa86SJames Wright 
79*a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
80*a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
81*a78efa86SJames Wright }
82*a78efa86SJames Wright 
83d1c51a42SJames Wright PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
84a515125bSLeila Ghaffari   WindType             wind_type;
85c51f031aSJames Wright   AdvectionICType      advectionic_type;
86a515125bSLeila Ghaffari   BubbleContinuityType bubble_continuity_type;
87a515125bSLeila Ghaffari   StabilizationType    stab;
8857272ee0SJames Wright   StabilizationTauType stab_tau;
893636f6a4SJames Wright   SetupContextAdv      setup_context;
90a515125bSLeila Ghaffari   User                 user = *(User *)ctx;
91b4c37c5cSJames Wright   MPI_Comm             comm = user->comm;
92b4c37c5cSJames Wright   Ceed                 ceed = user->ceed;
93a515125bSLeila Ghaffari   PetscBool            implicit;
9415a3537eSJed Brown   AdvectionContext     advection_ctx;
9515a3537eSJed Brown   CeedQFunctionContext advection_context;
969529d636SJames Wright   PetscInt             dim;
97a515125bSLeila Ghaffari 
9815a3537eSJed Brown   PetscFunctionBeginUser;
992b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
1002b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &advection_ctx));
1019529d636SJames Wright   PetscCall(DMGetDimension(dm, &dim));
102a515125bSLeila Ghaffari 
103a515125bSLeila Ghaffari   // ------------------------------------------------------
104a515125bSLeila Ghaffari   //               SET UP ADVECTION
105a515125bSLeila Ghaffari   // ------------------------------------------------------
1069529d636SJames Wright   switch (dim) {
1079529d636SJames Wright     case 2:
1089529d636SJames Wright       problem->dim                               = 2;
1099529d636SJames Wright       problem->q_data_size_vol                   = 5;
1109529d636SJames Wright       problem->q_data_size_sur                   = 3;
1119529d636SJames Wright       problem->setup_vol.qfunction               = Setup2d;
1129529d636SJames Wright       problem->setup_vol.qfunction_loc           = Setup2d_loc;
1139529d636SJames Wright       problem->setup_sur.qfunction               = SetupBoundary2d;
1149529d636SJames Wright       problem->setup_sur.qfunction_loc           = SetupBoundary2d_loc;
1159529d636SJames Wright       problem->ics.qfunction                     = ICsAdvection2d;
1169529d636SJames Wright       problem->ics.qfunction_loc                 = ICsAdvection2d_loc;
1179529d636SJames Wright       problem->apply_vol_rhs.qfunction           = RHS_Advection2d;
1189529d636SJames Wright       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection2d_loc;
1199529d636SJames Wright       problem->apply_vol_ifunction.qfunction     = IFunction_Advection2d;
1209529d636SJames Wright       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc;
1219529d636SJames Wright       problem->apply_inflow.qfunction            = Advection2d_InOutFlow;
1229529d636SJames Wright       problem->apply_inflow.qfunction_loc        = Advection2d_InOutFlow_loc;
1239529d636SJames Wright       problem->non_zero_time                     = PETSC_TRUE;
1249529d636SJames Wright       problem->print_info                        = PRINT_ADVECTION;
1259529d636SJames Wright       break;
1269529d636SJames Wright     case 3:
127a515125bSLeila Ghaffari       problem->dim                               = 3;
128a515125bSLeila Ghaffari       problem->q_data_size_vol                   = 10;
129493642f1SJames Wright       problem->q_data_size_sur                   = 10;
1309785fe93SJed Brown       problem->setup_vol.qfunction               = Setup;
1319785fe93SJed Brown       problem->setup_vol.qfunction_loc           = Setup_loc;
1329785fe93SJed Brown       problem->setup_sur.qfunction               = SetupBoundary;
1339785fe93SJed Brown       problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
1349785fe93SJed Brown       problem->ics.qfunction                     = ICsAdvection;
1359785fe93SJed Brown       problem->ics.qfunction_loc                 = ICsAdvection_loc;
1369529d636SJames Wright       problem->apply_vol_rhs.qfunction           = RHS_Advection;
1379529d636SJames Wright       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection_loc;
1389785fe93SJed Brown       problem->apply_vol_ifunction.qfunction     = IFunction_Advection;
1399785fe93SJed Brown       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc;
1409785fe93SJed Brown       problem->apply_inflow.qfunction            = Advection_InOutFlow;
1419785fe93SJed Brown       problem->apply_inflow.qfunction_loc        = Advection_InOutFlow_loc;
142a515125bSLeila Ghaffari       problem->non_zero_time                     = PETSC_FALSE;
143a515125bSLeila Ghaffari       problem->print_info                        = PRINT_ADVECTION;
1449529d636SJames Wright       break;
1459529d636SJames Wright   }
146a515125bSLeila Ghaffari 
147a515125bSLeila Ghaffari   // ------------------------------------------------------
148a515125bSLeila Ghaffari   //             Create the libCEED context
149a515125bSLeila Ghaffari   // ------------------------------------------------------
150a515125bSLeila Ghaffari   CeedScalar rc          = 1000.;  // m (Radius of bubble)
151a515125bSLeila Ghaffari   CeedScalar CtauS       = 0.;     // dimensionless
152059d1c40SJames Wright   PetscBool  strong_form = PETSC_FALSE;
153a515125bSLeila Ghaffari   CeedScalar E_wind      = 1.e6;  // J
15457272ee0SJames Wright   CeedScalar Ctau_a      = PetscPowScalarInt(user->app_ctx->degree, 2);
15557272ee0SJames Wright   CeedScalar Ctau_t      = 0.;
156a515125bSLeila Ghaffari   PetscReal  wind[3]     = {1., 0, 0};  // m/s
15705a512bdSLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
1582b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
1599529d636SJames Wright   for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
16005a512bdSLeila Ghaffari 
161a515125bSLeila Ghaffari   // ------------------------------------------------------
162a515125bSLeila Ghaffari   //             Create the PETSc context
163a515125bSLeila Ghaffari   // ------------------------------------------------------
164a515125bSLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
165a515125bSLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
166a515125bSLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
167a515125bSLeila Ghaffari   PetscScalar Joule;
168a515125bSLeila Ghaffari 
169a515125bSLeila Ghaffari   // ------------------------------------------------------
170a515125bSLeila Ghaffari   //              Command line Options
171a515125bSLeila Ghaffari   // ------------------------------------------------------
1721485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
173a515125bSLeila Ghaffari   // -- Physics
1742b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
175a515125bSLeila Ghaffari   PetscBool translation;
1762b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type,
1772b916ea7SJeremy L Thompson                              &translation));
178a515125bSLeila Ghaffari   PetscInt  n = problem->dim;
179a515125bSLeila Ghaffari   PetscBool user_wind;
1802b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
1812b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
182059d1c40SJames Wright   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
183059d1c40SJames Wright                              NULL));
1842b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
185c51f031aSJames Wright   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes,
186c51f031aSJames Wright                              (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
1879529d636SJames Wright   bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE;
1889529d636SJames Wright   PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type,
1899529d636SJames Wright                              (PetscEnum *)&bubble_continuity_type, NULL));
1902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
19157272ee0SJames Wright   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
19257272ee0SJames Wright                              (PetscEnum *)&stab_tau, NULL));
19357272ee0SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
19457272ee0SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL));
1952b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
196a515125bSLeila Ghaffari 
197a515125bSLeila Ghaffari   // -- Units
1982b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
199a515125bSLeila Ghaffari   meter = fabs(meter);
2002b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
201a515125bSLeila Ghaffari   kilogram = fabs(kilogram);
2022b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
203a515125bSLeila Ghaffari   second = fabs(second);
204a515125bSLeila Ghaffari 
205a515125bSLeila Ghaffari   // -- Warnings
206a515125bSLeila Ghaffari   if (wind_type == WIND_ROTATION && user_wind) {
2072b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
208a515125bSLeila Ghaffari   }
209c51f031aSJames Wright   if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) {
210a515125bSLeila Ghaffari     wind[2] = 0;
211c51f031aSJames Wright     PetscCall(
212c51f031aSJames Wright         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
213a515125bSLeila Ghaffari   }
214a515125bSLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
2152b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
216a515125bSLeila Ghaffari   }
2171485969bSJeremy L Thompson   PetscOptionsEnd();
218a515125bSLeila Ghaffari 
219*a78efa86SJames Wright   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
220*a78efa86SJames Wright 
221a515125bSLeila Ghaffari   // ------------------------------------------------------
222a515125bSLeila Ghaffari   //           Set up the PETSc context
223a515125bSLeila Ghaffari   // ------------------------------------------------------
224a515125bSLeila Ghaffari   // -- Define derived units
225a515125bSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
226a515125bSLeila Ghaffari 
227a515125bSLeila Ghaffari   user->units->meter    = meter;
228a515125bSLeila Ghaffari   user->units->kilogram = kilogram;
229a515125bSLeila Ghaffari   user->units->second   = second;
230a515125bSLeila Ghaffari   user->units->Joule    = Joule;
231a515125bSLeila Ghaffari 
232a515125bSLeila Ghaffari   // ------------------------------------------------------
233a515125bSLeila Ghaffari   //           Set up the libCEED context
234a515125bSLeila Ghaffari   // ------------------------------------------------------
235a515125bSLeila Ghaffari   // -- Scale variables to desired units
236a515125bSLeila Ghaffari   E_wind *= Joule;
237a515125bSLeila Ghaffari   rc = fabs(rc) * meter;
2389529d636SJames Wright   for (PetscInt i = 0; i < problem->dim; i++) {
23905a512bdSLeila Ghaffari     wind[i] *= (meter / second);
24005a512bdSLeila Ghaffari     domain_size[i] *= meter;
24105a512bdSLeila Ghaffari   }
24205a512bdSLeila Ghaffari   problem->dm_scale = meter;
243a515125bSLeila Ghaffari 
244a515125bSLeila Ghaffari   // -- Setup Context
245a515125bSLeila Ghaffari   setup_context->rc                     = rc;
24605a512bdSLeila Ghaffari   setup_context->lx                     = domain_size[0];
24705a512bdSLeila Ghaffari   setup_context->ly                     = domain_size[1];
2489529d636SJames Wright   setup_context->lz                     = problem->dim == 3 ? domain_size[2] : 0.;
249a515125bSLeila Ghaffari   setup_context->wind[0]                = wind[0];
250a515125bSLeila Ghaffari   setup_context->wind[1]                = wind[1];
2519529d636SJames Wright   setup_context->wind[2]                = problem->dim == 3 ? wind[2] : 0.;
252a515125bSLeila Ghaffari   setup_context->wind_type              = wind_type;
253c51f031aSJames Wright   setup_context->initial_condition_type = advectionic_type;
254a515125bSLeila Ghaffari   setup_context->bubble_continuity_type = bubble_continuity_type;
255a515125bSLeila Ghaffari   setup_context->time                   = 0;
256a515125bSLeila Ghaffari 
257a515125bSLeila Ghaffari   // -- QFunction Context
258a515125bSLeila Ghaffari   user->phys->implicit             = implicit;
25915a3537eSJed Brown   advection_ctx->CtauS             = CtauS;
26015a3537eSJed Brown   advection_ctx->E_wind            = E_wind;
26115a3537eSJed Brown   advection_ctx->implicit          = implicit;
26215a3537eSJed Brown   advection_ctx->strong_form       = strong_form;
26315a3537eSJed Brown   advection_ctx->stabilization     = stab;
26457272ee0SJames Wright   advection_ctx->stabilization_tau = stab_tau;
26557272ee0SJames Wright   advection_ctx->Ctau_a            = Ctau_a;
26657272ee0SJames Wright   advection_ctx->Ctau_t            = Ctau_t;
267a515125bSLeila Ghaffari 
268b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
269b4c37c5cSJames Wright   PetscCallCeed(ceed,
270b4c37c5cSJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
271b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
272a515125bSLeila Ghaffari 
273b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context));
274b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
275b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc));
27657272ee0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
27757272ee0SJames Wright                                                          "Size of timestep, delta t"));
27815a3537eSJed Brown   problem->apply_vol_rhs.qfunction_context = advection_context;
279b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context));
280b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context));
281d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
282a515125bSLeila Ghaffari }
283a515125bSLeila Ghaffari 
2842d49c0afSJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx) {
2852d49c0afSJames Wright   MPI_Comm         comm = user->comm;
286b4c37c5cSJames Wright   Ceed             ceed = user->ceed;
2873636f6a4SJames Wright   SetupContextAdv  setup_ctx;
28815a3537eSJed Brown   AdvectionContext advection_ctx;
289a515125bSLeila Ghaffari 
29015a3537eSJed Brown   PetscFunctionBeginUser;
291b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx));
292b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx));
2932b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
294a515125bSLeila Ghaffari                         "  Problem:\n"
295a515125bSLeila Ghaffari                         "    Problem Name                       : %s\n"
296a515125bSLeila Ghaffari                         "    Stabilization                      : %s\n"
297c51f031aSJames Wright                         "    Initial Condition Type             : %s\n"
298a515125bSLeila Ghaffari                         "    Bubble Continuity                  : %s\n"
299a515125bSLeila Ghaffari                         "    Wind Type                          : %s\n",
300c51f031aSJames Wright                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type],
301c51f031aSJames Wright                         BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type]));
302a515125bSLeila Ghaffari 
30315a3537eSJed Brown   if (setup_ctx->wind_type == WIND_TRANSLATION) {
3049529d636SJames Wright     switch (problem->dim) {
3059529d636SJames Wright       case 2:
3069529d636SJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1]));
3079529d636SJames Wright         break;
3089529d636SJames Wright       case 3:
3099529d636SJames Wright         PetscCall(
3109529d636SJames Wright             PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2]));
3119529d636SJames Wright         break;
3129529d636SJames Wright     }
313a515125bSLeila Ghaffari   }
314b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx));
315b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx));
316d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
317a515125bSLeila Ghaffari }
318