1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Utility functions for setting up Blasius Boundary Layer 10 11 #include "../navierstokes.h" 12 #include "../qfunctions/blasius.h" 13 #include "stg_shur14.h" 14 15 static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, 16 const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, 17 PetscInt *nynodes) { 18 PetscErrorCode ierr; 19 PetscInt ndims, dims[2]; 20 FILE *fp; 21 const PetscInt char_array_len = 512; 22 char line[char_array_len]; 23 char **array; 24 PetscReal *node_locs; 25 PetscFunctionBeginUser; 26 27 ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); 28 ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 29 ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 30 31 for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 32 if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 33 *nynodes = dims[0]; 34 ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr); 35 36 for (PetscInt i=0; i<dims[0]; i++) { 37 ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 38 ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 39 if (ndims < dims[1]) SETERRQ(comm, -1, 40 "Line %d of %s does not contain enough columns (%d instead of %d)", i, 41 path, ndims, dims[1]); 42 43 node_locs[i] = (PetscReal) atof(array[0]); 44 } 45 ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 46 *pynodes = node_locs; 47 PetscFunctionReturn(0); 48 } 49 50 /* \brief Modify the domain and mesh for blasius 51 * 52 * Modifies mesh such that `N` elements are within `refine_height` with a 53 * geometric growth ratio of `growth`. Excess elements are then distributed 54 * linearly in logspace to the top surface. 55 * 56 * The top surface is also angled downwards, so that it may be used as an 57 * outflow. It's angle is controlled by `top_angle` (in units of degrees). 58 * 59 * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 60 * locations. If it is NULL, then the modified coordinate values will be set in 61 * the array, along with `num_node_locs`. 62 */ 63 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 64 PetscReal growth, PetscInt N, 65 PetscReal refine_height, PetscReal top_angle, 66 PetscReal *node_locs[], PetscInt *num_node_locs) { 67 PetscInt ierr, narr, ncoords; 68 PetscReal domain_min[3], domain_max[3], domain_size[3]; 69 PetscScalar *arr_coords; 70 Vec vec_coords; 71 PetscFunctionBeginUser; 72 73 PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 74 75 // Get domain boundary information 76 ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 77 for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 78 79 // Get coords array from DM 80 ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 81 ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 82 ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 83 84 PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 85 ncoords = narr/dim; 86 87 // Get mesh information 88 PetscInt nmax = 3, faces[3]; 89 ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 90 NULL); CHKERRQ(ierr); 91 // Get element size of the box mesh, for indexing each node 92 const PetscReal dybox = domain_size[1]/faces[1]; 93 94 if (!*node_locs) { 95 // Calculate the first element height 96 PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 97 98 // Calculate log of sizing outside BL 99 PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 100 101 *num_node_locs = faces[1] + 1; 102 PetscReal *temp_node_locs; 103 ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr); 104 105 for (PetscInt i=0; i<ncoords; i++) { 106 PetscInt y_box_index = round(coords[i][1]/dybox); 107 if (y_box_index <= N) { 108 coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 109 * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 110 } else { 111 PetscInt j = y_box_index - N; 112 coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 113 * exp(log(refine_height) + logdy*j); 114 } 115 if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) 116 temp_node_locs[y_box_index] = coords[i][1]; 117 } 118 119 *node_locs = temp_node_locs; 120 } else { 121 // Error checking 122 if (*num_node_locs < faces[1] +1) 123 SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 124 "There are %d + 1 nodes, but only %d locations given", 125 faces[1]+1, *num_node_locs); 126 if (*num_node_locs > faces[1] +1) { 127 ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 128 "than the mesh has nodes (%d). This maybe unintended.", 129 *num_node_locs, faces[1]+1); CHKERRQ(ierr); 130 } 131 PetscScalar max_y = (*node_locs)[faces[1]]; 132 133 for (PetscInt i=0; i<ncoords; i++) { 134 // Determine which y-node we're at 135 PetscInt y_box_index = round(coords[i][1]/dybox); 136 coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y) 137 * (*node_locs)[y_box_index]; 138 } 139 } 140 141 ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 142 ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 143 144 PetscFunctionReturn(0); 145 } 146 147 PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 148 149 PetscInt ierr; 150 User user = *(User *)ctx; 151 MPI_Comm comm = PETSC_COMM_WORLD; 152 PetscBool use_stg = PETSC_FALSE; 153 BlasiusContext blasius_ctx; 154 NewtonianIdealGasContext newtonian_ig_ctx; 155 CeedQFunctionContext blasius_context; 156 157 PetscFunctionBeginUser; 158 ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 159 ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 160 161 // ------------------------------------------------------ 162 // SET UP Blasius 163 // ------------------------------------------------------ 164 problem->ics.qfunction = ICsBlasius; 165 problem->ics.qfunction_loc = ICsBlasius_loc; 166 167 CeedScalar Uinf = 40; // m/s 168 CeedScalar delta0 = 4.2e-4; // m 169 CeedScalar theta0 = 288.; // K 170 CeedScalar P0 = 1.01e5; // Pa 171 PetscBool weakT = PETSC_FALSE; // weak density or temperature 172 PetscReal mesh_refine_height = 5.9e-4; // m 173 PetscReal mesh_growth = 1.08; // [-] 174 PetscInt mesh_Ndelta = 45; // [-] 175 PetscReal mesh_top_angle = 5; // degrees 176 char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 177 178 PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 179 ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 180 NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 181 ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 182 NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 183 ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 184 NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 185 ierr = PetscOptionsScalar("-theta0", "Wall temperature", 186 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 187 ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 188 NULL, P0, &P0, NULL); CHKERRQ(ierr); 189 ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 190 "Velocity at boundary layer edge", 191 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 192 ierr = PetscOptionsScalar("-platemesh_refine_height", 193 "Height of boundary layer mesh refinement", 194 NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 195 ierr = PetscOptionsScalar("-platemesh_growth", 196 "Geometric growth rate of boundary layer mesh", 197 NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 198 ierr = PetscOptionsScalar("-platemesh_top_angle", 199 "Geometric top_angle rate of boundary layer mesh", 200 NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 201 ierr = PetscOptionsString("-platemesh_y_node_locs_path", 202 "Path to file with y node locations. " 203 "If empty, will use the algorithmic mesh warping.", NULL, 204 mesh_ynodes_path, mesh_ynodes_path, 205 sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 206 ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 207 NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 208 PetscOptionsEnd(); 209 210 PetscScalar meter = user->units->meter; 211 PetscScalar second = user->units->second; 212 PetscScalar Kelvin = user->units->Kelvin; 213 PetscScalar Pascal = user->units->Pascal; 214 215 theta0 *= Kelvin; 216 P0 *= Pascal; 217 Uinf *= meter / second; 218 delta0 *= meter; 219 220 PetscReal *mesh_ynodes = NULL; 221 PetscInt mesh_nynodes = 0; 222 if (strcmp(mesh_ynodes_path, "")) { 223 ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 224 CHKERRQ(ierr); 225 } 226 ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 227 mesh_refine_height, mesh_top_angle, &mesh_ynodes, 228 &mesh_nynodes); CHKERRQ(ierr); 229 230 // Some properties depend on parameters from NewtonianIdealGas 231 CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 232 CEED_MEM_HOST, &newtonian_ig_ctx); 233 234 blasius_ctx->weakT = weakT; 235 blasius_ctx->Uinf = Uinf; 236 blasius_ctx->delta0 = delta0; 237 blasius_ctx->theta0 = theta0; 238 blasius_ctx->P0 = P0; 239 newtonian_ig_ctx->P0 = P0; 240 blasius_ctx->implicit = user->phys->implicit; 241 blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 242 243 { 244 PetscReal domain_min[3]; 245 ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); 246 blasius_ctx->x_inflow = domain_min[0]; 247 } 248 249 CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 250 &newtonian_ig_ctx); 251 252 CeedQFunctionContextCreate(user->ceed, &blasius_context); 253 CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 254 CEED_USE_POINTER, 255 sizeof(*blasius_ctx), blasius_ctx); 256 CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 257 FreeContextPetsc); 258 259 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 260 problem->ics.qfunction_context = blasius_context; 261 if (use_stg) { 262 ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, 263 mesh_nynodes); CHKERRQ(ierr); 264 } else { 265 CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context); 266 problem->apply_inflow.qfunction = Blasius_Inflow; 267 problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 268 problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 269 problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 270 CeedQFunctionContextReferenceCopy(blasius_context, 271 &problem->apply_inflow.qfunction_context); 272 CeedQFunctionContextReferenceCopy(blasius_context, 273 &problem->apply_inflow_jacobian.qfunction_context); 274 } 275 ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278