1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2*ae2b091fSJames 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; 24a78efa86SJames Wright CeedQFunctionContext qf_ctx = 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 43a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 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 58a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 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 73d438a78dSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx)); 74a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 75a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 76a78efa86SJames Wright } 77a78efa86SJames Wright 78991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 79a515125bSLeila Ghaffari WindType wind_type; 80c51f031aSJames Wright AdvectionICType advectionic_type; 81a515125bSLeila Ghaffari BubbleContinuityType bubble_continuity_type; 82a515125bSLeila Ghaffari StabilizationType stab; 8357272ee0SJames Wright StabilizationTauType stab_tau; 843636f6a4SJames Wright SetupContextAdv setup_context; 85a515125bSLeila Ghaffari User user = *(User *)ctx; 86b4c37c5cSJames Wright MPI_Comm comm = user->comm; 87b4c37c5cSJames Wright Ceed ceed = user->ceed; 88a515125bSLeila Ghaffari PetscBool implicit; 8915a3537eSJed Brown AdvectionContext advection_ctx; 9015a3537eSJed Brown CeedQFunctionContext advection_context; 919529d636SJames Wright PetscInt dim; 92a515125bSLeila Ghaffari 9315a3537eSJed Brown PetscFunctionBeginUser; 942b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 952b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 969529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 97a515125bSLeila Ghaffari 98a515125bSLeila Ghaffari // ------------------------------------------------------ 99a515125bSLeila Ghaffari // SET UP ADVECTION 100a515125bSLeila Ghaffari // ------------------------------------------------------ 1019529d636SJames Wright switch (dim) { 1029529d636SJames Wright case 2: 1039529d636SJames Wright problem->dim = 2; 1049529d636SJames Wright problem->ics.qfunction = ICsAdvection2d; 1059529d636SJames Wright problem->ics.qfunction_loc = ICsAdvection2d_loc; 1069529d636SJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection2d; 1079529d636SJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection2d_loc; 1089529d636SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Advection2d; 1099529d636SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc; 1109529d636SJames Wright problem->apply_inflow.qfunction = Advection2d_InOutFlow; 1119529d636SJames Wright problem->apply_inflow.qfunction_loc = Advection2d_InOutFlow_loc; 11258ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 1139529d636SJames Wright problem->print_info = PRINT_ADVECTION; 1149529d636SJames Wright break; 1159529d636SJames Wright case 3: 116a515125bSLeila Ghaffari problem->dim = 3; 1179785fe93SJed Brown problem->ics.qfunction = ICsAdvection; 1189785fe93SJed Brown problem->ics.qfunction_loc = ICsAdvection_loc; 1199529d636SJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection; 1209529d636SJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection_loc; 1219785fe93SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Advection; 1229785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc; 1239785fe93SJed Brown problem->apply_inflow.qfunction = Advection_InOutFlow; 1249785fe93SJed Brown problem->apply_inflow.qfunction_loc = Advection_InOutFlow_loc; 12558ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 126a515125bSLeila Ghaffari problem->print_info = PRINT_ADVECTION; 1279529d636SJames Wright break; 1289529d636SJames Wright } 129a515125bSLeila Ghaffari 130a515125bSLeila Ghaffari // ------------------------------------------------------ 131a515125bSLeila Ghaffari // Create the libCEED context 132a515125bSLeila Ghaffari // ------------------------------------------------------ 133a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 134a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 135059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 136a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 13757272ee0SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 13857272ee0SJames Wright CeedScalar Ctau_t = 0.; 139a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 140c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 14187d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 1422b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 1439529d636SJames Wright for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 14405a512bdSLeila Ghaffari 145a515125bSLeila Ghaffari // ------------------------------------------------------ 146a515125bSLeila Ghaffari // Create the PETSc context 147a515125bSLeila Ghaffari // ------------------------------------------------------ 148a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 149a515125bSLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 150a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 151a515125bSLeila Ghaffari PetscScalar Joule; 152a515125bSLeila Ghaffari 153a515125bSLeila Ghaffari // ------------------------------------------------------ 154a515125bSLeila Ghaffari // Command line Options 155a515125bSLeila Ghaffari // ------------------------------------------------------ 1561485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 157a515125bSLeila Ghaffari // -- Physics 1582b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 159a515125bSLeila Ghaffari PetscBool translation; 1602b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type, 1612b916ea7SJeremy L Thompson &translation)); 162a515125bSLeila Ghaffari PetscInt n = problem->dim; 163a515125bSLeila Ghaffari PetscBool user_wind; 1642b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 165c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 1662b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 167059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 168059d1c40SJames Wright NULL)); 1692b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 170c51f031aSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes, 171c51f031aSJames Wright (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 1729529d636SJames Wright bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE; 1739529d636SJames Wright PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type, 1749529d636SJames Wright (PetscEnum *)&bubble_continuity_type, NULL)); 1752b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 17657272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 17757272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 17857272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 17957272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL)); 1802b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 181a515125bSLeila Ghaffari 182a515125bSLeila Ghaffari // -- Units 1832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 184a515125bSLeila Ghaffari meter = fabs(meter); 1852b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 186a515125bSLeila Ghaffari kilogram = fabs(kilogram); 1872b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 188a515125bSLeila Ghaffari second = fabs(second); 189a515125bSLeila Ghaffari 190a515125bSLeila Ghaffari // -- Warnings 191a515125bSLeila Ghaffari if (wind_type == WIND_ROTATION && user_wind) { 1922b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 193a515125bSLeila Ghaffari } 194c51f031aSJames Wright if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) { 195a515125bSLeila Ghaffari wind[2] = 0; 196c51f031aSJames Wright PetscCall( 197c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 198a515125bSLeila Ghaffari } 199a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 2002b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 201a515125bSLeila Ghaffari } 2021485969bSJeremy L Thompson PetscOptionsEnd(); 203a515125bSLeila Ghaffari 204a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 205a78efa86SJames Wright 206a515125bSLeila Ghaffari // ------------------------------------------------------ 207a515125bSLeila Ghaffari // Set up the PETSc context 208a515125bSLeila Ghaffari // ------------------------------------------------------ 209a515125bSLeila Ghaffari // -- Define derived units 210a515125bSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 211a515125bSLeila Ghaffari 212a515125bSLeila Ghaffari user->units->meter = meter; 213a515125bSLeila Ghaffari user->units->kilogram = kilogram; 214a515125bSLeila Ghaffari user->units->second = second; 215a515125bSLeila Ghaffari user->units->Joule = Joule; 216a515125bSLeila Ghaffari 217a515125bSLeila Ghaffari // ------------------------------------------------------ 218a515125bSLeila Ghaffari // Set up the libCEED context 219a515125bSLeila Ghaffari // ------------------------------------------------------ 220a515125bSLeila Ghaffari // -- Scale variables to desired units 221a515125bSLeila Ghaffari E_wind *= Joule; 222a515125bSLeila Ghaffari rc = fabs(rc) * meter; 2239529d636SJames Wright for (PetscInt i = 0; i < problem->dim; i++) { 22405a512bdSLeila Ghaffari wind[i] *= (meter / second); 22505a512bdSLeila Ghaffari domain_size[i] *= meter; 22605a512bdSLeila Ghaffari } 22705a512bdSLeila Ghaffari problem->dm_scale = meter; 228a515125bSLeila Ghaffari 229a515125bSLeila Ghaffari // -- Setup Context 230a515125bSLeila Ghaffari setup_context->rc = rc; 23105a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 23205a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 2339529d636SJames Wright setup_context->lz = problem->dim == 3 ? domain_size[2] : 0.; 234a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 235a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 2369529d636SJames Wright setup_context->wind[2] = problem->dim == 3 ? wind[2] : 0.; 237a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 238c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 239a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 240a515125bSLeila Ghaffari setup_context->time = 0; 241a515125bSLeila Ghaffari 242a515125bSLeila Ghaffari // -- QFunction Context 243a515125bSLeila Ghaffari user->phys->implicit = implicit; 24415a3537eSJed Brown advection_ctx->CtauS = CtauS; 24515a3537eSJed Brown advection_ctx->E_wind = E_wind; 24615a3537eSJed Brown advection_ctx->implicit = implicit; 24715a3537eSJed Brown advection_ctx->strong_form = strong_form; 24815a3537eSJed Brown advection_ctx->stabilization = stab; 24957272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 25057272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 25157272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 252c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 253a515125bSLeila Ghaffari 254b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 255b4c37c5cSJames Wright PetscCallCeed(ceed, 256b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 257b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 258a515125bSLeila Ghaffari 259b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context)); 260b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 261b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc)); 26257272ee0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 26357272ee0SJames Wright "Size of timestep, delta t")); 26415a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = advection_context; 265b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context)); 266b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context)); 267d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 268a515125bSLeila Ghaffari } 269a515125bSLeila Ghaffari 270991aef52SJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) { 2712d49c0afSJames Wright MPI_Comm comm = user->comm; 272b4c37c5cSJames Wright Ceed ceed = user->ceed; 2733636f6a4SJames Wright SetupContextAdv setup_ctx; 27415a3537eSJed Brown AdvectionContext advection_ctx; 275a515125bSLeila Ghaffari 27615a3537eSJed Brown PetscFunctionBeginUser; 277b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx)); 278b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx)); 2792b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 280a515125bSLeila Ghaffari " Problem:\n" 281a515125bSLeila Ghaffari " Problem Name : %s\n" 282a515125bSLeila Ghaffari " Stabilization : %s\n" 283c51f031aSJames Wright " Initial Condition Type : %s\n" 284a515125bSLeila Ghaffari " Bubble Continuity : %s\n" 285a515125bSLeila Ghaffari " Wind Type : %s\n", 286c51f031aSJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type], 287c51f031aSJames Wright BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type])); 288a515125bSLeila Ghaffari 28915a3537eSJed Brown if (setup_ctx->wind_type == WIND_TRANSLATION) { 2909529d636SJames Wright switch (problem->dim) { 2919529d636SJames Wright case 2: 2929529d636SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1])); 2939529d636SJames Wright break; 2949529d636SJames Wright case 3: 2959529d636SJames Wright PetscCall( 2969529d636SJames Wright PetscPrintf(comm, " Background Wind : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2])); 2979529d636SJames Wright break; 2989529d636SJames Wright } 299a515125bSLeila Ghaffari } 300b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx)); 301b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx)); 302d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 303a515125bSLeila Ghaffari } 304