xref: /libCEED/examples/fluids/problems/advection.c (revision d1d777230063fd92fc0f547b5443317655656d64)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Utility functions for setting up ADVECTION
1077841947SLeila Ghaffari 
112b730f8bSJeremy L Thompson #include "../qfunctions/advection.h"
122b730f8bSJeremy L Thompson 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
1677841947SLeila Ghaffari #include "../navierstokes.h"
1777841947SLeila Ghaffari #include "../qfunctions/setupgeo.h"
1893639ffbSJames Wright #include "../qfunctions/setupgeo2d.h"
1977841947SLeila Ghaffari 
20b4e9a8f8SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping
21b4e9a8f8SJames Wright //
22b4e9a8f8SJames Wright // Only used for SUPG stabilization
23b4e9a8f8SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) {
24b4e9a8f8SJames Wright   Ceed                 ceed = user->ceed;
25b4e9a8f8SJames Wright   CeedInt              num_comp_q, q_data_size;
26b4e9a8f8SJames Wright   CeedQFunction        qf_mass;
27b4e9a8f8SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
28b4e9a8f8SJames Wright   CeedBasis            basis_q;
29b4e9a8f8SJames Wright   CeedVector           q_data;
30b4e9a8f8SJames Wright   CeedQFunctionContext qf_ctx = NULL;
31b4e9a8f8SJames Wright   PetscInt             dim;
32b4e9a8f8SJames Wright 
33b4e9a8f8SJames Wright   PetscFunctionBeginUser;
34b4e9a8f8SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
35b4e9a8f8SJames Wright   {  // Get restriction and basis from the RHS function
36b4e9a8f8SJames Wright     CeedOperator     *sub_ops;
37b4e9a8f8SJames Wright     CeedOperatorField field;
38b4e9a8f8SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
39b4e9a8f8SJames Wright 
40b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
41b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
42b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
43b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
44b4e9a8f8SJames Wright 
45b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
46b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
47b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
48b4e9a8f8SJames Wright 
49b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx));
50b4e9a8f8SJames Wright   }
51b4e9a8f8SJames Wright 
52b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
53b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
54b4e9a8f8SJames Wright 
55b4e9a8f8SJames Wright   switch (dim) {
56b4e9a8f8SJames Wright     case 2:
57b4e9a8f8SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
58b4e9a8f8SJames Wright       break;
59b4e9a8f8SJames Wright     case 3:
60b4e9a8f8SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
61b4e9a8f8SJames Wright       break;
62b4e9a8f8SJames Wright   }
63b4e9a8f8SJames Wright 
64b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx));
65b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
66b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
67b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
68b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
69b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
70b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
71b4e9a8f8SJames Wright 
72b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
73b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
74b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
75b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
76b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
77b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
78b4e9a8f8SJames Wright 
79b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
80b4e9a8f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
81b4e9a8f8SJames Wright }
82b4e9a8f8SJames Wright 
8346de7363SJames Wright PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
8477841947SLeila Ghaffari   WindType             wind_type;
857b77ddfdSJames Wright   AdvectionICType      advectionic_type;
8677841947SLeila Ghaffari   BubbleContinuityType bubble_continuity_type;
8777841947SLeila Ghaffari   StabilizationType    stab;
88b18328c4SJames Wright   StabilizationTauType stab_tau;
8997baf651SJames Wright   SetupContextAdv      setup_context;
9077841947SLeila Ghaffari   User                 user = *(User *)ctx;
91a424bcd0SJames Wright   MPI_Comm             comm = user->comm;
92a424bcd0SJames Wright   Ceed                 ceed = user->ceed;
9377841947SLeila Ghaffari   PetscBool            implicit;
94841e4c73SJed Brown   AdvectionContext     advection_ctx;
95841e4c73SJed Brown   CeedQFunctionContext advection_context;
9693639ffbSJames Wright   PetscInt             dim;
9777841947SLeila Ghaffari 
98841e4c73SJed Brown   PetscFunctionBeginUser;
992b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
1002b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &advection_ctx));
10193639ffbSJames Wright   PetscCall(DMGetDimension(dm, &dim));
10277841947SLeila Ghaffari 
10377841947SLeila Ghaffari   // ------------------------------------------------------
10477841947SLeila Ghaffari   //               SET UP ADVECTION
10577841947SLeila Ghaffari   // ------------------------------------------------------
10693639ffbSJames Wright   switch (dim) {
10793639ffbSJames Wright     case 2:
10893639ffbSJames Wright       problem->dim                               = 2;
10993639ffbSJames Wright       problem->q_data_size_vol                   = 5;
11093639ffbSJames Wright       problem->q_data_size_sur                   = 3;
11193639ffbSJames Wright       problem->setup_vol.qfunction               = Setup2d;
11293639ffbSJames Wright       problem->setup_vol.qfunction_loc           = Setup2d_loc;
11393639ffbSJames Wright       problem->setup_sur.qfunction               = SetupBoundary2d;
11493639ffbSJames Wright       problem->setup_sur.qfunction_loc           = SetupBoundary2d_loc;
11593639ffbSJames Wright       problem->ics.qfunction                     = ICsAdvection2d;
11693639ffbSJames Wright       problem->ics.qfunction_loc                 = ICsAdvection2d_loc;
11793639ffbSJames Wright       problem->apply_vol_rhs.qfunction           = RHS_Advection2d;
11893639ffbSJames Wright       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection2d_loc;
11993639ffbSJames Wright       problem->apply_vol_ifunction.qfunction     = IFunction_Advection2d;
12093639ffbSJames Wright       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc;
12193639ffbSJames Wright       problem->apply_inflow.qfunction            = Advection2d_InOutFlow;
12293639ffbSJames Wright       problem->apply_inflow.qfunction_loc        = Advection2d_InOutFlow_loc;
12393639ffbSJames Wright       problem->non_zero_time                     = PETSC_TRUE;
12493639ffbSJames Wright       problem->print_info                        = PRINT_ADVECTION;
12593639ffbSJames Wright       break;
12693639ffbSJames Wright     case 3:
12777841947SLeila Ghaffari       problem->dim                               = 3;
12877841947SLeila Ghaffari       problem->q_data_size_vol                   = 10;
129ba6664aeSJames Wright       problem->q_data_size_sur                   = 10;
13091e5af17SJed Brown       problem->setup_vol.qfunction               = Setup;
13191e5af17SJed Brown       problem->setup_vol.qfunction_loc           = Setup_loc;
13291e5af17SJed Brown       problem->setup_sur.qfunction               = SetupBoundary;
13391e5af17SJed Brown       problem->setup_sur.qfunction_loc           = SetupBoundary_loc;
13491e5af17SJed Brown       problem->ics.qfunction                     = ICsAdvection;
13591e5af17SJed Brown       problem->ics.qfunction_loc                 = ICsAdvection_loc;
13693639ffbSJames Wright       problem->apply_vol_rhs.qfunction           = RHS_Advection;
13793639ffbSJames Wright       problem->apply_vol_rhs.qfunction_loc       = RHS_Advection_loc;
13891e5af17SJed Brown       problem->apply_vol_ifunction.qfunction     = IFunction_Advection;
13991e5af17SJed Brown       problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc;
14091e5af17SJed Brown       problem->apply_inflow.qfunction            = Advection_InOutFlow;
14191e5af17SJed Brown       problem->apply_inflow.qfunction_loc        = Advection_InOutFlow_loc;
14277841947SLeila Ghaffari       problem->non_zero_time                     = PETSC_FALSE;
14377841947SLeila Ghaffari       problem->print_info                        = PRINT_ADVECTION;
14493639ffbSJames Wright       break;
14593639ffbSJames Wright   }
14677841947SLeila Ghaffari 
14777841947SLeila Ghaffari   // ------------------------------------------------------
14877841947SLeila Ghaffari   //             Create the libCEED context
14977841947SLeila Ghaffari   // ------------------------------------------------------
15077841947SLeila Ghaffari   CeedScalar rc              = 1000.;  // m (Radius of bubble)
15177841947SLeila Ghaffari   CeedScalar CtauS           = 0.;     // dimensionless
152b767acfbSJames Wright   PetscBool  strong_form     = PETSC_FALSE;
15377841947SLeila Ghaffari   CeedScalar E_wind          = 1.e6;  // J
154b18328c4SJames Wright   CeedScalar Ctau_a          = PetscPowScalarInt(user->app_ctx->degree, 2);
155b18328c4SJames Wright   CeedScalar Ctau_t          = 0.;
15677841947SLeila Ghaffari   PetscReal  wind[3]         = {1., 0, 0};  // m/s
157*d1d77723SJames Wright   CeedScalar diffusion_coeff = 0.;
1581864f1c2SLeila Ghaffari   PetscReal  domain_min[3], domain_max[3], domain_size[3];
1592b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
16093639ffbSJames Wright   for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
1611864f1c2SLeila Ghaffari 
16277841947SLeila Ghaffari   // ------------------------------------------------------
16377841947SLeila Ghaffari   //             Create the PETSc context
16477841947SLeila Ghaffari   // ------------------------------------------------------
16577841947SLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
16677841947SLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
16777841947SLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
16877841947SLeila Ghaffari   PetscScalar Joule;
16977841947SLeila Ghaffari 
17077841947SLeila Ghaffari   // ------------------------------------------------------
17177841947SLeila Ghaffari   //              Command line Options
17277841947SLeila Ghaffari   // ------------------------------------------------------
17367490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
17477841947SLeila Ghaffari   // -- Physics
1752b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
17677841947SLeila Ghaffari   PetscBool translation;
1772b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type,
1782b730f8bSJeremy L Thompson                              &translation));
17977841947SLeila Ghaffari   PetscInt  n = problem->dim;
18077841947SLeila Ghaffari   PetscBool user_wind;
1812b730f8bSJeremy L Thompson   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
182*d1d77723SJames Wright   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
1832b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
184b767acfbSJames Wright   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
185b767acfbSJames Wright                              NULL));
1862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
1877b77ddfdSJames Wright   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes,
1887b77ddfdSJames Wright                              (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
18993639ffbSJames Wright   bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE;
19093639ffbSJames Wright   PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type,
19193639ffbSJames Wright                              (PetscEnum *)&bubble_continuity_type, NULL));
1922b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
193b18328c4SJames Wright   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
194b18328c4SJames Wright                              (PetscEnum *)&stab_tau, NULL));
195b18328c4SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
196b18328c4SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL));
1972b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
19877841947SLeila Ghaffari 
19977841947SLeila Ghaffari   // -- Units
2002b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
20177841947SLeila Ghaffari   meter = fabs(meter);
2022b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
20377841947SLeila Ghaffari   kilogram = fabs(kilogram);
2042b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
20577841947SLeila Ghaffari   second = fabs(second);
20677841947SLeila Ghaffari 
20777841947SLeila Ghaffari   // -- Warnings
20877841947SLeila Ghaffari   if (wind_type == WIND_ROTATION && user_wind) {
2092b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
21077841947SLeila Ghaffari   }
2117b77ddfdSJames Wright   if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) {
21277841947SLeila Ghaffari     wind[2] = 0;
2137b77ddfdSJames Wright     PetscCall(
2147b77ddfdSJames Wright         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
21577841947SLeila Ghaffari   }
21677841947SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
2172b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
21877841947SLeila Ghaffari   }
21967490bc6SJeremy L Thompson   PetscOptionsEnd();
22077841947SLeila Ghaffari 
221b4e9a8f8SJames Wright   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
222b4e9a8f8SJames Wright 
22377841947SLeila Ghaffari   // ------------------------------------------------------
22477841947SLeila Ghaffari   //           Set up the PETSc context
22577841947SLeila Ghaffari   // ------------------------------------------------------
22677841947SLeila Ghaffari   // -- Define derived units
22777841947SLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
22877841947SLeila Ghaffari 
22977841947SLeila Ghaffari   user->units->meter    = meter;
23077841947SLeila Ghaffari   user->units->kilogram = kilogram;
23177841947SLeila Ghaffari   user->units->second   = second;
23277841947SLeila Ghaffari   user->units->Joule    = Joule;
23377841947SLeila Ghaffari 
23477841947SLeila Ghaffari   // ------------------------------------------------------
23577841947SLeila Ghaffari   //           Set up the libCEED context
23677841947SLeila Ghaffari   // ------------------------------------------------------
23777841947SLeila Ghaffari   // -- Scale variables to desired units
23877841947SLeila Ghaffari   E_wind *= Joule;
23977841947SLeila Ghaffari   rc = fabs(rc) * meter;
24093639ffbSJames Wright   for (PetscInt i = 0; i < problem->dim; i++) {
2411864f1c2SLeila Ghaffari     wind[i] *= (meter / second);
2421864f1c2SLeila Ghaffari     domain_size[i] *= meter;
2431864f1c2SLeila Ghaffari   }
2441864f1c2SLeila Ghaffari   problem->dm_scale = meter;
24577841947SLeila Ghaffari 
24677841947SLeila Ghaffari   // -- Setup Context
24777841947SLeila Ghaffari   setup_context->rc                     = rc;
2481864f1c2SLeila Ghaffari   setup_context->lx                     = domain_size[0];
2491864f1c2SLeila Ghaffari   setup_context->ly                     = domain_size[1];
25093639ffbSJames Wright   setup_context->lz                     = problem->dim == 3 ? domain_size[2] : 0.;
25177841947SLeila Ghaffari   setup_context->wind[0]                = wind[0];
25277841947SLeila Ghaffari   setup_context->wind[1]                = wind[1];
25393639ffbSJames Wright   setup_context->wind[2]                = problem->dim == 3 ? wind[2] : 0.;
25477841947SLeila Ghaffari   setup_context->wind_type              = wind_type;
2557b77ddfdSJames Wright   setup_context->initial_condition_type = advectionic_type;
25677841947SLeila Ghaffari   setup_context->bubble_continuity_type = bubble_continuity_type;
25777841947SLeila Ghaffari   setup_context->time                   = 0;
25877841947SLeila Ghaffari 
25977841947SLeila Ghaffari   // -- QFunction Context
26077841947SLeila Ghaffari   user->phys->implicit             = implicit;
261841e4c73SJed Brown   advection_ctx->CtauS             = CtauS;
262841e4c73SJed Brown   advection_ctx->E_wind            = E_wind;
263841e4c73SJed Brown   advection_ctx->implicit          = implicit;
264841e4c73SJed Brown   advection_ctx->strong_form       = strong_form;
265841e4c73SJed Brown   advection_ctx->stabilization     = stab;
266b18328c4SJames Wright   advection_ctx->stabilization_tau = stab_tau;
267b18328c4SJames Wright   advection_ctx->Ctau_a            = Ctau_a;
268b18328c4SJames Wright   advection_ctx->Ctau_t            = Ctau_t;
269*d1d77723SJames Wright   advection_ctx->diffusion_coeff   = diffusion_coeff;
27077841947SLeila Ghaffari 
271a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context));
272a424bcd0SJames Wright   PetscCallCeed(ceed,
273a424bcd0SJames Wright                 CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
274a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc));
27577841947SLeila Ghaffari 
276a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context));
277a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
278a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc));
279b18328c4SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
280b18328c4SJames Wright                                                          "Size of timestep, delta t"));
281841e4c73SJed Brown   problem->apply_vol_rhs.qfunction_context = advection_context;
282a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context));
283a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context));
284ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
28577841947SLeila Ghaffari }
28677841947SLeila Ghaffari 
287367c849eSJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx) {
288367c849eSJames Wright   MPI_Comm         comm = user->comm;
289a424bcd0SJames Wright   Ceed             ceed = user->ceed;
29097baf651SJames Wright   SetupContextAdv  setup_ctx;
291841e4c73SJed Brown   AdvectionContext advection_ctx;
29277841947SLeila Ghaffari 
293841e4c73SJed Brown   PetscFunctionBeginUser;
294a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx));
295a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx));
2962b730f8bSJeremy L Thompson   PetscCall(PetscPrintf(comm,
29777841947SLeila Ghaffari                         "  Problem:\n"
29877841947SLeila Ghaffari                         "    Problem Name                       : %s\n"
29977841947SLeila Ghaffari                         "    Stabilization                      : %s\n"
3007b77ddfdSJames Wright                         "    Initial Condition Type             : %s\n"
30177841947SLeila Ghaffari                         "    Bubble Continuity                  : %s\n"
30277841947SLeila Ghaffari                         "    Wind Type                          : %s\n",
3037b77ddfdSJames Wright                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type],
3047b77ddfdSJames Wright                         BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type]));
30577841947SLeila Ghaffari 
306841e4c73SJed Brown   if (setup_ctx->wind_type == WIND_TRANSLATION) {
30793639ffbSJames Wright     switch (problem->dim) {
30893639ffbSJames Wright       case 2:
30993639ffbSJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1]));
31093639ffbSJames Wright         break;
31193639ffbSJames Wright       case 3:
31293639ffbSJames Wright         PetscCall(
31393639ffbSJames Wright             PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2]));
31493639ffbSJames Wright         break;
31593639ffbSJames Wright     }
31677841947SLeila Ghaffari   }
317a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx));
318a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx));
319ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32077841947SLeila Ghaffari }
321