xref: /honee/problems/blasius.c (revision 04e40bb60650195adcc92556a3eb81ec7887ccc8)
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"
122b916ea7SJeremy L Thompson 
132b916ea7SJeremy 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;
212b916ea7SJeremy 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);
40*04e40bb6SJeremy L Thompson     // mu and rho generally depend on h.
41*04e40bb6SJeremy L Thompson     // We naively assume constant mu.
420d850f2eSLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
430d850f2eSLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
440d850f2eSLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
450d850f2eSLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
460d850f2eSLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
470d850f2eSLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
480d850f2eSLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
490d850f2eSLeila Ghaffari     };
500d850f2eSLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
512b916ea7SJeremy 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]);
52e0d1a4dfSLeila Ghaffari   }
53e0d1a4dfSLeila Ghaffari 
54e0d1a4dfSLeila Ghaffari   // h - left end boundary condition
55e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
56aef1eb53SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->T_inf;
57e0d1a4dfSLeila Ghaffari 
58e0d1a4dfSLeila Ghaffari   // h - right end boundary condition
59e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
60e0d1a4dfSLeila Ghaffari   r[N + 1] = h[0] - 1.;
61e0d1a4dfSLeila Ghaffari 
62e0d1a4dfSLeila Ghaffari   // Restore vectors
63e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
64e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
65e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
66e0d1a4dfSLeila Ghaffari }
67e0d1a4dfSLeila Ghaffari 
68e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
69e0d1a4dfSLeila Ghaffari   SNES                snes;
70e0d1a4dfSLeila Ghaffari   Vec                 sol, res;
71e0d1a4dfSLeila Ghaffari   PetscReal          *w;
72e0d1a4dfSLeila Ghaffari   PetscInt            N = blasius->n_cheb;
730d850f2eSLeila Ghaffari   SNESConvergedReason reason;
74e0d1a4dfSLeila Ghaffari   const PetscScalar  *cheb_coefs;
75e0d1a4dfSLeila Ghaffari   PetscFunctionBegin;
760d850f2eSLeila Ghaffari 
770d850f2eSLeila Ghaffari   // Allocate memory
78e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
79e0d1a4dfSLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
800d850f2eSLeila Ghaffari 
810d850f2eSLeila Ghaffari   // Snes solve
82e0d1a4dfSLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
83e0d1a4dfSLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
84e0d1a4dfSLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
85e0d1a4dfSLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
860d850f2eSLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
870d850f2eSLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
88e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
89e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
900d850f2eSLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
91e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
92e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
930d850f2eSLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
942b916ea7SJeremy L Thompson   if (reason < 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
950d850f2eSLeila Ghaffari 
960d850f2eSLeila Ghaffari   // Assign Chebyshev coefficients
97e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
98e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
99e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
1000d850f2eSLeila Ghaffari 
1010d850f2eSLeila Ghaffari   // Destroy objects
102e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
103e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
104e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
105e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
106e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
107e0d1a4dfSLeila Ghaffari }
108e0d1a4dfSLeila Ghaffari 
1092b916ea7SJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
110f362fe62SJames Wright   PetscInt       ndims, dims[2];
111f362fe62SJames Wright   FILE          *fp;
112f362fe62SJames Wright   const PetscInt char_array_len = 512;
113f362fe62SJames Wright   char           line[char_array_len];
114f362fe62SJames Wright   char         **array;
115f362fe62SJames Wright   PetscReal     *node_locs;
116f362fe62SJames Wright   PetscFunctionBeginUser;
117f362fe62SJames Wright 
1182b916ea7SJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1192b916ea7SJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1202b916ea7SJeremy L Thompson   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
121f362fe62SJames Wright 
122f362fe62SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
123f362fe62SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
124f362fe62SJames Wright   *nynodes = dims[0];
1252b916ea7SJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
126f362fe62SJames Wright 
127f362fe62SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
1282b916ea7SJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1292b916ea7SJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1302b916ea7SJeremy L Thompson     if (ndims < dims[1])
1312b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1322b916ea7SJeremy L Thompson               ndims, dims[1]);
133f362fe62SJames Wright 
134f362fe62SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
135f362fe62SJames Wright   }
1362b916ea7SJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
137f362fe62SJames Wright   *pynodes = node_locs;
138f362fe62SJames Wright   PetscFunctionReturn(0);
139f362fe62SJames Wright }
140f362fe62SJames Wright 
141bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
142bb8a0c61SJames Wright  *
143*04e40bb6SJeremy L Thompson  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
144493642f1SJames Wright  * linearly in logspace to the top surface.
145bb8a0c61SJames Wright  *
146*04e40bb6SJeremy L Thompson  * The top surface is also angled downwards, so that it may be used as an outflow.
147*04e40bb6SJeremy L Thompson  * It's angle is controlled by `top_angle` (in units of degrees).
148f362fe62SJames Wright  *
149*04e40bb6SJeremy L Thompson  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
150*04e40bb6SJeremy L Thompson  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
151bb8a0c61SJames Wright  */
1522b916ea7SJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
1538e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
1542b916ea7SJeremy L Thompson   PetscInt     narr, ncoords;
155bb8a0c61SJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
156bb8a0c61SJames Wright   PetscScalar *arr_coords;
157bb8a0c61SJames Wright   Vec          vec_coords;
158bb8a0c61SJames Wright   PetscFunctionBeginUser;
159bb8a0c61SJames Wright 
160bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
161bb8a0c61SJames Wright 
162bb8a0c61SJames Wright   // Get domain boundary information
1632b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
164493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
165bb8a0c61SJames Wright 
166bb8a0c61SJames Wright   // Get coords array from DM
1672b916ea7SJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1682b916ea7SJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1692b916ea7SJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
170bb8a0c61SJames Wright 
171bb8a0c61SJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
172bb8a0c61SJames Wright   ncoords                   = narr / dim;
173bb8a0c61SJames Wright 
174bb8a0c61SJames Wright   // Get mesh information
175bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
1762b916ea7SJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
177f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
178f362fe62SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
179bb8a0c61SJames Wright 
1808e63b9cbSJames Wright   if (!*node_locs) {
181bb8a0c61SJames Wright     // Calculate the first element height
182bb8a0c61SJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
183bb8a0c61SJames Wright 
184bb8a0c61SJames Wright     // Calculate log of sizing outside BL
185bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
186bb8a0c61SJames Wright 
1878e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1888e63b9cbSJames Wright     PetscReal *temp_node_locs;
1892b916ea7SJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
1908e63b9cbSJames Wright 
191493642f1SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
192bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
193bb8a0c61SJames Wright       if (y_box_index <= N) {
1942b916ea7SJeremy L Thompson         coords[i][1] =
1952b916ea7SJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
196bb8a0c61SJames Wright       } else {
197bb8a0c61SJames Wright         PetscInt j   = y_box_index - N;
1982b916ea7SJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
199f362fe62SJames Wright       }
2002b916ea7SJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
201f362fe62SJames Wright     }
2028e63b9cbSJames Wright 
2038e63b9cbSJames Wright     *node_locs = temp_node_locs;
204f362fe62SJames Wright   } else {
205f362fe62SJames Wright     // Error checking
2062b916ea7SJeremy L Thompson     if (*num_node_locs < faces[1] + 1) {
2072b916ea7SJeremy L Thompson       SETERRQ(comm, -1,
2082b916ea7SJeremy L Thompson               "The y_node_locs_path has too few locations; "
209f362fe62SJames Wright               "There are %d + 1 nodes, but only %d locations given",
2108e63b9cbSJames Wright               faces[1] + 1, *num_node_locs);
2112b916ea7SJeremy L Thompson     }
2128e63b9cbSJames Wright     if (*num_node_locs > faces[1] + 1) {
2132b916ea7SJeremy L Thompson       PetscCall(PetscPrintf(comm,
2142b916ea7SJeremy L Thompson                             "WARNING: y_node_locs_path has more locations (%d) "
215e4677755SJames Wright                             "than the mesh has nodes (%d). This maybe unintended.\n",
2162b916ea7SJeremy L Thompson                             *num_node_locs, faces[1] + 1));
217f362fe62SJames Wright     }
2188e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
219f362fe62SJames Wright 
220f362fe62SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
221f362fe62SJames Wright       // Determine which y-node we're at
222f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2232b916ea7SJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
224bb8a0c61SJames Wright     }
225bb8a0c61SJames Wright   }
226bb8a0c61SJames Wright 
2272b916ea7SJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2282b916ea7SJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
229bb8a0c61SJames Wright 
230bb8a0c61SJames Wright   PetscFunctionReturn(0);
231bb8a0c61SJames Wright }
232bb8a0c61SJames Wright 
233d1c51a42SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
234bb8a0c61SJames Wright   User                     user    = *(User *)ctx;
235bb8a0c61SJames Wright   MPI_Comm                 comm    = PETSC_COMM_WORLD;
236493642f1SJames Wright   PetscBool                use_stg = PETSC_FALSE;
23715a3537eSJed Brown   BlasiusContext           blasius_ctx;
23815a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
23915a3537eSJed Brown   CeedQFunctionContext     blasius_context;
24015a3537eSJed Brown 
241bb8a0c61SJames Wright   PetscFunctionBeginUser;
242d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
2432b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &blasius_ctx));
244bb8a0c61SJames Wright 
245bb8a0c61SJames Wright   // ------------------------------------------------------
246bb8a0c61SJames Wright   //               SET UP Blasius
247bb8a0c61SJames Wright   // ------------------------------------------------------
2489785fe93SJed Brown   problem->ics.qfunction     = ICsBlasius;
2499785fe93SJed Brown   problem->ics.qfunction_loc = ICsBlasius_loc;
250bb8a0c61SJames Wright 
251aef1eb53SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
252aef1eb53SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
2530d850f2eSLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
254e0d1a4dfSLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
255bb8a0c61SJames Wright   CeedScalar P0                                   = 1.01e5;       // Pa
256e0d1a4dfSLeila Ghaffari   CeedInt    N                                    = 20;           // Number of Chebyshev terms
2572acc7cbcSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2587470235eSJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2597470235eSJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2607470235eSJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2617470235eSJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
262f362fe62SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
263bb8a0c61SJames Wright 
2643fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2652b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2662b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2672b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
2682b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2692b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2702b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
2712b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
2722b916ea7SJeremy 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,
2730d850f2eSLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2742b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
2752b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height,
2762b916ea7SJeremy L Thompson                                NULL));
2772b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
2782b916ea7SJeremy L Thompson   PetscCall(
2792b916ea7SJeremy L Thompson       PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
2802b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
281f362fe62SJames Wright                                "Path to file with y node locations. "
2822b916ea7SJeremy L Thompson                                "If empty, will use the algorithmic mesh warping.",
2832b916ea7SJeremy L Thompson                                NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
2842b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
285bb8a0c61SJames Wright   PetscOptionsEnd();
286bb8a0c61SJames Wright 
287bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
288bb8a0c61SJames Wright   PetscScalar second = user->units->second;
289bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
290bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
291bb8a0c61SJames Wright 
292aef1eb53SLeila Ghaffari   T_inf *= Kelvin;
293e0d1a4dfSLeila Ghaffari   T_wall *= Kelvin;
294bb8a0c61SJames Wright   P0 *= Pascal;
295aef1eb53SLeila Ghaffari   U_inf *= meter / second;
296bb8a0c61SJames Wright   delta0 *= meter;
297bb8a0c61SJames Wright 
298f362fe62SJames Wright   PetscReal *mesh_ynodes  = NULL;
299f362fe62SJames Wright   PetscInt   mesh_nynodes = 0;
300f362fe62SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
3012b916ea7SJeremy L Thompson     PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
302f362fe62SJames Wright   }
3032b916ea7SJeremy L Thompson   PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
304bb8a0c61SJames Wright 
30515a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
3062b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
307bb8a0c61SJames Wright 
308493642f1SJames Wright   blasius_ctx->weakT         = weakT;
309aef1eb53SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
310aef1eb53SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
311e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
31215a3537eSJed Brown   blasius_ctx->delta0        = delta0;
31315a3537eSJed Brown   blasius_ctx->P0            = P0;
314e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb        = N;
31504b9037bSJames Wright   newtonian_ig_ctx->P0       = P0;
31615a3537eSJed Brown   blasius_ctx->implicit      = user->phys->implicit;
31715a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
318493642f1SJames Wright 
319ef2c71fdSJames Wright   {
320e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3212b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
322ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
323e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
324ef2c71fdSJames Wright   }
3253a31d4d2SJames Wright   if (!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
326ef2c71fdSJames Wright 
3272b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
328bb8a0c61SJames Wright 
32915a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
3302b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx);
3312b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc);
332bb8a0c61SJames Wright 
33343bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
33415a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
335493642f1SJames Wright   if (use_stg) {
3362b916ea7SJeremy L Thompson     PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes));
3378085925cSJames Wright   } else {
3388085925cSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
3398085925cSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
3400aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
3410aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
3422b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context);
3432b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context);
344493642f1SJames Wright   }
3452b916ea7SJeremy L Thompson   PetscCall(PetscFree(mesh_ynodes));
346bb8a0c61SJames Wright   PetscFunctionReturn(0);
347bb8a0c61SJames Wright }
348