xref: /honee/problems/advection.c (revision c2d90829cabf1269f4b39c8bdf4d53aac009543c)
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
17a78efa86SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) {
18a78efa86SJames Wright   Ceed                 ceed = user->ceed;
19a78efa86SJames Wright   CeedInt              num_comp_q, q_data_size;
2087d3884fSJeremy L Thompson   CeedQFunction        qf_mass = NULL;
21a78efa86SJames Wright   CeedElemRestriction  elem_restr_q, elem_restr_qd_i;
22a78efa86SJames Wright   CeedBasis            basis_q;
23a78efa86SJames Wright   CeedVector           q_data;
24236a8ba6SJames Wright   CeedQFunctionContext qfctx = NULL;
25a78efa86SJames Wright   PetscInt             dim;
26a78efa86SJames Wright 
27a78efa86SJames Wright   PetscFunctionBeginUser;
28a78efa86SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
29a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
30a78efa86SJames Wright     CeedOperator     *sub_ops;
31a78efa86SJames Wright     CeedOperatorField field;
32a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
33a78efa86SJames Wright 
34a78efa86SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
35a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
36a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
37a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
38a78efa86SJames Wright 
39a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
40a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
41a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
42a78efa86SJames Wright 
43236a8ba6SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qfctx));
44a78efa86SJames Wright   }
45a78efa86SJames Wright 
46a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
47a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
48a78efa86SJames Wright 
49a78efa86SJames Wright   switch (dim) {
50a78efa86SJames Wright     case 2:
51a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass));
52a78efa86SJames Wright       break;
53a78efa86SJames Wright     case 3:
54a78efa86SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass));
55a78efa86SJames Wright       break;
56a78efa86SJames Wright   }
57a78efa86SJames Wright 
58236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qfctx));
59a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0));
60a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP));
61a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP));
62a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE));
63a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP));
64a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD));
65a78efa86SJames Wright 
66a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
67a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
68a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed));
69a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
70a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
71a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
72a78efa86SJames Wright 
73236a8ba6SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qfctx));
74a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
75a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
76a78efa86SJames Wright }
77a78efa86SJames Wright 
78*c2d90829SJames Wright /**
79*c2d90829SJames Wright   @brief Create RHS CeedOperator for direct projection of divergence of diffusive flux
80*c2d90829SJames Wright 
81*c2d90829SJames Wright   @param[in]  user           `User` object
82*c2d90829SJames Wright   @param[in]  ceed_data      `CeedData` object
83*c2d90829SJames Wright   @param[in]  diff_flux_proj `DivDiffFluxProjectionData` object
84*c2d90829SJames Wright   @param[out] op_rhs         Operator to calculate the RHS of the L^2 projection
85*c2d90829SJames Wright **/
86*c2d90829SJames Wright static PetscErrorCode DivDiffFluxProjectionCreateRHS_Direct_AdvDif(User user, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj,
87*c2d90829SJames Wright                                                                    CeedOperator *op_rhs) {
88*c2d90829SJames Wright   Ceed                 ceed       = user->ceed;
89*c2d90829SJames Wright   NodalProjectionData  projection = diff_flux_proj->projection;
90*c2d90829SJames Wright   CeedInt              num_comp_q;
91*c2d90829SJames Wright   PetscInt             dim, label_value = 0;
92*c2d90829SJames Wright   DMLabel              domain_label    = NULL;
93*c2d90829SJames Wright   CeedQFunctionContext advection_qfctx = NULL;
94*c2d90829SJames Wright 
95*c2d90829SJames Wright   PetscFunctionBeginUser;
96*c2d90829SJames Wright   // -- Get Pre-requisite things
97*c2d90829SJames Wright   PetscCall(DMGetDimension(projection->dm, &dim));
98*c2d90829SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
99*c2d90829SJames Wright 
100*c2d90829SJames Wright   {  // Get newtonian QF context
101*c2d90829SJames Wright     CeedOperator *sub_ops;
102*c2d90829SJames Wright     PetscInt      sub_op_index = 0;  // will be 0 for the volume op
103*c2d90829SJames Wright 
104*c2d90829SJames Wright     if (user->op_ifunction) PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_ifunction, &sub_ops));
105*c2d90829SJames Wright     else PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
106*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &advection_qfctx));
107*c2d90829SJames Wright   }
108*c2d90829SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_rhs));
109*c2d90829SJames Wright   {  // Add the volume integral CeedOperator
110*c2d90829SJames Wright     CeedQFunction       qf_rhs_volume = NULL;
111*c2d90829SJames Wright     CeedOperator        op_rhs_volume;
112*c2d90829SJames Wright     CeedVector          q_data;
113*c2d90829SJames Wright     CeedElemRestriction elem_restr_qd, elem_restr_diff_flux_volume = NULL;
114*c2d90829SJames Wright     CeedBasis           basis_diff_flux = NULL;
115*c2d90829SJames Wright     CeedInt             q_data_size;
116*c2d90829SJames Wright 
117*c2d90829SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(diff_flux_proj, &elem_restr_diff_flux_volume, &basis_diff_flux, NULL, NULL));
118*c2d90829SJames Wright     PetscCall(QDataGet(ceed, projection->dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
119*c2d90829SJames Wright                        &elem_restr_qd, &q_data, &q_data_size));
120*c2d90829SJames Wright     switch (dim) {
121*c2d90829SJames Wright       case 2:
122*c2d90829SJames Wright         PetscCallCeed(
123*c2d90829SJames Wright             ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_2D, DivDiffusiveFluxVolumeRHS_AdvDif_2D_loc, &qf_rhs_volume));
124*c2d90829SJames Wright         break;
125*c2d90829SJames Wright       case 3:
126*c2d90829SJames Wright         PetscCallCeed(
127*c2d90829SJames Wright             ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxVolumeRHS_AdvDif_3D, DivDiffusiveFluxVolumeRHS_AdvDif_3D_loc, &qf_rhs_volume));
128*c2d90829SJames Wright         break;
129*c2d90829SJames Wright     }
130*c2d90829SJames Wright     PetscCheck(qf_rhs_volume, user->comm, PETSC_ERR_SUP, "%s not valid for DM of dimension %" PetscInt_FMT, __func__, dim);
131*c2d90829SJames Wright 
132*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_volume, advection_qfctx));
133*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
134*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_volume, "qdata", q_data_size, CEED_EVAL_NONE));
135*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_volume, "diffusive flux RHS", projection->num_comp * dim, CEED_EVAL_GRAD));
136*c2d90829SJames Wright 
137*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_volume, NULL, NULL, &op_rhs_volume));
138*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
139*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
140*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_volume, "diffusive flux RHS", elem_restr_diff_flux_volume, basis_diff_flux, CEED_VECTOR_ACTIVE));
141*c2d90829SJames Wright 
142*c2d90829SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_volume));
143*c2d90829SJames Wright 
144*c2d90829SJames Wright     PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
145*c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
146*c2d90829SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_volume));
147*c2d90829SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
148*c2d90829SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_volume));
149*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_volume));
150*c2d90829SJames Wright   }
151*c2d90829SJames Wright 
152*c2d90829SJames Wright   {  // Add the boundary integral CeedOperator
153*c2d90829SJames Wright     CeedQFunction qf_rhs_boundary;
154*c2d90829SJames Wright     DMLabel       face_sets_label;
155*c2d90829SJames Wright     PetscInt      num_face_set_values, *face_set_values;
156*c2d90829SJames Wright     CeedInt       q_data_size;
157*c2d90829SJames Wright 
158*c2d90829SJames Wright     // -- Build RHS operator
159*c2d90829SJames Wright     switch (dim) {
160*c2d90829SJames Wright       case 2:
161*c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_2D, DivDiffusiveFluxBoundaryRHS_AdvDif_2D_loc,
162*c2d90829SJames Wright                                                         &qf_rhs_boundary));
163*c2d90829SJames Wright         break;
164*c2d90829SJames Wright       case 3:
165*c2d90829SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, DivDiffusiveFluxBoundaryRHS_AdvDif_3D, DivDiffusiveFluxBoundaryRHS_AdvDif_3D_loc,
166*c2d90829SJames Wright                                                         &qf_rhs_boundary));
167*c2d90829SJames Wright         break;
168*c2d90829SJames Wright     }
169*c2d90829SJames Wright 
170*c2d90829SJames Wright     PetscCall(QDataBoundaryGradientGetNumComponents(user->dm, &q_data_size));
171*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_boundary, advection_qfctx));
172*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
173*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_boundary, "qdata", q_data_size, CEED_EVAL_NONE));
174*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_boundary, "diffusive flux RHS", projection->num_comp, CEED_EVAL_INTERP));
175*c2d90829SJames Wright 
176*c2d90829SJames Wright     PetscCall(DMGetLabel(projection->dm, "Face Sets", &face_sets_label));
177*c2d90829SJames Wright     PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_sets_label, &num_face_set_values, &face_set_values));
178*c2d90829SJames Wright     for (PetscInt f = 0; f < num_face_set_values; f++) {
179*c2d90829SJames Wright       DMLabel  face_orientation_label;
180*c2d90829SJames Wright       PetscInt num_orientations_values, *orientation_values;
181*c2d90829SJames Wright 
182*c2d90829SJames Wright       {
183*c2d90829SJames Wright         char *face_orientation_label_name;
184*c2d90829SJames Wright 
185*c2d90829SJames Wright         PetscCall(DMPlexCreateFaceLabel(projection->dm, face_set_values[f], &face_orientation_label_name));
186*c2d90829SJames Wright         PetscCall(DMGetLabel(projection->dm, face_orientation_label_name, &face_orientation_label));
187*c2d90829SJames Wright         PetscCall(PetscFree(face_orientation_label_name));
188*c2d90829SJames Wright       }
189*c2d90829SJames Wright       PetscCall(DMLabelCreateGlobalValueArray(projection->dm, face_orientation_label, &num_orientations_values, &orientation_values));
190*c2d90829SJames Wright       for (PetscInt o = 0; o < num_orientations_values; o++) {
191*c2d90829SJames Wright         CeedOperator        op_rhs_boundary;
192*c2d90829SJames Wright         CeedBasis           basis_q, basis_diff_flux_boundary;
193*c2d90829SJames Wright         CeedElemRestriction elem_restr_qdata, elem_restr_q, elem_restr_diff_flux_boundary;
194*c2d90829SJames Wright         CeedVector          q_data;
195*c2d90829SJames Wright         CeedInt             q_data_size;
196*c2d90829SJames Wright         PetscInt            orientation = orientation_values[o], dm_field_q = 0, height_cell = 0, height_face = 1;
197*c2d90829SJames Wright 
198*c2d90829SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, user->dm, face_orientation_label, orientation, height_cell, dm_field_q, &elem_restr_q));
199*c2d90829SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCreate(ceed, user->dm, face_orientation_label, orientation, orientation, dm_field_q, &basis_q));
200*c2d90829SJames Wright         PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, face_orientation_label, orientation, height_face, 0,
201*c2d90829SJames Wright                                                   &elem_restr_diff_flux_boundary));
202*c2d90829SJames Wright         PetscCall(CreateBasisFromPlex(ceed, projection->dm, face_orientation_label, orientation, height_face, 0, &basis_diff_flux_boundary));
203*c2d90829SJames Wright         PetscCall(QDataBoundaryGradientGet(ceed, user->dm, face_orientation_label, orientation, ceed_data->x_coord, &elem_restr_qdata, &q_data,
204*c2d90829SJames Wright                                            &q_data_size));
205*c2d90829SJames Wright 
206*c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_boundary, NULL, NULL, &op_rhs_boundary));
207*c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "Grad_q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
208*c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "qdata", elem_restr_qdata, CEED_BASIS_NONE, q_data));
209*c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_boundary, "diffusive flux RHS", elem_restr_diff_flux_boundary, basis_diff_flux_boundary,
210*c2d90829SJames Wright                                                  CEED_VECTOR_ACTIVE));
211*c2d90829SJames Wright 
212*c2d90829SJames Wright         PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_rhs, op_rhs_boundary));
213*c2d90829SJames Wright 
214*c2d90829SJames Wright         PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_boundary));
215*c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qdata));
216*c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
217*c2d90829SJames Wright         PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux_boundary));
218*c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
219*c2d90829SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux_boundary));
220*c2d90829SJames Wright         PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
221*c2d90829SJames Wright       }
222*c2d90829SJames Wright       PetscCall(PetscFree(orientation_values));
223*c2d90829SJames Wright     }
224*c2d90829SJames Wright     PetscCall(PetscFree(face_set_values));
225*c2d90829SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_boundary));
226*c2d90829SJames Wright   }
227*c2d90829SJames Wright 
228*c2d90829SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&advection_qfctx));
229*c2d90829SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
230*c2d90829SJames Wright }
231*c2d90829SJames Wright 
232991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
2335f636aeaSJames Wright   AdvDifWindType             wind_type;
2345f636aeaSJames Wright   AdvDifICType               advectionic_type;
235108c6689SJames Wright   AdvDifBubbleContinuityType bubble_continuity_type = -1;
236a515125bSLeila Ghaffari   StabilizationType          stab;
23757272ee0SJames Wright   StabilizationTauType       stab_tau;
2383636f6a4SJames Wright   SetupContextAdv            setup_context;
239a515125bSLeila Ghaffari   User                       user = *(User *)ctx;
240b4c37c5cSJames Wright   MPI_Comm                   comm = user->comm;
241b4c37c5cSJames Wright   Ceed                       ceed = user->ceed;
242a515125bSLeila Ghaffari   PetscBool                  implicit;
24315a3537eSJed Brown   AdvectionContext           advection_ctx;
244e07531f7SJames Wright   CeedQFunctionContext       advection_qfctx;
2459529d636SJames Wright   PetscInt                   dim;
246a515125bSLeila Ghaffari 
24715a3537eSJed Brown   PetscFunctionBeginUser;
2482b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &setup_context));
2492b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &advection_ctx));
2509529d636SJames Wright   PetscCall(DMGetDimension(dm, &dim));
251a515125bSLeila Ghaffari 
252a515125bSLeila Ghaffari   // ------------------------------------------------------
253a515125bSLeila Ghaffari   //               SET UP ADVECTION
254a515125bSLeila Ghaffari   // ------------------------------------------------------
25528160fc2SJames Wright   problem->print_info        = PRINT_ADVECTION;
25628160fc2SJames Wright   problem->jac_data_size_vol = 0;
25728160fc2SJames Wright   problem->jac_data_size_sur = 0;
2589529d636SJames Wright   switch (dim) {
2599529d636SJames Wright     case 2:
260e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsAdvection2d;
261e07531f7SJames Wright       problem->ics.qf_loc                      = ICsAdvection2d_loc;
262e07531f7SJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHS_Advection2d;
263e07531f7SJames Wright       problem->apply_vol_rhs.qf_loc            = RHS_Advection2d_loc;
264e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection2d;
265e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Advection2d_loc;
266e07531f7SJames Wright       problem->apply_inflow.qf_func_ptr        = Advection2d_InOutFlow;
267e07531f7SJames Wright       problem->apply_inflow.qf_loc             = Advection2d_InOutFlow_loc;
26858ce1233SJames Wright       problem->compute_exact_solution_error    = PETSC_TRUE;
2699529d636SJames Wright       break;
2709529d636SJames Wright     case 3:
271e07531f7SJames Wright       problem->ics.qf_func_ptr                 = ICsAdvection;
272e07531f7SJames Wright       problem->ics.qf_loc                      = ICsAdvection_loc;
273e07531f7SJames Wright       problem->apply_vol_rhs.qf_func_ptr       = RHS_Advection;
274e07531f7SJames Wright       problem->apply_vol_rhs.qf_loc            = RHS_Advection_loc;
275e07531f7SJames Wright       problem->apply_vol_ifunction.qf_func_ptr = IFunction_Advection;
276e07531f7SJames Wright       problem->apply_vol_ifunction.qf_loc      = IFunction_Advection_loc;
277e07531f7SJames Wright       problem->apply_inflow.qf_func_ptr        = Advection_InOutFlow;
278e07531f7SJames Wright       problem->apply_inflow.qf_loc             = Advection_InOutFlow_loc;
27958ce1233SJames Wright       problem->compute_exact_solution_error    = PETSC_FALSE;
2809529d636SJames Wright       break;
2819529d636SJames Wright   }
282a515125bSLeila Ghaffari 
283*c2d90829SJames Wright   PetscCall(DivDiffFluxProjectionCreate(user, 1, &user->diff_flux_proj));
284*c2d90829SJames Wright   if (user->diff_flux_proj) {
285*c2d90829SJames Wright     DivDiffFluxProjectionData diff_flux_proj = user->diff_flux_proj;
286*c2d90829SJames Wright     NodalProjectionData       projection     = diff_flux_proj->projection;
287*c2d90829SJames Wright 
288*c2d90829SJames Wright     diff_flux_proj->CreateRHSOperator_Direct   = DivDiffFluxProjectionCreateRHS_Direct_AdvDif;
289*c2d90829SJames Wright     diff_flux_proj->CreateRHSOperator_Indirect = NULL;
290*c2d90829SJames Wright 
291*c2d90829SJames Wright     switch (user->diff_flux_proj->method) {
292*c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_DIRECT: {
293*c2d90829SJames Wright         PetscSection section;
294*c2d90829SJames Wright 
295*c2d90829SJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
296*c2d90829SJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
297*c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DivDiffusiveFlux_Scalar"));
298*c2d90829SJames Wright       } break;
299*c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_INDIRECT: {
300*c2d90829SJames Wright         PetscSection section;
301*c2d90829SJames Wright 
302*c2d90829SJames Wright         PetscCall(DMGetLocalSection(projection->dm, &section));
303*c2d90829SJames Wright         PetscCall(PetscSectionSetFieldName(section, 0, ""));
304*c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 0, "DiffusiveFlux_ScalarX"));
305*c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 1, "DiffusiveFlux_ScalarY"));
306*c2d90829SJames Wright         PetscCall(PetscSectionSetComponentName(section, 0, 2, "DiffusiveFlux_ScalarZ"));
307*c2d90829SJames Wright       } break;
308*c2d90829SJames Wright       case DIV_DIFF_FLUX_PROJ_NONE:
309*c2d90829SJames Wright         SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_method %s",
310*c2d90829SJames Wright                 DivDiffFluxProjectionMethods[user->app_ctx->divFdiffproj_method]);
311*c2d90829SJames Wright         break;
312*c2d90829SJames Wright     }
313*c2d90829SJames Wright   }
314*c2d90829SJames Wright 
315a515125bSLeila Ghaffari   // ------------------------------------------------------
316a515125bSLeila Ghaffari   //             Create the libCEED context
317a515125bSLeila Ghaffari   // ------------------------------------------------------
318a515125bSLeila Ghaffari   CeedScalar     rc              = 1000.;  // m (Radius of bubble)
319a515125bSLeila Ghaffari   CeedScalar     CtauS           = 0.;     // dimensionless
320059d1c40SJames Wright   PetscBool      strong_form     = PETSC_FALSE;
321a515125bSLeila Ghaffari   CeedScalar     E_wind          = 1.e6;  // J
32257272ee0SJames Wright   CeedScalar     Ctau_a          = PetscPowScalarInt(user->app_ctx->degree, 2);
323fbabb365SJames Wright   CeedScalar     Ctau_d          = PetscPowScalarInt(user->app_ctx->degree, 4);
32457272ee0SJames Wright   CeedScalar     Ctau_t          = 0.;
325a515125bSLeila Ghaffari   PetscReal      wind[3]         = {1., 0, 0};  // m/s
326c8d249deSJames Wright   CeedScalar     diffusion_coeff = 0.;
327a62be6baSJames Wright   CeedScalar     wave_frequency  = 2 * M_PI;
328a62be6baSJames Wright   CeedScalar     wave_phase      = 0;
329108c6689SJames Wright   AdvDifWaveType wave_type       = -1;
33087d3884fSJeremy L Thompson   PetscReal      domain_min[3], domain_max[3], domain_size[3] = {0.};
3312b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
33266d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) domain_size[i] = domain_max[i] - domain_min[i];
33305a512bdSLeila Ghaffari 
334a515125bSLeila Ghaffari   // ------------------------------------------------------
335a515125bSLeila Ghaffari   //             Create the PETSc context
336a515125bSLeila Ghaffari   // ------------------------------------------------------
337a515125bSLeila Ghaffari   PetscScalar meter    = 1e-2;  // 1 meter in scaled length units
338a515125bSLeila Ghaffari   PetscScalar kilogram = 1e-6;  // 1 kilogram in scaled mass units
339a515125bSLeila Ghaffari   PetscScalar second   = 1e-2;  // 1 second in scaled time units
340a515125bSLeila Ghaffari   PetscScalar Joule;
341a515125bSLeila Ghaffari 
342a515125bSLeila Ghaffari   // ------------------------------------------------------
343a515125bSLeila Ghaffari   //              Command line Options
344a515125bSLeila Ghaffari   // ------------------------------------------------------
3451485969bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL);
346a515125bSLeila Ghaffari   // -- Physics
3472b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL));
348a515125bSLeila Ghaffari   PetscBool translation;
3495f636aeaSJames Wright   PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, AdvDifWindTypes, (PetscEnum)(wind_type = ADVDIF_WIND_ROTATION),
3505f636aeaSJames Wright                              (PetscEnum *)&wind_type, &translation));
35166d54740SJames Wright   PetscInt  n = dim;
352a515125bSLeila Ghaffari   PetscBool user_wind;
3532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind));
354c8d249deSJames Wright   PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL));
3552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL));
356059d1c40SJames Wright   PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form,
357059d1c40SJames Wright                              NULL));
3582b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL));
3595f636aeaSJames Wright   PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvDifICTypes,
3605f636aeaSJames Wright                              (PetscEnum)(advectionic_type = ADVDIF_IC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL));
3612b916ea7SJeremy L Thompson   PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL));
36257272ee0SJames Wright   PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU),
36357272ee0SJames Wright                              (PetscEnum *)&stab_tau, NULL));
36457272ee0SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL));
365fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization, advection component", NULL, Ctau_a, &Ctau_a, NULL));
366fbabb365SJames Wright   PetscCall(PetscOptionsScalar("-Ctau_d", "Coefficient for the stabilization, diffusion component", NULL, Ctau_d, &Ctau_d, NULL));
3672b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL));
368a515125bSLeila Ghaffari 
369108c6689SJames Wright   if (advectionic_type == ADVDIF_IC_WAVE) {
370108c6689SJames Wright     PetscCall(PetscOptionsEnum("-wave_type", "Type of wave", NULL, AdvDifWaveTypes, (PetscEnum)(wave_type = ADVDIF_WAVE_SINE),
371108c6689SJames Wright                                (PetscEnum *)&wave_type, NULL));
372a62be6baSJames Wright     PetscCall(PetscOptionsScalar("-wave_frequency", "Frequency of sine wave", NULL, wave_frequency, &wave_frequency, NULL));
373a62be6baSJames Wright     PetscCall(PetscOptionsScalar("-wave_phase", "Length correction", NULL, wave_phase, &wave_phase, NULL));
374108c6689SJames Wright   }
375108c6689SJames Wright 
376108c6689SJames Wright   if (advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER || advectionic_type == ADVDIF_IC_BUBBLE_SPHERE) {
377108c6689SJames Wright     bubble_continuity_type = dim == 3 ? ADVDIF_BUBBLE_CONTINUITY_SMOOTH : ADVDIF_BUBBLE_CONTINUITY_COSINE;
378108c6689SJames Wright     PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, AdvDifBubbleContinuityTypes,
379108c6689SJames Wright                                (PetscEnum)bubble_continuity_type, (PetscEnum *)&bubble_continuity_type, NULL));
380108c6689SJames Wright   }
381a62be6baSJames Wright 
382a515125bSLeila Ghaffari   // -- Units
3832b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL));
384a515125bSLeila Ghaffari   meter = fabs(meter);
3852b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL));
386a515125bSLeila Ghaffari   kilogram = fabs(kilogram);
3872b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL));
388a515125bSLeila Ghaffari   second = fabs(second);
389a515125bSLeila Ghaffari 
390a515125bSLeila Ghaffari   // -- Warnings
3915f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_ROTATION && user_wind) {
3922b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n"));
393a515125bSLeila Ghaffari   }
3945f636aeaSJames Wright   if (wind_type == ADVDIF_WIND_TRANSLATION && advectionic_type == ADVDIF_IC_BUBBLE_CYLINDER && wind[2] != 0.) {
395a515125bSLeila Ghaffari     wind[2] = 0;
396c51f031aSJames Wright     PetscCall(
397c51f031aSJames Wright         PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n"));
398a515125bSLeila Ghaffari   }
399a515125bSLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
4002b916ea7SJeremy L Thompson     PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n"));
401a515125bSLeila Ghaffari   }
4021485969bSJeremy L Thompson   PetscOptionsEnd();
403a515125bSLeila Ghaffari 
404a78efa86SJames Wright   if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized;
405a78efa86SJames Wright 
406a515125bSLeila Ghaffari   // ------------------------------------------------------
407a515125bSLeila Ghaffari   //           Set up the PETSc context
408a515125bSLeila Ghaffari   // ------------------------------------------------------
409a515125bSLeila Ghaffari   // -- Define derived units
410a515125bSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
411a515125bSLeila Ghaffari 
412a515125bSLeila Ghaffari   user->units->meter    = meter;
413a515125bSLeila Ghaffari   user->units->kilogram = kilogram;
414a515125bSLeila Ghaffari   user->units->second   = second;
415a515125bSLeila Ghaffari   user->units->Joule    = Joule;
416a515125bSLeila Ghaffari 
417a515125bSLeila Ghaffari   // ------------------------------------------------------
418a515125bSLeila Ghaffari   //           Set up the libCEED context
419a515125bSLeila Ghaffari   // ------------------------------------------------------
420a515125bSLeila Ghaffari   // -- Scale variables to desired units
421a515125bSLeila Ghaffari   E_wind *= Joule;
422a515125bSLeila Ghaffari   rc = fabs(rc) * meter;
42366d54740SJames Wright   for (PetscInt i = 0; i < dim; i++) {
42405a512bdSLeila Ghaffari     wind[i] *= (meter / second);
42505a512bdSLeila Ghaffari     domain_size[i] *= meter;
42605a512bdSLeila Ghaffari   }
427a515125bSLeila Ghaffari 
428a515125bSLeila Ghaffari   // -- Setup Context
429a515125bSLeila Ghaffari   setup_context->rc                     = rc;
43005a512bdSLeila Ghaffari   setup_context->lx                     = domain_size[0];
43105a512bdSLeila Ghaffari   setup_context->ly                     = domain_size[1];
43266d54740SJames Wright   setup_context->lz                     = dim == 3 ? domain_size[2] : 0.;
433a515125bSLeila Ghaffari   setup_context->wind[0]                = wind[0];
434a515125bSLeila Ghaffari   setup_context->wind[1]                = wind[1];
43566d54740SJames Wright   setup_context->wind[2]                = dim == 3 ? wind[2] : 0.;
436a515125bSLeila Ghaffari   setup_context->wind_type              = wind_type;
437c51f031aSJames Wright   setup_context->initial_condition_type = advectionic_type;
438a515125bSLeila Ghaffari   setup_context->bubble_continuity_type = bubble_continuity_type;
439a515125bSLeila Ghaffari   setup_context->time                   = 0;
440a62be6baSJames Wright   setup_context->wave_frequency         = wave_frequency;
441a62be6baSJames Wright   setup_context->wave_phase             = wave_phase;
442a62be6baSJames Wright   setup_context->wave_type              = wave_type;
443a515125bSLeila Ghaffari 
444a515125bSLeila Ghaffari   // -- QFunction Context
445a515125bSLeila Ghaffari   user->phys->implicit             = implicit;
44615a3537eSJed Brown   advection_ctx->CtauS             = CtauS;
44715a3537eSJed Brown   advection_ctx->E_wind            = E_wind;
44815a3537eSJed Brown   advection_ctx->implicit          = implicit;
44915a3537eSJed Brown   advection_ctx->strong_form       = strong_form;
45015a3537eSJed Brown   advection_ctx->stabilization     = stab;
45157272ee0SJames Wright   advection_ctx->stabilization_tau = stab_tau;
45257272ee0SJames Wright   advection_ctx->Ctau_a            = Ctau_a;
453fbabb365SJames Wright   advection_ctx->Ctau_d            = Ctau_d;
45457272ee0SJames Wright   advection_ctx->Ctau_t            = Ctau_t;
455c8d249deSJames Wright   advection_ctx->diffusion_coeff   = diffusion_coeff;
456a515125bSLeila Ghaffari 
457e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfctx));
458e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(problem->ics.qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context));
459e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfctx, CEED_MEM_HOST, FreeContextPetsc));
460a515125bSLeila Ghaffari 
461e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_qfctx));
462e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx));
463e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_qfctx, CEED_MEM_HOST, FreeContextPetsc));
464e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_qfctx, "timestep size", offsetof(struct AdvectionContext_, dt), 1,
46557272ee0SJames Wright                                                          "Size of timestep, delta t"));
466e07531f7SJames Wright   problem->apply_vol_rhs.qfctx = advection_qfctx;
467e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_vol_ifunction.qfctx));
468e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_qfctx, &problem->apply_inflow.qfctx));
469d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
470a515125bSLeila Ghaffari }
471a515125bSLeila Ghaffari 
472991aef52SJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) {
4732d49c0afSJames Wright   MPI_Comm         comm = user->comm;
474b4c37c5cSJames Wright   Ceed             ceed = user->ceed;
4753636f6a4SJames Wright   SetupContextAdv  setup_ctx;
47615a3537eSJed Brown   AdvectionContext advection_ctx;
47766d54740SJames Wright   PetscInt         dim;
478a515125bSLeila Ghaffari 
47915a3537eSJed Brown   PetscFunctionBeginUser;
48066d54740SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
481e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx));
482e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
4832b916ea7SJeremy L Thompson   PetscCall(PetscPrintf(comm,
484a515125bSLeila Ghaffari                         "  Problem:\n"
485a515125bSLeila Ghaffari                         "    Problem Name                       : %s\n"
486a515125bSLeila Ghaffari                         "    Stabilization                      : %s\n"
4875bcc2007SJames Wright                         "    Stabilization Tau                  : %s\n"
488a515125bSLeila Ghaffari                         "    Wind Type                          : %s\n",
4895bcc2007SJames Wright                         app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization],
4905f636aeaSJames Wright                         StabilizationTauTypes[advection_ctx->stabilization_tau], AdvDifWindTypes[setup_ctx->wind_type]));
491a515125bSLeila Ghaffari 
4925f636aeaSJames Wright   if (setup_ctx->wind_type == ADVDIF_WIND_TRANSLATION) {
493c781aee8SJames Wright     CeedScalar *wind = setup_ctx->wind;
49466d54740SJames Wright     switch (dim) {
4959529d636SJames Wright       case 2:
496c781aee8SJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f\n", wind[0], wind[1]));
4979529d636SJames Wright         break;
4989529d636SJames Wright       case 3:
499c781aee8SJames Wright         PetscCall(PetscPrintf(comm, "    Background Wind                    : %f,%f,%f\n", wind[0], wind[1], wind[2]));
5009529d636SJames Wright         break;
5019529d636SJames Wright     }
502a515125bSLeila Ghaffari   }
503c781aee8SJames Wright 
5045f636aeaSJames Wright   PetscCall(PetscPrintf(comm, "    Initial Condition Type             : %s\n", AdvDifICTypes[setup_ctx->initial_condition_type]));
505c781aee8SJames Wright   switch (setup_ctx->initial_condition_type) {
5065f636aeaSJames Wright     case ADVDIF_IC_SKEW:
5075f636aeaSJames Wright     case ADVDIF_IC_COSINE_HILL:
508c781aee8SJames Wright       break;
5095f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_SPHERE:
5105f636aeaSJames Wright     case ADVDIF_IC_BUBBLE_CYLINDER:
5115f636aeaSJames Wright       PetscCall(PetscPrintf(comm, "    Bubble Continuity                  : %s\n", AdvDifBubbleContinuityTypes[setup_ctx->bubble_continuity_type]));
512c781aee8SJames Wright       break;
5135f636aeaSJames Wright     case ADVDIF_IC_WAVE:
514c781aee8SJames Wright       PetscCall(PetscPrintf(comm, "    Wave Type                          : %s\n", AdvDifWaveTypes[setup_ctx->wave_type]));
515c781aee8SJames Wright       break;
516c781aee8SJames Wright   }
517c781aee8SJames Wright 
518e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx));
519e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
520d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
521a515125bSLeila Ghaffari }
522