1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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 18b4e9a8f8SJames Wright // @brief Create CeedOperator for stabilized mass KSP for explicit timestepping 19b4e9a8f8SJames Wright // 20b4e9a8f8SJames Wright // Only used for SUPG stabilization 21b4e9a8f8SJames Wright PetscErrorCode CreateKSPMassOperator_AdvectionStabilized(User user, CeedOperator *op_mass) { 22b4e9a8f8SJames Wright Ceed ceed = user->ceed; 23b4e9a8f8SJames Wright CeedInt num_comp_q, q_data_size; 24b4e9a8f8SJames Wright CeedQFunction qf_mass; 25b4e9a8f8SJames Wright CeedElemRestriction elem_restr_q, elem_restr_qd_i; 26b4e9a8f8SJames Wright CeedBasis basis_q; 27b4e9a8f8SJames Wright CeedVector q_data; 28b4e9a8f8SJames Wright CeedQFunctionContext qf_ctx = NULL; 29b4e9a8f8SJames Wright PetscInt dim; 30b4e9a8f8SJames Wright 31b4e9a8f8SJames Wright PetscFunctionBeginUser; 32b4e9a8f8SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 33b4e9a8f8SJames Wright { // Get restriction and basis from the RHS function 34b4e9a8f8SJames Wright CeedOperator *sub_ops; 35b4e9a8f8SJames Wright CeedOperatorField field; 36b4e9a8f8SJames Wright PetscInt sub_op_index = 0; // will be 0 for the volume op 37b4e9a8f8SJames Wright 38*ed094490SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(user->op_rhs_ctx->op, &sub_ops)); 39b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field)); 40681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_q, &basis_q, NULL)); 41b4e9a8f8SJames Wright 42b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field)); 43681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_qd_i, NULL, &q_data)); 44b4e9a8f8SJames Wright 45b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorGetContext(sub_ops[sub_op_index], &qf_ctx)); 46b4e9a8f8SJames Wright } 47b4e9a8f8SJames Wright 48b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q)); 49b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size)); 50b4e9a8f8SJames Wright 51b4e9a8f8SJames Wright switch (dim) { 52b4e9a8f8SJames Wright case 2: 53b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection2D, MassFunction_Advection2D_loc, &qf_mass)); 54b4e9a8f8SJames Wright break; 55b4e9a8f8SJames Wright case 3: 56b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, MassFunction_Advection, MassFunction_Advection_loc, &qf_mass)); 57b4e9a8f8SJames Wright break; 58b4e9a8f8SJames Wright } 59b4e9a8f8SJames Wright 60b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_mass, qf_ctx)); 61b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_mass, 0)); 62b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q_dot", 5, CEED_EVAL_INTERP)); 63b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "q", 5, CEED_EVAL_INTERP)); 64b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_mass, "qdata", q_data_size, CEED_EVAL_NONE)); 65b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "v", 5, CEED_EVAL_INTERP)); 66b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_mass, "Grad_v", 5 * dim, CEED_EVAL_GRAD)); 67b4e9a8f8SJames Wright 68b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass)); 69b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q_dot", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 70b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "q", elem_restr_q, basis_q, user->q_ceed)); 71b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data)); 72b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 73b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE)); 74b4e9a8f8SJames Wright 75681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 76681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q)); 77681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i)); 78681d0ea7SJeremy L Thompson PetscCallCeed(ceed, CeedBasisDestroy(&basis_q)); 797de238d3SJeremy L Thompson PetscCallCeed(ceed, CeedQFunctionContextDestroy(&qf_ctx)); 80b4e9a8f8SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 81b4e9a8f8SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 82b4e9a8f8SJames Wright } 83b4e9a8f8SJames Wright 84731c13d7SJames Wright PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 8577841947SLeila Ghaffari WindType wind_type; 867b77ddfdSJames Wright AdvectionICType advectionic_type; 8777841947SLeila Ghaffari BubbleContinuityType bubble_continuity_type; 8877841947SLeila Ghaffari StabilizationType stab; 89b18328c4SJames Wright StabilizationTauType stab_tau; 9097baf651SJames Wright SetupContextAdv setup_context; 9177841947SLeila Ghaffari User user = *(User *)ctx; 92a424bcd0SJames Wright MPI_Comm comm = user->comm; 93a424bcd0SJames Wright Ceed ceed = user->ceed; 9477841947SLeila Ghaffari PetscBool implicit; 95841e4c73SJed Brown AdvectionContext advection_ctx; 96841e4c73SJed Brown CeedQFunctionContext advection_context; 9793639ffbSJames Wright PetscInt dim; 9877841947SLeila Ghaffari 99841e4c73SJed Brown PetscFunctionBeginUser; 1002b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &setup_context)); 1012b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &advection_ctx)); 10293639ffbSJames Wright PetscCall(DMGetDimension(dm, &dim)); 10377841947SLeila Ghaffari 10477841947SLeila Ghaffari // ------------------------------------------------------ 10577841947SLeila Ghaffari // SET UP ADVECTION 10677841947SLeila Ghaffari // ------------------------------------------------------ 10793639ffbSJames Wright switch (dim) { 10893639ffbSJames Wright case 2: 10993639ffbSJames Wright problem->dim = 2; 11093639ffbSJames Wright problem->ics.qfunction = ICsAdvection2d; 11193639ffbSJames Wright problem->ics.qfunction_loc = ICsAdvection2d_loc; 11293639ffbSJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection2d; 11393639ffbSJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection2d_loc; 11493639ffbSJames Wright problem->apply_vol_ifunction.qfunction = IFunction_Advection2d; 11593639ffbSJames Wright problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection2d_loc; 11693639ffbSJames Wright problem->apply_inflow.qfunction = Advection2d_InOutFlow; 11793639ffbSJames Wright problem->apply_inflow.qfunction_loc = Advection2d_InOutFlow_loc; 118de327db4SJames Wright problem->compute_exact_solution_error = PETSC_TRUE; 11993639ffbSJames Wright problem->print_info = PRINT_ADVECTION; 12093639ffbSJames Wright break; 12193639ffbSJames Wright case 3: 12277841947SLeila Ghaffari problem->dim = 3; 12391e5af17SJed Brown problem->ics.qfunction = ICsAdvection; 12491e5af17SJed Brown problem->ics.qfunction_loc = ICsAdvection_loc; 12593639ffbSJames Wright problem->apply_vol_rhs.qfunction = RHS_Advection; 12693639ffbSJames Wright problem->apply_vol_rhs.qfunction_loc = RHS_Advection_loc; 12791e5af17SJed Brown problem->apply_vol_ifunction.qfunction = IFunction_Advection; 12891e5af17SJed Brown problem->apply_vol_ifunction.qfunction_loc = IFunction_Advection_loc; 12991e5af17SJed Brown problem->apply_inflow.qfunction = Advection_InOutFlow; 13091e5af17SJed Brown problem->apply_inflow.qfunction_loc = Advection_InOutFlow_loc; 131de327db4SJames Wright problem->compute_exact_solution_error = PETSC_FALSE; 13277841947SLeila Ghaffari problem->print_info = PRINT_ADVECTION; 13393639ffbSJames Wright break; 13493639ffbSJames Wright } 13577841947SLeila Ghaffari 13677841947SLeila Ghaffari // ------------------------------------------------------ 13777841947SLeila Ghaffari // Create the libCEED context 13877841947SLeila Ghaffari // ------------------------------------------------------ 13977841947SLeila Ghaffari CeedScalar rc = 1000.; // m (Radius of bubble) 14077841947SLeila Ghaffari CeedScalar CtauS = 0.; // dimensionless 141b767acfbSJames Wright PetscBool strong_form = PETSC_FALSE; 14277841947SLeila Ghaffari CeedScalar E_wind = 1.e6; // J 143b18328c4SJames Wright CeedScalar Ctau_a = PetscPowScalarInt(user->app_ctx->degree, 2); 144b18328c4SJames Wright CeedScalar Ctau_t = 0.; 14577841947SLeila Ghaffari PetscReal wind[3] = {1., 0, 0}; // m/s 146d1d77723SJames Wright CeedScalar diffusion_coeff = 0.; 1471864f1c2SLeila Ghaffari PetscReal domain_min[3], domain_max[3], domain_size[3]; 1482b730f8bSJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 14993639ffbSJames Wright for (PetscInt i = 0; i < problem->dim; i++) domain_size[i] = domain_max[i] - domain_min[i]; 1501864f1c2SLeila Ghaffari 15177841947SLeila Ghaffari // ------------------------------------------------------ 15277841947SLeila Ghaffari // Create the PETSc context 15377841947SLeila Ghaffari // ------------------------------------------------------ 15477841947SLeila Ghaffari PetscScalar meter = 1e-2; // 1 meter in scaled length units 15577841947SLeila Ghaffari PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 15677841947SLeila Ghaffari PetscScalar second = 1e-2; // 1 second in scaled time units 15777841947SLeila Ghaffari PetscScalar Joule; 15877841947SLeila Ghaffari 15977841947SLeila Ghaffari // ------------------------------------------------------ 16077841947SLeila Ghaffari // Command line Options 16177841947SLeila Ghaffari // ------------------------------------------------------ 16267490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Options for ADVECTION problem", NULL); 16377841947SLeila Ghaffari // -- Physics 1642b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", NULL, rc, &rc, NULL)); 16577841947SLeila Ghaffari PetscBool translation; 1662b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-wind_type", "Wind type in Advection", NULL, WindTypes, (PetscEnum)(wind_type = WIND_ROTATION), (PetscEnum *)&wind_type, 1672b730f8bSJeremy L Thompson &translation)); 16877841947SLeila Ghaffari PetscInt n = problem->dim; 16977841947SLeila Ghaffari PetscBool user_wind; 1702b730f8bSJeremy L Thompson PetscCall(PetscOptionsRealArray("-wind_translation", "Constant wind vector", NULL, wind, &n, &user_wind)); 171d1d77723SJames Wright PetscCall(PetscOptionsScalar("-diffusion_coeff", "Diffusion coefficient", NULL, diffusion_coeff, &diffusion_coeff, NULL)); 1722b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-CtauS", "Scale coefficient for tau (nondimensional)", NULL, CtauS, &CtauS, NULL)); 173b767acfbSJames Wright PetscCall(PetscOptionsBool("-strong_form", "Strong (true) or weak/integrated by parts (false) advection residual", NULL, strong_form, &strong_form, 174b767acfbSJames Wright NULL)); 1752b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-E_wind", "Total energy of inflow wind", NULL, E_wind, &E_wind, NULL)); 1767b77ddfdSJames Wright PetscCall(PetscOptionsEnum("-advection_ic_type", "Initial condition for Advection problem", NULL, AdvectionICTypes, 1777b77ddfdSJames Wright (PetscEnum)(advectionic_type = ADVECTIONIC_BUBBLE_SPHERE), (PetscEnum *)&advectionic_type, NULL)); 17893639ffbSJames Wright bubble_continuity_type = problem->dim == 3 ? BUBBLE_CONTINUITY_SMOOTH : BUBBLE_CONTINUITY_COSINE; 17993639ffbSJames Wright PetscCall(PetscOptionsEnum("-bubble_continuity", "Smooth, back_sharp, or thick", NULL, BubbleContinuityTypes, (PetscEnum)bubble_continuity_type, 18093639ffbSJames Wright (PetscEnum *)&bubble_continuity_type, NULL)); 1812b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-stab", "Stabilization method", NULL, StabilizationTypes, (PetscEnum)(stab = STAB_NONE), (PetscEnum *)&stab, NULL)); 182b18328c4SJames Wright PetscCall(PetscOptionsEnum("-stab_tau", "Stabilization constant, tau", NULL, StabilizationTauTypes, (PetscEnum)(stab_tau = STAB_TAU_CTAU), 183b18328c4SJames Wright (PetscEnum *)&stab_tau, NULL)); 184b18328c4SJames Wright PetscCall(PetscOptionsScalar("-Ctau_t", "Stabilization time constant", NULL, Ctau_t, &Ctau_t, NULL)); 185b18328c4SJames Wright PetscCall(PetscOptionsScalar("-Ctau_a", "Coefficient for the stabilization ", NULL, Ctau_a, &Ctau_a, NULL)); 1862b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", NULL, implicit = PETSC_FALSE, &implicit, NULL)); 18777841947SLeila Ghaffari 18877841947SLeila Ghaffari // -- Units 1892b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_meter", "1 meter in scaled length units", NULL, meter, &meter, NULL)); 19077841947SLeila Ghaffari meter = fabs(meter); 1912b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_kilogram", "1 kilogram in scaled mass units", NULL, kilogram, &kilogram, NULL)); 19277841947SLeila Ghaffari kilogram = fabs(kilogram); 1932b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-units_second", "1 second in scaled time units", NULL, second, &second, NULL)); 19477841947SLeila Ghaffari second = fabs(second); 19577841947SLeila Ghaffari 19677841947SLeila Ghaffari // -- Warnings 19777841947SLeila Ghaffari if (wind_type == WIND_ROTATION && user_wind) { 1982b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -wind_translation only with -wind_type translation\n")); 19977841947SLeila Ghaffari } 2007b77ddfdSJames Wright if (wind_type == WIND_TRANSLATION && advectionic_type == ADVECTIONIC_BUBBLE_CYLINDER && wind[2] != 0.) { 20177841947SLeila Ghaffari wind[2] = 0; 2027b77ddfdSJames Wright PetscCall( 2037b77ddfdSJames Wright PetscPrintf(comm, "Warning! Background wind in the z direction should be zero (-wind_translation x,x,0) with -advection_ic_type cylinder\n")); 20477841947SLeila Ghaffari } 20577841947SLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 2062b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, "Warning! Use -CtauS only with -stab su or -stab supg\n")); 20777841947SLeila Ghaffari } 20867490bc6SJeremy L Thompson PetscOptionsEnd(); 20977841947SLeila Ghaffari 210b4e9a8f8SJames Wright if (stab == STAB_SUPG) problem->create_mass_operator = CreateKSPMassOperator_AdvectionStabilized; 211b4e9a8f8SJames Wright 21277841947SLeila Ghaffari // ------------------------------------------------------ 21377841947SLeila Ghaffari // Set up the PETSc context 21477841947SLeila Ghaffari // ------------------------------------------------------ 21577841947SLeila Ghaffari // -- Define derived units 21677841947SLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 21777841947SLeila Ghaffari 21877841947SLeila Ghaffari user->units->meter = meter; 21977841947SLeila Ghaffari user->units->kilogram = kilogram; 22077841947SLeila Ghaffari user->units->second = second; 22177841947SLeila Ghaffari user->units->Joule = Joule; 22277841947SLeila Ghaffari 22377841947SLeila Ghaffari // ------------------------------------------------------ 22477841947SLeila Ghaffari // Set up the libCEED context 22577841947SLeila Ghaffari // ------------------------------------------------------ 22677841947SLeila Ghaffari // -- Scale variables to desired units 22777841947SLeila Ghaffari E_wind *= Joule; 22877841947SLeila Ghaffari rc = fabs(rc) * meter; 22993639ffbSJames Wright for (PetscInt i = 0; i < problem->dim; i++) { 2301864f1c2SLeila Ghaffari wind[i] *= (meter / second); 2311864f1c2SLeila Ghaffari domain_size[i] *= meter; 2321864f1c2SLeila Ghaffari } 2331864f1c2SLeila Ghaffari problem->dm_scale = meter; 23477841947SLeila Ghaffari 23577841947SLeila Ghaffari // -- Setup Context 23677841947SLeila Ghaffari setup_context->rc = rc; 2371864f1c2SLeila Ghaffari setup_context->lx = domain_size[0]; 2381864f1c2SLeila Ghaffari setup_context->ly = domain_size[1]; 23993639ffbSJames Wright setup_context->lz = problem->dim == 3 ? domain_size[2] : 0.; 24077841947SLeila Ghaffari setup_context->wind[0] = wind[0]; 24177841947SLeila Ghaffari setup_context->wind[1] = wind[1]; 24293639ffbSJames Wright setup_context->wind[2] = problem->dim == 3 ? wind[2] : 0.; 24377841947SLeila Ghaffari setup_context->wind_type = wind_type; 2447b77ddfdSJames Wright setup_context->initial_condition_type = advectionic_type; 24577841947SLeila Ghaffari setup_context->bubble_continuity_type = bubble_continuity_type; 24677841947SLeila Ghaffari setup_context->time = 0; 24777841947SLeila Ghaffari 24877841947SLeila Ghaffari // -- QFunction Context 24977841947SLeila Ghaffari user->phys->implicit = implicit; 250841e4c73SJed Brown advection_ctx->CtauS = CtauS; 251841e4c73SJed Brown advection_ctx->E_wind = E_wind; 252841e4c73SJed Brown advection_ctx->implicit = implicit; 253841e4c73SJed Brown advection_ctx->strong_form = strong_form; 254841e4c73SJed Brown advection_ctx->stabilization = stab; 255b18328c4SJames Wright advection_ctx->stabilization_tau = stab_tau; 256b18328c4SJames Wright advection_ctx->Ctau_a = Ctau_a; 257b18328c4SJames Wright advection_ctx->Ctau_t = Ctau_t; 258d1d77723SJames Wright advection_ctx->diffusion_coeff = diffusion_coeff; 25977841947SLeila Ghaffari 260a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &problem->ics.qfunction_context)); 261a424bcd0SJames Wright PetscCallCeed(ceed, 262a424bcd0SJames Wright CeedQFunctionContextSetData(problem->ics.qfunction_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*setup_context), setup_context)); 263a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(problem->ics.qfunction_context, CEED_MEM_HOST, FreeContextPetsc)); 26477841947SLeila Ghaffari 265a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &advection_context)); 266a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(advection_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*advection_ctx), advection_ctx)); 267a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(advection_context, CEED_MEM_HOST, FreeContextPetsc)); 268b18328c4SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(advection_context, "timestep size", offsetof(struct AdvectionContext_, dt), 1, 269b18328c4SJames Wright "Size of timestep, delta t")); 270841e4c73SJed Brown problem->apply_vol_rhs.qfunction_context = advection_context; 271a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_vol_ifunction.qfunction_context)); 272a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(advection_context, &problem->apply_inflow.qfunction_context)); 273ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 27477841947SLeila Ghaffari } 27577841947SLeila Ghaffari 276731c13d7SJames Wright PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx) { 277367c849eSJames Wright MPI_Comm comm = user->comm; 278a424bcd0SJames Wright Ceed ceed = user->ceed; 27997baf651SJames Wright SetupContextAdv setup_ctx; 280841e4c73SJed Brown AdvectionContext advection_ctx; 28177841947SLeila Ghaffari 282841e4c73SJed Brown PetscFunctionBeginUser; 283a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfunction_context, CEED_MEM_HOST, &setup_ctx)); 284a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &advection_ctx)); 2852b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 28677841947SLeila Ghaffari " Problem:\n" 28777841947SLeila Ghaffari " Problem Name : %s\n" 28877841947SLeila Ghaffari " Stabilization : %s\n" 2897b77ddfdSJames Wright " Initial Condition Type : %s\n" 29077841947SLeila Ghaffari " Bubble Continuity : %s\n" 29177841947SLeila Ghaffari " Wind Type : %s\n", 2927b77ddfdSJames Wright app_ctx->problem_name, StabilizationTypes[advection_ctx->stabilization], AdvectionICTypes[setup_ctx->initial_condition_type], 2937b77ddfdSJames Wright BubbleContinuityTypes[setup_ctx->bubble_continuity_type], WindTypes[setup_ctx->wind_type])); 29477841947SLeila Ghaffari 295841e4c73SJed Brown if (setup_ctx->wind_type == WIND_TRANSLATION) { 29693639ffbSJames Wright switch (problem->dim) { 29793639ffbSJames Wright case 2: 29893639ffbSJames Wright PetscCall(PetscPrintf(comm, " Background Wind : %f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1])); 29993639ffbSJames Wright break; 30093639ffbSJames Wright case 3: 30193639ffbSJames Wright PetscCall( 30293639ffbSJames Wright PetscPrintf(comm, " Background Wind : %f,%f,%f\n", setup_ctx->wind[0], setup_ctx->wind[1], setup_ctx->wind[2])); 30393639ffbSJames Wright break; 30493639ffbSJames Wright } 30577841947SLeila Ghaffari } 306a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfunction_context, &setup_ctx)); 307a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &advection_ctx)); 308ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 30977841947SLeila Ghaffari } 310