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" 13493642f1SJames Wright #include "stg_shur14.h" 14bb8a0c61SJames Wright 15f362fe62SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, 16f362fe62SJames Wright const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, 17f362fe62SJames Wright PetscInt *nynodes) { 18f362fe62SJames Wright PetscErrorCode ierr; 19f362fe62SJames Wright PetscInt ndims, dims[2]; 20f362fe62SJames Wright FILE *fp; 21f362fe62SJames Wright const PetscInt char_array_len = 512; 22f362fe62SJames Wright char line[char_array_len]; 23f362fe62SJames Wright char **array; 24f362fe62SJames Wright PetscReal *node_locs; 25f362fe62SJames Wright PetscFunctionBeginUser; 26f362fe62SJames Wright 27f362fe62SJames Wright ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); 28f362fe62SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 29f362fe62SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 30f362fe62SJames Wright 31f362fe62SJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 32f362fe62SJames Wright if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 33f362fe62SJames Wright *nynodes = dims[0]; 34f362fe62SJames Wright ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr); 35f362fe62SJames Wright 36f362fe62SJames Wright for (PetscInt i=0; i<dims[0]; i++) { 37f362fe62SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 38f362fe62SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 39f362fe62SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 40f362fe62SJames Wright "Line %d of %s does not contain enough columns (%d instead of %d)", i, 41f362fe62SJames Wright path, ndims, dims[1]); 42f362fe62SJames Wright 43f362fe62SJames Wright node_locs[i] = (PetscReal) atof(array[0]); 44f362fe62SJames Wright } 45f362fe62SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 46f362fe62SJames Wright *pynodes = node_locs; 47f362fe62SJames Wright PetscFunctionReturn(0); 48f362fe62SJames Wright } 49f362fe62SJames Wright 50bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius 51bb8a0c61SJames Wright * 52493642f1SJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 53493642f1SJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 54493642f1SJames Wright * linearly in logspace to the top surface. 55bb8a0c61SJames Wright * 56bb8a0c61SJames Wright * The top surface is also angled downwards, so that it may be used as an 57493642f1SJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 58f362fe62SJames Wright * 59f362fe62SJames Wright * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 60f362fe62SJames Wright * locations. 61bb8a0c61SJames Wright */ 62f362fe62SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 63f362fe62SJames Wright PetscReal growth, PetscInt N, 64f362fe62SJames Wright PetscReal refine_height, PetscReal top_angle, 65f362fe62SJames Wright PetscReal node_locs[], PetscInt num_node_locs) { 66bb8a0c61SJames Wright PetscInt ierr, narr, ncoords; 67bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 68bb8a0c61SJames Wright PetscScalar *arr_coords; 69bb8a0c61SJames Wright Vec vec_coords; 70bb8a0c61SJames Wright PetscFunctionBeginUser; 71bb8a0c61SJames Wright 72bb8a0c61SJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 73bb8a0c61SJames Wright 74bb8a0c61SJames Wright // Get domain boundary information 75bb8a0c61SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 76493642f1SJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 77bb8a0c61SJames Wright 78bb8a0c61SJames Wright // Get coords array from DM 79bb8a0c61SJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 80bb8a0c61SJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 81bb8a0c61SJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 82bb8a0c61SJames Wright 83bb8a0c61SJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 84bb8a0c61SJames Wright ncoords = narr/dim; 85bb8a0c61SJames Wright 86bb8a0c61SJames Wright // Get mesh information 87bb8a0c61SJames Wright PetscInt nmax = 3, faces[3]; 88bb8a0c61SJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 89bb8a0c61SJames Wright NULL); CHKERRQ(ierr); 90f362fe62SJames Wright // Get element size of the box mesh, for indexing each node 91f362fe62SJames Wright const PetscReal dybox = domain_size[1]/faces[1]; 92bb8a0c61SJames Wright 93f362fe62SJames Wright if (!node_locs) { 94bb8a0c61SJames Wright // Calculate the first element height 95bb8a0c61SJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 96bb8a0c61SJames Wright 97bb8a0c61SJames Wright // Calculate log of sizing outside BL 98bb8a0c61SJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 99bb8a0c61SJames Wright 100493642f1SJames Wright for (PetscInt i=0; i<ncoords; i++) { 101bb8a0c61SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 102bb8a0c61SJames Wright if (y_box_index <= N) { 103029c6dfaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 104f362fe62SJames Wright * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 105bb8a0c61SJames Wright } else { 106bb8a0c61SJames Wright PetscInt j = y_box_index - N; 107029c6dfaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 108f362fe62SJames Wright * exp(log(refine_height) + logdy*j); 109f362fe62SJames Wright } 110f362fe62SJames Wright } 111f362fe62SJames Wright } else { 112f362fe62SJames Wright // Error checking 113f362fe62SJames Wright if (num_node_locs < faces[1] +1) 114f362fe62SJames Wright SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 115f362fe62SJames Wright "There are %d + 1 nodes, but only %d locations given", 116f362fe62SJames Wright faces[1]+1, num_node_locs); 117f362fe62SJames Wright if (num_node_locs > faces[1] +1) { 118f362fe62SJames Wright ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 119f362fe62SJames Wright "than the mesh has nodes (%d). This maybe unintended.", 120f362fe62SJames Wright num_node_locs, faces[1]+1); CHKERRQ(ierr); 121f362fe62SJames Wright } 122029c6dfaSJames Wright PetscScalar max_y = node_locs[faces[1]]; 123f362fe62SJames Wright 124f362fe62SJames Wright for (PetscInt i=0; i<ncoords; i++) { 125f362fe62SJames Wright // Determine which y-node we're at 126f362fe62SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 127029c6dfaSJames Wright // coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 128029c6dfaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y) 129f362fe62SJames Wright * node_locs[y_box_index]; 130bb8a0c61SJames Wright } 131bb8a0c61SJames Wright } 132bb8a0c61SJames Wright 133bb8a0c61SJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 134bb8a0c61SJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 135bb8a0c61SJames Wright 136bb8a0c61SJames Wright PetscFunctionReturn(0); 137bb8a0c61SJames Wright } 138bb8a0c61SJames Wright 139b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 140bb8a0c61SJames Wright 141bb8a0c61SJames Wright PetscInt ierr; 142bb8a0c61SJames Wright User user = *(User *)ctx; 143bb8a0c61SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 144493642f1SJames Wright PetscBool use_stg = PETSC_FALSE; 14515a3537eSJed Brown BlasiusContext blasius_ctx; 14615a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 14715a3537eSJed Brown CeedQFunctionContext blasius_context; 14815a3537eSJed Brown 149bb8a0c61SJames Wright PetscFunctionBeginUser; 150b7f03f12SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 15115a3537eSJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 152bb8a0c61SJames Wright 153bb8a0c61SJames Wright // ------------------------------------------------------ 154bb8a0c61SJames Wright // SET UP Blasius 155bb8a0c61SJames Wright // ------------------------------------------------------ 1569785fe93SJed Brown problem->ics.qfunction = ICsBlasius; 1579785fe93SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 158f0b65372SJed Brown problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 159f0b65372SJed Brown problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 160bb8a0c61SJames Wright 161bb8a0c61SJames Wright CeedScalar Uinf = 40; // m/s 162bb8a0c61SJames Wright CeedScalar delta0 = 4.2e-4; // m 163bb8a0c61SJames Wright CeedScalar theta0 = 288.; // K 164bb8a0c61SJames Wright CeedScalar P0 = 1.01e5; // Pa 1652acc7cbcSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 1667470235eSJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 1677470235eSJames Wright PetscReal mesh_growth = 1.08; // [-] 1687470235eSJames Wright PetscInt mesh_Ndelta = 45; // [-] 1697470235eSJames Wright PetscReal mesh_top_angle = 5; // degrees 170f362fe62SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 171bb8a0c61SJames Wright 172bb8a0c61SJames Wright PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 1732acc7cbcSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 1742acc7cbcSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 175bb8a0c61SJames Wright ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 176bb8a0c61SJames Wright NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 177bb8a0c61SJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 178bb8a0c61SJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 179bb8a0c61SJames Wright ierr = PetscOptionsScalar("-theta0", "Wall temperature", 180bb8a0c61SJames Wright NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 181bb8a0c61SJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 182bb8a0c61SJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 1837470235eSJames Wright ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 1847470235eSJames Wright "Velocity at boundary layer edge", 1857470235eSJames Wright NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 1867470235eSJames Wright ierr = PetscOptionsScalar("-platemesh_refine_height", 187bb8a0c61SJames Wright "Height of boundary layer mesh refinement", 1887470235eSJames Wright NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 1897470235eSJames Wright ierr = PetscOptionsScalar("-platemesh_growth", 190bb8a0c61SJames Wright "Geometric growth rate of boundary layer mesh", 1917470235eSJames Wright NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 1927470235eSJames Wright ierr = PetscOptionsScalar("-platemesh_top_angle", 193bb8a0c61SJames Wright "Geometric top_angle rate of boundary layer mesh", 1947470235eSJames Wright NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 195f362fe62SJames Wright ierr = PetscOptionsString("-platemesh_y_node_locs_path", 196f362fe62SJames Wright "Path to file with y node locations. " 197f362fe62SJames Wright "If empty, will use the algorithmic mesh warping.", NULL, 198f362fe62SJames Wright mesh_ynodes_path, mesh_ynodes_path, 199f362fe62SJames Wright sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 200493642f1SJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 201493642f1SJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 202bb8a0c61SJames Wright PetscOptionsEnd(); 203bb8a0c61SJames Wright 204bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 205bb8a0c61SJames Wright PetscScalar second = user->units->second; 206bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 207bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 208bb8a0c61SJames Wright 209bb8a0c61SJames Wright theta0 *= Kelvin; 210bb8a0c61SJames Wright P0 *= Pascal; 211bb8a0c61SJames Wright Uinf *= meter / second; 212bb8a0c61SJames Wright delta0 *= meter; 213bb8a0c61SJames Wright 214f362fe62SJames Wright PetscReal *mesh_ynodes = NULL; 215f362fe62SJames Wright PetscInt mesh_nynodes = 0; 216f362fe62SJames Wright if (strcmp(mesh_ynodes_path, "")) { 217f362fe62SJames Wright ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 218bb8a0c61SJames Wright CHKERRQ(ierr); 219f362fe62SJames Wright } 220f362fe62SJames Wright ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 221f362fe62SJames Wright mesh_refine_height, mesh_top_angle, mesh_ynodes, 222f362fe62SJames Wright mesh_nynodes); CHKERRQ(ierr); 223bb8a0c61SJames Wright 22415a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 22515a3537eSJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 22615a3537eSJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 227bb8a0c61SJames Wright 228493642f1SJames Wright blasius_ctx->weakT = weakT; 22915a3537eSJed Brown blasius_ctx->Uinf = Uinf; 23015a3537eSJed Brown blasius_ctx->delta0 = delta0; 23115a3537eSJed Brown blasius_ctx->theta0 = theta0; 23215a3537eSJed Brown blasius_ctx->P0 = P0; 23304b9037bSJames Wright newtonian_ig_ctx->P0 = P0; 23415a3537eSJed Brown blasius_ctx->implicit = user->phys->implicit; 23515a3537eSJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 236493642f1SJames Wright 237ef2c71fdSJames Wright { 238ef2c71fdSJames Wright PetscReal domain_min[3]; 239ef2c71fdSJames Wright ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); 240ef2c71fdSJames Wright blasius_ctx->x_inflow = domain_min[0]; 241ef2c71fdSJames Wright } 242ef2c71fdSJames Wright 24315a3537eSJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 24415a3537eSJed Brown &newtonian_ig_ctx); 245bb8a0c61SJames Wright 24615a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 24715a3537eSJed Brown CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 248bb8a0c61SJames Wright CEED_USE_POINTER, 24915a3537eSJed Brown sizeof(*blasius_ctx), blasius_ctx); 25015a3537eSJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 25115a3537eSJed Brown FreeContextPetsc); 252bb8a0c61SJames Wright 253*43bff553SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 25415a3537eSJed Brown problem->ics.qfunction_context = blasius_context; 25515a3537eSJed Brown CeedQFunctionContextReferenceCopy(blasius_context, 256f0b65372SJed Brown &problem->apply_inflow_jacobian.qfunction_context); 257493642f1SJames Wright if (use_stg) { 258e6098bcdSJames Wright ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, 259e6098bcdSJames Wright mesh_nynodes); CHKERRQ(ierr); 2608085925cSJames Wright } else { 2618085925cSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 2628085925cSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 2638085925cSJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 2648085925cSJames Wright &problem->apply_inflow.qfunction_context); 265493642f1SJames Wright } 266e6098bcdSJames Wright ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 267bb8a0c61SJames Wright PetscFunctionReturn(0); 268bb8a0c61SJames Wright } 269