xref: /libCEED/examples/fluids/problems/blasius.c (revision 2518f336a2bcc61653a1548f8a480d533f99d612)
188626eedSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
288626eedSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388626eedSJames Wright //
488626eedSJames Wright // SPDX-License-Identifier: BSD-2-Clause
588626eedSJames Wright //
688626eedSJames Wright // This file is part of CEED:  http://github.com/ceed
788626eedSJames Wright 
888626eedSJames Wright /// @file
988626eedSJames Wright /// Utility functions for setting up Blasius Boundary Layer
1088626eedSJames Wright 
1188626eedSJames Wright #include "../navierstokes.h"
1288626eedSJames Wright #include "../qfunctions/blasius.h"
13ba6664aeSJames Wright #include "stg_shur14.h"
1488626eedSJames Wright 
15*2518f336SLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
16*2518f336SLeila Ghaffari   const BlasiusContext blasius = (BlasiusContext)ctx;
17*2518f336SLeila Ghaffari   const PetscScalar *Tf, *Th;  // Chebyshev coefficients
18*2518f336SLeila Ghaffari   PetscScalar       *r, f[4], h[4];
19*2518f336SLeila Ghaffari   PetscInt          N = blasius->n_cheb;
20*2518f336SLeila Ghaffari   PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->Tinf, blasius->Uinf),
21*2518f336SLeila Ghaffari               Pr = Prandtl(&blasius->newtonian_ctx),
22*2518f336SLeila Ghaffari               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
23*2518f336SLeila Ghaffari   PetscFunctionBegin;
24*2518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
25*2518f336SLeila Ghaffari   Th = Tf + N;
26*2518f336SLeila Ghaffari   PetscCall(VecGetArray(R, &r));
27*2518f336SLeila Ghaffari 
28*2518f336SLeila Ghaffari   // Left boundary conditions f = f' = 0
29*2518f336SLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
30*2518f336SLeila Ghaffari   r[0] = f[0];
31*2518f336SLeila Ghaffari   r[1] = f[1];
32*2518f336SLeila Ghaffari 
33*2518f336SLeila Ghaffari   // f - right end boundary condition
34*2518f336SLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
35*2518f336SLeila Ghaffari   r[2] = f[1]  - 1.;
36*2518f336SLeila Ghaffari 
37*2518f336SLeila Ghaffari   for (int i=0; i<N-3; i++) {
38*2518f336SLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
39*2518f336SLeila Ghaffari     r[3+i] = 2*f[3] + f[2] * f[0];
40*2518f336SLeila Ghaffari     ChebyshevEval(N-1, Th, blasius->X[i], blasius->eta_max, h);
41*2518f336SLeila Ghaffari     r[N+2+i] = h[2] + Pr * f[0] * h[1] +
42*2518f336SLeila Ghaffari                Pr * (gamma - 1) * PetscSqr(Ma * f[2]);
43*2518f336SLeila Ghaffari   }
44*2518f336SLeila Ghaffari 
45*2518f336SLeila Ghaffari   // h - left end boundary condition
46*2518f336SLeila Ghaffari   ChebyshevEval(N-1, Th, -1., blasius->eta_max, h);
47*2518f336SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->Tinf;
48*2518f336SLeila Ghaffari 
49*2518f336SLeila Ghaffari   // h - right end boundary condition
50*2518f336SLeila Ghaffari   ChebyshevEval(N-1, Th, 1., blasius->eta_max, h);
51*2518f336SLeila Ghaffari   r[N+1] = h[0] - 1.;
52*2518f336SLeila Ghaffari 
53*2518f336SLeila Ghaffari   // Restore vectors
54*2518f336SLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
55*2518f336SLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
56*2518f336SLeila Ghaffari   PetscFunctionReturn(0);
57*2518f336SLeila Ghaffari }
58*2518f336SLeila Ghaffari 
59*2518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
60*2518f336SLeila Ghaffari   SNES      snes;
61*2518f336SLeila Ghaffari   Vec       sol, res;
62*2518f336SLeila Ghaffari   PetscReal *w;
63*2518f336SLeila Ghaffari   PetscInt  N = blasius->n_cheb;
64*2518f336SLeila Ghaffari   const PetscScalar *cheb_coefs;
65*2518f336SLeila Ghaffari   PetscFunctionBegin;
66*2518f336SLeila Ghaffari   PetscCall(PetscMalloc2(N-3, &blasius->X, N-3, &w));
67*2518f336SLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N-3, -1., 1., blasius->X, w));
68*2518f336SLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
69*2518f336SLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
70*2518f336SLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2*N-1));
71*2518f336SLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
72*2518f336SLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
73*2518f336SLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
74*2518f336SLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
75*2518f336SLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
76*2518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
77*2518f336SLeila Ghaffari   for (int i=0; i<N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
78*2518f336SLeila Ghaffari   for (int i=0; i<N-1; i++) blasius->Th_cheb[i] = cheb_coefs[i+N];
79*2518f336SLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
80*2518f336SLeila Ghaffari   PetscCall(VecDestroy(&sol));
81*2518f336SLeila Ghaffari   PetscCall(VecDestroy(&res));
82*2518f336SLeila Ghaffari   PetscCall(SNESDestroy(&snes));
83*2518f336SLeila Ghaffari   PetscFunctionReturn(0);
84*2518f336SLeila Ghaffari }
85*2518f336SLeila Ghaffari 
86d2714514SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
87d2714514SJames Wright                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
88d2714514SJames Wright                                    PetscInt *nynodes) {
89d2714514SJames Wright   PetscErrorCode ierr;
90d2714514SJames Wright   PetscInt ndims, dims[2];
91d2714514SJames Wright   FILE *fp;
92d2714514SJames Wright   const PetscInt char_array_len = 512;
93d2714514SJames Wright   char line[char_array_len];
94d2714514SJames Wright   char **array;
95d2714514SJames Wright   PetscReal *node_locs;
96d2714514SJames Wright   PetscFunctionBeginUser;
97d2714514SJames Wright 
98d2714514SJames Wright   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
99d2714514SJames Wright   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
100d2714514SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
101d2714514SJames Wright 
102d2714514SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
103d2714514SJames Wright   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
104d2714514SJames Wright   *nynodes = dims[0];
105d2714514SJames Wright   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
106d2714514SJames Wright 
107d2714514SJames Wright   for (PetscInt i=0; i<dims[0]; i++) {
108d2714514SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
109d2714514SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
110d2714514SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
111990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
112990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT ")", i,
113d2714514SJames Wright                                    path, ndims, dims[1]);
114d2714514SJames Wright 
115d2714514SJames Wright     node_locs[i] = (PetscReal) atof(array[0]);
116d2714514SJames Wright   }
117d2714514SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
118d2714514SJames Wright   *pynodes = node_locs;
119d2714514SJames Wright   PetscFunctionReturn(0);
120d2714514SJames Wright }
121d2714514SJames Wright 
12288626eedSJames Wright /* \brief Modify the domain and mesh for blasius
12388626eedSJames Wright  *
124ba6664aeSJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
125ba6664aeSJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
126ba6664aeSJames Wright  * linearly in logspace to the top surface.
12788626eedSJames Wright  *
12888626eedSJames Wright  * The top surface is also angled downwards, so that it may be used as an
129ba6664aeSJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
130d2714514SJames Wright  *
131d2714514SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
132798d42b8SJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
133798d42b8SJames Wright  * the array, along with `num_node_locs`.
13488626eedSJames Wright  */
135d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
136d2714514SJames Wright                                  PetscReal growth, PetscInt N,
137d2714514SJames Wright                                  PetscReal refine_height, PetscReal top_angle,
138798d42b8SJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
13988626eedSJames Wright   PetscInt ierr, narr, ncoords;
14088626eedSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
14188626eedSJames Wright   PetscScalar *arr_coords;
14288626eedSJames Wright   Vec vec_coords;
14388626eedSJames Wright   PetscFunctionBeginUser;
14488626eedSJames Wright 
14588626eedSJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
14688626eedSJames Wright 
14788626eedSJames Wright   // Get domain boundary information
14888626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
149ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
15088626eedSJames Wright 
15188626eedSJames Wright   // Get coords array from DM
15288626eedSJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
15388626eedSJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
15488626eedSJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
15588626eedSJames Wright 
15688626eedSJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
15788626eedSJames Wright   ncoords = narr/dim;
15888626eedSJames Wright 
15988626eedSJames Wright   // Get mesh information
16088626eedSJames Wright   PetscInt nmax = 3, faces[3];
16188626eedSJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
16288626eedSJames Wright                                  NULL); CHKERRQ(ierr);
163d2714514SJames Wright   // Get element size of the box mesh, for indexing each node
164d2714514SJames Wright   const PetscReal dybox = domain_size[1]/faces[1];
16588626eedSJames Wright 
166798d42b8SJames Wright   if (!*node_locs) {
16788626eedSJames Wright     // Calculate the first element height
16888626eedSJames Wright     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
16988626eedSJames Wright 
17088626eedSJames Wright     // Calculate log of sizing outside BL
17188626eedSJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
17288626eedSJames Wright 
173798d42b8SJames Wright     *num_node_locs = faces[1] + 1;
174798d42b8SJames Wright     PetscReal *temp_node_locs;
175798d42b8SJames Wright     ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr);
176798d42b8SJames Wright 
177ba6664aeSJames Wright     for (PetscInt i=0; i<ncoords; i++) {
17888626eedSJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
17988626eedSJames Wright       if (y_box_index <= N) {
180855641eaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
181d2714514SJames Wright                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
18288626eedSJames Wright       } else {
18388626eedSJames Wright         PetscInt j = y_box_index - N;
184855641eaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
185d2714514SJames Wright                        * exp(log(refine_height) + logdy*j);
186d2714514SJames Wright       }
187798d42b8SJames Wright       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2])
188798d42b8SJames Wright         temp_node_locs[y_box_index] = coords[i][1];
189d2714514SJames Wright     }
190798d42b8SJames Wright 
191798d42b8SJames Wright     *node_locs = temp_node_locs;
192d2714514SJames Wright   } else {
193d2714514SJames Wright     // Error checking
194798d42b8SJames Wright     if (*num_node_locs < faces[1] +1)
195d2714514SJames Wright       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
196d2714514SJames Wright               "There are %d + 1 nodes, but only %d locations given",
197798d42b8SJames Wright               faces[1]+1, *num_node_locs);
198798d42b8SJames Wright     if (*num_node_locs > faces[1] +1) {
199d2714514SJames Wright       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
200c10dcd27SJames Wright                          "than the mesh has nodes (%d). This maybe unintended.\n",
201798d42b8SJames Wright                          *num_node_locs, faces[1]+1); CHKERRQ(ierr);
202d2714514SJames Wright     }
203798d42b8SJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
204d2714514SJames Wright 
205d2714514SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
206d2714514SJames Wright       // Determine which y-node we're at
207d2714514SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
208855641eaSJames Wright       coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y)
209798d42b8SJames Wright                      * (*node_locs)[y_box_index];
21088626eedSJames Wright     }
21188626eedSJames Wright   }
21288626eedSJames Wright 
21388626eedSJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
21488626eedSJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
21588626eedSJames Wright 
21688626eedSJames Wright   PetscFunctionReturn(0);
21788626eedSJames Wright }
21888626eedSJames Wright 
219a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
22088626eedSJames Wright 
22188626eedSJames Wright   PetscInt ierr;
22288626eedSJames Wright   User      user    = *(User *)ctx;
22388626eedSJames Wright   MPI_Comm  comm    = PETSC_COMM_WORLD;
224ba6664aeSJames Wright   PetscBool use_stg = PETSC_FALSE;
225841e4c73SJed Brown   BlasiusContext blasius_ctx;
226841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
227841e4c73SJed Brown   CeedQFunctionContext blasius_context;
228841e4c73SJed Brown 
22988626eedSJames Wright   PetscFunctionBeginUser;
230a0add3c9SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
231841e4c73SJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
23288626eedSJames Wright 
23388626eedSJames Wright   // ------------------------------------------------------
23488626eedSJames Wright   //               SET UP Blasius
23588626eedSJames Wright   // ------------------------------------------------------
23691e5af17SJed Brown   problem->ics.qfunction        = ICsBlasius;
23791e5af17SJed Brown   problem->ics.qfunction_loc    = ICsBlasius_loc;
23888626eedSJames Wright 
23988626eedSJames Wright   CeedScalar Uinf               = 40;          // m/s
240*2518f336SLeila Ghaffari   CeedScalar Tinf               = 288.;        // K
241*2518f336SLeila Ghaffari   CeedScalar T_wall             = 400.;        // K
242*2518f336SLeila Ghaffari   CeedScalar delta0             = 4.2e-3;      // m
24388626eedSJames Wright   CeedScalar theta0             = 288.;        // K
24488626eedSJames Wright   CeedScalar P0                 = 1.01e5;      // Pa
245*2518f336SLeila Ghaffari   CeedInt    N                  = 20;          // Number of Chebyshev terms
246871db79fSKenneth E. Jansen   PetscBool  weakT              = PETSC_FALSE; // weak density or temperature
2477b8d3891SJames Wright   PetscReal  mesh_refine_height = 5.9e-4;      // m
2487b8d3891SJames Wright   PetscReal  mesh_growth        = 1.08;        // [-]
2497b8d3891SJames Wright   PetscInt   mesh_Ndelta        = 45;          // [-]
2507b8d3891SJames Wright   PetscReal  mesh_top_angle     = 5;           // degrees
251d2714514SJames Wright   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
25288626eedSJames Wright 
2537b37518fSJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
254871db79fSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
255871db79fSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
25688626eedSJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
25788626eedSJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
258*2518f336SLeila Ghaffari   ierr = PetscOptionsScalar("-Tinf", "Temperature at boundary layer edge",
259*2518f336SLeila Ghaffari                             NULL, Tinf, &Tinf, NULL); CHKERRQ(ierr);
260*2518f336SLeila Ghaffari   ierr = PetscOptionsScalar("-T_wall", "Temperature at wall",
261*2518f336SLeila Ghaffari                             NULL, T_wall, &T_wall, NULL); CHKERRQ(ierr);
26288626eedSJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
26388626eedSJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
26488626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
26588626eedSJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
26688626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
26788626eedSJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
268*2518f336SLeila Ghaffari   ierr = PetscOptionsInt("-N_Chebyshev", "Number of Chebyshev terms",
269*2518f336SLeila Ghaffari                          NULL, N, &N, NULL); CHKERRQ(ierr);
2707b8d3891SJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
2717b8d3891SJames Wright                                 "Velocity at boundary layer edge",
2727b8d3891SJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
2737b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
27488626eedSJames Wright                             "Height of boundary layer mesh refinement",
2757b8d3891SJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
2767b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
27788626eedSJames Wright                             "Geometric growth rate of boundary layer mesh",
2787b8d3891SJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
2797b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
28088626eedSJames Wright                             "Geometric top_angle rate of boundary layer mesh",
2817b8d3891SJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
282d2714514SJames Wright   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
283d2714514SJames Wright                             "Path to file with y node locations. "
284d2714514SJames Wright                             "If empty, will use the algorithmic mesh warping.", NULL,
285d2714514SJames Wright                             mesh_ynodes_path, mesh_ynodes_path,
286d2714514SJames Wright                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
287ba6664aeSJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
288ba6664aeSJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
28988626eedSJames Wright   PetscOptionsEnd();
29088626eedSJames Wright 
29188626eedSJames Wright   PetscScalar meter  = user->units->meter;
29288626eedSJames Wright   PetscScalar second = user->units->second;
29388626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
29488626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
29588626eedSJames Wright 
29688626eedSJames Wright   theta0 *= Kelvin;
297*2518f336SLeila Ghaffari   Tinf   *= Kelvin;
298*2518f336SLeila Ghaffari   T_wall *= Kelvin;
29988626eedSJames Wright   P0     *= Pascal;
30088626eedSJames Wright   Uinf   *= meter / second;
30188626eedSJames Wright   delta0 *= meter;
30288626eedSJames Wright 
303d2714514SJames Wright   PetscReal *mesh_ynodes = NULL;
304d2714514SJames Wright   PetscInt  mesh_nynodes = 0;
305d2714514SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
306d2714514SJames Wright     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
30788626eedSJames Wright     CHKERRQ(ierr);
308d2714514SJames Wright   }
309d2714514SJames Wright   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
310798d42b8SJames Wright                     mesh_refine_height, mesh_top_angle, &mesh_ynodes,
311798d42b8SJames Wright                     &mesh_nynodes); CHKERRQ(ierr);
31288626eedSJames Wright 
313841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
314841e4c73SJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
315841e4c73SJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
31688626eedSJames Wright 
317ba6664aeSJames Wright   blasius_ctx->weakT         = weakT;
318841e4c73SJed Brown   blasius_ctx->Uinf          = Uinf;
319*2518f336SLeila Ghaffari   blasius_ctx->Tinf          = Tinf;
320*2518f336SLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
321841e4c73SJed Brown   blasius_ctx->delta0        = delta0;
322841e4c73SJed Brown   blasius_ctx->theta0        = theta0;
323841e4c73SJed Brown   blasius_ctx->P0            = P0;
324*2518f336SLeila Ghaffari   blasius_ctx->n_cheb        = N;
32530e9fa81SJames Wright   newtonian_ig_ctx->P0       = P0;
326841e4c73SJed Brown   blasius_ctx->implicit      = user->phys->implicit;
327841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
328ba6664aeSJames Wright 
329f1122ed0SJames Wright   {
330*2518f336SLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
331*2518f336SLeila Ghaffari     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
332f1122ed0SJames Wright     blasius_ctx->x_inflow = domain_min[0];
333*2518f336SLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
334f1122ed0SJames Wright   }
335*2518f336SLeila Ghaffari   PetscCall(PetscMalloc2(blasius_ctx->n_cheb, &blasius_ctx->Tf_cheb,
336*2518f336SLeila Ghaffari                          blasius_ctx->n_cheb-1, &blasius_ctx->Th_cheb));
337*2518f336SLeila Ghaffari   PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
338f1122ed0SJames Wright 
339841e4c73SJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
340841e4c73SJed Brown                                   &newtonian_ig_ctx);
34188626eedSJames Wright 
342841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
343*2518f336SLeila Ghaffari   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER,
344841e4c73SJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
345841e4c73SJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
346841e4c73SJed Brown                                      FreeContextPetsc);
34788626eedSJames Wright 
348b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
349841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
350ba6664aeSJames Wright   if (use_stg) {
3514e139266SJames Wright     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes,
3524e139266SJames Wright                     mesh_nynodes); CHKERRQ(ierr);
35365dd5cafSJames Wright   } else {
35465dd5cafSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
35565dd5cafSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
356fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
357fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
35865dd5cafSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
35965dd5cafSJames Wright                                       &problem->apply_inflow.qfunction_context);
360fc02c281SJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
361fc02c281SJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
362ba6664aeSJames Wright   }
3634e139266SJames Wright   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
36488626eedSJames Wright   PetscFunctionReturn(0);
36588626eedSJames Wright }
366