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" 13*ba6664aeSJames Wright #include "stg_shur14.h" 1488626eedSJames Wright 1588626eedSJames Wright /* \brief Modify the domain and mesh for blasius 1688626eedSJames Wright * 17*ba6664aeSJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 18*ba6664aeSJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 19*ba6664aeSJames Wright * linearly in logspace to the top surface. 2088626eedSJames Wright * 2188626eedSJames Wright * The top surface is also angled downwards, so that it may be used as an 22*ba6664aeSJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 2388626eedSJames Wright */ 2488626eedSJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N, 2588626eedSJames Wright PetscReal refine_height, PetscReal top_angle) { 2688626eedSJames Wright 2788626eedSJames Wright PetscInt ierr, narr, ncoords; 2888626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 2988626eedSJames Wright PetscScalar *arr_coords; 3088626eedSJames Wright Vec vec_coords; 3188626eedSJames Wright PetscFunctionBeginUser; 3288626eedSJames Wright 3388626eedSJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 3488626eedSJames Wright 3588626eedSJames Wright // Get domain boundary information 3688626eedSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 37*ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 3888626eedSJames Wright 3988626eedSJames Wright // Get coords array from DM 4088626eedSJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 4188626eedSJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 4288626eedSJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 4388626eedSJames Wright 4488626eedSJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 4588626eedSJames Wright ncoords = narr/dim; 4688626eedSJames Wright 4788626eedSJames Wright // Get mesh information 4888626eedSJames Wright PetscInt nmax = 3, faces[3]; 4988626eedSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 5088626eedSJames Wright NULL); CHKERRQ(ierr); 5188626eedSJames Wright 5288626eedSJames Wright // Calculate the first element height 5388626eedSJames Wright PetscReal dybox = domain_size[1]/faces[1]; 5488626eedSJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 5588626eedSJames Wright 5688626eedSJames Wright // Calculate log of sizing outside BL 5788626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 5888626eedSJames Wright 59*ba6664aeSJames Wright for(PetscInt i=0; i<ncoords; i++) { 6088626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 6188626eedSJames Wright if(y_box_index <= N) { 6288626eedSJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 6388626eedSJames Wright dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1); 6488626eedSJames Wright } else { 6588626eedSJames Wright PetscInt j = y_box_index - N; 6688626eedSJames Wright coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) * 6788626eedSJames Wright exp(log(refine_height) + logdy*j); 6888626eedSJames Wright } 6988626eedSJames Wright } 7088626eedSJames Wright 7188626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 7288626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 7388626eedSJames Wright 7488626eedSJames Wright PetscFunctionReturn(0); 7588626eedSJames Wright } 7688626eedSJames Wright 77a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 7888626eedSJames Wright 7988626eedSJames Wright PetscInt ierr; 8088626eedSJames Wright User user = *(User *)ctx; 8188626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 82*ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 83841e4c73SJed Brown BlasiusContext blasius_ctx; 84841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 85841e4c73SJed Brown CeedQFunctionContext blasius_context; 86841e4c73SJed Brown 8788626eedSJames Wright PetscFunctionBeginUser; 88a0add3c9SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 89841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 9088626eedSJames Wright 9188626eedSJames Wright // ------------------------------------------------------ 9288626eedSJames Wright // SET UP Blasius 9388626eedSJames Wright // ------------------------------------------------------ 94841e4c73SJed Brown CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 9591e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 9691e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 9791e5af17SJed Brown problem->apply_outflow.qfunction = Blasius_Outflow; 9891e5af17SJed Brown problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc; 99*ba6664aeSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 100*ba6664aeSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 10188626eedSJames Wright 10288626eedSJames Wright // CeedScalar mu = .04; // Pa s, dynamic viscosity 10388626eedSJames Wright CeedScalar Uinf = 40; // m/s 10488626eedSJames Wright CeedScalar delta0 = 4.2e-4; // m 10588626eedSJames Wright PetscReal refine_height = 5.9e-4; // m 10688626eedSJames Wright PetscReal growth = 1.08; // [-] 10788626eedSJames Wright PetscInt Ndelta = 45; // [-] 10888626eedSJames Wright PetscReal top_angle = 5; // degrees 10988626eedSJames Wright CeedScalar theta0 = 288.; // K 11088626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 111871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 11288626eedSJames Wright 11388626eedSJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 114871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 115871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 11688626eedSJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 11788626eedSJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 11888626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 11988626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 12088626eedSJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 12188626eedSJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 12288626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 12388626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 12488626eedSJames Wright ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge", 12588626eedSJames Wright NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr); 12688626eedSJames Wright ierr = PetscOptionsScalar("-refine_height", 12788626eedSJames Wright "Height of boundary layer mesh refinement", 12888626eedSJames Wright NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr); 12988626eedSJames Wright ierr = PetscOptionsScalar("-growth", 13088626eedSJames Wright "Geometric growth rate of boundary layer mesh", 13188626eedSJames Wright NULL, growth, &growth, NULL); CHKERRQ(ierr); 13288626eedSJames Wright ierr = PetscOptionsScalar("-top_angle", 13388626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 13488626eedSJames Wright NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr); 135*ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 136*ba6664aeSJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 13788626eedSJames Wright PetscOptionsEnd(); 13888626eedSJames Wright 13988626eedSJames Wright PetscScalar meter = user->units->meter; 14088626eedSJames Wright PetscScalar second = user->units->second; 14188626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 14288626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 14388626eedSJames Wright 14488626eedSJames Wright theta0 *= Kelvin; 14588626eedSJames Wright P0 *= Pascal; 14688626eedSJames Wright Uinf *= meter / second; 14788626eedSJames Wright delta0 *= meter; 14888626eedSJames Wright 14988626eedSJames Wright ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle); 15088626eedSJames Wright CHKERRQ(ierr); 15188626eedSJames Wright 152841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 153841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 154841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 15588626eedSJames Wright 156*ba6664aeSJames Wright blasius_ctx->weakT = weakT; 157841e4c73SJed Brown blasius_ctx->Uinf = Uinf; 158841e4c73SJed Brown blasius_ctx->delta0 = delta0; 159841e4c73SJed Brown blasius_ctx->theta0 = theta0; 160841e4c73SJed Brown blasius_ctx->P0 = P0; 161841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 162841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 163*ba6664aeSJames Wright 164841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 165841e4c73SJed Brown &newtonian_ig_ctx); 16688626eedSJames Wright 167841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 168841e4c73SJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 16988626eedSJames Wright CEED_USE_POINTER, 170841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 171841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 172841e4c73SJed Brown FreeContextPetsc); 17388626eedSJames Wright 174841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 175841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 176841e4c73SJed Brown &problem->apply_inflow.qfunction_context); 177841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 178841e4c73SJed Brown &problem->apply_outflow.qfunction_context); 179*ba6664aeSJames Wright if (use_stg) { 180*ba6664aeSJames Wright ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0); CHKERRQ(ierr); 181*ba6664aeSJames Wright } 18288626eedSJames Wright PetscFunctionReturn(0); 18388626eedSJames Wright } 184