xref: /libCEED/examples/fluids/problems/blasius.c (revision 2249ac91e41bcb38055c29b6427350f9214a7158)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson #include <petscdt.h>
1649aac155SJeremy 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;
25ff9b3c0eSJames Wright   State                S_infty = blasius->S_infty;
26ff9b3c0eSJames Wright   CeedScalar           U_infty = sqrt(Dot3(S_infty.Y.velocity, S_infty.Y.velocity));
27f17d818dSJames Wright 
28f17d818dSJames Wright   PetscFunctionBeginUser;
29ff9b3c0eSJames Wright   PetscScalar Ma = Mach(&blasius->newtonian_ctx, S_infty.Y.temperature, U_infty), Pr = Prandtl(&blasius->newtonian_ctx),
302518f336SLeila Ghaffari               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
31f17d818dSJames Wright 
322518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
332518f336SLeila Ghaffari   Th = Tf + N;
342518f336SLeila Ghaffari   PetscCall(VecGetArray(R, &r));
352518f336SLeila Ghaffari 
362518f336SLeila Ghaffari   // Left boundary conditions f = f' = 0
372518f336SLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
382518f336SLeila Ghaffari   r[0] = f[0];
392518f336SLeila Ghaffari   r[1] = f[1];
402518f336SLeila Ghaffari 
412518f336SLeila Ghaffari   // f - right end boundary condition
422518f336SLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
432518f336SLeila Ghaffari   r[2] = f[1] - 1.;
442518f336SLeila Ghaffari 
452518f336SLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
462518f336SLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
472518f336SLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
48ea61e9acSJeremy L Thompson     // mu and rho generally depend on h.
49ea61e9acSJeremy L Thompson     // We naively assume constant mu.
5007d14e58SLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
5107d14e58SLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
5207d14e58SLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
5307d14e58SLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
5407d14e58SLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
5507d14e58SLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
5607d14e58SLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
5707d14e58SLeila Ghaffari     };
5807d14e58SLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
592b730f8bSJeremy 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]);
602518f336SLeila Ghaffari   }
612518f336SLeila Ghaffari 
622518f336SLeila Ghaffari   // h - left end boundary condition
632518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
64ff9b3c0eSJames Wright   r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature;
652518f336SLeila Ghaffari 
662518f336SLeila Ghaffari   // h - right end boundary condition
672518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
682518f336SLeila Ghaffari   r[N + 1] = h[0] - 1.;
692518f336SLeila Ghaffari 
702518f336SLeila Ghaffari   // Restore vectors
712518f336SLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
722518f336SLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
73ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
742518f336SLeila Ghaffari }
752518f336SLeila Ghaffari 
762518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
772518f336SLeila Ghaffari   SNES                snes;
782518f336SLeila Ghaffari   Vec                 sol, res;
792518f336SLeila Ghaffari   PetscReal          *w;
802518f336SLeila Ghaffari   PetscInt            N = blasius->n_cheb;
8107d14e58SLeila Ghaffari   SNESConvergedReason reason;
822518f336SLeila Ghaffari   const PetscScalar  *cheb_coefs;
8307d14e58SLeila Ghaffari 
84f17d818dSJames Wright   PetscFunctionBeginUser;
8507d14e58SLeila Ghaffari   // Allocate memory
862518f336SLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
872518f336SLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
8807d14e58SLeila Ghaffari 
8907d14e58SLeila Ghaffari   // Snes solve
902518f336SLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
912518f336SLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
922518f336SLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
932518f336SLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
9407d14e58SLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
9507d14e58SLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
962518f336SLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
972518f336SLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
9807d14e58SLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
992518f336SLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
1002518f336SLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
10107d14e58SLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
1020e654f56SJames Wright   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
10307d14e58SLeila Ghaffari 
10407d14e58SLeila Ghaffari   // Assign Chebyshev coefficients
1052518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
1062518f336SLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
1072518f336SLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
10807d14e58SLeila Ghaffari 
10907d14e58SLeila Ghaffari   // Destroy objects
1102518f336SLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
1112518f336SLeila Ghaffari   PetscCall(VecDestroy(&sol));
1122518f336SLeila Ghaffari   PetscCall(VecDestroy(&res));
1132518f336SLeila Ghaffari   PetscCall(SNESDestroy(&snes));
114ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1152518f336SLeila Ghaffari }
1162518f336SLeila Ghaffari 
1172b730f8bSJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
118457e73b2SJames Wright   int            ndims, dims[2];
119d2714514SJames Wright   FILE          *fp;
120d2714514SJames Wright   const PetscInt char_array_len = 512;
121d2714514SJames Wright   char           line[char_array_len];
122d2714514SJames Wright   char         **array;
123d2714514SJames Wright   PetscReal     *node_locs;
124d2714514SJames Wright 
125f17d818dSJames Wright   PetscFunctionBeginUser;
1262b730f8bSJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1272b730f8bSJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1282b730f8bSJeremy L Thompson   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
129d2714514SJames Wright 
130d2714514SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
131d2714514SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
132d2714514SJames Wright   *nynodes = dims[0];
1332b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
134d2714514SJames Wright 
135d2714514SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
1362b730f8bSJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1372b730f8bSJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1380e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
139457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]);
140d2714514SJames Wright 
141d2714514SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
142d2714514SJames Wright   }
1432b730f8bSJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
144d2714514SJames Wright   *pynodes = node_locs;
145ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
146d2714514SJames Wright }
147d2714514SJames Wright 
14888626eedSJames Wright /* \brief Modify the domain and mesh for blasius
14988626eedSJames Wright  *
150ea61e9acSJeremy L Thompson  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
151ba6664aeSJames Wright  * linearly in logspace to the top surface.
15288626eedSJames Wright  *
153ea61e9acSJeremy L Thompson  * The top surface is also angled downwards, so that it may be used as an outflow.
154ea61e9acSJeremy L Thompson  * It's angle is controlled by `top_angle` (in units of degrees).
155d2714514SJames Wright  *
156ea61e9acSJeremy L Thompson  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
157ea61e9acSJeremy L Thompson  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
15888626eedSJames Wright  */
1592b730f8bSJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
160798d42b8SJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
1612b730f8bSJeremy L Thompson   PetscInt     narr, ncoords;
16288626eedSJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
16388626eedSJames Wright   PetscScalar *arr_coords;
16488626eedSJames Wright   Vec          vec_coords;
165f17d818dSJames Wright 
16688626eedSJames Wright   PetscFunctionBeginUser;
16788626eedSJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
16888626eedSJames Wright   // Get domain boundary information
1692b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
170ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
17188626eedSJames Wright 
17288626eedSJames Wright   // Get coords array from DM
1732b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1742b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1752b730f8bSJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
17688626eedSJames Wright 
17788626eedSJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
17888626eedSJames Wright   ncoords                   = narr / dim;
17988626eedSJames Wright 
18088626eedSJames Wright   // Get mesh information
18188626eedSJames Wright   PetscInt nmax = 3, faces[3];
1822b730f8bSJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
183d2714514SJames Wright   // Get element size of the box mesh, for indexing each node
184d2714514SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
18588626eedSJames Wright 
186798d42b8SJames Wright   if (!*node_locs) {
18788626eedSJames Wright     // Calculate the first element height
18888626eedSJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
18988626eedSJames Wright 
19088626eedSJames Wright     // Calculate log of sizing outside BL
19188626eedSJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
19288626eedSJames Wright 
193798d42b8SJames Wright     *num_node_locs = faces[1] + 1;
194798d42b8SJames Wright     PetscReal *temp_node_locs;
1952b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
196798d42b8SJames Wright 
197ba6664aeSJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
19888626eedSJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
19988626eedSJames Wright       if (y_box_index <= N) {
2002b730f8bSJeremy L Thompson         coords[i][1] =
2012b730f8bSJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
20288626eedSJames Wright       } else {
20388626eedSJames Wright         PetscInt j   = y_box_index - N;
2042b730f8bSJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
205d2714514SJames Wright       }
2062b730f8bSJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
207d2714514SJames Wright     }
208798d42b8SJames Wright 
209798d42b8SJames Wright     *node_locs = temp_node_locs;
210d2714514SJames Wright   } else {
2110e654f56SJames Wright     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
212457e73b2SJames Wright                "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given",
213457e73b2SJames Wright                faces[1] + 1, *num_node_locs);
214798d42b8SJames Wright     if (*num_node_locs > faces[1] + 1) {
2152b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
216457e73b2SJames Wright                             "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") "
217457e73b2SJames Wright                             "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n",
2182b730f8bSJeremy L Thompson                             *num_node_locs, faces[1] + 1));
219d2714514SJames Wright     }
220798d42b8SJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
221d2714514SJames Wright 
222d2714514SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
223d2714514SJames Wright       // Determine which y-node we're at
224d2714514SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2252b730f8bSJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
22688626eedSJames Wright     }
22788626eedSJames Wright   }
22888626eedSJames Wright 
2292b730f8bSJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2302b730f8bSJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
231ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
23288626eedSJames Wright }
23388626eedSJames Wright 
234731c13d7SJames Wright PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
23588626eedSJames Wright   User                     user    = *(User *)ctx;
2360a37e1beSJames Wright   MPI_Comm                 comm    = user->comm;
237a424bcd0SJames Wright   Ceed                     ceed    = user->ceed;
238ba6664aeSJames Wright   PetscBool                use_stg = PETSC_FALSE;
239841e4c73SJed Brown   BlasiusContext           blasius_ctx;
240841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
241841e4c73SJed Brown   CeedQFunctionContext     blasius_context;
242841e4c73SJed Brown 
24388626eedSJames Wright   PetscFunctionBeginUser;
24446de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
2452b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &blasius_ctx));
24688626eedSJames Wright 
24788626eedSJames Wright   // ------------------------------------------------------
24888626eedSJames Wright   //               SET UP Blasius
24988626eedSJames Wright   // ------------------------------------------------------
25091e5af17SJed Brown   problem->ics.qfunction     = ICsBlasius;
25191e5af17SJed Brown   problem->ics.qfunction_loc = ICsBlasius_loc;
25288626eedSJames Wright 
253fb455ff0SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
254fb455ff0SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
25507d14e58SLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
2562518f336SLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
257ff9b3c0eSJames Wright   CeedScalar P_inf                                = 1.01e5;       // Pa
258457e73b2SJames Wright   PetscInt   N                                    = 20;           // Number of Chebyshev terms
259871db79fSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2607b8d3891SJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2617b8d3891SJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2627b8d3891SJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2637b8d3891SJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
264d2714514SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
265*2249ac91SJames Wright   PetscBool  P0_set;
26688626eedSJames Wright 
2677b37518fSJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2682b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2692b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2702b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
271*2249ac91SJames Wright   PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set));  // For maintaining behavior of -P0 flag (which is deprecated)
272*2249ac91SJames Wright   PetscCall(
273*2249ac91SJames Wright       PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0",
274*2249ac91SJames Wright                              "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure"));
275*2249ac91SJames Wright   PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL));
2762b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2772b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2782b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
2792b730f8bSJeremy 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,
28007d14e58SLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2812526956eSJames Wright   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
2822b730f8bSJeremy L Thompson     PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
2832526956eSJames Wright     PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height,
2842526956eSJames Wright                                  &mesh_refine_height, NULL));
2852b730f8bSJeremy L Thompson     PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
2862b730f8bSJeremy L Thompson     PetscCall(
2872b730f8bSJeremy L Thompson         PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
2882b730f8bSJeremy L Thompson     PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
289d2714514SJames Wright                                  "Path to file with y node locations. "
2902b730f8bSJeremy L Thompson                                  "If empty, will use the algorithmic mesh warping.",
2912b730f8bSJeremy L Thompson                                  NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
2922526956eSJames Wright   }
2932b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
29488626eedSJames Wright   PetscOptionsEnd();
29588626eedSJames Wright 
29688626eedSJames Wright   PetscScalar meter  = user->units->meter;
29788626eedSJames Wright   PetscScalar second = user->units->second;
29888626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
29988626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
30088626eedSJames Wright 
301fb455ff0SLeila Ghaffari   T_inf *= Kelvin;
3022518f336SLeila Ghaffari   T_wall *= Kelvin;
303ff9b3c0eSJames Wright   P_inf *= Pascal;
304fb455ff0SLeila Ghaffari   U_inf *= meter / second;
30588626eedSJames Wright   delta0 *= meter;
30688626eedSJames Wright 
3072526956eSJames Wright   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
308d2714514SJames Wright     PetscReal *mesh_ynodes  = NULL;
309d2714514SJames Wright     PetscInt   mesh_nynodes = 0;
3109309e21cSJames Wright     if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
3112b730f8bSJeremy L Thompson     PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
312f8839eb4SJames Wright     PetscCall(PetscFree(mesh_ynodes));
3139309e21cSJames Wright   }
31488626eedSJames Wright 
315841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
316a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
31788626eedSJames Wright 
318ff9b3c0eSJames Wright   StatePrimitive Y_inf = {
319ff9b3c0eSJames Wright       .pressure = P_inf, .velocity = {U_inf, 0, 0},
320ff9b3c0eSJames Wright            .temperature = T_inf
321ff9b3c0eSJames Wright   };
322ff9b3c0eSJames Wright   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
323ff9b3c0eSJames Wright 
324ba6664aeSJames Wright   blasius_ctx->weakT    = weakT;
3252518f336SLeila Ghaffari   blasius_ctx->T_wall   = T_wall;
326841e4c73SJed Brown   blasius_ctx->delta0   = delta0;
327ff9b3c0eSJames Wright   blasius_ctx->S_infty  = S_infty;
3282518f336SLeila Ghaffari   blasius_ctx->n_cheb   = N;
329841e4c73SJed Brown   blasius_ctx->implicit = user->phys->implicit;
330*2249ac91SJames Wright   if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf;  // For maintaining behavior of -P0 flag (which is deprecated)
331841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
332ba6664aeSJames Wright 
333f1122ed0SJames Wright   {
3342518f336SLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3352b730f8bSJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
336f1122ed0SJames Wright     blasius_ctx->x_inflow = domain_min[0];
3372518f336SLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
338f1122ed0SJames Wright   }
3392029a677SJames Wright   PetscBool diff_filter_mms = PETSC_FALSE;
3402029a677SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
3412029a677SJames Wright   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
342f1122ed0SJames Wright 
343a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
34488626eedSJames Wright 
345a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &blasius_context));
346a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
347a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc));
34888626eedSJames Wright 
349a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
350841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
351ba6664aeSJames Wright   if (use_stg) {
352ff9b3c0eSJames Wright     PetscCall(SetupStg(comm, dm, problem, user, weakT, S_infty.Y.temperature, S_infty.Y.pressure));
3532029a677SJames Wright   } else if (diff_filter_mms) {
354cbef7084SJames Wright     PetscCall(DifferentialFilterMmsICSetup(problem));
35565dd5cafSJames Wright   } else {
35665dd5cafSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
35765dd5cafSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
358fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
359fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
360a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context));
361a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context));
362ba6664aeSJames Wright   }
363ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
36488626eedSJames Wright }
365