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. 61 */ 62 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 63 PetscReal growth, PetscInt N, 64 PetscReal refine_height, PetscReal top_angle, 65 PetscReal node_locs[], PetscInt num_node_locs) { 66 PetscInt ierr, narr, ncoords; 67 PetscReal domain_min[3], domain_max[3], domain_size[3]; 68 PetscScalar *arr_coords; 69 Vec vec_coords; 70 PetscFunctionBeginUser; 71 72 PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 73 74 // Get domain boundary information 75 ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 76 for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 77 78 // Get coords array from DM 79 ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 80 ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 81 ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 82 83 PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 84 ncoords = narr/dim; 85 86 // Get mesh information 87 PetscInt nmax = 3, faces[3]; 88 ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 89 NULL); CHKERRQ(ierr); 90 // Get element size of the box mesh, for indexing each node 91 const PetscReal dybox = domain_size[1]/faces[1]; 92 93 if (!node_locs) { 94 // Calculate the first element height 95 PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 96 97 // Calculate log of sizing outside BL 98 PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 99 100 for (PetscInt i=0; i<ncoords; i++) { 101 PetscInt y_box_index = round(coords[i][1]/dybox); 102 if (y_box_index <= N) { 103 coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 104 * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 105 } else { 106 PetscInt j = y_box_index - N; 107 coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 108 * exp(log(refine_height) + logdy*j); 109 } 110 } 111 } else { 112 // Error checking 113 if (num_node_locs < faces[1] +1) 114 SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 115 "There are %d + 1 nodes, but only %d locations given", 116 faces[1]+1, num_node_locs); 117 if (num_node_locs > faces[1] +1) { 118 ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 119 "than the mesh has nodes (%d). This maybe unintended.", 120 num_node_locs, faces[1]+1); CHKERRQ(ierr); 121 } 122 123 for (PetscInt i=0; i<ncoords; i++) { 124 // Determine which y-node we're at 125 PetscInt y_box_index = round(coords[i][1]/dybox); 126 coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff) 127 * node_locs[y_box_index]; 128 } 129 } 130 131 ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 132 ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 133 134 PetscFunctionReturn(0); 135 } 136 137 PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 138 139 PetscInt ierr; 140 User user = *(User *)ctx; 141 MPI_Comm comm = PETSC_COMM_WORLD; 142 PetscBool use_stg = PETSC_FALSE; 143 BlasiusContext blasius_ctx; 144 NewtonianIdealGasContext newtonian_ig_ctx; 145 CeedQFunctionContext blasius_context; 146 147 PetscFunctionBeginUser; 148 ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 149 ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 150 151 // ------------------------------------------------------ 152 // SET UP Blasius 153 // ------------------------------------------------------ 154 CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 155 problem->ics.qfunction = ICsBlasius; 156 problem->ics.qfunction_loc = ICsBlasius_loc; 157 problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 158 problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 159 160 CeedScalar Uinf = 40; // m/s 161 CeedScalar delta0 = 4.2e-4; // m 162 CeedScalar theta0 = 288.; // K 163 CeedScalar P0 = 1.01e5; // Pa 164 PetscBool weakT = PETSC_FALSE; // weak density or temperature 165 PetscReal mesh_refine_height = 5.9e-4; // m 166 PetscReal mesh_growth = 1.08; // [-] 167 PetscInt mesh_Ndelta = 45; // [-] 168 PetscReal mesh_top_angle = 5; // degrees 169 char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 170 171 PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); 172 ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 173 NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 174 ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", 175 NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); 176 ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 177 NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 178 ierr = PetscOptionsScalar("-theta0", "Wall temperature", 179 NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 180 ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 181 NULL, P0, &P0, NULL); CHKERRQ(ierr); 182 ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 183 "Velocity at boundary layer edge", 184 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 185 ierr = PetscOptionsScalar("-platemesh_refine_height", 186 "Height of boundary layer mesh refinement", 187 NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 188 ierr = PetscOptionsScalar("-platemesh_growth", 189 "Geometric growth rate of boundary layer mesh", 190 NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 191 ierr = PetscOptionsScalar("-platemesh_top_angle", 192 "Geometric top_angle rate of boundary layer mesh", 193 NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 194 ierr = PetscOptionsString("-platemesh_y_node_locs_path", 195 "Path to file with y node locations. " 196 "If empty, will use the algorithmic mesh warping.", NULL, 197 mesh_ynodes_path, mesh_ynodes_path, 198 sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 199 ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 200 NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 201 PetscOptionsEnd(); 202 203 PetscScalar meter = user->units->meter; 204 PetscScalar second = user->units->second; 205 PetscScalar Kelvin = user->units->Kelvin; 206 PetscScalar Pascal = user->units->Pascal; 207 208 theta0 *= Kelvin; 209 P0 *= Pascal; 210 Uinf *= meter / second; 211 delta0 *= meter; 212 213 PetscReal *mesh_ynodes = NULL; 214 PetscInt mesh_nynodes = 0; 215 if (strcmp(mesh_ynodes_path, "")) { 216 ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 217 CHKERRQ(ierr); 218 } 219 ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 220 mesh_refine_height, mesh_top_angle, mesh_ynodes, 221 mesh_nynodes); CHKERRQ(ierr); 222 223 // Some properties depend on parameters from NewtonianIdealGas 224 CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 225 CEED_MEM_HOST, &newtonian_ig_ctx); 226 227 blasius_ctx->weakT = weakT; 228 blasius_ctx->Uinf = Uinf; 229 blasius_ctx->delta0 = delta0; 230 blasius_ctx->theta0 = theta0; 231 blasius_ctx->P0 = P0; 232 newtonian_ig_ctx->P0 = P0; 233 blasius_ctx->implicit = user->phys->implicit; 234 blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 235 236 { 237 PetscReal domain_min[3]; 238 ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); 239 blasius_ctx->x_inflow = domain_min[0]; 240 } 241 242 CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 243 &newtonian_ig_ctx); 244 245 CeedQFunctionContextCreate(user->ceed, &blasius_context); 246 CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, 247 CEED_USE_POINTER, 248 sizeof(*blasius_ctx), blasius_ctx); 249 CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 250 FreeContextPetsc); 251 252 problem->ics.qfunction_context = blasius_context; 253 CeedQFunctionContextReferenceCopy(blasius_context, 254 &problem->apply_inflow_jacobian.qfunction_context); 255 if (use_stg) { 256 ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, 257 mesh_nynodes); CHKERRQ(ierr); 258 } else { 259 problem->apply_inflow.qfunction = Blasius_Inflow; 260 problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 261 CeedQFunctionContextReferenceCopy(blasius_context, 262 &problem->apply_inflow.qfunction_context); 263 } 264 ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 265 PetscFunctionReturn(0); 266 } 267