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 16*149fb536SJames Wright #include <navierstokes.h> 17a515125bSLeila Ghaffari 18a78efa86SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 19a78efa86SJames Wright // 20a78efa86SJames Wright // Only used for SUPG stabilization 21a78efa86SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) { 22a78efa86SJames Wright Ceed ceed = user->ceed; 23a78efa86SJames Wright CeedInt num_comp_q, q_data_size; 2487d3884fSJeremy L Thompson CeedQFunction qf_mass = NULL; 25a78efa86SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 26a78efa86SJames Wright CeedBasis basis_q; 27a78efa86SJames Wright CeedVector q_data; 28a78efa86SJames Wright CeedQFunctionContext qf_ctx = NULL; 29a78efa86SJames Wright PetscInt dim; 30a78efa86SJames Wright 31a78efa86SJames Wright PetscFunctionBeginUser; 32a78efa86SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 33a78efa86SJames Wright { // Get restriction and basis from the RHS function 34a78efa86SJames Wright CeedOperator *sub_ops; 35a78efa86SJames Wright CeedOperatorField field; 36a78efa86SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 37a78efa86SJames Wright 38a78efa86SJames Wright PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops)); 39a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 40a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q)); 41a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q)); 42a78efa86SJames Wright 43a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 44a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i)); 45a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data)); 46a78efa86SJames Wright 47a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 48a78efa86SJames Wright } 49a78efa86SJames Wright 50a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 51a78efa86SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 52a78efa86SJames Wright 53a78efa86SJames Wright switch (dim) { 54a78efa86SJames Wright case 2: 55a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 56a78efa86SJames Wright break; 57a78efa86SJames Wright case 3: 58a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 59a78efa86SJames Wright break; 60a78efa86SJames Wright } 61a78efa86SJames Wright 62a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 63a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 64a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 65a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 66a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 67a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 68a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 69a78efa86SJames Wright 70a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 71a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 72a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 73a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 74a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 75a78efa86SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 76a78efa86SJames Wright 77a78efa86SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 78a78efa86SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 79a78efa86SJames Wright } 80a78efa86SJames Wright 81991aef52SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 82a515125bSLeila Ghaffari WindType wind_type; 83c51f031aSJames Wright AdvectionICType advectionic_type; 84a515125bSLeila Ghaffari BubbleContinuityType bubble_continuity_type; 85a515125bSLeila Ghaffari StabilizationType stab; 8657272ee0SJames Wright StabilizationTauType stab_tau; 873636f6a4SJames Wright SetupContextAdv setup_context; 88a515125bSLeila Ghaffari User user = *(User *)ctx; 89b4c37c5cSJames Wright MPI_Comm comm = user->comm; 90b4c37c5cSJames Wright Ceed ceed = user->ceed; 91a515125bSLeila Ghaffari PetscBool implicit; 9215a3537eSJed Brown AdvectionContext advection_ctx; 9315a3537eSJed Brown CeedQFunctionContext advection_context; 949529d636SJames Wright PetscInt dim; 95a515125bSLeila Ghaffari 9615a3537eSJed Brown PetscFunctionBeginUser; 972b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 982b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 999529d636SJames Wright PetscCall(DMGetDimension(dm, &dim)); 100a515125bSLeila Ghaffari 101a515125bSLeila Ghaffari // ------------------------------------------------------ 102a515125bSLeila Ghaffari // SET UP ADVECTION 103a515125bSLeila Ghaffari // ------------------------------------------------------ 1049529d636SJames Wright switch (dim) { 1059529d636SJames Wright case 2: 1069529d636SJames Wright problem->dim = 2; 1079529d636SJames Wright problem->ics.qfunction = ICsAdvection2d; 1089529d636SJames Wright problem->ics.qfunction_loc = ICsAdvection2d_loc; 1099529d636SJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection2d; 1109529d636SJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection2d_loc; 1119529d636SJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Advection2d; 1129529d636SJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc; 1139529d636SJames Wright problem->apply_inflow.qfunction = Advection2d_InOutFlow; 1149529d636SJames Wright problem->apply_inflow.qfunction_loc = Advection2d_InOutFlow_loc; 11558ce1233SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 1169529d636SJames Wright problem->print_info = PRINT_ADVECTION; 1179529d636SJames Wright break; 1189529d636SJames Wright case 3: 119a515125bSLeila Ghaffari problem->dim = 3; 1209785fe93SJed Brown problem->ics.qfunction = ICsAdvection; 1219785fe93SJed Brown problem->ics.qfunction_loc = ICsAdvection_loc; 1229529d636SJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection; 1239529d636SJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection_loc; 1249785fe93SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Advection; 1259785fe93SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc; 1269785fe93SJed Brown problem->apply_inflow.qfunction = Advection_InOutFlow; 1279785fe93SJed Brown problem->apply_inflow.qfunction_loc = Advection_InOutFlow_loc; 12858ce1233SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 129a515125bSLeila Ghaffari problem->print_info = PRINT_ADVECTION; 1309529d636SJames Wright break; 1319529d636SJames Wright } 132a515125bSLeila Ghaffari 133a515125bSLeila Ghaffari // ------------------------------------------------------ 134a515125bSLeila Ghaffari // Create the libCEED context 135a515125bSLeila Ghaffari // ------------------------------------------------------ 136a515125bSLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 137a515125bSLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 138059d1c40SJames Wright PetscBool strong_form = PETSC_FALSE; 139a515125bSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 14057272ee0SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 14157272ee0SJames Wright CeedScalar Ctau_t = 0.; 142a515125bSLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 143c8d249deSJames Wright CeedScalar diffusion_coeff = 0.; 14487d3884fSJeremy L Thompson PetscReal domain_min[3], domain_max[3], domain_size[3] = {0.}; 1452b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 1469529d636SJames Wright for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 14705a512bdSLeila Ghaffari 148a515125bSLeila Ghaffari // ------------------------------------------------------ 149a515125bSLeila Ghaffari // Create the PETSc context 150a515125bSLeila Ghaffari // ------------------------------------------------------ 151a515125bSLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 152a515125bSLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 153a515125bSLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 154a515125bSLeila Ghaffari PetscScalar Joule; 155a515125bSLeila Ghaffari 156a515125bSLeila Ghaffari // ------------------------------------------------------ 157a515125bSLeila Ghaffari // Command line Options 158a515125bSLeila Ghaffari // ------------------------------------------------------ 1591485969bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 160a515125bSLeila Ghaffari // -- Physics 1612b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 162a515125bSLeila Ghaffari PetscBool translation; 1632b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type, 1642b916ea7SJeremy L Thompson &translation)); 165a515125bSLeila Ghaffari PetscInt n = problem->dim; 166a515125bSLeila Ghaffari PetscBool user_wind; 1672b916ea7SJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 168c8d249deSJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 1692b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 170059d1c40SJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 171059d1c40SJames Wright NULL)); 1722b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 173c51f031aSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes, 174c51f031aSJames Wright (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 1759529d636SJames Wright bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE; 1769529d636SJames Wright PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type, 1779529d636SJames Wright (PetscEnum *)&bubble_continuity_type, NULL)); 1782b916ea7SJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 17957272ee0SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 18057272ee0SJames Wright (PetscEnum *)&stab_tau, NULL)); 18157272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 18257272ee0SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL)); 1832b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 184a515125bSLeila Ghaffari 185a515125bSLeila Ghaffari // -- Units 1862b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 187a515125bSLeila Ghaffari meter = fabs(meter); 1882b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 189a515125bSLeila Ghaffari kilogram = fabs(kilogram); 1902b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 191a515125bSLeila Ghaffari second = fabs(second); 192a515125bSLeila Ghaffari 193a515125bSLeila Ghaffari // -- Warnings 194a515125bSLeila Ghaffari if (wind_type == WIND_ROTATION && user_wind) { 1952b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 196a515125bSLeila Ghaffari } 197c51f031aSJames Wright if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) { 198a515125bSLeila Ghaffari wind[2] = 0; 199c51f031aSJames Wright PetscCall( 200c51f031aSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 201a515125bSLeila Ghaffari } 202a515125bSLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 2032b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 204a515125bSLeila Ghaffari } 2051485969bSJeremy L Thompson PetscOptionsEnd(); 206a515125bSLeila Ghaffari 207a78efa86SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 208a78efa86SJames Wright 209a515125bSLeila Ghaffari // ------------------------------------------------------ 210a515125bSLeila Ghaffari // Set up the PETSc context 211a515125bSLeila Ghaffari // ------------------------------------------------------ 212a515125bSLeila Ghaffari // -- Define derived units 213a515125bSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 214a515125bSLeila Ghaffari 215a515125bSLeila Ghaffari user->units->meter = meter; 216a515125bSLeila Ghaffari user->units->kilogram = kilogram; 217a515125bSLeila Ghaffari user->units->second = second; 218a515125bSLeila Ghaffari user->units->Joule = Joule; 219a515125bSLeila Ghaffari 220a515125bSLeila Ghaffari // ------------------------------------------------------ 221a515125bSLeila Ghaffari // Set up the libCEED context 222a515125bSLeila Ghaffari // ------------------------------------------------------ 223a515125bSLeila Ghaffari // -- Scale variables to desired units 224a515125bSLeila Ghaffari E_wind *= Joule; 225a515125bSLeila Ghaffari rc = fabs(rc) * meter; 2269529d636SJames Wright for (PetscInt i = 0; i < problem->dim; i++) { 22705a512bdSLeila Ghaffari wind[i] *= (meter / second); 22805a512bdSLeila Ghaffari domain_size[i] *= meter; 22905a512bdSLeila Ghaffari } 23005a512bdSLeila Ghaffari problem->dm_scale = meter; 231a515125bSLeila Ghaffari 232a515125bSLeila Ghaffari // -- Setup Context 233a515125bSLeila Ghaffari setup_context->rc = rc; 23405a512bdSLeila Ghaffari setup_context->lx = domain_size[0]; 23505a512bdSLeila Ghaffari setup_context->ly = domain_size[1]; 2369529d636SJames Wright setup_context->lz = problem->dim == 3 ? domain_size[2] : 0.; 237a515125bSLeila Ghaffari setup_context->wind[0] = wind[0]; 238a515125bSLeila Ghaffari setup_context->wind[1] = wind[1]; 2399529d636SJames Wright setup_context->wind[2] = problem->dim == 3 ? wind[2] : 0.; 240a515125bSLeila Ghaffari setup_context->wind_type = wind_type; 241c51f031aSJames Wright setup_context->initial_condition_type = advectionic_type; 242a515125bSLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 243a515125bSLeila Ghaffari setup_context->time = 0; 244a515125bSLeila Ghaffari 245a515125bSLeila Ghaffari // -- QFunction Context 246a515125bSLeila Ghaffari user->phys->implicit = implicit; 24715a3537eSJed Brown advection_ctx->CtauS = CtauS; 24815a3537eSJed Brown advection_ctx->E_wind = E_wind; 24915a3537eSJed Brown advection_ctx->implicit = implicit; 25015a3537eSJed Brown advection_ctx->strong_form = strong_form; 25115a3537eSJed Brown advection_ctx->stabilization = stab; 25257272ee0SJames Wright advection_ctx->stabilization_tau = stab_tau; 25357272ee0SJames Wright advection_ctx->Ctau_a = Ctau_a; 25457272ee0SJames Wright advection_ctx->Ctau_t = Ctau_t; 255c8d249deSJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 256a515125bSLeila Ghaffari 257b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 258b4c37c5cSJames Wright PetscCallCeed(ceed, 259b4c37c5cSJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 260b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 261a515125bSLeila Ghaffari 262b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context)); 263b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 264b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc)); 26557272ee0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 26657272ee0SJames Wright "Size of timestep, delta t")); 26715a3537eSJed Brown problem->apply_vol_rhs.qfunction_context = advection_context; 268b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context)); 269b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context)); 270d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 271a515125bSLeila Ghaffari } 272a515125bSLeila Ghaffari 273991aef52SJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) { 2742d49c0afSJames Wright MPI_Comm comm = user->comm; 275b4c37c5cSJames Wright Ceed ceed = user->ceed; 2763636f6a4SJames Wright SetupContextAdv setup_ctx; 27715a3537eSJed Brown AdvectionContext advection_ctx; 278a515125bSLeila Ghaffari 27915a3537eSJed Brown PetscFunctionBeginUser; 280b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx)); 281b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx)); 2822b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 283a515125bSLeila Ghaffari " Problem:\n" 284a515125bSLeila Ghaffari " Problem Name : %s\n" 285a515125bSLeila Ghaffari " Stabilization : %s\n" 286c51f031aSJames Wright " Initial Condition Type : %s\n" 287a515125bSLeila Ghaffari " Bubble Continuity : %s\n" 288a515125bSLeila Ghaffari " Wind Type : %s\n", 289c51f031aSJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type], 290c51f031aSJames Wright BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type])); 291a515125bSLeila Ghaffari 29215a3537eSJed Brown if (setup_ctx->wind_type == WIND_TRANSLATION) { 2939529d636SJames Wright switch (problem->dim) { 2949529d636SJames Wright case 2: 2959529d636SJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1])); 2969529d636SJames Wright break; 2979529d636SJames Wright case 3: 2989529d636SJames Wright PetscCall( 2999529d636SJames Wright PetscPrintf(comm, " Background Wind : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2])); 3009529d636SJames Wright break; 3019529d636SJames Wright } 302a515125bSLeila Ghaffari } 303b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx)); 304b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx)); 305d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 306a515125bSLeila Ghaffari } 307