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` 60d2714514SJames Wright * locations. 6188626eedSJames Wright */ 62d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 63d2714514SJames Wright PetscReal growth, PetscInt N, 64d2714514SJames Wright PetscReal refine_height, PetscReal top_angle, 65d2714514SJames Wright PetscReal node_locs[], PetscInt num_node_locs) { 6688626eedSJames Wright 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 94d2714514SJames 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 101ba6664aeSJames Wright for (PetscInt i=0; i<ncoords; i++) { 10288626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 10388626eedSJames Wright if (y_box_index <= N) { 104d2714514SJames Wright coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 105d2714514SJames Wright * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 10688626eedSJames Wright } else { 10788626eedSJames Wright PetscInt j = y_box_index - N; 108d2714514SJames Wright coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 109d2714514SJames Wright * exp(log(refine_height) + logdy*j); 110d2714514SJames Wright } 111d2714514SJames Wright } 112d2714514SJames Wright } else { 113d2714514SJames Wright // Error checking 114d2714514SJames Wright if (num_node_locs < faces[1] +1) 115d2714514SJames Wright SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 116d2714514SJames Wright "There are %d + 1 nodes, but only %d locations given", 117d2714514SJames Wright faces[1]+1, num_node_locs); 118d2714514SJames Wright if (num_node_locs > faces[1] +1) { 119d2714514SJames Wright ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 120d2714514SJames Wright "than the mesh has nodes (%d). This maybe unintended.", 121d2714514SJames Wright num_node_locs, faces[1]+1); CHKERRQ(ierr); 122d2714514SJames Wright } 123d2714514SJames Wright 124d2714514SJames Wright for (PetscInt i=0; i<ncoords; i++) { 125d2714514SJames Wright // Determine which y-node we're at 126d2714514SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 127d2714514SJames Wright coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 128d2714514SJames Wright * node_locs[y_box_index]; 12988626eedSJames Wright } 13088626eedSJames Wright } 13188626eedSJames Wright 13288626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 13388626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 13488626eedSJames Wright 13588626eedSJames Wright PetscFunctionReturn(0); 13688626eedSJames Wright } 13788626eedSJames Wright 138a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 13988626eedSJames Wright 14088626eedSJames Wright PetscInt ierr; 14188626eedSJames Wright User user = *(User *)ctx; 14288626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 143ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 144841e4c73SJed Brown BlasiusContext blasius_ctx; 145841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 146841e4c73SJed Brown CeedQFunctionContext blasius_context; 147841e4c73SJed Brown 14888626eedSJames Wright PetscFunctionBeginUser; 149a0add3c9SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 150841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 15188626eedSJames Wright 15288626eedSJames Wright // ------------------------------------------------------ 15388626eedSJames Wright // SET UP Blasius 15488626eedSJames Wright // ------------------------------------------------------ 155841e4c73SJed Brown CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 15691e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 15791e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 15891e5af17SJed Brown problem->apply_outflow.qfunction = Blasius_Outflow; 15991e5af17SJed Brown problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc; 160ba6664aeSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 161ba6664aeSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 16288626eedSJames Wright 16388626eedSJames Wright CeedScalar Uinf = 40; // m/s 16488626eedSJames Wright CeedScalar delta0 = 4.2e-4; // m 16588626eedSJames Wright CeedScalar theta0 = 288.; // K 16688626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 167871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 1687b8d3891SJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 1697b8d3891SJames Wright PetscReal mesh_growth = 1.08; // [-] 1707b8d3891SJames Wright PetscInt mesh_Ndelta = 45; // [-] 1717b8d3891SJames Wright PetscReal mesh_top_angle = 5; // degrees 172d2714514SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 17388626eedSJames Wright 17488626eedSJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 175871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 176871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 17788626eedSJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 17888626eedSJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 17988626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 18088626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 18188626eedSJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 18288626eedSJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 18388626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 18488626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 1857b8d3891SJames Wright ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 1867b8d3891SJames Wright "Velocity at boundary layer edge", 1877b8d3891SJames Wright NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 1887b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_refine_height", 18988626eedSJames Wright "Height of boundary layer mesh refinement", 1907b8d3891SJames Wright NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 1917b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_growth", 19288626eedSJames Wright "Geometric growth rate of boundary layer mesh", 1937b8d3891SJames Wright NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 1947b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_top_angle", 19588626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 1967b8d3891SJames Wright NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 197d2714514SJames Wright ierr = PetscOptionsString("-platemesh_y_node_locs_path", 198d2714514SJames Wright "Path to file with y node locations. " 199d2714514SJames Wright "If empty, will use the algorithmic mesh warping.", NULL, 200d2714514SJames Wright mesh_ynodes_path, mesh_ynodes_path, 201d2714514SJames Wright sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 202ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 203ba6664aeSJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 20488626eedSJames Wright PetscOptionsEnd(); 20588626eedSJames Wright 20688626eedSJames Wright PetscScalar meter = user->units->meter; 20788626eedSJames Wright PetscScalar second = user->units->second; 20888626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 20988626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 21088626eedSJames Wright 21188626eedSJames Wright theta0 *= Kelvin; 21288626eedSJames Wright P0 *= Pascal; 21388626eedSJames Wright Uinf *= meter / second; 21488626eedSJames Wright delta0 *= meter; 21588626eedSJames Wright 216d2714514SJames Wright PetscReal *mesh_ynodes = NULL; 217d2714514SJames Wright PetscInt mesh_nynodes = 0; 218d2714514SJames Wright if (strcmp(mesh_ynodes_path, "")) { 219d2714514SJames Wright ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 22088626eedSJames Wright CHKERRQ(ierr); 221d2714514SJames Wright } 222d2714514SJames Wright ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 223d2714514SJames Wright mesh_refine_height, mesh_top_angle, mesh_ynodes, 224d2714514SJames Wright mesh_nynodes); CHKERRQ(ierr); 22588626eedSJames Wright 226841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 227841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 228841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 22988626eedSJames Wright 230ba6664aeSJames Wright blasius_ctx->weakT = weakT; 231841e4c73SJed Brown blasius_ctx->Uinf = Uinf; 232841e4c73SJed Brown blasius_ctx->delta0 = delta0; 233841e4c73SJed Brown blasius_ctx->theta0 = theta0; 234841e4c73SJed Brown blasius_ctx->P0 = P0; 235841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 236841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 237ba6664aeSJames Wright 238*f1122ed0SJames Wright { 239*f1122ed0SJames Wright PetscReal domain_min[3]; 240*f1122ed0SJames Wright ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); 241*f1122ed0SJames Wright blasius_ctx->x_inflow = domain_min[0]; 242*f1122ed0SJames Wright } 243*f1122ed0SJames Wright 244841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 245841e4c73SJed Brown &newtonian_ig_ctx); 24688626eedSJames Wright 247841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 248841e4c73SJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 24988626eedSJames Wright CEED_USE_POINTER, 250841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 251841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 252841e4c73SJed Brown FreeContextPetsc); 25388626eedSJames Wright 254841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 255841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 256841e4c73SJed Brown &problem->apply_inflow.qfunction_context); 257841e4c73SJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 258841e4c73SJed Brown &problem->apply_outflow.qfunction_context); 259ba6664aeSJames Wright if (use_stg) { 2604e139266SJames Wright ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, 2614e139266SJames Wright mesh_nynodes); CHKERRQ(ierr); 262ba6664aeSJames Wright } 2634e139266SJames Wright ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 26488626eedSJames Wright PetscFunctionReturn(0); 26588626eedSJames Wright } 266