1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up Channel flow 10 11 #include "../navierstokes.h" 12 #include "../qfunctions/channel.h" 13 14 PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx) { 15 16 PetscInt ierr; 17 User user = *(User *)ctx; 18 MPI_Comm comm = PETSC_COMM_WORLD; 19 ChannelContext channel_ctx; 20 NewtonianIdealGasContext newtonian_ig_ctx; 21 CeedQFunctionContext channel_context; 22 23 PetscFunctionBeginUser; 24 ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 25 ierr = PetscCalloc1(1, &channel_ctx); CHKERRQ(ierr); 26 27 // ------------------------------------------------------ 28 // SET UP Channel 29 // ------------------------------------------------------ 30 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 31 problem->ics.qfunction = ICsChannel; 32 problem->ics.qfunction_loc = ICsChannel_loc; 33 problem->apply_inflow.qfunction = Channel_Inflow; 34 problem->apply_inflow.qfunction_loc = Channel_Inflow_loc; 35 problem->apply_outflow.qfunction = Channel_Outflow; 36 problem->apply_outflow.qfunction_loc = Channel_Outflow_loc; 37 38 // -- Command Line Options 39 CeedScalar umax = 10.; // m/s 40 CeedScalar theta0 = 300.; // K 41 CeedScalar P0 = 1.e5; // Pa 42 PetscReal body_force_scale = 1.; 43 PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 44 ierr = PetscOptionsScalar("-umax", "Centerline velocity of the Channel", 45 NULL, umax, &umax, NULL); CHKERRQ(ierr); 46 ierr = PetscOptionsScalar("-theta0", "Wall temperature", 47 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 48 ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 49 NULL, P0, &P0, NULL); CHKERRQ(ierr); 50 ierr = PetscOptionsReal("-body_force_scale", "Multiplier for body force", 51 NULL, body_force_scale=1, &body_force_scale, NULL); CHKERRQ(ierr); 52 PetscOptionsEnd(); 53 54 PetscScalar meter = user->units->meter; 55 PetscScalar second = user->units->second; 56 PetscScalar Kelvin = user->units->Kelvin; 57 PetscScalar Pascal = user->units->Pascal; 58 59 theta0 *= Kelvin; 60 P0 *= Pascal; 61 umax *= meter / second; 62 63 //-- Setup Problem information 64 CeedScalar H, center; 65 { 66 PetscReal domain_min[3], domain_max[3], domain_size[3]; 67 ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 68 for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 69 70 H = 0.5*domain_size[1]*meter; 71 center = H + domain_min[1]*meter; 72 } 73 74 // Some properties depend on parameters from NewtonianIdealGas 75 CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 76 CEED_MEM_HOST, &newtonian_ig_ctx); 77 78 channel_ctx->center = center; 79 channel_ctx->H = H; 80 channel_ctx->theta0 = theta0; 81 channel_ctx->P0 = P0; 82 channel_ctx->umax = umax; 83 channel_ctx->implicit = user->phys->implicit; 84 channel_ctx->B = body_force_scale * 2 * umax*newtonian_ig_ctx->mu / (H*H); 85 86 { 87 // Calculate Body force 88 CeedScalar cv = newtonian_ig_ctx->cv, 89 cp = newtonian_ig_ctx->cp; 90 CeedScalar Rd = cp - cv; 91 CeedScalar rho = P0 / (Rd*theta0); 92 CeedScalar g[] = {channel_ctx->B / rho, 0., 0.}; 93 ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 94 } 95 channel_ctx->newtonian_ctx = *newtonian_ig_ctx; 96 CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 97 &newtonian_ig_ctx); 98 99 CeedQFunctionContextCreate(user->ceed, &channel_context); 100 CeedQFunctionContextSetData(channel_context, CEED_MEM_HOST, 101 CEED_USE_POINTER, 102 sizeof(*channel_ctx), channel_ctx); 103 CeedQFunctionContextSetDataDestroy(channel_context, CEED_MEM_HOST, 104 FreeContextPetsc); 105 106 problem->ics.qfunction_context = channel_context; 107 CeedQFunctionContextReferenceCopy(channel_context, 108 &problem->apply_inflow.qfunction_context); 109 CeedQFunctionContextReferenceCopy(channel_context, 110 &problem->apply_outflow.qfunction_context); 111 PetscFunctionReturn(0); 112 } 113