188626eedSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 288626eedSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388626eedSJames Wright // 488626eedSJames Wright // SPDX-License-Identifier: BSD-2-Clause 588626eedSJames Wright // 688626eedSJames Wright // This file is part of CEED: http://github.com/ceed 788626eedSJames Wright 888626eedSJames Wright /// @file 988626eedSJames Wright /// Utility functions for setting up Blasius Boundary Layer 1088626eedSJames Wright 1188626eedSJames Wright #include "../navierstokes.h" 1288626eedSJames Wright #include "../qfunctions/blasius.h" 1388626eedSJames Wright 1488626eedSJames Wright /* \brief Modify the domain and mesh for blasius 1588626eedSJames Wright * 1688626eedSJames Wright * Modifies mesh such that `N` elements are within 1.2*`delta0` with a geometric 1788626eedSJames Wright * growth ratio of `growth`. Excess elements are then geometrically distributed 1888626eedSJames Wright * to the top surface. 1988626eedSJames Wright * 2088626eedSJames Wright * The top surface is also angled downwards, so that it may be used as an 2188626eedSJames Wright * outflow. It's angle is controlled by top_angle (in units of degrees). 2288626eedSJames Wright */ 2388626eedSJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N, 2488626eedSJames Wright PetscReal refine_height, PetscReal top_angle) { 2588626eedSJames Wright 2688626eedSJames Wright PetscInt ierr, narr, ncoords; 2788626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 2888626eedSJames Wright PetscScalar *arr_coords; 2988626eedSJames Wright Vec vec_coords; 3088626eedSJames Wright PetscFunctionBeginUser; 3188626eedSJames Wright 3288626eedSJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 3388626eedSJames Wright 3488626eedSJames Wright // Get domain boundary information 3588626eedSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 3688626eedSJames Wright for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 3788626eedSJames Wright 3888626eedSJames Wright // Get coords array from DM 3988626eedSJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 4088626eedSJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 4188626eedSJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 4288626eedSJames Wright 4388626eedSJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 4488626eedSJames Wright ncoords = narr/dim; 4588626eedSJames Wright 4688626eedSJames Wright // Get mesh information 4788626eedSJames Wright PetscInt nmax = 3, faces[3]; 4888626eedSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 4988626eedSJames Wright NULL); CHKERRQ(ierr); 5088626eedSJames Wright 5188626eedSJames Wright // Calculate the first element height 5288626eedSJames Wright PetscReal dybox = domain_size[1]/faces[1]; 5388626eedSJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 5488626eedSJames Wright 5588626eedSJames Wright // Calculate log of sizing outside BL 5688626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 5788626eedSJames Wright 5888626eedSJames Wright for(int i=0; i<ncoords; i++) { 5988626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 6088626eedSJames Wright if(y_box_index <= N) { 6188626eedSJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 6288626eedSJames Wright dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1); 6388626eedSJames Wright } else { 6488626eedSJames Wright PetscInt j = y_box_index - N; 6588626eedSJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 6688626eedSJames Wright exp(log(refine_height) + logdy*j); 6788626eedSJames Wright } 6888626eedSJames Wright } 6988626eedSJames Wright 7088626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 7188626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 7288626eedSJames Wright 7388626eedSJames Wright PetscFunctionReturn(0); 7488626eedSJames Wright } 7588626eedSJames Wright 7688626eedSJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *setup_ctx, 7788626eedSJames Wright void *ctx) { 7888626eedSJames Wright 7988626eedSJames Wright PetscInt ierr; 8088626eedSJames Wright User user = *(User *)ctx; 8188626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 82*841e4c73SJed Brown BlasiusContext blasius_ctx; 83*841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 84*841e4c73SJed Brown CeedQFunctionContext blasius_context; 85*841e4c73SJed Brown 8688626eedSJames Wright PetscFunctionBeginUser; 87*841e4c73SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr); 88*841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 8988626eedSJames Wright 9088626eedSJames Wright // ------------------------------------------------------ 9188626eedSJames Wright // SET UP Blasius 9288626eedSJames Wright // ------------------------------------------------------ 93*841e4c73SJed Brown CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 9491e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 9591e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 9691e5af17SJed Brown problem->apply_inflow.qfunction = Blasius_Inflow; 9791e5af17SJed Brown problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 9891e5af17SJed Brown problem->apply_outflow.qfunction = Blasius_Outflow; 9991e5af17SJed Brown problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc; 10088626eedSJames Wright 10188626eedSJames Wright // CeedScalar mu = .04; // Pa s, dynamic viscosity 10288626eedSJames Wright CeedScalar Uinf = 40; // m/s 10388626eedSJames Wright CeedScalar delta0 = 4.2e-4; // m 10488626eedSJames Wright PetscReal refine_height = 5.9e-4; // m 10588626eedSJames Wright PetscReal growth = 1.08; // [-] 10688626eedSJames Wright PetscInt Ndelta = 45; // [-] 10788626eedSJames Wright PetscReal top_angle = 5; // degrees 10888626eedSJames Wright CeedScalar theta0 = 288.; // K 10988626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 110871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 11188626eedSJames Wright 11288626eedSJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 113871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 114871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 11588626eedSJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 11688626eedSJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 11788626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 11888626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 11988626eedSJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 12088626eedSJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 12188626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 12288626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 12388626eedSJames Wright ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge", 12488626eedSJames Wright NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr); 12588626eedSJames Wright ierr = PetscOptionsScalar("-refine_height", 12688626eedSJames Wright "Height of boundary layer mesh refinement", 12788626eedSJames Wright NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr); 12888626eedSJames Wright ierr = PetscOptionsScalar("-growth", 12988626eedSJames Wright "Geometric growth rate of boundary layer mesh", 13088626eedSJames Wright NULL, growth, &growth, NULL); CHKERRQ(ierr); 13188626eedSJames Wright ierr = PetscOptionsScalar("-top_angle", 13288626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 13388626eedSJames Wright NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr); 13488626eedSJames Wright PetscOptionsEnd(); 13588626eedSJames Wright 13688626eedSJames Wright PetscScalar meter = user->units->meter; 13788626eedSJames Wright PetscScalar second = user->units->second; 13888626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 13988626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 14088626eedSJames Wright 14188626eedSJames Wright theta0 *= Kelvin; 14288626eedSJames Wright P0 *= Pascal; 14388626eedSJames Wright Uinf *= meter / second; 14488626eedSJames Wright delta0 *= meter; 14588626eedSJames Wright 14688626eedSJames Wright ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle); 14788626eedSJames Wright CHKERRQ(ierr); 14888626eedSJames Wright 149*841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 150*841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 151*841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 15288626eedSJames Wright 153*841e4c73SJed Brown blasius_ctx->weakT = !!weakT; 154*841e4c73SJed Brown blasius_ctx->Uinf = Uinf; 155*841e4c73SJed Brown blasius_ctx->delta0 = delta0; 156*841e4c73SJed Brown blasius_ctx->theta0 = theta0; 157*841e4c73SJed Brown blasius_ctx->P0 = P0; 158*841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 159*841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 160*841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 161*841e4c73SJed Brown &newtonian_ig_ctx); 16288626eedSJames Wright 163*841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 164*841e4c73SJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 16588626eedSJames Wright CEED_USE_POINTER, 166*841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 167*841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 168*841e4c73SJed Brown FreeContextPetsc); 16988626eedSJames Wright 170*841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 171*841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 172*841e4c73SJed Brown &problem->apply_inflow.qfunction_context); 173*841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 174*841e4c73SJed Brown &problem->apply_outflow.qfunction_context); 17588626eedSJames Wright PetscFunctionReturn(0); 17688626eedSJames Wright } 177