1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3bb8a0c61SJames Wright 4bb8a0c61SJames Wright /// @file 5bb8a0c61SJames Wright /// Utility functions for setting up Channel flow 6bb8a0c61SJames Wright 7bb8a0c61SJames Wright #include "../qfunctions/channel.h" 8bb8a0c61SJames Wright 9e419654dSJeremy L Thompson #include <ceed.h> 10e419654dSJeremy L Thompson #include <petscdm.h> 11e419654dSJeremy L Thompson 12149fb536SJames Wright #include <navierstokes.h> 13bb8a0c61SJames Wright 14991aef52SJames Wright PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 15bb8a0c61SJames Wright User user = *(User *)ctx; 16b4c37c5cSJames Wright MPI_Comm comm = user->comm; 17b4c37c5cSJames Wright Ceed ceed = user->ceed; 1815a3537eSJed Brown ChannelContext channel_ctx; 1915a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 20*e07531f7SJames Wright CeedQFunctionContext channel_qfctx; 2115a3537eSJed Brown 22bb8a0c61SJames Wright PetscFunctionBeginUser; 23d1c51a42SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 242b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &channel_ctx)); 25bb8a0c61SJames Wright 26bb8a0c61SJames Wright // ------------------------------------------------------ 27bb8a0c61SJames Wright // SET UP Channel 28bb8a0c61SJames Wright // ------------------------------------------------------ 29*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 30*e07531f7SJames Wright problem->ics.qf_func_ptr = ICsChannel; 31*e07531f7SJames Wright problem->ics.qf_loc = ICsChannel_loc; 323636f6a4SJames Wright if (user->phys->state_var == STATEVAR_CONSERVATIVE) { 33*e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Channel_Inflow; 34*e07531f7SJames Wright problem->apply_inflow.qf_loc = Channel_Inflow_loc; 35*e07531f7SJames Wright problem->apply_outflow.qf_func_ptr = Channel_Outflow; 36*e07531f7SJames Wright problem->apply_outflow.qf_loc = Channel_Outflow_loc; 37cbe60e31SLeila Ghaffari } 38bb8a0c61SJames Wright 39bb8a0c61SJames Wright // -- Command Line Options 40bb8a0c61SJames Wright CeedScalar umax = 10.; // m/s 41bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K 42bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa 4310786903SJed Brown PetscReal body_force_scale = 1.; 44bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 452b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-umax", "Centerline velocity of the Channel", NULL, umax, &umax, NULL)); 462b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-theta0", "Wall temperature", NULL, theta0, &theta0, NULL)); 472b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL)); 482b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-body_force_scale", "Multiplier for body force", NULL, body_force_scale = 1, &body_force_scale, NULL)); 49bb8a0c61SJames Wright PetscOptionsEnd(); 50bb8a0c61SJames Wright 51bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 52bb8a0c61SJames Wright PetscScalar second = user->units->second; 53bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 54bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 55bb8a0c61SJames Wright 56bb8a0c61SJames Wright theta0 *= Kelvin; 57bb8a0c61SJames Wright P0 *= Pascal; 58bb8a0c61SJames Wright umax *= meter / second; 59bb8a0c61SJames Wright 60bb8a0c61SJames Wright //-- Setup Problem information 61bb8a0c61SJames Wright CeedScalar H, center; 62bb8a0c61SJames Wright { 63bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 642b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 65493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 66bb8a0c61SJames Wright 67bb8a0c61SJames Wright H = 0.5 * domain_size[1] * meter; 68bb8a0c61SJames Wright center = H + domain_min[1] * meter; 69bb8a0c61SJames Wright } 70bb8a0c61SJames Wright 7115a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 72*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 7315a3537eSJed Brown 7415a3537eSJed Brown channel_ctx->center = center; 7515a3537eSJed Brown channel_ctx->H = H; 7615a3537eSJed Brown channel_ctx->theta0 = theta0; 7715a3537eSJed Brown channel_ctx->P0 = P0; 7815a3537eSJed Brown channel_ctx->umax = umax; 7915a3537eSJed Brown channel_ctx->implicit = user->phys->implicit; 8010786903SJed Brown channel_ctx->B = body_force_scale * 2 * umax * newtonian_ig_ctx->mu / (H * H); 81bb8a0c61SJames Wright 82bb8a0c61SJames Wright { 83bb8a0c61SJames Wright // Calculate Body force 842b916ea7SJeremy L Thompson CeedScalar cv = newtonian_ig_ctx->cv, cp = newtonian_ig_ctx->cp; 85bb8a0c61SJames Wright CeedScalar Rd = cp - cv; 86bb8a0c61SJames Wright CeedScalar rho = P0 / (Rd * theta0); 8715a3537eSJed Brown CeedScalar g[] = {channel_ctx->B / rho, 0., 0.}; 882b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 89bb8a0c61SJames Wright } 9015a3537eSJed Brown channel_ctx->newtonian_ctx = *newtonian_ig_ctx; 91*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 92bb8a0c61SJames Wright 93*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &channel_qfctx)); 94*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(channel_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*channel_ctx), channel_ctx)); 95*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(channel_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 9615a3537eSJed Brown 97*e07531f7SJames Wright problem->ics.qfctx = channel_qfctx; 98*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &problem->apply_inflow.qfctx)); 99*e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &problem->apply_outflow.qfctx)); 100d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 101bb8a0c61SJames Wright } 102