xref: /honee/problems/blasius.c (revision aef1eb5380903bd55ea4d009b07c78f59013fb38)
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;
20*aef1eb53SLeila 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     r[3+i] = 2*f[3] + f[2] * f[0];
40e0d1a4dfSLeila Ghaffari     ChebyshevEval(N-1, Th, blasius->X[i], blasius->eta_max, h);
41e0d1a4dfSLeila Ghaffari     r[N+2+i] = h[2] + Pr * f[0] * h[1] +
42e0d1a4dfSLeila Ghaffari                Pr * (gamma - 1) * PetscSqr(Ma * f[2]);
43e0d1a4dfSLeila Ghaffari   }
44e0d1a4dfSLeila Ghaffari 
45e0d1a4dfSLeila Ghaffari   // h - left end boundary condition
46e0d1a4dfSLeila Ghaffari   ChebyshevEval(N-1, Th, -1., blasius->eta_max, h);
47*aef1eb53SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->T_inf;
48e0d1a4dfSLeila Ghaffari 
49e0d1a4dfSLeila Ghaffari   // h - right end boundary condition
50e0d1a4dfSLeila Ghaffari   ChebyshevEval(N-1, Th, 1., blasius->eta_max, h);
51e0d1a4dfSLeila Ghaffari   r[N+1] = h[0] - 1.;
52e0d1a4dfSLeila Ghaffari 
53e0d1a4dfSLeila Ghaffari   // Restore vectors
54e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
55e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
56e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
57e0d1a4dfSLeila Ghaffari }
58e0d1a4dfSLeila Ghaffari 
59e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
60e0d1a4dfSLeila Ghaffari   SNES      snes;
61e0d1a4dfSLeila Ghaffari   Vec       sol, res;
62e0d1a4dfSLeila Ghaffari   PetscReal *w;
63e0d1a4dfSLeila Ghaffari   PetscInt  N = blasius->n_cheb;
64e0d1a4dfSLeila Ghaffari   const PetscScalar *cheb_coefs;
65e0d1a4dfSLeila Ghaffari   PetscFunctionBegin;
66e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(N-3, &blasius->X, N-3, &w));
67e0d1a4dfSLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N-3, -1., 1., blasius->X, w));
68e0d1a4dfSLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
69e0d1a4dfSLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
70e0d1a4dfSLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2*N-1));
71e0d1a4dfSLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
72e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
73e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
74e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
75e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
76e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
77e0d1a4dfSLeila Ghaffari   for (int i=0; i<N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
78e0d1a4dfSLeila Ghaffari   for (int i=0; i<N-1; i++) blasius->Th_cheb[i] = cheb_coefs[i+N];
79e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
80e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
81e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
82e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
83e0d1a4dfSLeila Ghaffari   PetscFunctionReturn(0);
84e0d1a4dfSLeila Ghaffari }
85e0d1a4dfSLeila Ghaffari 
86f362fe62SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
87f362fe62SJames Wright                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
88f362fe62SJames Wright                                    PetscInt *nynodes) {
89f362fe62SJames Wright   PetscErrorCode ierr;
90f362fe62SJames Wright   PetscInt ndims, dims[2];
91f362fe62SJames Wright   FILE *fp;
92f362fe62SJames Wright   const PetscInt char_array_len = 512;
93f362fe62SJames Wright   char line[char_array_len];
94f362fe62SJames Wright   char **array;
95f362fe62SJames Wright   PetscReal *node_locs;
96f362fe62SJames Wright   PetscFunctionBeginUser;
97f362fe62SJames Wright 
98f362fe62SJames Wright   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
99f362fe62SJames Wright   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
100f362fe62SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
101f362fe62SJames Wright 
102f362fe62SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
103f362fe62SJames Wright   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
104f362fe62SJames Wright   *nynodes = dims[0];
105f362fe62SJames Wright   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
106f362fe62SJames Wright 
107f362fe62SJames Wright   for (PetscInt i=0; i<dims[0]; i++) {
108f362fe62SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
109f362fe62SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
110f362fe62SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
11145aa3cadSJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
11245aa3cadSJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT ")", i,
113f362fe62SJames Wright                                    path, ndims, dims[1]);
114f362fe62SJames Wright 
115f362fe62SJames Wright     node_locs[i] = (PetscReal) atof(array[0]);
116f362fe62SJames Wright   }
117f362fe62SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
118f362fe62SJames Wright   *pynodes = node_locs;
119f362fe62SJames Wright   PetscFunctionReturn(0);
120f362fe62SJames Wright }
121f362fe62SJames Wright 
122bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
123bb8a0c61SJames Wright  *
124493642f1SJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
125493642f1SJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
126493642f1SJames Wright  * linearly in logspace to the top surface.
127bb8a0c61SJames Wright  *
128bb8a0c61SJames Wright  * The top surface is also angled downwards, so that it may be used as an
129493642f1SJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
130f362fe62SJames Wright  *
131f362fe62SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
1328e63b9cbSJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
1338e63b9cbSJames Wright  * the array, along with `num_node_locs`.
134bb8a0c61SJames Wright  */
135f362fe62SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
136f362fe62SJames Wright                                  PetscReal growth, PetscInt N,
137f362fe62SJames Wright                                  PetscReal refine_height, PetscReal top_angle,
1388e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
139bb8a0c61SJames Wright   PetscInt ierr, narr, ncoords;
140bb8a0c61SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
141bb8a0c61SJames Wright   PetscScalar *arr_coords;
142bb8a0c61SJames Wright   Vec vec_coords;
143bb8a0c61SJames Wright   PetscFunctionBeginUser;
144bb8a0c61SJames Wright 
145bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
146bb8a0c61SJames Wright 
147bb8a0c61SJames Wright   // Get domain boundary information
148bb8a0c61SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
149493642f1SJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
150bb8a0c61SJames Wright 
151bb8a0c61SJames Wright   // Get coords array from DM
152bb8a0c61SJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
153bb8a0c61SJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
154bb8a0c61SJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
155bb8a0c61SJames Wright 
156bb8a0c61SJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
157bb8a0c61SJames Wright   ncoords = narr/dim;
158bb8a0c61SJames Wright 
159bb8a0c61SJames Wright   // Get mesh information
160bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
161bb8a0c61SJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
162bb8a0c61SJames Wright                                  NULL); CHKERRQ(ierr);
163f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
164f362fe62SJames Wright   const PetscReal dybox = domain_size[1]/faces[1];
165bb8a0c61SJames Wright 
1668e63b9cbSJames Wright   if (!*node_locs) {
167bb8a0c61SJames Wright     // Calculate the first element height
168bb8a0c61SJames Wright     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
169bb8a0c61SJames Wright 
170bb8a0c61SJames Wright     // Calculate log of sizing outside BL
171bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
172bb8a0c61SJames Wright 
1738e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1748e63b9cbSJames Wright     PetscReal *temp_node_locs;
1758e63b9cbSJames Wright     ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr);
1768e63b9cbSJames Wright 
177493642f1SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
178bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
179bb8a0c61SJames Wright       if (y_box_index <= N) {
180029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
181f362fe62SJames Wright                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
182bb8a0c61SJames Wright       } else {
183bb8a0c61SJames Wright         PetscInt j = y_box_index - N;
184029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
185f362fe62SJames Wright                        * exp(log(refine_height) + logdy*j);
186f362fe62SJames Wright       }
1878e63b9cbSJames Wright       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2])
1888e63b9cbSJames Wright         temp_node_locs[y_box_index] = coords[i][1];
189f362fe62SJames Wright     }
1908e63b9cbSJames Wright 
1918e63b9cbSJames Wright     *node_locs = temp_node_locs;
192f362fe62SJames Wright   } else {
193f362fe62SJames Wright     // Error checking
1948e63b9cbSJames Wright     if (*num_node_locs < faces[1] +1)
195f362fe62SJames Wright       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
196f362fe62SJames Wright               "There are %d + 1 nodes, but only %d locations given",
1978e63b9cbSJames Wright               faces[1]+1, *num_node_locs);
1988e63b9cbSJames Wright     if (*num_node_locs > faces[1] +1) {
199f362fe62SJames Wright       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
200e4677755SJames Wright                          "than the mesh has nodes (%d). This maybe unintended.\n",
2018e63b9cbSJames Wright                          *num_node_locs, faces[1]+1); CHKERRQ(ierr);
202f362fe62SJames Wright     }
2038e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
204f362fe62SJames Wright 
205f362fe62SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
206f362fe62SJames Wright       // Determine which y-node we're at
207f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
208029c6dfaSJames Wright       coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y)
2098e63b9cbSJames Wright                      * (*node_locs)[y_box_index];
210bb8a0c61SJames Wright     }
211bb8a0c61SJames Wright   }
212bb8a0c61SJames Wright 
213bb8a0c61SJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
214bb8a0c61SJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
215bb8a0c61SJames Wright 
216bb8a0c61SJames Wright   PetscFunctionReturn(0);
217bb8a0c61SJames Wright }
218bb8a0c61SJames Wright 
219b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
220bb8a0c61SJames Wright 
221bb8a0c61SJames Wright   PetscInt ierr;
222bb8a0c61SJames Wright   User      user    = *(User *)ctx;
223bb8a0c61SJames Wright   MPI_Comm  comm    = PETSC_COMM_WORLD;
224493642f1SJames Wright   PetscBool use_stg = PETSC_FALSE;
22515a3537eSJed Brown   BlasiusContext blasius_ctx;
22615a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
22715a3537eSJed Brown   CeedQFunctionContext blasius_context;
22815a3537eSJed Brown 
229bb8a0c61SJames Wright   PetscFunctionBeginUser;
230b7f03f12SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
23115a3537eSJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
232bb8a0c61SJames Wright 
233bb8a0c61SJames Wright   // ------------------------------------------------------
234bb8a0c61SJames Wright   //               SET UP Blasius
235bb8a0c61SJames Wright   // ------------------------------------------------------
2369785fe93SJed Brown   problem->ics.qfunction        = ICsBlasius;
2379785fe93SJed Brown   problem->ics.qfunction_loc    = ICsBlasius_loc;
238bb8a0c61SJames Wright 
239*aef1eb53SLeila Ghaffari   CeedScalar U_inf              = 40;          // m/s
240*aef1eb53SLeila Ghaffari   CeedScalar T_inf              = 288.;        // K
241e0d1a4dfSLeila Ghaffari   CeedScalar T_wall             = 400.;        // K
242e0d1a4dfSLeila Ghaffari   CeedScalar delta0             = 4.2e-3;      // m
243bb8a0c61SJames Wright   CeedScalar P0                 = 1.01e5;      // Pa
244e0d1a4dfSLeila Ghaffari   CeedInt    N                  = 20;          // Number of Chebyshev terms
2452acc7cbcSKenneth E. Jansen   PetscBool  weakT              = PETSC_FALSE; // weak density or temperature
2467470235eSJames Wright   PetscReal  mesh_refine_height = 5.9e-4;      // m
2477470235eSJames Wright   PetscReal  mesh_growth        = 1.08;        // [-]
2487470235eSJames Wright   PetscInt   mesh_Ndelta        = 45;          // [-]
2497470235eSJames Wright   PetscReal  mesh_top_angle     = 5;           // degrees
250f362fe62SJames Wright   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
251bb8a0c61SJames Wright 
2523fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2532acc7cbcSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
2542acc7cbcSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
255*aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-velocity_infinity",
256*aef1eb53SLeila Ghaffari                             "Velocity at boundary layer edge",
257*aef1eb53SLeila Ghaffari                             NULL, U_inf, &U_inf, NULL); CHKERRQ(ierr);
258*aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-temperature_infinity",
259*aef1eb53SLeila Ghaffari                             "Temperature at boundary layer edge",
260*aef1eb53SLeila Ghaffari                             NULL, T_inf, &T_inf, NULL); CHKERRQ(ierr);
261*aef1eb53SLeila Ghaffari   ierr = PetscOptionsScalar("-temperature_wall", "Temperature at wall",
262e0d1a4dfSLeila Ghaffari                             NULL, T_wall, &T_wall, NULL); CHKERRQ(ierr);
263bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
264bb8a0c61SJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
265bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
266bb8a0c61SJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
267e0d1a4dfSLeila Ghaffari   ierr = PetscOptionsInt("-N_Chebyshev", "Number of Chebyshev terms",
268e0d1a4dfSLeila Ghaffari                          NULL, N, &N, NULL); CHKERRQ(ierr);
2697470235eSJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
2707470235eSJames Wright                                 "Velocity at boundary layer edge",
2717470235eSJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
2727470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
273bb8a0c61SJames Wright                             "Height of boundary layer mesh refinement",
2747470235eSJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
2757470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
276bb8a0c61SJames Wright                             "Geometric growth rate of boundary layer mesh",
2777470235eSJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
2787470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
279bb8a0c61SJames Wright                             "Geometric top_angle rate of boundary layer mesh",
2807470235eSJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
281f362fe62SJames Wright   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
282f362fe62SJames Wright                             "Path to file with y node locations. "
283f362fe62SJames Wright                             "If empty, will use the algorithmic mesh warping.", NULL,
284f362fe62SJames Wright                             mesh_ynodes_path, mesh_ynodes_path,
285f362fe62SJames Wright                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
286493642f1SJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
287493642f1SJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
288bb8a0c61SJames Wright   PetscOptionsEnd();
289bb8a0c61SJames Wright 
290bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
291bb8a0c61SJames Wright   PetscScalar second = user->units->second;
292bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
293bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
294bb8a0c61SJames Wright 
295*aef1eb53SLeila Ghaffari   T_inf   *= Kelvin;
296e0d1a4dfSLeila Ghaffari   T_wall *= Kelvin;
297bb8a0c61SJames Wright   P0     *= Pascal;
298*aef1eb53SLeila Ghaffari   U_inf   *= meter / second;
299bb8a0c61SJames Wright   delta0 *= meter;
300bb8a0c61SJames Wright 
301f362fe62SJames Wright   PetscReal *mesh_ynodes = NULL;
302f362fe62SJames Wright   PetscInt  mesh_nynodes = 0;
303f362fe62SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
304f362fe62SJames Wright     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
305bb8a0c61SJames Wright     CHKERRQ(ierr);
306f362fe62SJames Wright   }
307f362fe62SJames Wright   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
3088e63b9cbSJames Wright                     mesh_refine_height, mesh_top_angle, &mesh_ynodes,
3098e63b9cbSJames Wright                     &mesh_nynodes); CHKERRQ(ierr);
310bb8a0c61SJames Wright 
31115a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
31215a3537eSJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
31315a3537eSJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
314bb8a0c61SJames Wright 
315493642f1SJames Wright   blasius_ctx->weakT         = weakT;
316*aef1eb53SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
317*aef1eb53SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
318e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
31915a3537eSJed Brown   blasius_ctx->delta0        = delta0;
32015a3537eSJed Brown   blasius_ctx->P0            = P0;
321e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb        = N;
32204b9037bSJames Wright   newtonian_ig_ctx->P0       = P0;
32315a3537eSJed Brown   blasius_ctx->implicit      = user->phys->implicit;
32415a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
325493642f1SJames Wright 
326ef2c71fdSJames Wright   {
327e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
328e0d1a4dfSLeila Ghaffari     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
329ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
330e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
331ef2c71fdSJames Wright   }
332e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(blasius_ctx->n_cheb, &blasius_ctx->Tf_cheb,
333e0d1a4dfSLeila Ghaffari                          blasius_ctx->n_cheb-1, &blasius_ctx->Th_cheb));
334e0d1a4dfSLeila Ghaffari   PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
335ef2c71fdSJames Wright 
33615a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
33715a3537eSJed Brown                                   &newtonian_ig_ctx);
338bb8a0c61SJames Wright 
33915a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
340e0d1a4dfSLeila Ghaffari   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER,
34115a3537eSJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
34215a3537eSJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
34315a3537eSJed Brown                                      FreeContextPetsc);
344bb8a0c61SJames Wright 
34543bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
34615a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
347493642f1SJames Wright   if (use_stg) {
348*aef1eb53SLeila Ghaffari     ierr = SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes,
349e6098bcdSJames Wright                     mesh_nynodes); CHKERRQ(ierr);
3508085925cSJames Wright   } else {
3518085925cSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
3528085925cSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
3530aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
3540aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
3558085925cSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
3568085925cSJames Wright                                       &problem->apply_inflow.qfunction_context);
3570aa41abfSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
3580aa41abfSJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
359493642f1SJames Wright   }
360e6098bcdSJames Wright   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
361bb8a0c61SJames Wright   PetscFunctionReturn(0);
362bb8a0c61SJames Wright }
363