xref: /honee/problems/blasius.c (revision 0d850f2ecb3e7cca04d509aa253ec815011355d0)
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 "../navierstokes.h"
12bb8a0c61SJames Wright #include "../qfunctions/blasius.h"
13493642f1SJames Wright #include "stg_shur14.h"
14bb8a0c61SJames Wright 
15e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
16e0d1a4dfSLeila Ghaffari   const BlasiusContext blasius = (BlasiusContext)ctx;
17e0d1a4dfSLeila Ghaffari   const PetscScalar *Tf, *Th;  // Chebyshev coefficients
18e0d1a4dfSLeila Ghaffari   PetscScalar       *r, f[4], h[4];
19e0d1a4dfSLeila Ghaffari   PetscInt          N = blasius->n_cheb;
20aef1eb53SLeila Ghaffari   PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf),
21e0d1a4dfSLeila Ghaffari               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*0d850f2eSLeila Ghaffari     // mu and rho generally depend on h. We naively assume constant mu.
41*0d850f2eSLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
42*0d850f2eSLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
43*0d850f2eSLeila Ghaffari     const PetscScalar mu_tilde[2] = {1, 0};
44*0d850f2eSLeila Ghaffari     const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])};
45*0d850f2eSLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
46*0d850f2eSLeila Ghaffari       mu_tilde[0] *rho_tilde[0],
47*0d850f2eSLeila Ghaffari       mu_tilde[1] *rho_tilde[0] + mu_tilde[0] *rho_tilde[1],
48*0d850f2eSLeila Ghaffari     };
49*0d850f2eSLeila Ghaffari     r[3+i] = 2*(mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
50*0d850f2eSLeila Ghaffari     r[N+2+i] = (mu_rho_tilde[0] * h[2] + mu_rho_tilde[1] * h[1]) + Pr * f[0] * h[1]
51*0d850f2eSLeila Ghaffari                + 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;
73*0d850f2eSLeila Ghaffari   SNESConvergedReason reason;
74e0d1a4dfSLeila Ghaffari   const PetscScalar *cheb_coefs;
75e0d1a4dfSLeila Ghaffari   PetscFunctionBegin;
76*0d850f2eSLeila Ghaffari 
77*0d850f2eSLeila 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));
80*0d850f2eSLeila Ghaffari 
81*0d850f2eSLeila 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));
86*0d850f2eSLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
87*0d850f2eSLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
88e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
89e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
90*0d850f2eSLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
91e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
92e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
93*0d850f2eSLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
94*0d850f2eSLeila Ghaffari   if (reason < 0)
95*0d850f2eSLeila Ghaffari     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED,
96*0d850f2eSLeila Ghaffari             "The Chebyshev solve failed.\n");
97*0d850f2eSLeila Ghaffari 
98*0d850f2eSLeila Ghaffari   // Assign Chebyshev coefficients
99e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
100e0d1a4dfSLeila Ghaffari   for (int i=0; i<N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
101e0d1a4dfSLeila Ghaffari   for (int i=0; i<N-1; i++) blasius->Th_cheb[i] = cheb_coefs[i+N];
102*0d850f2eSLeila Ghaffari 
103*0d850f2eSLeila Ghaffari   // Destroy objects
104e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
105e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
106e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
107e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
108e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
109e0d1a4dfSLeila Ghaffari }
110e0d1a4dfSLeila Ghaffari 
111f362fe62SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
112f362fe62SJames Wright                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
113f362fe62SJames Wright                                    PetscInt *nynodes) {
114f362fe62SJames Wright   PetscErrorCode ierr;
115f362fe62SJames Wright   PetscInt ndims, dims[2];
116f362fe62SJames Wright   FILE *fp;
117f362fe62SJames Wright   const PetscInt char_array_len = 512;
118f362fe62SJames Wright   char line[char_array_len];
119f362fe62SJames Wright   char **array;
120f362fe62SJames Wright   PetscReal *node_locs;
121f362fe62SJames Wright   PetscFunctionBeginUser;
122f362fe62SJames Wright 
123f362fe62SJames Wright   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
124f362fe62SJames Wright   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
125f362fe62SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
126f362fe62SJames Wright 
127f362fe62SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
128f362fe62SJames Wright   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
129f362fe62SJames Wright   *nynodes = dims[0];
130f362fe62SJames Wright   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
131f362fe62SJames Wright 
132f362fe62SJames Wright   for (PetscInt i=0; i<dims[0]; i++) {
133f362fe62SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
134f362fe62SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
135f362fe62SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
13645aa3cadSJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
13745aa3cadSJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT ")", i,
138f362fe62SJames Wright                                    path, ndims, dims[1]);
139f362fe62SJames Wright 
140f362fe62SJames Wright     node_locs[i] = (PetscReal) atof(array[0]);
141f362fe62SJames Wright   }
142f362fe62SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
143f362fe62SJames Wright   *pynodes = node_locs;
144f362fe62SJames Wright   PetscFunctionReturn(0);
145f362fe62SJames Wright }
146f362fe62SJames Wright 
147bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
148bb8a0c61SJames Wright  *
149493642f1SJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
150493642f1SJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
151493642f1SJames Wright  * linearly in logspace to the top surface.
152bb8a0c61SJames Wright  *
153bb8a0c61SJames Wright  * The top surface is also angled downwards, so that it may be used as an
154493642f1SJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
155f362fe62SJames Wright  *
156f362fe62SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
1578e63b9cbSJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
1588e63b9cbSJames Wright  * the array, along with `num_node_locs`.
159bb8a0c61SJames Wright  */
160f362fe62SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
161f362fe62SJames Wright                                  PetscReal growth, PetscInt N,
162f362fe62SJames Wright                                  PetscReal refine_height, PetscReal top_angle,
1638e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
164bb8a0c61SJames Wright   PetscInt ierr, narr, ncoords;
165bb8a0c61SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
166bb8a0c61SJames Wright   PetscScalar *arr_coords;
167bb8a0c61SJames Wright   Vec vec_coords;
168bb8a0c61SJames Wright   PetscFunctionBeginUser;
169bb8a0c61SJames Wright 
170bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
171bb8a0c61SJames Wright 
172bb8a0c61SJames Wright   // Get domain boundary information
173bb8a0c61SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
174493642f1SJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
175bb8a0c61SJames Wright 
176bb8a0c61SJames Wright   // Get coords array from DM
177bb8a0c61SJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
178bb8a0c61SJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
179bb8a0c61SJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
180bb8a0c61SJames Wright 
181bb8a0c61SJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
182bb8a0c61SJames Wright   ncoords = narr/dim;
183bb8a0c61SJames Wright 
184bb8a0c61SJames Wright   // Get mesh information
185bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
186bb8a0c61SJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
187bb8a0c61SJames Wright                                  NULL); CHKERRQ(ierr);
188f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
189f362fe62SJames Wright   const PetscReal dybox = domain_size[1]/faces[1];
190bb8a0c61SJames Wright 
1918e63b9cbSJames Wright   if (!*node_locs) {
192bb8a0c61SJames Wright     // Calculate the first element height
193bb8a0c61SJames Wright     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
194bb8a0c61SJames Wright 
195bb8a0c61SJames Wright     // Calculate log of sizing outside BL
196bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
197bb8a0c61SJames Wright 
1988e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1998e63b9cbSJames Wright     PetscReal *temp_node_locs;
2008e63b9cbSJames Wright     ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr);
2018e63b9cbSJames Wright 
202493642f1SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
203bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
204bb8a0c61SJames Wright       if (y_box_index <= N) {
205029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
206f362fe62SJames Wright                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
207bb8a0c61SJames Wright       } else {
208bb8a0c61SJames Wright         PetscInt j = y_box_index - N;
209029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
210f362fe62SJames Wright                        * exp(log(refine_height) + logdy*j);
211f362fe62SJames Wright       }
2128e63b9cbSJames Wright       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2])
2138e63b9cbSJames Wright         temp_node_locs[y_box_index] = coords[i][1];
214f362fe62SJames Wright     }
2158e63b9cbSJames Wright 
2168e63b9cbSJames Wright     *node_locs = temp_node_locs;
217f362fe62SJames Wright   } else {
218f362fe62SJames Wright     // Error checking
2198e63b9cbSJames Wright     if (*num_node_locs < faces[1] +1)
220f362fe62SJames Wright       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
221f362fe62SJames Wright               "There are %d + 1 nodes, but only %d locations given",
2228e63b9cbSJames Wright               faces[1]+1, *num_node_locs);
2238e63b9cbSJames Wright     if (*num_node_locs > faces[1] +1) {
224f362fe62SJames Wright       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
225e4677755SJames Wright                          "than the mesh has nodes (%d). This maybe unintended.\n",
2268e63b9cbSJames Wright                          *num_node_locs, faces[1]+1); CHKERRQ(ierr);
227f362fe62SJames Wright     }
2288e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
229f362fe62SJames Wright 
230f362fe62SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
231f362fe62SJames Wright       // Determine which y-node we're at
232f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
233029c6dfaSJames Wright       coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y)
2348e63b9cbSJames Wright                      * (*node_locs)[y_box_index];
235bb8a0c61SJames Wright     }
236bb8a0c61SJames Wright   }
237bb8a0c61SJames Wright 
238bb8a0c61SJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
239bb8a0c61SJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
240bb8a0c61SJames Wright 
241bb8a0c61SJames Wright   PetscFunctionReturn(0);
242bb8a0c61SJames Wright }
243bb8a0c61SJames Wright 
244b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
245bb8a0c61SJames Wright 
246bb8a0c61SJames Wright   PetscInt ierr;
247bb8a0c61SJames Wright   User      user    = *(User *)ctx;
248bb8a0c61SJames Wright   MPI_Comm  comm    = PETSC_COMM_WORLD;
249493642f1SJames Wright   PetscBool use_stg = PETSC_FALSE;
25015a3537eSJed Brown   BlasiusContext blasius_ctx;
25115a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
25215a3537eSJed Brown   CeedQFunctionContext blasius_context;
25315a3537eSJed Brown 
254bb8a0c61SJames Wright   PetscFunctionBeginUser;
255b7f03f12SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
25615a3537eSJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
257bb8a0c61SJames Wright 
258bb8a0c61SJames Wright   // ------------------------------------------------------
259bb8a0c61SJames Wright   //               SET UP Blasius
260bb8a0c61SJames Wright   // ------------------------------------------------------
2619785fe93SJed Brown   problem->ics.qfunction        = ICsBlasius;
2629785fe93SJed Brown   problem->ics.qfunction_loc    = ICsBlasius_loc;
263bb8a0c61SJames Wright 
264aef1eb53SLeila Ghaffari   CeedScalar U_inf              = 40;          // m/s
265aef1eb53SLeila Ghaffari   CeedScalar T_inf              = 288.;        // K
266*0d850f2eSLeila Ghaffari   CeedScalar T_wall             = 288.;        // K
267e0d1a4dfSLeila Ghaffari   CeedScalar delta0             = 4.2e-3;      // m
268bb8a0c61SJames Wright   CeedScalar P0                 = 1.01e5;      // Pa
269e0d1a4dfSLeila Ghaffari   CeedInt    N                  = 20;          // Number of Chebyshev terms
2702acc7cbcSKenneth E. Jansen   PetscBool  weakT              = PETSC_FALSE; // weak density or temperature
2717470235eSJames Wright   PetscReal  mesh_refine_height = 5.9e-4;      // m
2727470235eSJames Wright   PetscReal  mesh_growth        = 1.08;        // [-]
2737470235eSJames Wright   PetscInt   mesh_Ndelta        = 45;          // [-]
2747470235eSJames Wright   PetscReal  mesh_top_angle     = 5;           // degrees
275f362fe62SJames Wright   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
276bb8a0c61SJames Wright 
2773fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2782acc7cbcSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
2792acc7cbcSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
280aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-velocity_infinity",
281aef1eb53SLeila Ghaffari                             "Velocity at boundary layer edge",
282aef1eb53SLeila Ghaffari                             NULL, U_inf, &U_inf, NULL); CHKERRQ(ierr);
283aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-temperature_infinity",
284aef1eb53SLeila Ghaffari                             "Temperature at boundary layer edge",
285aef1eb53SLeila Ghaffari                             NULL, T_inf, &T_inf, NULL); CHKERRQ(ierr);
286aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-temperature_wall", "Temperature at wall",
287e0d1a4dfSLeila Ghaffari                             NULL, T_wall, &T_wall, NULL); CHKERRQ(ierr);
288bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
289bb8a0c61SJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
290bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
291bb8a0c61SJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
292*0d850f2eSLeila Ghaffari   ierr = PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms",
293e0d1a4dfSLeila Ghaffari                          NULL, N, &N, NULL); CHKERRQ(ierr);
294*0d850f2eSLeila Ghaffari   PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV,
295*0d850f2eSLeila Ghaffari              comm, PETSC_ERR_ARG_OUTOFRANGE,
296*0d850f2eSLeila Ghaffari              "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N,
297*0d850f2eSLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2987470235eSJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
2997470235eSJames Wright                                 "Velocity at boundary layer edge",
3007470235eSJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
3017470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
302bb8a0c61SJames Wright                             "Height of boundary layer mesh refinement",
3037470235eSJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
3047470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
305bb8a0c61SJames Wright                             "Geometric growth rate of boundary layer mesh",
3067470235eSJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
3077470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
308bb8a0c61SJames Wright                             "Geometric top_angle rate of boundary layer mesh",
3097470235eSJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
310f362fe62SJames Wright   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
311f362fe62SJames Wright                             "Path to file with y node locations. "
312f362fe62SJames Wright                             "If empty, will use the algorithmic mesh warping.", NULL,
313f362fe62SJames Wright                             mesh_ynodes_path, mesh_ynodes_path,
314f362fe62SJames Wright                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
315493642f1SJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
316493642f1SJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
317bb8a0c61SJames Wright   PetscOptionsEnd();
318bb8a0c61SJames Wright 
319bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
320bb8a0c61SJames Wright   PetscScalar second = user->units->second;
321bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
322bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
323bb8a0c61SJames Wright 
324aef1eb53SLeila Ghaffari   T_inf   *= Kelvin;
325e0d1a4dfSLeila Ghaffari   T_wall *= Kelvin;
326bb8a0c61SJames Wright   P0     *= Pascal;
327aef1eb53SLeila Ghaffari   U_inf   *= meter / second;
328bb8a0c61SJames Wright   delta0 *= meter;
329bb8a0c61SJames Wright 
330f362fe62SJames Wright   PetscReal *mesh_ynodes = NULL;
331f362fe62SJames Wright   PetscInt  mesh_nynodes = 0;
332f362fe62SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
333f362fe62SJames Wright     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
334bb8a0c61SJames Wright     CHKERRQ(ierr);
335f362fe62SJames Wright   }
336f362fe62SJames Wright   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
3378e63b9cbSJames Wright                     mesh_refine_height, mesh_top_angle, &mesh_ynodes,
3388e63b9cbSJames Wright                     &mesh_nynodes); CHKERRQ(ierr);
339bb8a0c61SJames Wright 
34015a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
34115a3537eSJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
34215a3537eSJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
343bb8a0c61SJames Wright 
344493642f1SJames Wright   blasius_ctx->weakT         = weakT;
345aef1eb53SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
346aef1eb53SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
347e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
34815a3537eSJed Brown   blasius_ctx->delta0        = delta0;
34915a3537eSJed Brown   blasius_ctx->P0            = P0;
350e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb        = N;
35104b9037bSJames Wright   newtonian_ig_ctx->P0       = P0;
35215a3537eSJed Brown   blasius_ctx->implicit      = user->phys->implicit;
35315a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
354493642f1SJames Wright 
355ef2c71fdSJames Wright   {
356e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
357e0d1a4dfSLeila Ghaffari     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
358ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
359e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
360ef2c71fdSJames Wright   }
361e0d1a4dfSLeila Ghaffari   PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
362ef2c71fdSJames Wright 
36315a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
36415a3537eSJed Brown                                   &newtonian_ig_ctx);
365bb8a0c61SJames Wright 
36615a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
367e0d1a4dfSLeila Ghaffari   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER,
36815a3537eSJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
36915a3537eSJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
37015a3537eSJed Brown                                      FreeContextPetsc);
371bb8a0c61SJames Wright 
37243bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
37315a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
374493642f1SJames Wright   if (use_stg) {
375aef1eb53SLeila Ghaffari     ierr = SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes,
376e6098bcdSJames Wright                     mesh_nynodes); CHKERRQ(ierr);
3778085925cSJames Wright   } else {
3788085925cSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
3798085925cSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
3800aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
3810aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
3828085925cSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
3838085925cSJames Wright                                       &problem->apply_inflow.qfunction_context);
3840aa41abfSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
3850aa41abfSJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
386493642f1SJames Wright   }
387e6098bcdSJames Wright   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
388bb8a0c61SJames Wright   PetscFunctionReturn(0);
389bb8a0c61SJames Wright }
390