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 "../qfunctions/blasius.h" 12*2b916ea7SJeremy L Thompson 13*2b916ea7SJeremy L Thompson #include "../navierstokes.h" 14493642f1SJames Wright #include "stg_shur14.h" 15bb8a0c61SJames Wright 16e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) { 17e0d1a4dfSLeila Ghaffari const BlasiusContext blasius = (BlasiusContext)ctx; 18e0d1a4dfSLeila Ghaffari const PetscScalar *Tf, *Th; // Chebyshev coefficients 19e0d1a4dfSLeila Ghaffari PetscScalar *r, f[4], h[4]; 20e0d1a4dfSLeila Ghaffari PetscInt N = blasius->n_cheb; 21*2b916ea7SJeremy L Thompson PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx), 22e0d1a4dfSLeila Ghaffari gamma = HeatCapacityRatio(&blasius->newtonian_ctx); 23e0d1a4dfSLeila Ghaffari PetscFunctionBegin; 24e0d1a4dfSLeila Ghaffari PetscCall(VecGetArrayRead(X, &Tf)); 25e0d1a4dfSLeila Ghaffari Th = Tf + N; 26e0d1a4dfSLeila Ghaffari PetscCall(VecGetArray(R, &r)); 27e0d1a4dfSLeila Ghaffari 28e0d1a4dfSLeila Ghaffari // Left boundary conditions f = f' = 0 29e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, -1., blasius->eta_max, f); 30e0d1a4dfSLeila Ghaffari r[0] = f[0]; 31e0d1a4dfSLeila Ghaffari r[1] = f[1]; 32e0d1a4dfSLeila Ghaffari 33e0d1a4dfSLeila Ghaffari // f - right end boundary condition 34e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, 1., blasius->eta_max, f); 35e0d1a4dfSLeila Ghaffari r[2] = f[1] - 1.; 36e0d1a4dfSLeila Ghaffari 37e0d1a4dfSLeila Ghaffari for (int i = 0; i < N - 3; i++) { 38e0d1a4dfSLeila Ghaffari ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f); 39e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h); 400d850f2eSLeila Ghaffari // mu and rho generally depend on h. We naively assume constant mu. 410d850f2eSLeila Ghaffari // For an ideal gas at constant pressure, density is inversely proportional to enthalpy. 420d850f2eSLeila Ghaffari // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here. 430d850f2eSLeila Ghaffari const PetscScalar mu_tilde[2] = {1, 0}; 440d850f2eSLeila Ghaffari const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])}; 450d850f2eSLeila Ghaffari const PetscScalar mu_rho_tilde[2] = { 460d850f2eSLeila Ghaffari mu_tilde[0] * rho_tilde[0], 470d850f2eSLeila Ghaffari mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1], 480d850f2eSLeila Ghaffari }; 490d850f2eSLeila Ghaffari r[3 + i] = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0]; 50*2b916ea7SJeremy 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]); 51e0d1a4dfSLeila Ghaffari } 52e0d1a4dfSLeila Ghaffari 53e0d1a4dfSLeila Ghaffari // h - left end boundary condition 54e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h); 55aef1eb53SLeila Ghaffari r[N] = h[0] - blasius->T_wall / blasius->T_inf; 56e0d1a4dfSLeila Ghaffari 57e0d1a4dfSLeila Ghaffari // h - right end boundary condition 58e0d1a4dfSLeila Ghaffari ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h); 59e0d1a4dfSLeila Ghaffari r[N + 1] = h[0] - 1.; 60e0d1a4dfSLeila Ghaffari 61e0d1a4dfSLeila Ghaffari // Restore vectors 62e0d1a4dfSLeila Ghaffari PetscCall(VecRestoreArrayRead(X, &Tf)); 63e0d1a4dfSLeila Ghaffari PetscCall(VecRestoreArray(R, &r)); 64e0d1a4dfSLeila Ghaffari PetscFunctionReturn(0); 65e0d1a4dfSLeila Ghaffari } 66e0d1a4dfSLeila Ghaffari 67e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 68e0d1a4dfSLeila Ghaffari SNES snes; 69e0d1a4dfSLeila Ghaffari Vec sol, res; 70e0d1a4dfSLeila Ghaffari PetscReal *w; 71e0d1a4dfSLeila Ghaffari PetscInt N = blasius->n_cheb; 720d850f2eSLeila Ghaffari SNESConvergedReason reason; 73e0d1a4dfSLeila Ghaffari const PetscScalar *cheb_coefs; 74e0d1a4dfSLeila Ghaffari PetscFunctionBegin; 750d850f2eSLeila Ghaffari 760d850f2eSLeila Ghaffari // Allocate memory 77e0d1a4dfSLeila Ghaffari PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w)); 78e0d1a4dfSLeila Ghaffari PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w)); 790d850f2eSLeila Ghaffari 800d850f2eSLeila Ghaffari // Snes solve 81e0d1a4dfSLeila Ghaffari PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 82e0d1a4dfSLeila Ghaffari PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 83e0d1a4dfSLeila Ghaffari PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1)); 84e0d1a4dfSLeila Ghaffari PetscCall(VecSetFromOptions(sol)); 850d850f2eSLeila Ghaffari // Constant relative enthalpy 1 as initial guess 860d850f2eSLeila Ghaffari PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES)); 87e0d1a4dfSLeila Ghaffari PetscCall(VecDuplicate(sol, &res)); 88e0d1a4dfSLeila Ghaffari PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 890d850f2eSLeila Ghaffari PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_")); 90e0d1a4dfSLeila Ghaffari PetscCall(SNESSetFromOptions(snes)); 91e0d1a4dfSLeila Ghaffari PetscCall(SNESSolve(snes, NULL, sol)); 920d850f2eSLeila Ghaffari PetscCall(SNESGetConvergedReason(snes, &reason)); 93*2b916ea7SJeremy L Thompson if (reason < 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n"); 940d850f2eSLeila Ghaffari 950d850f2eSLeila Ghaffari // Assign Chebyshev coefficients 96e0d1a4dfSLeila Ghaffari PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 97e0d1a4dfSLeila Ghaffari for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 98e0d1a4dfSLeila Ghaffari for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N]; 990d850f2eSLeila Ghaffari 1000d850f2eSLeila Ghaffari // Destroy objects 101e0d1a4dfSLeila Ghaffari PetscCall(PetscFree2(blasius->X, w)); 102e0d1a4dfSLeila Ghaffari PetscCall(VecDestroy(&sol)); 103e0d1a4dfSLeila Ghaffari PetscCall(VecDestroy(&res)); 104e0d1a4dfSLeila Ghaffari PetscCall(SNESDestroy(&snes)); 105e0d1a4dfSLeila Ghaffari PetscFunctionReturn(0); 106e0d1a4dfSLeila Ghaffari } 107e0d1a4dfSLeila Ghaffari 108*2b916ea7SJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) { 109f362fe62SJames Wright PetscInt ndims, dims[2]; 110f362fe62SJames Wright FILE *fp; 111f362fe62SJames Wright const PetscInt char_array_len = 512; 112f362fe62SJames Wright char line[char_array_len]; 113f362fe62SJames Wright char **array; 114f362fe62SJames Wright PetscReal *node_locs; 115f362fe62SJames Wright PetscFunctionBeginUser; 116f362fe62SJames Wright 117*2b916ea7SJeremy L Thompson PetscCall(PetscFOpen(comm, path, "r", &fp)); 118*2b916ea7SJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 119*2b916ea7SJeremy L Thompson PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 120f362fe62SJames Wright 121f362fe62SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 122f362fe62SJames Wright if (ndims < 2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 123f362fe62SJames Wright *nynodes = dims[0]; 124*2b916ea7SJeremy L Thompson PetscCall(PetscMalloc1(*nynodes, &node_locs)); 125f362fe62SJames Wright 126f362fe62SJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 127*2b916ea7SJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 128*2b916ea7SJeremy L Thompson PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 129*2b916ea7SJeremy L Thompson if (ndims < dims[1]) 130*2b916ea7SJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 131*2b916ea7SJeremy L Thompson ndims, dims[1]); 132f362fe62SJames Wright 133f362fe62SJames Wright node_locs[i] = (PetscReal)atof(array[0]); 134f362fe62SJames Wright } 135*2b916ea7SJeremy L Thompson PetscCall(PetscFClose(comm, fp)); 136f362fe62SJames Wright *pynodes = node_locs; 137f362fe62SJames Wright PetscFunctionReturn(0); 138f362fe62SJames Wright } 139f362fe62SJames Wright 140bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius 141bb8a0c61SJames Wright * 142493642f1SJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 143493642f1SJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 144493642f1SJames Wright * linearly in logspace to the top surface. 145bb8a0c61SJames Wright * 146bb8a0c61SJames Wright * The top surface is also angled downwards, so that it may be used as an 147493642f1SJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 148f362fe62SJames Wright * 149f362fe62SJames Wright * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 1508e63b9cbSJames Wright * locations. If it is NULL, then the modified coordinate values will be set in 1518e63b9cbSJames Wright * the array, along with `num_node_locs`. 152bb8a0c61SJames Wright */ 153*2b916ea7SJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle, 1548e63b9cbSJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 155*2b916ea7SJeremy L Thompson PetscInt narr, ncoords; 156bb8a0c61SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 157bb8a0c61SJames Wright PetscScalar *arr_coords; 158bb8a0c61SJames Wright Vec vec_coords; 159bb8a0c61SJames Wright PetscFunctionBeginUser; 160bb8a0c61SJames Wright 161bb8a0c61SJames Wright PetscReal angle_coeff = tan(top_angle * (M_PI / 180)); 162bb8a0c61SJames Wright 163bb8a0c61SJames Wright // Get domain boundary information 164*2b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 165493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 166bb8a0c61SJames Wright 167bb8a0c61SJames Wright // Get coords array from DM 168*2b916ea7SJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &vec_coords)); 169*2b916ea7SJeremy L Thompson PetscCall(VecGetLocalSize(vec_coords, &narr)); 170*2b916ea7SJeremy L Thompson PetscCall(VecGetArray(vec_coords, &arr_coords)); 171bb8a0c61SJames Wright 172bb8a0c61SJames Wright PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords; 173bb8a0c61SJames Wright ncoords = narr / dim; 174bb8a0c61SJames Wright 175bb8a0c61SJames Wright // Get mesh information 176bb8a0c61SJames Wright PetscInt nmax = 3, faces[3]; 177*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 178f362fe62SJames Wright // Get element size of the box mesh, for indexing each node 179f362fe62SJames Wright const PetscReal dybox = domain_size[1] / faces[1]; 180bb8a0c61SJames Wright 1818e63b9cbSJames Wright if (!*node_locs) { 182bb8a0c61SJames Wright // Calculate the first element height 183bb8a0c61SJames Wright PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1); 184bb8a0c61SJames Wright 185bb8a0c61SJames Wright // Calculate log of sizing outside BL 186bb8a0c61SJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 187bb8a0c61SJames Wright 1888e63b9cbSJames Wright *num_node_locs = faces[1] + 1; 1898e63b9cbSJames Wright PetscReal *temp_node_locs; 190*2b916ea7SJeremy L Thompson PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs)); 1918e63b9cbSJames Wright 192493642f1SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 193bb8a0c61SJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 194bb8a0c61SJames Wright if (y_box_index <= N) { 195*2b916ea7SJeremy L Thompson coords[i][1] = 196*2b916ea7SJeremy L Thompson (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1); 197bb8a0c61SJames Wright } else { 198bb8a0c61SJames Wright PetscInt j = y_box_index - N; 199*2b916ea7SJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j); 200f362fe62SJames Wright } 201*2b916ea7SJeremy L Thompson if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1]; 202f362fe62SJames Wright } 2038e63b9cbSJames Wright 2048e63b9cbSJames Wright *node_locs = temp_node_locs; 205f362fe62SJames Wright } else { 206f362fe62SJames Wright // Error checking 207*2b916ea7SJeremy L Thompson if (*num_node_locs < faces[1] + 1) { 208*2b916ea7SJeremy L Thompson SETERRQ(comm, -1, 209*2b916ea7SJeremy L Thompson "The y_node_locs_path has too few locations; " 210f362fe62SJames Wright "There are %d + 1 nodes, but only %d locations given", 2118e63b9cbSJames Wright faces[1] + 1, *num_node_locs); 212*2b916ea7SJeremy L Thompson } 2138e63b9cbSJames Wright if (*num_node_locs > faces[1] + 1) { 214*2b916ea7SJeremy L Thompson PetscCall(PetscPrintf(comm, 215*2b916ea7SJeremy L Thompson "WARNING: y_node_locs_path has more locations (%d) " 216e4677755SJames Wright "than the mesh has nodes (%d). This maybe unintended.\n", 217*2b916ea7SJeremy L Thompson *num_node_locs, faces[1] + 1)); 218f362fe62SJames Wright } 2198e63b9cbSJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 220f362fe62SJames Wright 221f362fe62SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 222f362fe62SJames Wright // Determine which y-node we're at 223f362fe62SJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 224*2b916ea7SJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index]; 225bb8a0c61SJames Wright } 226bb8a0c61SJames Wright } 227bb8a0c61SJames Wright 228*2b916ea7SJeremy L Thompson PetscCall(VecRestoreArray(vec_coords, &arr_coords)); 229*2b916ea7SJeremy L Thompson PetscCall(DMSetCoordinatesLocal(dm, vec_coords)); 230bb8a0c61SJames Wright 231bb8a0c61SJames Wright PetscFunctionReturn(0); 232bb8a0c61SJames Wright } 233bb8a0c61SJames Wright 234b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 235bb8a0c61SJames Wright User user = *(User *)ctx; 236bb8a0c61SJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 237493642f1SJames Wright PetscBool use_stg = PETSC_FALSE; 23815a3537eSJed Brown BlasiusContext blasius_ctx; 23915a3537eSJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 24015a3537eSJed Brown CeedQFunctionContext blasius_context; 24115a3537eSJed Brown 242bb8a0c61SJames Wright PetscFunctionBeginUser; 243*2b916ea7SJeremy L Thompson PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx)); 244*2b916ea7SJeremy L Thompson PetscCall(PetscCalloc1(1, &blasius_ctx)); 245bb8a0c61SJames Wright 246bb8a0c61SJames Wright // ------------------------------------------------------ 247bb8a0c61SJames Wright // SET UP Blasius 248bb8a0c61SJames Wright // ------------------------------------------------------ 2499785fe93SJed Brown problem->ics.qfunction = ICsBlasius; 2509785fe93SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 251bb8a0c61SJames Wright 252aef1eb53SLeila Ghaffari CeedScalar U_inf = 40; // m/s 253aef1eb53SLeila Ghaffari CeedScalar T_inf = 288.; // K 2540d850f2eSLeila Ghaffari CeedScalar T_wall = 288.; // K 255e0d1a4dfSLeila Ghaffari CeedScalar delta0 = 4.2e-3; // m 256bb8a0c61SJames Wright CeedScalar P0 = 1.01e5; // Pa 257e0d1a4dfSLeila Ghaffari CeedInt N = 20; // Number of Chebyshev terms 2582acc7cbcSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 2597470235eSJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 2607470235eSJames Wright PetscReal mesh_growth = 1.08; // [-] 2617470235eSJames Wright PetscInt mesh_Ndelta = 45; // [-] 2627470235eSJames Wright PetscReal mesh_top_angle = 5; // degrees 263f362fe62SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 264bb8a0c61SJames Wright 2653fd71269SJames Wright PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 266*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL)); 267*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL)); 268*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL)); 269*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL)); 270*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL)); 271*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL)); 272*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL)); 273*2b916ea7SJeremy 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, 2740d850f2eSLeila Ghaffari BLASIUS_MAX_N_CHEBYSHEV); 275*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1)); 276*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height, 277*2b916ea7SJeremy L Thompson NULL)); 278*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL)); 279*2b916ea7SJeremy L Thompson PetscCall( 280*2b916ea7SJeremy L Thompson PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL)); 281*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-platemesh_y_node_locs_path", 282f362fe62SJames Wright "Path to file with y node locations. " 283*2b916ea7SJeremy L Thompson "If empty, will use the algorithmic mesh warping.", 284*2b916ea7SJeremy L Thompson NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL)); 285*2b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL)); 286bb8a0c61SJames Wright PetscOptionsEnd(); 287bb8a0c61SJames Wright 288bb8a0c61SJames Wright PetscScalar meter = user->units->meter; 289bb8a0c61SJames Wright PetscScalar second = user->units->second; 290bb8a0c61SJames Wright PetscScalar Kelvin = user->units->Kelvin; 291bb8a0c61SJames Wright PetscScalar Pascal = user->units->Pascal; 292bb8a0c61SJames Wright 293aef1eb53SLeila Ghaffari T_inf *= Kelvin; 294e0d1a4dfSLeila Ghaffari T_wall *= Kelvin; 295bb8a0c61SJames Wright P0 *= Pascal; 296aef1eb53SLeila Ghaffari U_inf *= meter / second; 297bb8a0c61SJames Wright delta0 *= meter; 298bb8a0c61SJames Wright 299f362fe62SJames Wright PetscReal *mesh_ynodes = NULL; 300f362fe62SJames Wright PetscInt mesh_nynodes = 0; 301f362fe62SJames Wright if (strcmp(mesh_ynodes_path, "")) { 302*2b916ea7SJeremy L Thompson PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes)); 303f362fe62SJames Wright } 304*2b916ea7SJeremy L Thompson PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes)); 305bb8a0c61SJames Wright 30615a3537eSJed Brown // Some properties depend on parameters from NewtonianIdealGas 307*2b916ea7SJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 308bb8a0c61SJames Wright 309493642f1SJames Wright blasius_ctx->weakT = weakT; 310aef1eb53SLeila Ghaffari blasius_ctx->U_inf = U_inf; 311aef1eb53SLeila Ghaffari blasius_ctx->T_inf = T_inf; 312e0d1a4dfSLeila Ghaffari blasius_ctx->T_wall = T_wall; 31315a3537eSJed Brown blasius_ctx->delta0 = delta0; 31415a3537eSJed Brown blasius_ctx->P0 = P0; 315e0d1a4dfSLeila Ghaffari blasius_ctx->n_cheb = N; 31604b9037bSJames Wright newtonian_ig_ctx->P0 = P0; 31715a3537eSJed Brown blasius_ctx->implicit = user->phys->implicit; 31815a3537eSJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 319493642f1SJames Wright 320ef2c71fdSJames Wright { 321e0d1a4dfSLeila Ghaffari PetscReal domain_min[3], domain_max[3]; 322*2b916ea7SJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 323ef2c71fdSJames Wright blasius_ctx->x_inflow = domain_min[0]; 324e0d1a4dfSLeila Ghaffari blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 325ef2c71fdSJames Wright } 3263a31d4d2SJames Wright if (!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 327ef2c71fdSJames Wright 328*2b916ea7SJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 329bb8a0c61SJames Wright 33015a3537eSJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 331*2b916ea7SJeremy L Thompson CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx); 332*2b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc); 333bb8a0c61SJames Wright 33443bff553SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 33515a3537eSJed Brown problem->ics.qfunction_context = blasius_context; 336493642f1SJames Wright if (use_stg) { 337*2b916ea7SJeremy L Thompson PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes)); 3388085925cSJames Wright } else { 3398085925cSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 3408085925cSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 3410aa41abfSJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 3420aa41abfSJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 343*2b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context); 344*2b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context); 345493642f1SJames Wright } 346*2b916ea7SJeremy L Thompson PetscCall(PetscFree(mesh_ynodes)); 347bb8a0c61SJames Wright PetscFunctionReturn(0); 348bb8a0c61SJames Wright } 349