1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 "../qfunctions/blasius.h" 122b916ea7SJeremy L Thompson 13e419654dSJeremy L Thompson #include <ceed.h> 14e419654dSJeremy L Thompson #include <petscdm.h> 15e419654dSJeremy L Thompson #include <petscdt.h> 16e419654dSJeremy L Thompson 17149fb536SJames Wright #include <navierstokes.h> 18493642f1SJames Wright #include "stg_shur14.h" 19bb8a0c61SJames Wright 20e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) { 21e0d1a4dfSLeila Ghaffari const BlasiusContext blasius = (BlasiusContext)ctx; 22e0d1a4dfSLeila Ghaffari const PetscScalar *Tf, *Th; // Chebyshev coefficients 23e0d1a4dfSLeila Ghaffari PetscScalar *r, f[4], h[4]; 24e0d1a4dfSLeila Ghaffari PetscInt N = blasius->n_cheb; 25fcb2c22aSJames Wright State S_infty = blasius->S_infty; 26*64667825SJames Wright CeedScalar U_infty = Norm3(S_infty.Y.velocity); 2706f41313SJames Wright 2806f41313SJames Wright PetscFunctionBeginUser; 29fcb2c22aSJames Wright PetscScalar Ma = Mach(&blasius->newtonian_ctx, S_infty.Y.temperature, U_infty), Pr = Prandtl(&blasius->newtonian_ctx), 30e0d1a4dfSLeila Ghaffari gamma = HeatCapacityRatio(&blasius->newtonian_ctx); 3106f41313SJames Wright 32e0d1a4dfSLeila Ghaffari PetscCall(VecGetArrayRead(X, &Tf)); 33e0d1a4dfSLeila Ghaffari Th = Tf + N; 34e0d1a4dfSLeila Ghaffari PetscCall(VecGetArray(R, &r)); 35e0d1a4dfSLeila Ghaffari 36e0d1a4dfSLeila Ghaffari // Left boundary conditions f = f' = 0 37e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, -1., blasius->eta_max, f); 38e0d1a4dfSLeila Ghaffari r[0] = f[0]; 39e0d1a4dfSLeila Ghaffari r[1] = f[1]; 40e0d1a4dfSLeila Ghaffari 41e0d1a4dfSLeila Ghaffari // f - right end boundary condition 42e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, 1., blasius->eta_max, f); 43e0d1a4dfSLeila Ghaffari r[2] = f[1] - 1.; 44e0d1a4dfSLeila Ghaffari 45e0d1a4dfSLeila Ghaffari for (int i = 0; i < N - 3; i++) { 46e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f); 47e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h); 4804e40bb6SJeremy L Thompson // mu and rho generally depend on h. 4904e40bb6SJeremy L Thompson // We naively assume constant mu. 500d850f2eSLeila Ghaffari // For an ideal gas at constant pressure, density is inversely proportional to enthalpy. 510d850f2eSLeila Ghaffari // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here. 520d850f2eSLeila Ghaffari const PetscScalar mu_tilde[2] = {1, 0}; 530d850f2eSLeila Ghaffari const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])}; 540d850f2eSLeila Ghaffari const PetscScalar mu_rho_tilde[2] = { 550d850f2eSLeila Ghaffari mu_tilde[0] * rho_tilde[0], 560d850f2eSLeila Ghaffari mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1], 570d850f2eSLeila Ghaffari }; 580d850f2eSLeila Ghaffari r[3 + i] = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0]; 592b916ea7SJeremy L Thompson 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]); 60e0d1a4dfSLeila Ghaffari } 61e0d1a4dfSLeila Ghaffari 62e0d1a4dfSLeila Ghaffari // h - left end boundary condition 63e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h); 64fcb2c22aSJames Wright r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature; 65e0d1a4dfSLeila Ghaffari 66e0d1a4dfSLeila Ghaffari // h - right end boundary condition 67e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h); 68e0d1a4dfSLeila Ghaffari r[N + 1] = h[0] - 1.; 69e0d1a4dfSLeila Ghaffari 70e0d1a4dfSLeila Ghaffari // Restore vectors 71e0d1a4dfSLeila Ghaffari PetscCall(VecRestoreArrayRead(X, &Tf)); 72e0d1a4dfSLeila Ghaffari PetscCall(VecRestoreArray(R, &r)); 73d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 74e0d1a4dfSLeila Ghaffari } 75e0d1a4dfSLeila Ghaffari 76e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 77e0d1a4dfSLeila Ghaffari SNES snes; 78e0d1a4dfSLeila Ghaffari Vec sol, res; 79e0d1a4dfSLeila Ghaffari PetscReal *w; 80e0d1a4dfSLeila Ghaffari PetscInt N = blasius->n_cheb; 810d850f2eSLeila Ghaffari SNESConvergedReason reason; 82e0d1a4dfSLeila Ghaffari const PetscScalar *cheb_coefs; 830d850f2eSLeila Ghaffari 8406f41313SJames Wright PetscFunctionBeginUser; 850d850f2eSLeila Ghaffari // Allocate memory 86e0d1a4dfSLeila Ghaffari PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w)); 87e0d1a4dfSLeila Ghaffari PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w)); 880d850f2eSLeila Ghaffari 890d850f2eSLeila Ghaffari // Snes solve 90e0d1a4dfSLeila Ghaffari PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 91e0d1a4dfSLeila Ghaffari PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 92e0d1a4dfSLeila Ghaffari PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1)); 93e0d1a4dfSLeila Ghaffari PetscCall(VecSetFromOptions(sol)); 940d850f2eSLeila Ghaffari // Constant relative enthalpy 1 as initial guess 950d850f2eSLeila Ghaffari PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES)); 96e0d1a4dfSLeila Ghaffari PetscCall(VecDuplicate(sol, &res)); 97e0d1a4dfSLeila Ghaffari PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 980d850f2eSLeila Ghaffari PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_")); 99e0d1a4dfSLeila Ghaffari PetscCall(SNESSetFromOptions(snes)); 100e0d1a4dfSLeila Ghaffari PetscCall(SNESSolve(snes, NULL, sol)); 1010d850f2eSLeila Ghaffari PetscCall(SNESGetConvergedReason(snes, &reason)); 1025d28dccaSJames Wright PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n"); 1030d850f2eSLeila Ghaffari 1040d850f2eSLeila Ghaffari // Assign Chebyshev coefficients 105e0d1a4dfSLeila Ghaffari PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 106e0d1a4dfSLeila Ghaffari for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 107e0d1a4dfSLeila Ghaffari for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N]; 1080d850f2eSLeila Ghaffari 1090d850f2eSLeila Ghaffari // Destroy objects 110e0d1a4dfSLeila Ghaffari PetscCall(PetscFree2(blasius->X, w)); 111e0d1a4dfSLeila Ghaffari PetscCall(VecDestroy(&sol)); 112e0d1a4dfSLeila Ghaffari PetscCall(VecDestroy(&res)); 113e0d1a4dfSLeila Ghaffari PetscCall(SNESDestroy(&snes)); 114d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 115e0d1a4dfSLeila Ghaffari } 116e0d1a4dfSLeila Ghaffari 1172b916ea7SJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) { 11887d3884fSJeremy L Thompson int ndims, dims[2] = {0.}; 119f362fe62SJames Wright FILE *fp; 120f362fe62SJames Wright const PetscInt char_array_len = 512; 121f362fe62SJames Wright char line[char_array_len]; 122f362fe62SJames Wright PetscReal *node_locs; 123f362fe62SJames Wright 12406f41313SJames Wright PetscFunctionBeginUser; 1252b916ea7SJeremy L Thompson PetscCall(PetscFOpen(comm, path, "r", &fp)); 1262b916ea7SJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 127039caf0dSJames Wright { 128039caf0dSJames Wright char **array; 129f362fe62SJames Wright 130039caf0dSJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 131f362fe62SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 132039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 133039caf0dSJames Wright } 134f362fe62SJames Wright if (ndims < 2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 135f362fe62SJames Wright *nynodes = dims[0]; 1362b916ea7SJeremy L Thompson PetscCall(PetscMalloc1(*nynodes, &node_locs)); 137f362fe62SJames Wright 138f362fe62SJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 139039caf0dSJames Wright char **array; 140039caf0dSJames Wright 1412b916ea7SJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 1422b916ea7SJeremy L Thompson PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1435d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 144defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]); 145f362fe62SJames Wright 146f362fe62SJames Wright node_locs[i] = (PetscReal)atof(array[0]); 147039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 148f362fe62SJames Wright } 1492b916ea7SJeremy L Thompson PetscCall(PetscFClose(comm, fp)); 150f362fe62SJames Wright *pynodes = node_locs; 151d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 152f362fe62SJames Wright } 153f362fe62SJames Wright 154bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius 155bb8a0c61SJames Wright * 15604e40bb6SJeremy L Thompson * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed 157493642f1SJames Wright * linearly in logspace to the top surface. 158bb8a0c61SJames Wright * 15904e40bb6SJeremy L Thompson * The top surface is also angled downwards, so that it may be used as an outflow. 16004e40bb6SJeremy L Thompson * It's angle is controlled by `top_angle` (in units of degrees). 161f362fe62SJames Wright * 16204e40bb6SJeremy L Thompson * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations. 16304e40bb6SJeremy L Thompson * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`. 164bb8a0c61SJames Wright */ 1652b916ea7SJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle, 1668e63b9cbSJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 1672b916ea7SJeremy L Thompson PetscInt narr, ncoords; 168bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 169bb8a0c61SJames Wright PetscScalar *arr_coords; 170bb8a0c61SJames Wright Vec vec_coords; 17106f41313SJames Wright 172bb8a0c61SJames Wright PetscFunctionBeginUser; 173bb8a0c61SJames Wright PetscReal angle_coeff = tan(top_angle * (M_PI / 180)); 174bb8a0c61SJames Wright // Get domain boundary information 1752b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 176493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 177bb8a0c61SJames Wright 178bb8a0c61SJames Wright // Get coords array from DM 1792b916ea7SJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &vec_coords)); 1802b916ea7SJeremy L Thompson PetscCall(VecGetLocalSize(vec_coords, &narr)); 1812b916ea7SJeremy L Thompson PetscCall(VecGetArray(vec_coords, &arr_coords)); 182bb8a0c61SJames Wright 183bb8a0c61SJames Wright PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords; 184bb8a0c61SJames Wright ncoords = narr / dim; 185bb8a0c61SJames Wright 186bb8a0c61SJames Wright // Get mesh information 187bb8a0c61SJames Wright PetscInt nmax = 3, faces[3]; 1882b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 189f362fe62SJames Wright // Get element size of the box mesh, for indexing each node 190f362fe62SJames Wright const PetscReal dybox = domain_size[1] / faces[1]; 191bb8a0c61SJames Wright 1928e63b9cbSJames Wright if (!*node_locs) { 193bb8a0c61SJames Wright // Calculate the first element height 194bb8a0c61SJames Wright PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1); 195bb8a0c61SJames Wright 196bb8a0c61SJames Wright // Calculate log of sizing outside BL 197bb8a0c61SJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 198bb8a0c61SJames Wright 1998e63b9cbSJames Wright *num_node_locs = faces[1] + 1; 2008e63b9cbSJames Wright PetscReal *temp_node_locs; 2012b916ea7SJeremy L Thompson PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs)); 2028e63b9cbSJames Wright 203493642f1SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 204bb8a0c61SJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 205bb8a0c61SJames Wright if (y_box_index <= N) { 2062b916ea7SJeremy L Thompson coords[i][1] = 2072b916ea7SJeremy L Thompson (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1); 208bb8a0c61SJames Wright } else { 209bb8a0c61SJames Wright PetscInt j = y_box_index - N; 2102b916ea7SJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j); 211f362fe62SJames Wright } 2122b916ea7SJeremy L Thompson if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1]; 213f362fe62SJames Wright } 2148e63b9cbSJames Wright 2158e63b9cbSJames Wright *node_locs = temp_node_locs; 216f362fe62SJames Wright } else { 2175d28dccaSJames Wright PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED, 218defe8520SJames Wright "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given", 219defe8520SJames Wright faces[1] + 1, *num_node_locs); 2208e63b9cbSJames Wright if (*num_node_locs > faces[1] + 1) { 2212b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 222defe8520SJames Wright "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") " 223defe8520SJames Wright "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n", 2242b916ea7SJeremy L Thompson *num_node_locs, faces[1] + 1)); 225f362fe62SJames Wright } 2268e63b9cbSJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 227f362fe62SJames Wright 228f362fe62SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 229f362fe62SJames Wright // Determine which y-node we're at 230f362fe62SJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 2312b916ea7SJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index]; 232bb8a0c61SJames Wright } 233bb8a0c61SJames Wright } 234bb8a0c61SJames Wright 2352b916ea7SJeremy L Thompson PetscCall(VecRestoreArray(vec_coords, &arr_coords)); 2362b916ea7SJeremy L Thompson PetscCall(DMSetCoordinatesLocal(dm, vec_coords)); 237d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 238bb8a0c61SJames Wright } 239bb8a0c61SJames Wright 240991aef52SJames Wright PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 241bb8a0c61SJames Wright User user = *(User *)ctx; 2424b0f6111SJames Wright MPI_Comm comm = user->comm; 243b4c37c5cSJames Wright Ceed ceed = user->ceed; 244493642f1SJames Wright PetscBool use_stg = PETSC_FALSE; 24515a3537eSJed Brown BlasiusContext blasius_ctx; 24615a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 24715a3537eSJed Brown CeedQFunctionContext blasius_context; 24815a3537eSJed Brown 249bb8a0c61SJames Wright PetscFunctionBeginUser; 250d1c51a42SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 2512b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &blasius_ctx)); 252bb8a0c61SJames Wright 253bb8a0c61SJames Wright // ------------------------------------------------------ 254bb8a0c61SJames Wright // SET UP Blasius 255bb8a0c61SJames Wright // ------------------------------------------------------ 2569785fe93SJed Brown problem->ics.qfunction = ICsBlasius; 2579785fe93SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 258bb8a0c61SJames Wright 259aef1eb53SLeila Ghaffari CeedScalar U_inf = 40; // m/s 260aef1eb53SLeila Ghaffari CeedScalar T_inf = 288.; // K 2610d850f2eSLeila Ghaffari CeedScalar T_wall = 288.; // K 262e0d1a4dfSLeila Ghaffari CeedScalar delta0 = 4.2e-3; // m 263fcb2c22aSJames Wright CeedScalar P_inf = 1.01e5; // Pa 264defe8520SJames Wright PetscInt N = 20; // Number of Chebyshev terms 2652acc7cbcSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 2667470235eSJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 2677470235eSJames Wright PetscReal mesh_growth = 1.08; // [-] 2687470235eSJames Wright PetscInt mesh_Ndelta = 45; // [-] 2697470235eSJames Wright PetscReal mesh_top_angle = 5; // degrees 270f362fe62SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 271a213b8aaSJames Wright PetscBool P0_set; 272bb8a0c61SJames Wright 2733fd71269SJames Wright PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 2742b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL)); 2752b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL)); 2762b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL)); 277a213b8aaSJames Wright PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set)); // For maintaining behavior of -P0 flag (which is deprecated) 278a213b8aaSJames Wright PetscCall( 279a213b8aaSJames Wright PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0", 280a213b8aaSJames Wright "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure")); 281a213b8aaSJames Wright PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL)); 2822b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL)); 2832b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL)); 2842b916ea7SJeremy L Thompson PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL)); 2852b916ea7SJeremy L Thompson PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N, 2860d850f2eSLeila Ghaffari BLASIUS_MAX_N_CHEBYSHEV); 287f31f4833SJames Wright if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) { 2882b916ea7SJeremy L Thompson PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1)); 289f31f4833SJames Wright PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, 290f31f4833SJames Wright &mesh_refine_height, NULL)); 2912b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL)); 2922b916ea7SJeremy L Thompson PetscCall( 2932b916ea7SJeremy L Thompson PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL)); 2942b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-platemesh_y_node_locs_path", 295f362fe62SJames Wright "Path to file with y node locations. " 2962b916ea7SJeremy L Thompson "If empty, will use the algorithmic mesh warping.", 2972b916ea7SJeremy L Thompson NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL)); 298f31f4833SJames Wright } 2992b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL)); 300bb8a0c61SJames Wright PetscOptionsEnd(); 301bb8a0c61SJames Wright 302bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 303bb8a0c61SJames Wright PetscScalar second = user->units->second; 304bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 305bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 306bb8a0c61SJames Wright 307aef1eb53SLeila Ghaffari T_inf *= Kelvin; 308e0d1a4dfSLeila Ghaffari T_wall *= Kelvin; 309fcb2c22aSJames Wright P_inf *= Pascal; 310aef1eb53SLeila Ghaffari U_inf *= meter / second; 311bb8a0c61SJames Wright delta0 *= meter; 312bb8a0c61SJames Wright 313f31f4833SJames Wright if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) { 314f362fe62SJames Wright PetscReal *mesh_ynodes = NULL; 315f362fe62SJames Wright PetscInt mesh_nynodes = 0; 316c029f0c5SJames Wright if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes)); 3172b916ea7SJeremy L Thompson PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes)); 3189ef62cddSJames Wright PetscCall(PetscFree(mesh_ynodes)); 319c029f0c5SJames Wright } 320bb8a0c61SJames Wright 32115a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 322b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 323bb8a0c61SJames Wright 324fcb2c22aSJames Wright StatePrimitive Y_inf = { 325fcb2c22aSJames Wright .pressure = P_inf, .velocity = {U_inf, 0, 0}, 326fcb2c22aSJames Wright .temperature = T_inf 327fcb2c22aSJames Wright }; 328fcb2c22aSJames Wright State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf); 329fcb2c22aSJames Wright 330493642f1SJames Wright blasius_ctx->weakT = weakT; 331e0d1a4dfSLeila Ghaffari blasius_ctx->T_wall = T_wall; 33215a3537eSJed Brown blasius_ctx->delta0 = delta0; 333fcb2c22aSJames Wright blasius_ctx->S_infty = S_infty; 334e0d1a4dfSLeila Ghaffari blasius_ctx->n_cheb = N; 33515a3537eSJed Brown blasius_ctx->implicit = user->phys->implicit; 336a213b8aaSJames Wright if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf; // For maintaining behavior of -P0 flag (which is deprecated) 33715a3537eSJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 338493642f1SJames Wright 339ef2c71fdSJames Wright { 340e0d1a4dfSLeila Ghaffari PetscReal domain_min[3], domain_max[3]; 3412b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 342ef2c71fdSJames Wright blasius_ctx->x_inflow = domain_min[0]; 343e0d1a4dfSLeila Ghaffari blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 344ef2c71fdSJames Wright } 34525c92e8fSJames Wright PetscBool diff_filter_mms = PETSC_FALSE; 34625c92e8fSJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL)); 34725c92e8fSJames Wright if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 348ef2c71fdSJames Wright 349b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 350bb8a0c61SJames Wright 351b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &blasius_context)); 352b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx)); 353b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc)); 354bb8a0c61SJames Wright 355b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 35615a3537eSJed Brown problem->ics.qfunction_context = blasius_context; 357493642f1SJames Wright if (use_stg) { 358fcb2c22aSJames Wright PetscCall(SetupStg(comm, dm, problem, user, weakT, S_infty.Y.temperature, S_infty.Y.pressure)); 35925c92e8fSJames Wright } else if (diff_filter_mms) { 36042454adaSJames Wright PetscCall(DifferentialFilterMmsICSetup(problem)); 3618085925cSJames Wright } else { 362727ec889SJames Wright PetscCheck((user->phys->state_var == STATEVAR_CONSERVATIVE) || (user->app_ctx->test_type == TESTTYPE_DIFF_FILTER), user->comm, 363727ec889SJames Wright PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow"); 3648085925cSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 3658085925cSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 3660aa41abfSJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 3670aa41abfSJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 368b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context)); 369b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context)); 370493642f1SJames Wright } 371d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 372bb8a0c61SJames Wright } 373