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" 13ba6664aeSJames Wright #include "stg_shur14.h" 1488626eedSJames Wright 15d2714514SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, 16d2714514SJames Wright const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, 17d2714514SJames Wright PetscInt *nynodes) { 18d2714514SJames Wright PetscErrorCode ierr; 19d2714514SJames Wright PetscInt ndims, dims[2]; 20d2714514SJames Wright FILE *fp; 21d2714514SJames Wright const PetscInt char_array_len = 512; 22d2714514SJames Wright char line[char_array_len]; 23d2714514SJames Wright char **array; 24d2714514SJames Wright PetscReal *node_locs; 25d2714514SJames Wright PetscFunctionBeginUser; 26d2714514SJames Wright 27d2714514SJames Wright ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); 28d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 29d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 30d2714514SJames Wright 31d2714514SJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 32d2714514SJames Wright if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 33d2714514SJames Wright *nynodes = dims[0]; 34d2714514SJames Wright ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr); 35d2714514SJames Wright 36d2714514SJames Wright for (PetscInt i=0; i<dims[0]; i++) { 37d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 38d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 39d2714514SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 40d2714514SJames Wright "Line %d of %s does not contain enough columns (%d instead of %d)", i, 41d2714514SJames Wright path, ndims, dims[1]); 42d2714514SJames Wright 43d2714514SJames Wright node_locs[i] = (PetscReal) atof(array[0]); 44d2714514SJames Wright } 45d2714514SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 46d2714514SJames Wright *pynodes = node_locs; 47d2714514SJames Wright PetscFunctionReturn(0); 48d2714514SJames Wright } 49d2714514SJames Wright 5088626eedSJames Wright /* \brief Modify the domain and mesh for blasius 5188626eedSJames Wright * 52ba6664aeSJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 53ba6664aeSJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 54ba6664aeSJames Wright * linearly in logspace to the top surface. 5588626eedSJames Wright * 5688626eedSJames Wright * The top surface is also angled downwards, so that it may be used as an 57ba6664aeSJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 58d2714514SJames Wright * 59d2714514SJames Wright * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 60*798d42b8SJames Wright * locations. If it is NULL, then the modified coordinate values will be set in 61*798d42b8SJames Wright * the array, along with `num_node_locs`. 6288626eedSJames Wright */ 63d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 64d2714514SJames Wright PetscReal growth, PetscInt N, 65d2714514SJames Wright PetscReal refine_height, PetscReal top_angle, 66*798d42b8SJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 6788626eedSJames Wright PetscInt ierr, narr, ncoords; 6888626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 6988626eedSJames Wright PetscScalar *arr_coords; 7088626eedSJames Wright Vec vec_coords; 7188626eedSJames Wright PetscFunctionBeginUser; 7288626eedSJames Wright 7388626eedSJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 7488626eedSJames Wright 7588626eedSJames Wright // Get domain boundary information 7688626eedSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 77ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 7888626eedSJames Wright 7988626eedSJames Wright // Get coords array from DM 8088626eedSJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 8188626eedSJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 8288626eedSJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 8388626eedSJames Wright 8488626eedSJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 8588626eedSJames Wright ncoords = narr/dim; 8688626eedSJames Wright 8788626eedSJames Wright // Get mesh information 8888626eedSJames Wright PetscInt nmax = 3, faces[3]; 8988626eedSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 9088626eedSJames Wright NULL); CHKERRQ(ierr); 91d2714514SJames Wright // Get element size of the box mesh, for indexing each node 92d2714514SJames Wright const PetscReal dybox = domain_size[1]/faces[1]; 9388626eedSJames Wright 94*798d42b8SJames Wright if (!*node_locs) { 9588626eedSJames Wright // Calculate the first element height 9688626eedSJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 9788626eedSJames Wright 9888626eedSJames Wright // Calculate log of sizing outside BL 9988626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 10088626eedSJames Wright 101*798d42b8SJames Wright *num_node_locs = faces[1] + 1; 102*798d42b8SJames Wright PetscReal *temp_node_locs; 103*798d42b8SJames Wright ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr); 104*798d42b8SJames Wright 105ba6664aeSJames Wright for (PetscInt i=0; i<ncoords; i++) { 10688626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 10788626eedSJames Wright if (y_box_index <= N) { 108855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 109d2714514SJames Wright * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 11088626eedSJames Wright } else { 11188626eedSJames Wright PetscInt j = y_box_index - N; 112855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 113d2714514SJames Wright * exp(log(refine_height) + logdy*j); 114d2714514SJames Wright } 115*798d42b8SJames Wright if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) 116*798d42b8SJames Wright temp_node_locs[y_box_index] = coords[i][1]; 117d2714514SJames Wright } 118*798d42b8SJames Wright 119*798d42b8SJames Wright *node_locs = temp_node_locs; 120d2714514SJames Wright } else { 121d2714514SJames Wright // Error checking 122*798d42b8SJames Wright if (*num_node_locs < faces[1] +1) 123d2714514SJames Wright SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 124d2714514SJames Wright "There are %d + 1 nodes, but only %d locations given", 125*798d42b8SJames Wright faces[1]+1, *num_node_locs); 126*798d42b8SJames Wright if (*num_node_locs > faces[1] +1) { 127d2714514SJames Wright ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 128d2714514SJames Wright "than the mesh has nodes (%d). This maybe unintended.", 129*798d42b8SJames Wright *num_node_locs, faces[1]+1); CHKERRQ(ierr); 130d2714514SJames Wright } 131*798d42b8SJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 132d2714514SJames Wright 133d2714514SJames Wright for (PetscInt i=0; i<ncoords; i++) { 134d2714514SJames Wright // Determine which y-node we're at 135d2714514SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 136855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y) 137*798d42b8SJames Wright * (*node_locs)[y_box_index]; 13888626eedSJames Wright } 13988626eedSJames Wright } 14088626eedSJames Wright 14188626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 14288626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 14388626eedSJames Wright 14488626eedSJames Wright PetscFunctionReturn(0); 14588626eedSJames Wright } 14688626eedSJames Wright 147a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 14888626eedSJames Wright 14988626eedSJames Wright PetscInt ierr; 15088626eedSJames Wright User user = *(User *)ctx; 15188626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 152ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 153841e4c73SJed Brown BlasiusContext blasius_ctx; 154841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 155841e4c73SJed Brown CeedQFunctionContext blasius_context; 156841e4c73SJed Brown 15788626eedSJames Wright PetscFunctionBeginUser; 158a0add3c9SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 159841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 16088626eedSJames Wright 16188626eedSJames Wright // ------------------------------------------------------ 16288626eedSJames Wright // SET UP Blasius 16388626eedSJames Wright // ------------------------------------------------------ 16491e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 16591e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 16688626eedSJames Wright 16788626eedSJames Wright CeedScalar Uinf = 40; // m/s 16888626eedSJames Wright CeedScalar delta0 = 4.2e-4; // m 16988626eedSJames Wright CeedScalar theta0 = 288.; // K 17088626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 171871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 1727b8d3891SJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 1737b8d3891SJames Wright PetscReal mesh_growth = 1.08; // [-] 1747b8d3891SJames Wright PetscInt mesh_Ndelta = 45; // [-] 1757b8d3891SJames Wright PetscReal mesh_top_angle = 5; // degrees 176d2714514SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 17788626eedSJames Wright 17888626eedSJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 179871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 180871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 18188626eedSJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 18288626eedSJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 18388626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 18488626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 18588626eedSJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 18688626eedSJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 18788626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 18888626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 1897b8d3891SJames Wright ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 1907b8d3891SJames Wright "Velocity at boundary layer edge", 1917b8d3891SJames Wright NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 1927b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_refine_height", 19388626eedSJames Wright "Height of boundary layer mesh refinement", 1947b8d3891SJames Wright NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 1957b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_growth", 19688626eedSJames Wright "Geometric growth rate of boundary layer mesh", 1977b8d3891SJames Wright NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 1987b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_top_angle", 19988626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 2007b8d3891SJames Wright NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 201d2714514SJames Wright ierr = PetscOptionsString("-platemesh_y_node_locs_path", 202d2714514SJames Wright "Path to file with y node locations. " 203d2714514SJames Wright "If empty, will use the algorithmic mesh warping.", NULL, 204d2714514SJames Wright mesh_ynodes_path, mesh_ynodes_path, 205d2714514SJames Wright sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 206ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 207ba6664aeSJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 20888626eedSJames Wright PetscOptionsEnd(); 20988626eedSJames Wright 21088626eedSJames Wright PetscScalar meter = user->units->meter; 21188626eedSJames Wright PetscScalar second = user->units->second; 21288626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 21388626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 21488626eedSJames Wright 21588626eedSJames Wright theta0 *= Kelvin; 21688626eedSJames Wright P0 *= Pascal; 21788626eedSJames Wright Uinf *= meter / second; 21888626eedSJames Wright delta0 *= meter; 21988626eedSJames Wright 220d2714514SJames Wright PetscReal *mesh_ynodes = NULL; 221d2714514SJames Wright PetscInt mesh_nynodes = 0; 222d2714514SJames Wright if (strcmp(mesh_ynodes_path, "")) { 223d2714514SJames Wright ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 22488626eedSJames Wright CHKERRQ(ierr); 225d2714514SJames Wright } 226d2714514SJames Wright ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 227*798d42b8SJames Wright mesh_refine_height, mesh_top_angle, &mesh_ynodes, 228*798d42b8SJames Wright &mesh_nynodes); CHKERRQ(ierr); 22988626eedSJames Wright 230841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 231841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 232841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 23388626eedSJames Wright 234ba6664aeSJames Wright blasius_ctx->weakT = weakT; 235841e4c73SJed Brown blasius_ctx->Uinf = Uinf; 236841e4c73SJed Brown blasius_ctx->delta0 = delta0; 237841e4c73SJed Brown blasius_ctx->theta0 = theta0; 238841e4c73SJed Brown blasius_ctx->P0 = P0; 23930e9fa81SJames Wright newtonian_ig_ctx->P0 = P0; 240841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 241841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 242ba6664aeSJames Wright 243f1122ed0SJames Wright { 244f1122ed0SJames Wright PetscReal domain_min[3]; 245f1122ed0SJames Wright ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); 246f1122ed0SJames Wright blasius_ctx->x_inflow = domain_min[0]; 247f1122ed0SJames Wright } 248f1122ed0SJames Wright 249841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 250841e4c73SJed Brown &newtonian_ig_ctx); 25188626eedSJames Wright 252841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 253841e4c73SJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 25488626eedSJames Wright CEED_USE_POINTER, 255841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 256841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 257841e4c73SJed Brown FreeContextPetsc); 25888626eedSJames Wright 259b77c53c9SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 260841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 261ba6664aeSJames Wright if (use_stg) { 2624e139266SJames Wright ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, 2634e139266SJames Wright mesh_nynodes); CHKERRQ(ierr); 26465dd5cafSJames Wright } else { 265fc02c281SJames Wright CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context); 26665dd5cafSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 26765dd5cafSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 268fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 269fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 27065dd5cafSJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 27165dd5cafSJames Wright &problem->apply_inflow.qfunction_context); 272fc02c281SJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 273fc02c281SJames Wright &problem->apply_inflow_jacobian.qfunction_context); 274ba6664aeSJames Wright } 2754e139266SJames Wright ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 27688626eedSJames Wright PetscFunctionReturn(0); 27788626eedSJames Wright } 278