xref: /honee/problems/blasius.c (revision 5d28dccaccae4dbbdfc8aa7c8439b84e1bbb591b)
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 
13e419654dSJeremy L Thompson #include <ceed.h>
14e419654dSJeremy L Thompson #include <petscdm.h>
15e419654dSJeremy L Thompson #include <petscdt.h>
16e419654dSJeremy L Thompson 
172b916ea7SJeremy L Thompson #include "../navierstokes.h"
18493642f1SJames Wright #include "stg_shur14.h"
19bb8a0c61SJames Wright 
20e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
21e0d1a4dfSLeila Ghaffari   const BlasiusContext blasius = (BlasiusContext)ctx;
22e0d1a4dfSLeila Ghaffari   const PetscScalar   *Tf, *Th;  // Chebyshev coefficients
23e0d1a4dfSLeila Ghaffari   PetscScalar         *r, f[4], h[4];
24e0d1a4dfSLeila Ghaffari   PetscInt             N  = blasius->n_cheb;
252b916ea7SJeremy L Thompson   PetscScalar          Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx),
26e0d1a4dfSLeila Ghaffari               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
27e0d1a4dfSLeila Ghaffari   PetscFunctionBegin;
28e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
29e0d1a4dfSLeila Ghaffari   Th = Tf + N;
30e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArray(R, &r));
31e0d1a4dfSLeila Ghaffari 
32e0d1a4dfSLeila Ghaffari   // Left boundary conditions f = f' = 0
33e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
34e0d1a4dfSLeila Ghaffari   r[0] = f[0];
35e0d1a4dfSLeila Ghaffari   r[1] = f[1];
36e0d1a4dfSLeila Ghaffari 
37e0d1a4dfSLeila Ghaffari   // f - right end boundary condition
38e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
39e0d1a4dfSLeila Ghaffari   r[2] = f[1] - 1.;
40e0d1a4dfSLeila Ghaffari 
41e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
42e0d1a4dfSLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
43e0d1a4dfSLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
4404e40bb6SJeremy L Thompson     // mu and rho generally depend on h.
4504e40bb6SJeremy L Thompson     // We naively assume constant mu.
460d850f2eSLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
470d850f2eSLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
480d850f2eSLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
490d850f2eSLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
500d850f2eSLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
510d850f2eSLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
520d850f2eSLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
530d850f2eSLeila Ghaffari     };
540d850f2eSLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
552b916ea7SJeremy 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]);
56e0d1a4dfSLeila Ghaffari   }
57e0d1a4dfSLeila Ghaffari 
58e0d1a4dfSLeila Ghaffari   // h - left end boundary condition
59e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
60aef1eb53SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->T_inf;
61e0d1a4dfSLeila Ghaffari 
62e0d1a4dfSLeila Ghaffari   // h - right end boundary condition
63e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
64e0d1a4dfSLeila Ghaffari   r[N + 1] = h[0] - 1.;
65e0d1a4dfSLeila Ghaffari 
66e0d1a4dfSLeila Ghaffari   // Restore vectors
67e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
68e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
69e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
70e0d1a4dfSLeila Ghaffari }
71e0d1a4dfSLeila Ghaffari 
72e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
73e0d1a4dfSLeila Ghaffari   SNES                snes;
74e0d1a4dfSLeila Ghaffari   Vec                 sol, res;
75e0d1a4dfSLeila Ghaffari   PetscReal          *w;
76e0d1a4dfSLeila Ghaffari   PetscInt            N = blasius->n_cheb;
770d850f2eSLeila Ghaffari   SNESConvergedReason reason;
78e0d1a4dfSLeila Ghaffari   const PetscScalar  *cheb_coefs;
79e0d1a4dfSLeila Ghaffari   PetscFunctionBegin;
800d850f2eSLeila Ghaffari 
810d850f2eSLeila Ghaffari   // Allocate memory
82e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
83e0d1a4dfSLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
840d850f2eSLeila Ghaffari 
850d850f2eSLeila Ghaffari   // Snes solve
86e0d1a4dfSLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
87e0d1a4dfSLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
88e0d1a4dfSLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
89e0d1a4dfSLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
900d850f2eSLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
910d850f2eSLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
92e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
93e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
940d850f2eSLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
95e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
96e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
970d850f2eSLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
98*5d28dccaSJames Wright   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
990d850f2eSLeila Ghaffari 
1000d850f2eSLeila Ghaffari   // Assign Chebyshev coefficients
101e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
102e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
103e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
1040d850f2eSLeila Ghaffari 
1050d850f2eSLeila Ghaffari   // Destroy objects
106e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
107e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
108e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
109e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
110e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
111e0d1a4dfSLeila Ghaffari }
112e0d1a4dfSLeila Ghaffari 
1132b916ea7SJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
114f362fe62SJames Wright   PetscInt       ndims, dims[2];
115f362fe62SJames Wright   FILE          *fp;
116f362fe62SJames Wright   const PetscInt char_array_len = 512;
117f362fe62SJames Wright   char           line[char_array_len];
118f362fe62SJames Wright   char         **array;
119f362fe62SJames Wright   PetscReal     *node_locs;
120f362fe62SJames Wright   PetscFunctionBeginUser;
121f362fe62SJames Wright 
1222b916ea7SJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1232b916ea7SJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1242b916ea7SJeremy L Thompson   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
125f362fe62SJames Wright 
126f362fe62SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
127f362fe62SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
128f362fe62SJames Wright   *nynodes = dims[0];
1292b916ea7SJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
130f362fe62SJames Wright 
131f362fe62SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
1322b916ea7SJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1332b916ea7SJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
134*5d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
135*5d28dccaSJames Wright                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1362b916ea7SJeremy L Thompson                ndims, dims[1]);
137f362fe62SJames Wright 
138f362fe62SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
139f362fe62SJames Wright   }
1402b916ea7SJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
141f362fe62SJames Wright   *pynodes = node_locs;
142f362fe62SJames Wright   PetscFunctionReturn(0);
143f362fe62SJames Wright }
144f362fe62SJames Wright 
145bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
146bb8a0c61SJames Wright  *
14704e40bb6SJeremy L Thompson  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
148493642f1SJames Wright  * linearly in logspace to the top surface.
149bb8a0c61SJames Wright  *
15004e40bb6SJeremy L Thompson  * The top surface is also angled downwards, so that it may be used as an outflow.
15104e40bb6SJeremy L Thompson  * It's angle is controlled by `top_angle` (in units of degrees).
152f362fe62SJames Wright  *
15304e40bb6SJeremy L Thompson  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
15404e40bb6SJeremy L Thompson  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
155bb8a0c61SJames Wright  */
1562b916ea7SJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
1578e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
1582b916ea7SJeremy L Thompson   PetscInt     narr, ncoords;
159bb8a0c61SJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
160bb8a0c61SJames Wright   PetscScalar *arr_coords;
161bb8a0c61SJames Wright   Vec          vec_coords;
162bb8a0c61SJames Wright   PetscFunctionBeginUser;
163bb8a0c61SJames Wright 
164bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
165bb8a0c61SJames Wright 
166bb8a0c61SJames Wright   // Get domain boundary information
1672b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
168493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
169bb8a0c61SJames Wright 
170bb8a0c61SJames Wright   // Get coords array from DM
1712b916ea7SJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1722b916ea7SJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1732b916ea7SJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
174bb8a0c61SJames Wright 
175bb8a0c61SJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
176bb8a0c61SJames Wright   ncoords                   = narr / dim;
177bb8a0c61SJames Wright 
178bb8a0c61SJames Wright   // Get mesh information
179bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
1802b916ea7SJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
181f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
182f362fe62SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
183bb8a0c61SJames Wright 
1848e63b9cbSJames Wright   if (!*node_locs) {
185bb8a0c61SJames Wright     // Calculate the first element height
186bb8a0c61SJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
187bb8a0c61SJames Wright 
188bb8a0c61SJames Wright     // Calculate log of sizing outside BL
189bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
190bb8a0c61SJames Wright 
1918e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1928e63b9cbSJames Wright     PetscReal *temp_node_locs;
1932b916ea7SJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
1948e63b9cbSJames Wright 
195493642f1SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
196bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
197bb8a0c61SJames Wright       if (y_box_index <= N) {
1982b916ea7SJeremy L Thompson         coords[i][1] =
1992b916ea7SJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
200bb8a0c61SJames Wright       } else {
201bb8a0c61SJames Wright         PetscInt j   = y_box_index - N;
2022b916ea7SJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
203f362fe62SJames Wright       }
2042b916ea7SJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
205f362fe62SJames Wright     }
2068e63b9cbSJames Wright 
2078e63b9cbSJames Wright     *node_locs = temp_node_locs;
208f362fe62SJames Wright   } else {
209*5d28dccaSJames Wright     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
210*5d28dccaSJames Wright                "The y_node_locs_path has too few locations; There are %d + 1 nodes, but only %d locations given", faces[1] + 1, *num_node_locs);
2118e63b9cbSJames Wright     if (*num_node_locs > faces[1] + 1) {
2122b916ea7SJeremy L Thompson       PetscCall(PetscPrintf(comm,
2132b916ea7SJeremy L Thompson                             "WARNING: y_node_locs_path has more locations (%d) "
214e4677755SJames Wright                             "than the mesh has nodes (%d). This maybe unintended.\n",
2152b916ea7SJeremy L Thompson                             *num_node_locs, faces[1] + 1));
216f362fe62SJames Wright     }
2178e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
218f362fe62SJames Wright 
219f362fe62SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
220f362fe62SJames Wright       // Determine which y-node we're at
221f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2222b916ea7SJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
223bb8a0c61SJames Wright     }
224bb8a0c61SJames Wright   }
225bb8a0c61SJames Wright 
2262b916ea7SJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2272b916ea7SJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
228bb8a0c61SJames Wright 
229bb8a0c61SJames Wright   PetscFunctionReturn(0);
230bb8a0c61SJames Wright }
231bb8a0c61SJames Wright 
232d1c51a42SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
233bb8a0c61SJames Wright   User                     user    = *(User *)ctx;
234bb8a0c61SJames Wright   MPI_Comm                 comm    = PETSC_COMM_WORLD;
235493642f1SJames Wright   PetscBool                use_stg = PETSC_FALSE;
23615a3537eSJed Brown   BlasiusContext           blasius_ctx;
23715a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
23815a3537eSJed Brown   CeedQFunctionContext     blasius_context;
23915a3537eSJed Brown 
240bb8a0c61SJames Wright   PetscFunctionBeginUser;
241d1c51a42SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
2422b916ea7SJeremy L Thompson   PetscCall(PetscCalloc1(1, &blasius_ctx));
243bb8a0c61SJames Wright 
244bb8a0c61SJames Wright   // ------------------------------------------------------
245bb8a0c61SJames Wright   //               SET UP Blasius
246bb8a0c61SJames Wright   // ------------------------------------------------------
2479785fe93SJed Brown   problem->ics.qfunction     = ICsBlasius;
2489785fe93SJed Brown   problem->ics.qfunction_loc = ICsBlasius_loc;
249bb8a0c61SJames Wright 
250aef1eb53SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
251aef1eb53SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
2520d850f2eSLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
253e0d1a4dfSLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
254bb8a0c61SJames Wright   CeedScalar P0                                   = 1.01e5;       // Pa
255e0d1a4dfSLeila Ghaffari   CeedInt    N                                    = 20;           // Number of Chebyshev terms
2562acc7cbcSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2577470235eSJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2587470235eSJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2597470235eSJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2607470235eSJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
261f362fe62SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
262bb8a0c61SJames Wright 
2633fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2642b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2652b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2662b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
2672b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2682b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2692b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
2702b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
2712b916ea7SJeremy 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,
2720d850f2eSLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2732b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
2742b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height,
2752b916ea7SJeremy L Thompson                                NULL));
2762b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
2772b916ea7SJeremy L Thompson   PetscCall(
2782b916ea7SJeremy L Thompson       PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
2792b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
280f362fe62SJames Wright                                "Path to file with y node locations. "
2812b916ea7SJeremy L Thompson                                "If empty, will use the algorithmic mesh warping.",
2822b916ea7SJeremy L Thompson                                NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
2832b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
284bb8a0c61SJames Wright   PetscOptionsEnd();
285bb8a0c61SJames Wright 
286bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
287bb8a0c61SJames Wright   PetscScalar second = user->units->second;
288bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
289bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
290bb8a0c61SJames Wright 
291aef1eb53SLeila Ghaffari   T_inf *= Kelvin;
292e0d1a4dfSLeila Ghaffari   T_wall *= Kelvin;
293bb8a0c61SJames Wright   P0 *= Pascal;
294aef1eb53SLeila Ghaffari   U_inf *= meter / second;
295bb8a0c61SJames Wright   delta0 *= meter;
296bb8a0c61SJames Wright 
297f362fe62SJames Wright   PetscReal *mesh_ynodes  = NULL;
298f362fe62SJames Wright   PetscInt   mesh_nynodes = 0;
299f362fe62SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
3002b916ea7SJeremy L Thompson     PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
301f362fe62SJames Wright   }
3022b916ea7SJeremy L Thompson   PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
303bb8a0c61SJames Wright 
30415a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
3052b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
306bb8a0c61SJames Wright 
307493642f1SJames Wright   blasius_ctx->weakT         = weakT;
308aef1eb53SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
309aef1eb53SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
310e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
31115a3537eSJed Brown   blasius_ctx->delta0        = delta0;
31215a3537eSJed Brown   blasius_ctx->P0            = P0;
313e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb        = N;
31404b9037bSJames Wright   newtonian_ig_ctx->P0       = P0;
31515a3537eSJed Brown   blasius_ctx->implicit      = user->phys->implicit;
31615a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
317493642f1SJames Wright 
318ef2c71fdSJames Wright   {
319e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3202b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
321ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
322e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
323ef2c71fdSJames Wright   }
3243a31d4d2SJames Wright   if (!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
325ef2c71fdSJames Wright 
3262b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
327bb8a0c61SJames Wright 
32815a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
3292b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx);
3302b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc);
331bb8a0c61SJames Wright 
33243bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
33315a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
334493642f1SJames Wright   if (use_stg) {
3352b916ea7SJeremy L Thompson     PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes));
3368085925cSJames Wright   } else {
3378085925cSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
3388085925cSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
3390aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
3400aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
3412b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context);
3422b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context);
343493642f1SJames Wright   }
3442b916ea7SJeremy L Thompson   PetscCall(PetscFree(mesh_ynodes));
345bb8a0c61SJames Wright   PetscFunctionReturn(0);
346bb8a0c61SJames Wright }
347