188626eedSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 288626eedSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 388626eedSJames Wright // 488626eedSJames Wright // SPDX-License-Identifier: BSD-2-Clause 588626eedSJames Wright // 688626eedSJames Wright // This file is part of CEED: http://github.com/ceed 788626eedSJames Wright 888626eedSJames Wright /// @file 988626eedSJames Wright /// Utility functions for setting up Blasius Boundary Layer 1088626eedSJames Wright 1188626eedSJames Wright #include "../qfunctions/blasius.h" 122b730f8bSJeremy L Thompson 13*49aac155SJeremy L Thompson #include <ceed.h> 14*49aac155SJeremy L Thompson #include <petscdm.h> 15*49aac155SJeremy L Thompson #include <petscdt.h> 16*49aac155SJeremy L Thompson 172b730f8bSJeremy L Thompson #include "../navierstokes.h" 18ba6664aeSJames Wright #include "stg_shur14.h" 1988626eedSJames Wright 202518f336SLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) { 212518f336SLeila Ghaffari const BlasiusContext blasius = (BlasiusContext)ctx; 222518f336SLeila Ghaffari const PetscScalar *Tf, *Th; // Chebyshev coefficients 232518f336SLeila Ghaffari PetscScalar *r, f[4], h[4]; 242518f336SLeila Ghaffari PetscInt N = blasius->n_cheb; 252b730f8bSJeremy L Thompson PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx), 262518f336SLeila Ghaffari gamma = HeatCapacityRatio(&blasius->newtonian_ctx); 272518f336SLeila Ghaffari PetscFunctionBegin; 282518f336SLeila Ghaffari PetscCall(VecGetArrayRead(X, &Tf)); 292518f336SLeila Ghaffari Th = Tf + N; 302518f336SLeila Ghaffari PetscCall(VecGetArray(R, &r)); 312518f336SLeila Ghaffari 322518f336SLeila Ghaffari // Left boundary conditions f = f' = 0 332518f336SLeila Ghaffari ChebyshevEval(N, Tf, -1., blasius->eta_max, f); 342518f336SLeila Ghaffari r[0] = f[0]; 352518f336SLeila Ghaffari r[1] = f[1]; 362518f336SLeila Ghaffari 372518f336SLeila Ghaffari // f - right end boundary condition 382518f336SLeila Ghaffari ChebyshevEval(N, Tf, 1., blasius->eta_max, f); 392518f336SLeila Ghaffari r[2] = f[1] - 1.; 402518f336SLeila Ghaffari 412518f336SLeila Ghaffari for (int i = 0; i < N - 3; i++) { 422518f336SLeila Ghaffari ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f); 432518f336SLeila Ghaffari ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h); 44ea61e9acSJeremy L Thompson // mu and rho generally depend on h. 45ea61e9acSJeremy L Thompson // We naively assume constant mu. 4607d14e58SLeila Ghaffari // For an ideal gas at constant pressure, density is inversely proportional to enthalpy. 4707d14e58SLeila Ghaffari // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here. 4807d14e58SLeila Ghaffari const PetscScalar mu_tilde[2] = {1, 0}; 4907d14e58SLeila Ghaffari const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])}; 5007d14e58SLeila Ghaffari const PetscScalar mu_rho_tilde[2] = { 5107d14e58SLeila Ghaffari mu_tilde[0] * rho_tilde[0], 5207d14e58SLeila Ghaffari mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1], 5307d14e58SLeila Ghaffari }; 5407d14e58SLeila Ghaffari r[3 + i] = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0]; 552b730f8bSJeremy 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]); 562518f336SLeila Ghaffari } 572518f336SLeila Ghaffari 582518f336SLeila Ghaffari // h - left end boundary condition 592518f336SLeila Ghaffari ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h); 60fb455ff0SLeila Ghaffari r[N] = h[0] - blasius->T_wall / blasius->T_inf; 612518f336SLeila Ghaffari 622518f336SLeila Ghaffari // h - right end boundary condition 632518f336SLeila Ghaffari ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h); 642518f336SLeila Ghaffari r[N + 1] = h[0] - 1.; 652518f336SLeila Ghaffari 662518f336SLeila Ghaffari // Restore vectors 672518f336SLeila Ghaffari PetscCall(VecRestoreArrayRead(X, &Tf)); 682518f336SLeila Ghaffari PetscCall(VecRestoreArray(R, &r)); 692518f336SLeila Ghaffari PetscFunctionReturn(0); 702518f336SLeila Ghaffari } 712518f336SLeila Ghaffari 722518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 732518f336SLeila Ghaffari SNES snes; 742518f336SLeila Ghaffari Vec sol, res; 752518f336SLeila Ghaffari PetscReal *w; 762518f336SLeila Ghaffari PetscInt N = blasius->n_cheb; 7707d14e58SLeila Ghaffari SNESConvergedReason reason; 782518f336SLeila Ghaffari const PetscScalar *cheb_coefs; 792518f336SLeila Ghaffari PetscFunctionBegin; 8007d14e58SLeila Ghaffari 8107d14e58SLeila Ghaffari // Allocate memory 822518f336SLeila Ghaffari PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w)); 832518f336SLeila Ghaffari PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w)); 8407d14e58SLeila Ghaffari 8507d14e58SLeila Ghaffari // Snes solve 862518f336SLeila Ghaffari PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 872518f336SLeila Ghaffari PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 882518f336SLeila Ghaffari PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1)); 892518f336SLeila Ghaffari PetscCall(VecSetFromOptions(sol)); 9007d14e58SLeila Ghaffari // Constant relative enthalpy 1 as initial guess 9107d14e58SLeila Ghaffari PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES)); 922518f336SLeila Ghaffari PetscCall(VecDuplicate(sol, &res)); 932518f336SLeila Ghaffari PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 9407d14e58SLeila Ghaffari PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_")); 952518f336SLeila Ghaffari PetscCall(SNESSetFromOptions(snes)); 962518f336SLeila Ghaffari PetscCall(SNESSolve(snes, NULL, sol)); 9707d14e58SLeila Ghaffari PetscCall(SNESGetConvergedReason(snes, &reason)); 982b730f8bSJeremy L Thompson if (reason < 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n"); 9907d14e58SLeila Ghaffari 10007d14e58SLeila Ghaffari // Assign Chebyshev coefficients 1012518f336SLeila Ghaffari PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 1022518f336SLeila Ghaffari for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 1032518f336SLeila Ghaffari for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N]; 10407d14e58SLeila Ghaffari 10507d14e58SLeila Ghaffari // Destroy objects 1062518f336SLeila Ghaffari PetscCall(PetscFree2(blasius->X, w)); 1072518f336SLeila Ghaffari PetscCall(VecDestroy(&sol)); 1082518f336SLeila Ghaffari PetscCall(VecDestroy(&res)); 1092518f336SLeila Ghaffari PetscCall(SNESDestroy(&snes)); 1102518f336SLeila Ghaffari PetscFunctionReturn(0); 1112518f336SLeila Ghaffari } 1122518f336SLeila Ghaffari 1132b730f8bSJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) { 114d2714514SJames Wright PetscInt ndims, dims[2]; 115d2714514SJames Wright FILE *fp; 116d2714514SJames Wright const PetscInt char_array_len = 512; 117d2714514SJames Wright char line[char_array_len]; 118d2714514SJames Wright char **array; 119d2714514SJames Wright PetscReal *node_locs; 120d2714514SJames Wright PetscFunctionBeginUser; 121d2714514SJames Wright 1222b730f8bSJeremy L Thompson PetscCall(PetscFOpen(comm, path, "r", &fp)); 1232b730f8bSJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 1242b730f8bSJeremy L Thompson PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 125d2714514SJames Wright 126d2714514SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 127d2714514SJames Wright if (ndims < 2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 128d2714514SJames Wright *nynodes = dims[0]; 1292b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(*nynodes, &node_locs)); 130d2714514SJames Wright 131d2714514SJames Wright for (PetscInt i = 0; i < dims[0]; i++) { 1322b730f8bSJeremy L Thompson PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 1332b730f8bSJeremy L Thompson PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 134*49aac155SJeremy L Thompson if (ndims < dims[1]) { 1352b730f8bSJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 1362b730f8bSJeremy L Thompson ndims, dims[1]); 137*49aac155SJeremy L Thompson } 138d2714514SJames Wright 139d2714514SJames Wright node_locs[i] = (PetscReal)atof(array[0]); 140d2714514SJames Wright } 1412b730f8bSJeremy L Thompson PetscCall(PetscFClose(comm, fp)); 142d2714514SJames Wright *pynodes = node_locs; 143d2714514SJames Wright PetscFunctionReturn(0); 144d2714514SJames Wright } 145d2714514SJames Wright 14688626eedSJames Wright /* \brief Modify the domain and mesh for blasius 14788626eedSJames Wright * 148ea61e9acSJeremy L Thompson * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed 149ba6664aeSJames Wright * linearly in logspace to the top surface. 15088626eedSJames Wright * 151ea61e9acSJeremy L Thompson * The top surface is also angled downwards, so that it may be used as an outflow. 152ea61e9acSJeremy L Thompson * It's angle is controlled by `top_angle` (in units of degrees). 153d2714514SJames Wright * 154ea61e9acSJeremy L Thompson * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations. 155ea61e9acSJeremy L Thompson * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`. 15688626eedSJames Wright */ 1572b730f8bSJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle, 158798d42b8SJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 1592b730f8bSJeremy L Thompson PetscInt narr, ncoords; 16088626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 16188626eedSJames Wright PetscScalar *arr_coords; 16288626eedSJames Wright Vec vec_coords; 16388626eedSJames Wright PetscFunctionBeginUser; 16488626eedSJames Wright 16588626eedSJames Wright PetscReal angle_coeff = tan(top_angle * (M_PI / 180)); 16688626eedSJames Wright 16788626eedSJames Wright // Get domain boundary information 1682b730f8bSJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 169ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 17088626eedSJames Wright 17188626eedSJames Wright // Get coords array from DM 1722b730f8bSJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &vec_coords)); 1732b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(vec_coords, &narr)); 1742b730f8bSJeremy L Thompson PetscCall(VecGetArray(vec_coords, &arr_coords)); 17588626eedSJames Wright 17688626eedSJames Wright PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords; 17788626eedSJames Wright ncoords = narr / dim; 17888626eedSJames Wright 17988626eedSJames Wright // Get mesh information 18088626eedSJames Wright PetscInt nmax = 3, faces[3]; 1812b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 182d2714514SJames Wright // Get element size of the box mesh, for indexing each node 183d2714514SJames Wright const PetscReal dybox = domain_size[1] / faces[1]; 18488626eedSJames Wright 185798d42b8SJames Wright if (!*node_locs) { 18688626eedSJames Wright // Calculate the first element height 18788626eedSJames Wright PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1); 18888626eedSJames Wright 18988626eedSJames Wright // Calculate log of sizing outside BL 19088626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 19188626eedSJames Wright 192798d42b8SJames Wright *num_node_locs = faces[1] + 1; 193798d42b8SJames Wright PetscReal *temp_node_locs; 1942b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs)); 195798d42b8SJames Wright 196ba6664aeSJames Wright for (PetscInt i = 0; i < ncoords; i++) { 19788626eedSJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 19888626eedSJames Wright if (y_box_index <= N) { 1992b730f8bSJeremy L Thompson coords[i][1] = 2002b730f8bSJeremy L Thompson (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1); 20188626eedSJames Wright } else { 20288626eedSJames Wright PetscInt j = y_box_index - N; 2032b730f8bSJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j); 204d2714514SJames Wright } 2052b730f8bSJeremy L Thompson if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1]; 206d2714514SJames Wright } 207798d42b8SJames Wright 208798d42b8SJames Wright *node_locs = temp_node_locs; 209d2714514SJames Wright } else { 210d2714514SJames Wright // Error checking 2112b730f8bSJeremy L Thompson if (*num_node_locs < faces[1] + 1) { 2122b730f8bSJeremy L Thompson SETERRQ(comm, -1, 2132b730f8bSJeremy L Thompson "The y_node_locs_path has too few locations; " 214d2714514SJames Wright "There are %d + 1 nodes, but only %d locations given", 215798d42b8SJames Wright faces[1] + 1, *num_node_locs); 2162b730f8bSJeremy L Thompson } 217798d42b8SJames Wright if (*num_node_locs > faces[1] + 1) { 2182b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 2192b730f8bSJeremy L Thompson "WARNING: y_node_locs_path has more locations (%d) " 220c10dcd27SJames Wright "than the mesh has nodes (%d). This maybe unintended.\n", 2212b730f8bSJeremy L Thompson *num_node_locs, faces[1] + 1)); 222d2714514SJames Wright } 223798d42b8SJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 224d2714514SJames Wright 225d2714514SJames Wright for (PetscInt i = 0; i < ncoords; i++) { 226d2714514SJames Wright // Determine which y-node we're at 227d2714514SJames Wright PetscInt y_box_index = round(coords[i][1] / dybox); 2282b730f8bSJeremy L Thompson coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index]; 22988626eedSJames Wright } 23088626eedSJames Wright } 23188626eedSJames Wright 2322b730f8bSJeremy L Thompson PetscCall(VecRestoreArray(vec_coords, &arr_coords)); 2332b730f8bSJeremy L Thompson PetscCall(DMSetCoordinatesLocal(dm, vec_coords)); 23488626eedSJames Wright 23588626eedSJames Wright PetscFunctionReturn(0); 23688626eedSJames Wright } 23788626eedSJames Wright 23846de7363SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) { 23988626eedSJames Wright User user = *(User *)ctx; 24088626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 241ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 242841e4c73SJed Brown BlasiusContext blasius_ctx; 243841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 244841e4c73SJed Brown CeedQFunctionContext blasius_context; 245841e4c73SJed Brown 24688626eedSJames Wright PetscFunctionBeginUser; 24746de7363SJames Wright PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 2482b730f8bSJeremy L Thompson PetscCall(PetscCalloc1(1, &blasius_ctx)); 24988626eedSJames Wright 25088626eedSJames Wright // ------------------------------------------------------ 25188626eedSJames Wright // SET UP Blasius 25288626eedSJames Wright // ------------------------------------------------------ 25391e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 25491e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 25588626eedSJames Wright 256fb455ff0SLeila Ghaffari CeedScalar U_inf = 40; // m/s 257fb455ff0SLeila Ghaffari CeedScalar T_inf = 288.; // K 25807d14e58SLeila Ghaffari CeedScalar T_wall = 288.; // K 2592518f336SLeila Ghaffari CeedScalar delta0 = 4.2e-3; // m 26088626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 2612518f336SLeila Ghaffari CeedInt N = 20; // Number of Chebyshev terms 262871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 2637b8d3891SJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 2647b8d3891SJames Wright PetscReal mesh_growth = 1.08; // [-] 2657b8d3891SJames Wright PetscInt mesh_Ndelta = 45; // [-] 2667b8d3891SJames Wright PetscReal mesh_top_angle = 5; // degrees 267d2714514SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 26888626eedSJames Wright 2697b37518fSJames Wright PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 2702b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL)); 2712b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL)); 2722b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL)); 2732b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL)); 2742b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL)); 2752b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL)); 2762b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL)); 2772b730f8bSJeremy 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, 27807d14e58SLeila Ghaffari BLASIUS_MAX_N_CHEBYSHEV); 2792b730f8bSJeremy L Thompson PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1)); 2802b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height, 2812b730f8bSJeremy L Thompson NULL)); 2822b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL)); 2832b730f8bSJeremy L Thompson PetscCall( 2842b730f8bSJeremy L Thompson PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL)); 2852b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-platemesh_y_node_locs_path", 286d2714514SJames Wright "Path to file with y node locations. " 2872b730f8bSJeremy L Thompson "If empty, will use the algorithmic mesh warping.", 2882b730f8bSJeremy L Thompson NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL)); 2892b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL)); 29088626eedSJames Wright PetscOptionsEnd(); 29188626eedSJames Wright 29288626eedSJames Wright PetscScalar meter = user->units->meter; 29388626eedSJames Wright PetscScalar second = user->units->second; 29488626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 29588626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 29688626eedSJames Wright 297fb455ff0SLeila Ghaffari T_inf *= Kelvin; 2982518f336SLeila Ghaffari T_wall *= Kelvin; 29988626eedSJames Wright P0 *= Pascal; 300fb455ff0SLeila Ghaffari U_inf *= meter / second; 30188626eedSJames Wright delta0 *= meter; 30288626eedSJames Wright 303d2714514SJames Wright PetscReal *mesh_ynodes = NULL; 304d2714514SJames Wright PetscInt mesh_nynodes = 0; 305d2714514SJames Wright if (strcmp(mesh_ynodes_path, "")) { 3062b730f8bSJeremy L Thompson PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes)); 307d2714514SJames Wright } 3082b730f8bSJeremy L Thompson PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes)); 30988626eedSJames Wright 310841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 3112b730f8bSJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 31288626eedSJames Wright 313ba6664aeSJames Wright blasius_ctx->weakT = weakT; 314fb455ff0SLeila Ghaffari blasius_ctx->U_inf = U_inf; 315fb455ff0SLeila Ghaffari blasius_ctx->T_inf = T_inf; 3162518f336SLeila Ghaffari blasius_ctx->T_wall = T_wall; 317841e4c73SJed Brown blasius_ctx->delta0 = delta0; 318841e4c73SJed Brown blasius_ctx->P0 = P0; 3192518f336SLeila Ghaffari blasius_ctx->n_cheb = N; 32030e9fa81SJames Wright newtonian_ig_ctx->P0 = P0; 321841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 322841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 323ba6664aeSJames Wright 324f1122ed0SJames Wright { 3252518f336SLeila Ghaffari PetscReal domain_min[3], domain_max[3]; 3262b730f8bSJeremy L Thompson PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 327f1122ed0SJames Wright blasius_ctx->x_inflow = domain_min[0]; 3282518f336SLeila Ghaffari blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 329f1122ed0SJames Wright } 3303ce1835eSJames Wright if (!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 331f1122ed0SJames Wright 3322b730f8bSJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 33388626eedSJames Wright 334841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 3352b730f8bSJeremy L Thompson CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx); 3362b730f8bSJeremy L Thompson CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc); 33788626eedSJames Wright 338b77c53c9SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 339841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 340ba6664aeSJames Wright if (use_stg) { 3412b730f8bSJeremy L Thompson PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes)); 34265dd5cafSJames Wright } else { 34365dd5cafSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 34465dd5cafSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 345fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 346fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 3472b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context); 3482b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context); 349ba6664aeSJames Wright } 3502b730f8bSJeremy L Thompson PetscCall(PetscFree(mesh_ynodes)); 35188626eedSJames Wright PetscFunctionReturn(0); 35288626eedSJames Wright } 353