1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up Blasius Boundary Layer 6 7 #include "../qfunctions/blasius.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 #include <petscdt.h> 12 13 #include <navierstokes.h> 14 #include "stg_shur14.h" 15 16 PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) { 17 const BlasiusContext blasius = (BlasiusContext)ctx; 18 const PetscScalar *Tf, *Th; // Chebyshev coefficients 19 PetscScalar *r, f[4], h[4]; 20 PetscInt N = blasius->n_cheb; 21 State S_infty = blasius->S_infty; 22 CeedScalar U_infty = Norm3(S_infty.Y.velocity); 23 24 PetscFunctionBeginUser; 25 PetscScalar Ma = Mach(&blasius->newtonian_ctx, S_infty.Y.temperature, U_infty), Pr = Prandtl(&blasius->newtonian_ctx), 26 gamma = HeatCapacityRatio(&blasius->newtonian_ctx); 27 28 PetscCall(VecGetArrayRead(X, &Tf)); 29 Th = Tf + N; 30 PetscCall(VecGetArray(R, &r)); 31 32 // Left boundary conditions f = f' = 0 33 ChebyshevEval(N, Tf, -1., blasius->eta_max, f); 34 r[0] = f[0]; 35 r[1] = f[1]; 36 37 // f - right end boundary condition 38 ChebyshevEval(N, Tf, 1., blasius->eta_max, f); 39 r[2] = f[1] - 1.; 40 41 for (int i = 0; i < N - 3; i++) { 42 ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f); 43 ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h); 44 // mu and rho generally depend on h. 45 // We naively assume constant mu. 46 // For an ideal gas at constant pressure, density is inversely proportional to enthalpy. 47 // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here. 48 const PetscScalar mu_tilde[2] = {1, 0}; 49 const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])}; 50 const PetscScalar mu_rho_tilde[2] = { 51 mu_tilde[0] * rho_tilde[0], 52 mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1], 53 }; 54 r[3 + i] = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0]; 55 r[N + 2 + i] = (mu_rho_tilde[0] * h[2] + mu_rho_tilde[1] * h[1]) + Pr * f[0] * h[1] + Pr * (gamma - 1) * mu_rho_tilde[0] * PetscSqr(Ma * f[2]); 56 } 57 58 // h - left end boundary condition 59 ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h); 60 r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature; 61 62 // h - right end boundary condition 63 ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h); 64 r[N + 1] = h[0] - 1.; 65 66 // Restore vectors 67 PetscCall(VecRestoreArrayRead(X, &Tf)); 68 PetscCall(VecRestoreArray(R, &r)); 69 PetscFunctionReturn(PETSC_SUCCESS); 70 } 71 72 PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 73 SNES snes; 74 Vec sol, res; 75 PetscReal *w; 76 PetscInt N = blasius->n_cheb; 77 SNESConvergedReason reason; 78 const PetscScalar *cheb_coefs; 79 80 PetscFunctionBeginUser; 81 // Allocate memory 82 PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w)); 83 PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w)); 84 85 // Snes solve 86 PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 87 PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 88 PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1)); 89 PetscCall(VecSetFromOptions(sol)); 90 // Constant relative enthalpy 1 as initial guess 91 PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES)); 92 PetscCall(VecDuplicate(sol, &res)); 93 PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 94 PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_")); 95 PetscCall(SNESSetFromOptions(snes)); 96 PetscCall(SNESSolve(snes, NULL, sol)); 97 PetscCall(SNESGetConvergedReason(snes, &reason)); 98 PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n"); 99 100 // Assign Chebyshev coefficients 101 PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 102 for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 103 for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N]; 104 105 // Destroy objects 106 PetscCall(PetscFree2(blasius->X, w)); 107 PetscCall(VecDestroy(&sol)); 108 PetscCall(VecDestroy(&res)); 109 PetscCall(SNESDestroy(&snes)); 110 PetscFunctionReturn(PETSC_SUCCESS); 111 } 112 113 static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) { 114 int ndims, dims[2] = {0.}; 115 FILE *fp; 116 const PetscInt char_array_len = 512; 117 char line[char_array_len]; 118 PetscReal *node_locs; 119 120 PetscFunctionBeginUser; 121 PetscCall(PetscFOpen(comm, path, "r", &fp)); 122 PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 123 { 124 char **array; 125 126 PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 127 for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 128 PetscCall(PetscStrToArrayDestroy(ndims, array)); 129 } 130 if (ndims < 2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 131 *nynodes = dims[0]; 132 PetscCall(PetscMalloc1(*nynodes, &node_locs)); 133 134 for (PetscInt i = 0; i < dims[0]; i++) { 135 char **array; 136 137 PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 138 PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 139 PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 140 "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]); 141 142 node_locs[i] = (PetscReal)atof(array[0]); 143 PetscCall(PetscStrToArrayDestroy(ndims, array)); 144 } 145 PetscCall(PetscFClose(comm, fp)); 146 *pynodes = node_locs; 147 PetscFunctionReturn(PETSC_SUCCESS); 148 } 149 150 /* \brief Modify the domain and mesh for blasius 151 * 152 * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed 153 * linearly in logspace to the top surface. 154 * 155 * The top surface is also angled downwards, so that it may be used as an outflow. 156 * It's angle is controlled by `top_angle` (in units of degrees). 157 * 158 * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations. 159 * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`. 160 */ 161 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle, 162 PetscReal *node_locs[], PetscInt *num_node_locs) { 163 PetscInt narr, ncoords, dim; 164 PetscReal domain_min[3], domain_max[3], domain_size[3]; 165 PetscScalar *arr_coords; 166 Vec vec_coords; 167 168 PetscFunctionBeginUser; 169 PetscCall(DMGetDimension(dm, &dim)); 170 PetscReal angle_coeff = tan(top_angle * (M_PI / 180)); 171 // Get domain boundary information 172 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 173 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 174 175 // Get coords array from DM 176 PetscCall(DMGetCoordinatesLocal(dm, &vec_coords)); 177 PetscCall(VecGetLocalSize(vec_coords, &narr)); 178 PetscCall(VecGetArray(vec_coords, &arr_coords)); 179 180 PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords; 181 ncoords = narr / dim; 182 183 // Get mesh information 184 PetscInt nmax = 3, faces[3]; 185 PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 186 // Get element size of the box mesh, for indexing each node 187 const PetscReal dybox = domain_size[1] / faces[1]; 188 189 if (!*node_locs) { 190 // Calculate the first element height 191 PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1); 192 193 // Calculate log of sizing outside BL 194 PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 195 196 *num_node_locs = faces[1] + 1; 197 PetscReal *temp_node_locs; 198 PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs)); 199 200 for (PetscInt i = 0; i < ncoords; i++) { 201 PetscInt y_box_index = round(coords[i][1] / dybox); 202 if (y_box_index <= N) { 203 coords[i][1] = 204 (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1); 205 } else { 206 PetscInt j = y_box_index - N; 207 coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j); 208 } 209 if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1]; 210 } 211 212 *node_locs = temp_node_locs; 213 } else { 214 PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED, 215 "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given", 216 faces[1] + 1, *num_node_locs); 217 if (*num_node_locs > faces[1] + 1) { 218 PetscCall(PetscPrintf(comm, 219 "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") " 220 "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n", 221 *num_node_locs, faces[1] + 1)); 222 } 223 PetscScalar max_y = (*node_locs)[faces[1]]; 224 225 for (PetscInt i = 0; i < ncoords; i++) { 226 // Determine which y-node we're at 227 PetscInt y_box_index = round(coords[i][1] / dybox); 228 coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index]; 229 } 230 } 231 232 PetscCall(VecRestoreArray(vec_coords, &arr_coords)); 233 PetscCall(DMSetCoordinatesLocal(dm, vec_coords)); 234 PetscFunctionReturn(PETSC_SUCCESS); 235 } 236 237 PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 238 Honee honee = *(Honee *)ctx; 239 MPI_Comm comm = honee->comm; 240 Ceed ceed = honee->ceed; 241 PetscBool use_stg = PETSC_FALSE; 242 BlasiusContext blasius_ctx; 243 NewtonianIdealGasContext newtonian_ig_ctx; 244 CeedQFunctionContext blasius_qfctx; 245 246 PetscFunctionBeginUser; 247 PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 248 PetscCall(PetscCalloc1(1, &blasius_ctx)); 249 250 // ------------------------------------------------------ 251 // SET UP Blasius 252 // ------------------------------------------------------ 253 problem->ics.qf_func_ptr = ICsBlasius; 254 problem->ics.qf_loc = ICsBlasius_loc; 255 256 CeedScalar U_inf = 40; // m/s 257 CeedScalar T_inf = 288.; // K 258 CeedScalar T_wall = 288.; // K 259 CeedScalar delta0 = 4.2e-3; // m 260 CeedScalar P_inf = 1.01e5; // Pa 261 PetscInt N = 20; // Number of Chebyshev terms 262 PetscBool weakT = PETSC_FALSE; // weak density or temperature 263 PetscReal mesh_refine_height = 5.9e-4; // m 264 PetscReal mesh_growth = 1.08; // [-] 265 PetscInt mesh_Ndelta = 45; // [-] 266 PetscReal mesh_top_angle = 5; // degrees 267 char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 268 PetscBool P0_set; 269 270 PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 271 PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL)); 272 PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL)); 273 PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL)); 274 PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set)); // For maintaining behavior of -P0 flag (which is deprecated) 275 PetscCall( 276 PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0", 277 "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure")); 278 PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL)); 279 PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL)); 280 PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL)); 281 PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL)); 282 PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N, 283 BLASIUS_MAX_N_CHEBYSHEV); 284 if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) { 285 PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1)); 286 PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, 287 &mesh_refine_height, NULL)); 288 PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL)); 289 PetscCall( 290 PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL)); 291 PetscCall(PetscOptionsString("-platemesh_y_node_locs_path", 292 "Path to file with y node locations. " 293 "If empty, will use the algorithmic mesh warping.", 294 NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL)); 295 } 296 PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL)); 297 PetscOptionsEnd(); 298 299 PetscScalar meter = honee->units->meter; 300 PetscScalar second = honee->units->second; 301 PetscScalar Kelvin = honee->units->Kelvin; 302 PetscScalar Pascal = honee->units->Pascal; 303 304 T_inf *= Kelvin; 305 T_wall *= Kelvin; 306 P_inf *= Pascal; 307 U_inf *= meter / second; 308 delta0 *= meter; 309 310 if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) { 311 PetscReal *mesh_ynodes = NULL; 312 PetscInt mesh_nynodes = 0; 313 if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes)); 314 PetscCall(ModifyMesh(comm, dm, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes)); 315 PetscCall(PetscFree(mesh_ynodes)); 316 } 317 318 // Some properties depend on parameters from NewtonianIdealGas 319 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 320 321 StatePrimitive Y_inf = { 322 .pressure = P_inf, .velocity = {U_inf, 0, 0}, 323 .temperature = T_inf 324 }; 325 State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf); 326 327 blasius_ctx->weakT = weakT; 328 blasius_ctx->T_wall = T_wall; 329 blasius_ctx->delta0 = delta0; 330 blasius_ctx->S_infty = S_infty; 331 blasius_ctx->n_cheb = N; 332 blasius_ctx->implicit = honee->phys->implicit; 333 if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf; // For maintaining behavior of -P0 flag (which is deprecated) 334 blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 335 336 { 337 PetscReal domain_min[3], domain_max[3]; 338 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 339 blasius_ctx->x_inflow = domain_min[0]; 340 blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 341 } 342 PetscBool diff_filter_mms = PETSC_FALSE; 343 PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL)); 344 if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 345 346 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 347 348 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &blasius_qfctx)); 349 PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx)); 350 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 351 352 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 353 problem->ics.qfctx = blasius_qfctx; 354 if (use_stg) { 355 PetscCall(SetupStg(comm, dm, problem, honee, weakT, S_infty.Y.temperature, S_infty.Y.pressure)); 356 } else if (diff_filter_mms) { 357 PetscCall(DifferentialFilterMmsICSetup(problem)); 358 } else { 359 PetscCheck((honee->phys->state_var == STATEVAR_CONSERVATIVE) || (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER), honee->comm, 360 PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow"); 361 problem->apply_inflow.qf_func_ptr = Blasius_Inflow; 362 problem->apply_inflow.qf_loc = Blasius_Inflow_loc; 363 problem->apply_inflow_jacobian.qf_func_ptr = Blasius_Inflow_Jacobian; 364 problem->apply_inflow_jacobian.qf_loc = Blasius_Inflow_Jacobian_loc; 365 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_qfctx, &problem->apply_inflow.qfctx)); 366 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_qfctx, &problem->apply_inflow_jacobian.qfctx)); 367 } 368 PetscFunctionReturn(PETSC_SUCCESS); 369 } 370