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