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, 15 void *ctx) { 16 17 PetscInt ierr; 18 User user = *(User *)ctx; 19 MPI_Comm comm = PETSC_COMM_WORLD; 20 ChannelContext channel_ctx; 21 NewtonianIdealGasContext newtonian_ig_ctx; 22 CeedQFunctionContext channel_context; 23 24 PetscFunctionBeginUser; 25 ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 26 ierr = PetscCalloc1(1, &channel_ctx); CHKERRQ(ierr); 27 28 // ------------------------------------------------------ 29 // SET UP Channel 30 // ------------------------------------------------------ 31 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 32 problem->ics.qfunction = ICsChannel; 33 problem->ics.qfunction_loc = ICsChannel_loc; 34 problem->apply_inflow.qfunction = Channel_Inflow; 35 problem->apply_inflow.qfunction_loc = Channel_Inflow_loc; 36 problem->apply_outflow.qfunction = Channel_Outflow; 37 problem->apply_outflow.qfunction_loc = Channel_Outflow_loc; 38 39 // -- Command Line Options 40 CeedScalar umax = 10.; // m/s 41 CeedScalar theta0 = 300.; // K 42 CeedScalar P0 = 1.e5; // Pa 43 PetscReal body_force_scale = 1.; 44 PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 45 ierr = PetscOptionsScalar("-umax", "Centerline velocity of the Channel", 46 NULL, umax, &umax, NULL); CHKERRQ(ierr); 47 ierr = PetscOptionsScalar("-theta0", "Wall temperature", 48 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 49 ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 50 NULL, P0, &P0, NULL); CHKERRQ(ierr); 51 ierr = PetscOptionsReal("-body_force_scale", "Multiplier for body force", 52 NULL, body_force_scale=1, &body_force_scale, NULL); CHKERRQ(ierr); 53 PetscOptionsEnd(); 54 55 PetscScalar meter = user->units->meter; 56 PetscScalar second = user->units->second; 57 PetscScalar Kelvin = user->units->Kelvin; 58 PetscScalar Pascal = user->units->Pascal; 59 60 theta0 *= Kelvin; 61 P0 *= Pascal; 62 umax *= meter / second; 63 64 //-- Setup Problem information 65 CeedScalar H, center; 66 { 67 PetscReal domain_min[3], domain_max[3], domain_size[3]; 68 ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 69 for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 70 71 H = 0.5*domain_size[1]*meter; 72 center = H + domain_min[1]*meter; 73 } 74 75 // Some properties depend on parameters from NewtonianIdealGas 76 CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 77 CEED_MEM_HOST, &newtonian_ig_ctx); 78 79 channel_ctx->center = center; 80 channel_ctx->H = H; 81 channel_ctx->theta0 = theta0; 82 channel_ctx->P0 = P0; 83 channel_ctx->umax = umax; 84 channel_ctx->implicit = user->phys->implicit; 85 channel_ctx->B = body_force_scale * 2 * umax*newtonian_ig_ctx->mu / (H*H); 86 87 { 88 // Calculate Body force 89 CeedScalar cv = newtonian_ig_ctx->cv, 90 cp = newtonian_ig_ctx->cp; 91 CeedScalar Rd = cp - cv; 92 CeedScalar rho = P0 / (Rd*theta0); 93 CeedScalar g[] = {channel_ctx->B / rho, 0., 0.}; 94 ierr = PetscArraycpy(newtonian_ig_ctx->g, g, 3); CHKERRQ(ierr); 95 } 96 channel_ctx->newtonian_ctx = *newtonian_ig_ctx; 97 CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 98 &newtonian_ig_ctx); 99 100 CeedQFunctionContextCreate(user->ceed, &channel_context); 101 CeedQFunctionContextSetData(channel_context, CEED_MEM_HOST, 102 CEED_USE_POINTER, 103 sizeof(*channel_ctx), channel_ctx); 104 CeedQFunctionContextSetDataDestroy(channel_context, CEED_MEM_HOST, 105 FreeContextPetsc); 106 107 problem->ics.qfunction_context = channel_context; 108 CeedQFunctionContextReferenceCopy(channel_context, 109 &problem->apply_inflow.qfunction_context); 110 CeedQFunctionContextReferenceCopy(channel_context, 111 &problem->apply_outflow.qfunction_context); 112 PetscFunctionReturn(0); 113 } 114