1*bb8a0c61SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*bb8a0c61SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*bb8a0c61SJames Wright // 4*bb8a0c61SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5*bb8a0c61SJames Wright // 6*bb8a0c61SJames Wright // This file is part of CEED: http://github.com/ceed 7*bb8a0c61SJames Wright 8*bb8a0c61SJames Wright /// @file 9*bb8a0c61SJames Wright /// Utility functions for setting up Blasius Boundary Layer 10*bb8a0c61SJames Wright 11*bb8a0c61SJames Wright #include "../navierstokes.h" 12*bb8a0c61SJames Wright #include "../qfunctions/blasius.h" 13*bb8a0c61SJames Wright 14*bb8a0c61SJames Wright #ifndef blasius_context_struct 15*bb8a0c61SJames Wright #define blasius_context_struct 16*bb8a0c61SJames Wright typedef struct BlasiusContext_ *BlasiusContext; 17*bb8a0c61SJames Wright struct BlasiusContext_ { 18*bb8a0c61SJames Wright bool implicit; // !< Using implicit timesteping or not 19*bb8a0c61SJames Wright CeedScalar delta0; // !< Boundary layer height at inflow 20*bb8a0c61SJames Wright CeedScalar Uinf; // !< Velocity at boundary layer edge 21*bb8a0c61SJames Wright CeedScalar P0; // !< Pressure at outflow 22*bb8a0c61SJames Wright CeedScalar theta0; // !< Temperature at inflow 23*bb8a0c61SJames Wright struct NewtonianIdealGasContext_ newtonian_ctx; 24*bb8a0c61SJames Wright }; 25*bb8a0c61SJames Wright #endif 26*bb8a0c61SJames Wright 27*bb8a0c61SJames Wright #ifndef M_PI 28*bb8a0c61SJames Wright #define M_PI 3.14159265358979323846 29*bb8a0c61SJames Wright #endif 30*bb8a0c61SJames Wright 31*bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius 32*bb8a0c61SJames Wright * 33*bb8a0c61SJames Wright * Modifies mesh such that `N` elements are within 1.2*`delta0` with a geometric 34*bb8a0c61SJames Wright * growth ratio of `growth`. Excess elements are then geometrically distributed 35*bb8a0c61SJames Wright * to the top surface. 36*bb8a0c61SJames Wright * 37*bb8a0c61SJames Wright * The top surface is also angled downwards, so that it may be used as an 38*bb8a0c61SJames Wright * outflow. It's angle is controlled by top_angle (in units of degrees). 39*bb8a0c61SJames Wright */ 40*bb8a0c61SJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N, 41*bb8a0c61SJames Wright PetscReal refine_height, PetscReal top_angle) { 42*bb8a0c61SJames Wright 43*bb8a0c61SJames Wright PetscInt ierr, narr, ncoords; 44*bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 45*bb8a0c61SJames Wright PetscScalar *arr_coords; 46*bb8a0c61SJames Wright Vec vec_coords; 47*bb8a0c61SJames Wright PetscFunctionBeginUser; 48*bb8a0c61SJames Wright 49*bb8a0c61SJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 50*bb8a0c61SJames Wright 51*bb8a0c61SJames Wright // Get domain boundary information 52*bb8a0c61SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 53*bb8a0c61SJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 54*bb8a0c61SJames Wright 55*bb8a0c61SJames Wright // Get coords array from DM 56*bb8a0c61SJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 57*bb8a0c61SJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 58*bb8a0c61SJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 59*bb8a0c61SJames Wright 60*bb8a0c61SJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 61*bb8a0c61SJames Wright ncoords = narr/dim; 62*bb8a0c61SJames Wright 63*bb8a0c61SJames Wright // Get mesh information 64*bb8a0c61SJames Wright PetscInt nmax = 3, faces[3]; 65*bb8a0c61SJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 66*bb8a0c61SJames Wright NULL); CHKERRQ(ierr); 67*bb8a0c61SJames Wright 68*bb8a0c61SJames Wright // Calculate the first element height 69*bb8a0c61SJames Wright PetscReal dybox = domain_size[1]/faces[1]; 70*bb8a0c61SJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 71*bb8a0c61SJames Wright 72*bb8a0c61SJames Wright // Calculate log of sizing outside BL 73*bb8a0c61SJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 74*bb8a0c61SJames Wright 75*bb8a0c61SJames Wright for(int i=0; i<ncoords; i++) { 76*bb8a0c61SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 77*bb8a0c61SJames Wright if(y_box_index <= N) { 78*bb8a0c61SJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 79*bb8a0c61SJames Wright dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1); 80*bb8a0c61SJames Wright } else { 81*bb8a0c61SJames Wright PetscInt j = y_box_index - N; 82*bb8a0c61SJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 83*bb8a0c61SJames Wright exp(log(refine_height) + logdy*j); 84*bb8a0c61SJames Wright } 85*bb8a0c61SJames Wright } 86*bb8a0c61SJames Wright 87*bb8a0c61SJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 88*bb8a0c61SJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 89*bb8a0c61SJames Wright 90*bb8a0c61SJames Wright PetscFunctionReturn(0); 91*bb8a0c61SJames Wright } 92*bb8a0c61SJames Wright 93*bb8a0c61SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *setup_ctx, 94*bb8a0c61SJames Wright void *ctx) { 95*bb8a0c61SJames Wright 96*bb8a0c61SJames Wright PetscInt ierr; 97*bb8a0c61SJames Wright ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr); 98*bb8a0c61SJames Wright User user = *(User *)ctx; 99*bb8a0c61SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 100*bb8a0c61SJames Wright PetscFunctionBeginUser; 101*bb8a0c61SJames Wright ierr = PetscCalloc1(1, &user->phys->blasius_ctx); CHKERRQ(ierr); 102*bb8a0c61SJames Wright 103*bb8a0c61SJames Wright // ------------------------------------------------------ 104*bb8a0c61SJames Wright // SET UP Blasius 105*bb8a0c61SJames Wright // ------------------------------------------------------ 106*bb8a0c61SJames Wright problem->ics = ICsBlasius; 107*bb8a0c61SJames Wright problem->ics_loc = ICsBlasius_loc; 108*bb8a0c61SJames Wright problem->apply_inflow = Blasius_Inflow; 109*bb8a0c61SJames Wright problem->apply_inflow_loc = Blasius_Inflow_loc; 110*bb8a0c61SJames Wright problem->apply_outflow = Blasius_Outflow; 111*bb8a0c61SJames Wright problem->apply_outflow_loc = Blasius_Outflow_loc; 112*bb8a0c61SJames Wright problem->setup_ctx = SetupContext_BLASIUS; 113*bb8a0c61SJames Wright 114*bb8a0c61SJames Wright // CeedScalar mu = .04; // Pa s, dynamic viscosity 115*bb8a0c61SJames Wright CeedScalar mu = 1.8e-5; // Pa s, dynamic viscosity 116*bb8a0c61SJames Wright CeedScalar Uinf = 40; // m/s 117*bb8a0c61SJames Wright CeedScalar delta0 = 4.2e-4; // m 118*bb8a0c61SJames Wright PetscReal refine_height = 5.9e-4; // m 119*bb8a0c61SJames Wright PetscReal growth = 1.08; // [-] 120*bb8a0c61SJames Wright PetscInt Ndelta = 45; // [-] 121*bb8a0c61SJames Wright PetscReal top_angle = 5; // degrees 122*bb8a0c61SJames Wright CeedScalar theta0 = 288.; // K 123*bb8a0c61SJames Wright CeedScalar P0 = 1.01e5; // Pa 124*bb8a0c61SJames Wright 125*bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 126*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 127*bb8a0c61SJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 128*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 129*bb8a0c61SJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 130*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 131*bb8a0c61SJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 132*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 133*bb8a0c61SJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 134*bb8a0c61SJames Wright ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge", 135*bb8a0c61SJames Wright NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr); 136*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-refine_height", 137*bb8a0c61SJames Wright "Height of boundary layer mesh refinement", 138*bb8a0c61SJames Wright NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr); 139*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-growth", 140*bb8a0c61SJames Wright "Geometric growth rate of boundary layer mesh", 141*bb8a0c61SJames Wright NULL, growth, &growth, NULL); CHKERRQ(ierr); 142*bb8a0c61SJames Wright ierr = PetscOptionsScalar("-top_angle", 143*bb8a0c61SJames Wright "Geometric top_angle rate of boundary layer mesh", 144*bb8a0c61SJames Wright NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr); 145*bb8a0c61SJames Wright PetscOptionsEnd(); 146*bb8a0c61SJames Wright 147*bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 148*bb8a0c61SJames Wright PetscScalar second = user->units->second; 149*bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 150*bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 151*bb8a0c61SJames Wright 152*bb8a0c61SJames Wright mu *= Pascal * second; 153*bb8a0c61SJames Wright theta0 *= Kelvin; 154*bb8a0c61SJames Wright P0 *= Pascal; 155*bb8a0c61SJames Wright Uinf *= meter / second; 156*bb8a0c61SJames Wright delta0 *= meter; 157*bb8a0c61SJames Wright 158*bb8a0c61SJames Wright ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle); 159*bb8a0c61SJames Wright CHKERRQ(ierr); 160*bb8a0c61SJames Wright 161*bb8a0c61SJames Wright user->phys->blasius_ctx->Uinf = Uinf; 162*bb8a0c61SJames Wright user->phys->blasius_ctx->delta0 = delta0; 163*bb8a0c61SJames Wright user->phys->blasius_ctx->theta0 = theta0; 164*bb8a0c61SJames Wright user->phys->blasius_ctx->P0 = P0; 165*bb8a0c61SJames Wright user->phys->blasius_ctx->implicit = user->phys->implicit; 166*bb8a0c61SJames Wright 167*bb8a0c61SJames Wright user->phys->newtonian_ig_ctx->mu = mu; 168*bb8a0c61SJames Wright user->phys->blasius_ctx->newtonian_ctx = *user->phys->newtonian_ig_ctx; 169*bb8a0c61SJames Wright PetscFunctionReturn(0); 170*bb8a0c61SJames Wright } 171*bb8a0c61SJames Wright 172*bb8a0c61SJames Wright PetscErrorCode SetupContext_BLASIUS(Ceed ceed, CeedData ceed_data, 173*bb8a0c61SJames Wright AppCtx app_ctx, SetupContext setup_ctx, Physics phys) { 174*bb8a0c61SJames Wright PetscFunctionBeginUser; 175*bb8a0c61SJames Wright PetscInt ierr; 176*bb8a0c61SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->setup_context); 177*bb8a0c61SJames Wright CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST, 178*bb8a0c61SJames Wright CEED_USE_POINTER, 179*bb8a0c61SJames Wright sizeof(*setup_ctx), setup_ctx); 180*bb8a0c61SJames Wright ierr = SetupContext_NEWTONIAN_IG(ceed, ceed_data, app_ctx, setup_ctx, phys); 181*bb8a0c61SJames Wright CHKERRQ(ierr); 182*bb8a0c61SJames Wright 183*bb8a0c61SJames Wright CeedQFunctionContextCreate(ceed, &ceed_data->blasius_context); 184*bb8a0c61SJames Wright CeedQFunctionContextSetData(ceed_data->blasius_context, CEED_MEM_HOST, 185*bb8a0c61SJames Wright CEED_USE_POINTER, sizeof(*phys->blasius_ctx), phys->blasius_ctx); 186*bb8a0c61SJames Wright phys->has_neumann = PETSC_TRUE; 187*bb8a0c61SJames Wright if (ceed_data->qf_ics) 188*bb8a0c61SJames Wright CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->blasius_context); 189*bb8a0c61SJames Wright if (ceed_data->qf_apply_inflow) 190*bb8a0c61SJames Wright CeedQFunctionSetContext(ceed_data->qf_apply_inflow, ceed_data->blasius_context); 191*bb8a0c61SJames Wright if (ceed_data->qf_apply_outflow) 192*bb8a0c61SJames Wright CeedQFunctionSetContext(ceed_data->qf_apply_outflow, 193*bb8a0c61SJames Wright ceed_data->blasius_context); 194*bb8a0c61SJames Wright PetscFunctionReturn(0); 195*bb8a0c61SJames Wright } 196*bb8a0c61SJames Wright 197