1bb8a0c61SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2bb8a0c61SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3bb8a0c61SJames Wright // 4bb8a0c61SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5bb8a0c61SJames Wright // 6bb8a0c61SJames Wright // This file is part of CEED: http://github.com/ceed 7bb8a0c61SJames Wright 8bb8a0c61SJames Wright /// @file 9bb8a0c61SJames Wright /// Utility functions for setting up Blasius Boundary Layer 10bb8a0c61SJames Wright 11bb8a0c61SJames Wright #include "../navierstokes.h" 12bb8a0c61SJames Wright #include "../qfunctions/blasius.h" 13*493642f1SJames Wright #include "stg_shur14.h" 14bb8a0c61SJames Wright 15bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius 16bb8a0c61SJames Wright * 17*493642f1SJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 18*493642f1SJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 19*493642f1SJames Wright * linearly in logspace to the top surface. 20bb8a0c61SJames Wright * 21bb8a0c61SJames Wright * The top surface is also angled downwards, so that it may be used as an 22*493642f1SJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 23bb8a0c61SJames Wright */ 24bb8a0c61SJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N, 25bb8a0c61SJames Wright PetscReal refine_height, PetscReal top_angle) { 26bb8a0c61SJames Wright 27bb8a0c61SJames Wright PetscInt ierr, narr, ncoords; 28bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 29bb8a0c61SJames Wright PetscScalar *arr_coords; 30bb8a0c61SJames Wright Vec vec_coords; 31bb8a0c61SJames Wright PetscFunctionBeginUser; 32bb8a0c61SJames Wright 33bb8a0c61SJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 34bb8a0c61SJames Wright 35bb8a0c61SJames Wright // Get domain boundary information 36bb8a0c61SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 37*493642f1SJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 38bb8a0c61SJames Wright 39bb8a0c61SJames Wright // Get coords array from DM 40bb8a0c61SJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 41bb8a0c61SJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 42bb8a0c61SJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 43bb8a0c61SJames Wright 44bb8a0c61SJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 45bb8a0c61SJames Wright ncoords = narr/dim; 46bb8a0c61SJames Wright 47bb8a0c61SJames Wright // Get mesh information 48bb8a0c61SJames Wright PetscInt nmax = 3, faces[3]; 49bb8a0c61SJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 50bb8a0c61SJames Wright NULL); CHKERRQ(ierr); 51bb8a0c61SJames Wright 52bb8a0c61SJames Wright // Calculate the first element height 53bb8a0c61SJames Wright PetscReal dybox = domain_size[1]/faces[1]; 54bb8a0c61SJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 55bb8a0c61SJames Wright 56bb8a0c61SJames Wright // Calculate log of sizing outside BL 57bb8a0c61SJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 58bb8a0c61SJames Wright 59*493642f1SJames Wright for(PetscInt i=0; i<ncoords; i++) { 60bb8a0c61SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 61bb8a0c61SJames Wright if(y_box_index <= N) { 62bb8a0c61SJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 63bb8a0c61SJames Wright dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1); 64bb8a0c61SJames Wright } else { 65bb8a0c61SJames Wright PetscInt j = y_box_index - N; 66bb8a0c61SJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 67bb8a0c61SJames Wright exp(log(refine_height) + logdy*j); 68bb8a0c61SJames Wright } 69bb8a0c61SJames Wright } 70bb8a0c61SJames Wright 71bb8a0c61SJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 72bb8a0c61SJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 73bb8a0c61SJames Wright 74bb8a0c61SJames Wright PetscFunctionReturn(0); 75bb8a0c61SJames Wright } 76bb8a0c61SJames Wright 77b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 78bb8a0c61SJames Wright 79bb8a0c61SJames Wright PetscInt ierr; 80bb8a0c61SJames Wright User user = *(User *)ctx; 81bb8a0c61SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 82*493642f1SJames Wright PetscBool use_stg = PETSC_FALSE; 8315a3537eSJed Brown BlasiusContext blasius_ctx; 8415a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 8515a3537eSJed Brown CeedQFunctionContext blasius_context; 8615a3537eSJed Brown 87bb8a0c61SJames Wright PetscFunctionBeginUser; 88b7f03f12SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 8915a3537eSJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 90bb8a0c61SJames Wright 91bb8a0c61SJames Wright // ------------------------------------------------------ 92bb8a0c61SJames Wright // SET UP Blasius 93bb8a0c61SJames Wright // ------------------------------------------------------ 9415a3537eSJed Brown CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 959785fe93SJed Brown problem->ics.qfunction = ICsBlasius; 969785fe93SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 979785fe93SJed Brown problem->apply_outflow.qfunction = Blasius_Outflow; 989785fe93SJed Brown problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc; 99*493642f1SJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 100*493642f1SJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 101bb8a0c61SJames Wright 102bb8a0c61SJames Wright // CeedScalar mu = .04; // Pa s, dynamic viscosity 103bb8a0c61SJames Wright CeedScalar Uinf = 40; // m/s 104bb8a0c61SJames Wright CeedScalar delta0 = 4.2e-4; // m 105bb8a0c61SJames Wright PetscReal refine_height = 5.9e-4; // m 106bb8a0c61SJames Wright PetscReal growth = 1.08; // [-] 107bb8a0c61SJames Wright PetscInt Ndelta = 45; // [-] 108bb8a0c61SJames Wright PetscReal top_angle = 5; // degrees 109bb8a0c61SJames Wright CeedScalar theta0 = 288.; // K 110bb8a0c61SJames Wright CeedScalar P0 = 1.01e5; // Pa 1112acc7cbcSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 112bb8a0c61SJames Wright 113bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 1142acc7cbcSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 1152acc7cbcSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 116bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 117bb8a0c61SJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 118bb8a0c61SJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 119bb8a0c61SJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 120bb8a0c61SJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 121bb8a0c61SJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 122bb8a0c61SJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 123bb8a0c61SJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 124bb8a0c61SJames Wright ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge", 125bb8a0c61SJames Wright NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr); 126bb8a0c61SJames Wright ierr = PetscOptionsScalar("-refine_height", 127bb8a0c61SJames Wright "Height of boundary layer mesh refinement", 128bb8a0c61SJames Wright NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr); 129bb8a0c61SJames Wright ierr = PetscOptionsScalar("-growth", 130bb8a0c61SJames Wright "Geometric growth rate of boundary layer mesh", 131bb8a0c61SJames Wright NULL, growth, &growth, NULL); CHKERRQ(ierr); 132bb8a0c61SJames Wright ierr = PetscOptionsScalar("-top_angle", 133bb8a0c61SJames Wright "Geometric top_angle rate of boundary layer mesh", 134bb8a0c61SJames Wright NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr); 135*493642f1SJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 136*493642f1SJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 137bb8a0c61SJames Wright PetscOptionsEnd(); 138bb8a0c61SJames Wright 139bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 140bb8a0c61SJames Wright PetscScalar second = user->units->second; 141bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 142bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 143bb8a0c61SJames Wright 144bb8a0c61SJames Wright theta0 *= Kelvin; 145bb8a0c61SJames Wright P0 *= Pascal; 146bb8a0c61SJames Wright Uinf *= meter / second; 147bb8a0c61SJames Wright delta0 *= meter; 148bb8a0c61SJames Wright 149bb8a0c61SJames Wright ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle); 150bb8a0c61SJames Wright CHKERRQ(ierr); 151bb8a0c61SJames Wright 15215a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 15315a3537eSJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 15415a3537eSJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 155bb8a0c61SJames Wright 156*493642f1SJames Wright blasius_ctx->weakT = weakT; 15715a3537eSJed Brown blasius_ctx->Uinf = Uinf; 15815a3537eSJed Brown blasius_ctx->delta0 = delta0; 15915a3537eSJed Brown blasius_ctx->theta0 = theta0; 16015a3537eSJed Brown blasius_ctx->P0 = P0; 16115a3537eSJed Brown blasius_ctx->implicit = user->phys->implicit; 16215a3537eSJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 163*493642f1SJames Wright 16415a3537eSJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 16515a3537eSJed Brown &newtonian_ig_ctx); 166bb8a0c61SJames Wright 16715a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 16815a3537eSJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 169bb8a0c61SJames Wright CEED_USE_POINTER, 17015a3537eSJed Brown sizeof(*blasius_ctx), blasius_ctx); 17115a3537eSJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 17215a3537eSJed Brown FreeContextPetsc); 173bb8a0c61SJames Wright 17415a3537eSJed Brown problem->ics.qfunction_context = blasius_context; 17515a3537eSJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 17615a3537eSJed Brown &problem->apply_inflow.qfunction_context); 17715a3537eSJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 17815a3537eSJed Brown &problem->apply_outflow.qfunction_context); 179*493642f1SJames Wright if (use_stg) { 180*493642f1SJames Wright ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0); CHKERRQ(ierr); 181*493642f1SJames Wright } 182bb8a0c61SJames Wright PetscFunctionReturn(0); 183bb8a0c61SJames Wright } 184