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 14*b3b24828SJames Wright static PetscErrorCode DivDiffFluxVerifyMesh(DM dm); 15*b3b24828SJames Wright 16991aef52SJames Wright PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 170c373b74SJames Wright Honee honee = *(Honee *)ctx; 180c373b74SJames Wright MPI_Comm comm = honee->comm; 190c373b74SJames Wright Ceed ceed = honee->ceed; 2015a3537eSJed Brown ChannelContext channel_ctx; 2115a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 22e07531f7SJames Wright CeedQFunctionContext channel_qfctx; 23*b3b24828SJames Wright PetscBool use_divdiff_verify_mesh = PETSC_FALSE; 2415a3537eSJed Brown 25bb8a0c61SJames Wright PetscFunctionBeginUser; 26d1c51a42SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 272b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &channel_ctx)); 28bb8a0c61SJames Wright 29*b3b24828SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-mesh_transform_channel_div_diff_projection_verify", &use_divdiff_verify_mesh, NULL)); 30*b3b24828SJames Wright if (use_divdiff_verify_mesh) PetscCall(DivDiffFluxVerifyMesh(dm)); 31*b3b24828SJames Wright 32bb8a0c61SJames Wright // ------------------------------------------------------ 33bb8a0c61SJames Wright // SET UP Channel 34bb8a0c61SJames Wright // ------------------------------------------------------ 35e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 36e07531f7SJames Wright problem->ics.qf_func_ptr = ICsChannel; 37e07531f7SJames Wright problem->ics.qf_loc = ICsChannel_loc; 380c373b74SJames Wright if (honee->phys->state_var == STATEVAR_CONSERVATIVE) { 39e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = Channel_Inflow; 40e07531f7SJames Wright problem->apply_inflow.qf_loc = Channel_Inflow_loc; 41e07531f7SJames Wright problem->apply_outflow.qf_func_ptr = Channel_Outflow; 42e07531f7SJames Wright problem->apply_outflow.qf_loc = Channel_Outflow_loc; 43cbe60e31SLeila Ghaffari } 44bb8a0c61SJames Wright 45bb8a0c61SJames Wright // -- Command Line Options 46bb8a0c61SJames Wright CeedScalar umax = 10.; // m/s 47bb8a0c61SJames Wright CeedScalar theta0 = 300.; // K 48bb8a0c61SJames Wright CeedScalar P0 = 1.e5; // Pa 4910786903SJed Brown PetscReal body_force_scale = 1.; 50bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 512b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-umax", "Centerline velocity of the Channel", NULL, umax, &umax, NULL)); 522b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-theta0", "Wall temperature", NULL, theta0, &theta0, NULL)); 532b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL)); 542b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-body_force_scale", "Multiplier for body force", NULL, body_force_scale = 1, &body_force_scale, NULL)); 55bb8a0c61SJames Wright PetscOptionsEnd(); 56bb8a0c61SJames Wright 570c373b74SJames Wright PetscScalar meter = honee->units->meter; 580c373b74SJames Wright PetscScalar second = honee->units->second; 590c373b74SJames Wright PetscScalar Kelvin = honee->units->Kelvin; 600c373b74SJames Wright PetscScalar Pascal = honee->units->Pascal; 61bb8a0c61SJames Wright 62bb8a0c61SJames Wright theta0 *= Kelvin; 63bb8a0c61SJames Wright P0 *= Pascal; 64bb8a0c61SJames Wright umax *= meter / second; 65bb8a0c61SJames Wright 66bb8a0c61SJames Wright //-- Setup Problem information 67bb8a0c61SJames Wright CeedScalar H, center; 68bb8a0c61SJames Wright { 69bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 702b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 71493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 72bb8a0c61SJames Wright 73bb8a0c61SJames Wright H = 0.5 * domain_size[1] * meter; 74bb8a0c61SJames Wright center = H + domain_min[1] * meter; 75bb8a0c61SJames Wright } 76bb8a0c61SJames Wright 7715a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 78e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 7915a3537eSJed Brown 8015a3537eSJed Brown channel_ctx->center = center; 8115a3537eSJed Brown channel_ctx->H = H; 8215a3537eSJed Brown channel_ctx->theta0 = theta0; 8315a3537eSJed Brown channel_ctx->P0 = P0; 8415a3537eSJed Brown channel_ctx->umax = umax; 850c373b74SJames Wright channel_ctx->implicit = honee->phys->implicit; 8610786903SJed Brown channel_ctx->B = body_force_scale * 2 * umax * newtonian_ig_ctx->mu / (H * H); 87bb8a0c61SJames Wright 88bb8a0c61SJames Wright { 89bb8a0c61SJames Wright // Calculate Body force 902b916ea7SJeremy L Thompson CeedScalar cv = newtonian_ig_ctx->cv, cp = newtonian_ig_ctx->cp; 91bb8a0c61SJames Wright CeedScalar Rd = cp - cv; 92bb8a0c61SJames Wright CeedScalar rho = P0 / (Rd * theta0); 9315a3537eSJed Brown CeedScalar g[] = {channel_ctx->B / rho, 0., 0.}; 942b916ea7SJeremy L Thompson PetscCall(PetscArraycpy(newtonian_ig_ctx->g, g, 3)); 95bb8a0c61SJames Wright } 9615a3537eSJed Brown channel_ctx->newtonian_ctx = *newtonian_ig_ctx; 97e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 98bb8a0c61SJames Wright 990c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &channel_qfctx)); 100e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(channel_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*channel_ctx), channel_ctx)); 101e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(channel_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 10215a3537eSJed Brown 103e07531f7SJames Wright problem->ics.qfctx = channel_qfctx; 104e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &problem->apply_inflow.qfctx)); 105e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(channel_qfctx, &problem->apply_outflow.qfctx)); 106d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 107bb8a0c61SJames Wright } 108*b3b24828SJames Wright 109*b3b24828SJames Wright // This function transforms the mesh coordinates to mimic the mesh used in 110*b3b24828SJames Wright // *A better consistency for low-order stabilized finite element methods* Jansen et. al. 1999 111*b3b24828SJames Wright // which is used to verify the projection of divergence of diffusive flux. See !27 and !94 for more details. 112*b3b24828SJames Wright static PetscErrorCode DivDiffFluxVerifyMesh(DM dm) { 113*b3b24828SJames Wright PetscInt narr, ncoords, dim; 114*b3b24828SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 115*b3b24828SJames Wright PetscScalar *arr_coords; 116*b3b24828SJames Wright Vec vec_coords; 117*b3b24828SJames Wright 118*b3b24828SJames Wright PetscFunctionBeginUser; 119*b3b24828SJames Wright PetscCall(DMGetDimension(dm, &dim)); 120*b3b24828SJames Wright // Get domain boundary information 121*b3b24828SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 122*b3b24828SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 123*b3b24828SJames Wright 124*b3b24828SJames Wright // Get coords array from DM 125*b3b24828SJames Wright PetscCall(DMGetCoordinatesLocal(dm, &vec_coords)); 126*b3b24828SJames Wright PetscCall(VecGetLocalSize(vec_coords, &narr)); 127*b3b24828SJames Wright PetscCall(VecGetArray(vec_coords, &arr_coords)); 128*b3b24828SJames Wright 129*b3b24828SJames Wright PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords; 130*b3b24828SJames Wright ncoords = narr / dim; 131*b3b24828SJames Wright 132*b3b24828SJames Wright // Get mesh information 133*b3b24828SJames Wright PetscInt nmax = 3, faces[3]; 134*b3b24828SJames Wright PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 135*b3b24828SJames Wright // Get element size of the box mesh, for indexing each node 136*b3b24828SJames Wright const PetscReal dxbox = domain_size[0] / (faces[0]); 137*b3b24828SJames Wright 138*b3b24828SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 139*b3b24828SJames Wright PetscInt x_box_index = round(coords[i][0] / dxbox); 140*b3b24828SJames Wright if (x_box_index % 2) { 141*b3b24828SJames Wright coords[i][0] = (x_box_index - 1) * dxbox + 0.5 * dxbox; 142*b3b24828SJames Wright } 143*b3b24828SJames Wright } 144*b3b24828SJames Wright 145*b3b24828SJames Wright PetscCall(VecRestoreArray(vec_coords, &arr_coords)); 146*b3b24828SJames Wright PetscCall(DMSetCoordinatesLocal(dm, vec_coords)); 147*b3b24828SJames Wright PetscFunctionReturn(0); 148*b3b24828SJames Wright } 149