xref: /libCEED/examples/fluids/problems/blasius.c (revision 46de7363edb07f0ea682e5e2fa0ce97210b98e34)
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 "../qfunctions/blasius.h"
122b730f8bSJeremy L Thompson 
132b730f8bSJeremy L Thompson #include "../navierstokes.h"
14ba6664aeSJames Wright #include "stg_shur14.h"
1588626eedSJames Wright 
162518f336SLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
172518f336SLeila Ghaffari   const BlasiusContext blasius = (BlasiusContext)ctx;
182518f336SLeila Ghaffari   const PetscScalar   *Tf, *Th;  // Chebyshev coefficients
192518f336SLeila Ghaffari   PetscScalar         *r, f[4], h[4];
202518f336SLeila Ghaffari   PetscInt             N  = blasius->n_cheb;
212b730f8bSJeremy L Thompson   PetscScalar          Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx),
222518f336SLeila Ghaffari               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
232518f336SLeila Ghaffari   PetscFunctionBegin;
242518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
252518f336SLeila Ghaffari   Th = Tf + N;
262518f336SLeila Ghaffari   PetscCall(VecGetArray(R, &r));
272518f336SLeila Ghaffari 
282518f336SLeila Ghaffari   // Left boundary conditions f = f' = 0
292518f336SLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
302518f336SLeila Ghaffari   r[0] = f[0];
312518f336SLeila Ghaffari   r[1] = f[1];
322518f336SLeila Ghaffari 
332518f336SLeila Ghaffari   // f - right end boundary condition
342518f336SLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
352518f336SLeila Ghaffari   r[2] = f[1] - 1.;
362518f336SLeila Ghaffari 
372518f336SLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
382518f336SLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
392518f336SLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
4007d14e58SLeila Ghaffari     // mu and rho generally depend on h. We naively assume constant mu.
4107d14e58SLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
4207d14e58SLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
4307d14e58SLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
4407d14e58SLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
4507d14e58SLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
4607d14e58SLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
4707d14e58SLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
4807d14e58SLeila Ghaffari     };
4907d14e58SLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
502b730f8bSJeremy 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]);
512518f336SLeila Ghaffari   }
522518f336SLeila Ghaffari 
532518f336SLeila Ghaffari   // h - left end boundary condition
542518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
55fb455ff0SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->T_inf;
562518f336SLeila Ghaffari 
572518f336SLeila Ghaffari   // h - right end boundary condition
582518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
592518f336SLeila Ghaffari   r[N + 1] = h[0] - 1.;
602518f336SLeila Ghaffari 
612518f336SLeila Ghaffari   // Restore vectors
622518f336SLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
632518f336SLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
642518f336SLeila Ghaffari   PetscFunctionReturn(0);
652518f336SLeila Ghaffari }
662518f336SLeila Ghaffari 
672518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
682518f336SLeila Ghaffari   SNES                snes;
692518f336SLeila Ghaffari   Vec                 sol, res;
702518f336SLeila Ghaffari   PetscReal          *w;
712518f336SLeila Ghaffari   PetscInt            N = blasius->n_cheb;
7207d14e58SLeila Ghaffari   SNESConvergedReason reason;
732518f336SLeila Ghaffari   const PetscScalar  *cheb_coefs;
742518f336SLeila Ghaffari   PetscFunctionBegin;
7507d14e58SLeila Ghaffari 
7607d14e58SLeila Ghaffari   // Allocate memory
772518f336SLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
782518f336SLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
7907d14e58SLeila Ghaffari 
8007d14e58SLeila Ghaffari   // Snes solve
812518f336SLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
822518f336SLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
832518f336SLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
842518f336SLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
8507d14e58SLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
8607d14e58SLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
872518f336SLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
882518f336SLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
8907d14e58SLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
902518f336SLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
912518f336SLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
9207d14e58SLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
932b730f8bSJeremy L Thompson   if (reason < 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
9407d14e58SLeila Ghaffari 
9507d14e58SLeila Ghaffari   // Assign Chebyshev coefficients
962518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
972518f336SLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
982518f336SLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
9907d14e58SLeila Ghaffari 
10007d14e58SLeila Ghaffari   // Destroy objects
1012518f336SLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
1022518f336SLeila Ghaffari   PetscCall(VecDestroy(&sol));
1032518f336SLeila Ghaffari   PetscCall(VecDestroy(&res));
1042518f336SLeila Ghaffari   PetscCall(SNESDestroy(&snes));
1052518f336SLeila Ghaffari   PetscFunctionReturn(0);
1062518f336SLeila Ghaffari }
1072518f336SLeila Ghaffari 
1082b730f8bSJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
109d2714514SJames Wright   PetscInt       ndims, dims[2];
110d2714514SJames Wright   FILE          *fp;
111d2714514SJames Wright   const PetscInt char_array_len = 512;
112d2714514SJames Wright   char           line[char_array_len];
113d2714514SJames Wright   char         **array;
114d2714514SJames Wright   PetscReal     *node_locs;
115d2714514SJames Wright   PetscFunctionBeginUser;
116d2714514SJames Wright 
1172b730f8bSJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1182b730f8bSJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1192b730f8bSJeremy L Thompson   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
120d2714514SJames Wright 
121d2714514SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
122d2714514SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
123d2714514SJames Wright   *nynodes = dims[0];
1242b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
125d2714514SJames Wright 
126d2714514SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
1272b730f8bSJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1282b730f8bSJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1292b730f8bSJeremy L Thompson     if (ndims < dims[1])
1302b730f8bSJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1312b730f8bSJeremy L Thompson               ndims, dims[1]);
132d2714514SJames Wright 
133d2714514SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
134d2714514SJames Wright   }
1352b730f8bSJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
136d2714514SJames Wright   *pynodes = node_locs;
137d2714514SJames Wright   PetscFunctionReturn(0);
138d2714514SJames Wright }
139d2714514SJames Wright 
14088626eedSJames Wright /* \brief Modify the domain and mesh for blasius
14188626eedSJames Wright  *
142ba6664aeSJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
143ba6664aeSJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
144ba6664aeSJames Wright  * linearly in logspace to the top surface.
14588626eedSJames Wright  *
14688626eedSJames Wright  * The top surface is also angled downwards, so that it may be used as an
147ba6664aeSJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
148d2714514SJames Wright  *
149d2714514SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
150798d42b8SJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
151798d42b8SJames Wright  * the array, along with `num_node_locs`.
15288626eedSJames Wright  */
1532b730f8bSJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
154798d42b8SJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
1552b730f8bSJeremy L Thompson   PetscInt     narr, ncoords;
15688626eedSJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
15788626eedSJames Wright   PetscScalar *arr_coords;
15888626eedSJames Wright   Vec          vec_coords;
15988626eedSJames Wright   PetscFunctionBeginUser;
16088626eedSJames Wright 
16188626eedSJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
16288626eedSJames Wright 
16388626eedSJames Wright   // Get domain boundary information
1642b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
165ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
16688626eedSJames Wright 
16788626eedSJames Wright   // Get coords array from DM
1682b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1692b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1702b730f8bSJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
17188626eedSJames Wright 
17288626eedSJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
17388626eedSJames Wright   ncoords                   = narr / dim;
17488626eedSJames Wright 
17588626eedSJames Wright   // Get mesh information
17688626eedSJames Wright   PetscInt nmax = 3, faces[3];
1772b730f8bSJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
178d2714514SJames Wright   // Get element size of the box mesh, for indexing each node
179d2714514SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
18088626eedSJames Wright 
181798d42b8SJames Wright   if (!*node_locs) {
18288626eedSJames Wright     // Calculate the first element height
18388626eedSJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
18488626eedSJames Wright 
18588626eedSJames Wright     // Calculate log of sizing outside BL
18688626eedSJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
18788626eedSJames Wright 
188798d42b8SJames Wright     *num_node_locs = faces[1] + 1;
189798d42b8SJames Wright     PetscReal *temp_node_locs;
1902b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
191798d42b8SJames Wright 
192ba6664aeSJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
19388626eedSJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
19488626eedSJames Wright       if (y_box_index <= N) {
1952b730f8bSJeremy L Thompson         coords[i][1] =
1962b730f8bSJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
19788626eedSJames Wright       } else {
19888626eedSJames Wright         PetscInt j   = y_box_index - N;
1992b730f8bSJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
200d2714514SJames Wright       }
2012b730f8bSJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
202d2714514SJames Wright     }
203798d42b8SJames Wright 
204798d42b8SJames Wright     *node_locs = temp_node_locs;
205d2714514SJames Wright   } else {
206d2714514SJames Wright     // Error checking
2072b730f8bSJeremy L Thompson     if (*num_node_locs < faces[1] + 1) {
2082b730f8bSJeremy L Thompson       SETERRQ(comm, -1,
2092b730f8bSJeremy L Thompson               "The y_node_locs_path has too few locations; "
210d2714514SJames Wright               "There are %d + 1 nodes, but only %d locations given",
211798d42b8SJames Wright               faces[1] + 1, *num_node_locs);
2122b730f8bSJeremy L Thompson     }
213798d42b8SJames Wright     if (*num_node_locs > faces[1] + 1) {
2142b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
2152b730f8bSJeremy L Thompson                             "WARNING: y_node_locs_path has more locations (%d) "
216c10dcd27SJames Wright                             "than the mesh has nodes (%d). This maybe unintended.\n",
2172b730f8bSJeremy L Thompson                             *num_node_locs, faces[1] + 1));
218d2714514SJames Wright     }
219798d42b8SJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
220d2714514SJames Wright 
221d2714514SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
222d2714514SJames Wright       // Determine which y-node we're at
223d2714514SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2242b730f8bSJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
22588626eedSJames Wright     }
22688626eedSJames Wright   }
22788626eedSJames Wright 
2282b730f8bSJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2292b730f8bSJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
23088626eedSJames Wright 
23188626eedSJames Wright   PetscFunctionReturn(0);
23288626eedSJames Wright }
23388626eedSJames Wright 
234*46de7363SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
23588626eedSJames Wright   User                     user    = *(User *)ctx;
23688626eedSJames Wright   MPI_Comm                 comm    = PETSC_COMM_WORLD;
237ba6664aeSJames Wright   PetscBool                use_stg = PETSC_FALSE;
238841e4c73SJed Brown   BlasiusContext           blasius_ctx;
239841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
240841e4c73SJed Brown   CeedQFunctionContext     blasius_context;
241841e4c73SJed Brown 
24288626eedSJames Wright   PetscFunctionBeginUser;
243*46de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
2442b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &blasius_ctx));
24588626eedSJames Wright 
24688626eedSJames Wright   // ------------------------------------------------------
24788626eedSJames Wright   //               SET UP Blasius
24888626eedSJames Wright   // ------------------------------------------------------
24991e5af17SJed Brown   problem->ics.qfunction     = ICsBlasius;
25091e5af17SJed Brown   problem->ics.qfunction_loc = ICsBlasius_loc;
25188626eedSJames Wright 
252fb455ff0SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
253fb455ff0SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
25407d14e58SLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
2552518f336SLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
25688626eedSJames Wright   CeedScalar P0                                   = 1.01e5;       // Pa
2572518f336SLeila Ghaffari   CeedInt    N                                    = 20;           // Number of Chebyshev terms
258871db79fSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2597b8d3891SJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2607b8d3891SJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2617b8d3891SJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2627b8d3891SJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
263d2714514SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
26488626eedSJames Wright 
2657b37518fSJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2662b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2672b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2682b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
2692b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2702b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2712b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
2722b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
2732b730f8bSJeremy 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,
27407d14e58SLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2752b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
2762b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height,
2772b730f8bSJeremy L Thompson                                NULL));
2782b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
2792b730f8bSJeremy L Thompson   PetscCall(
2802b730f8bSJeremy L Thompson       PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
2812b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
282d2714514SJames Wright                                "Path to file with y node locations. "
2832b730f8bSJeremy L Thompson                                "If empty, will use the algorithmic mesh warping.",
2842b730f8bSJeremy L Thompson                                NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
2852b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
28688626eedSJames Wright   PetscOptionsEnd();
28788626eedSJames Wright 
28888626eedSJames Wright   PetscScalar meter  = user->units->meter;
28988626eedSJames Wright   PetscScalar second = user->units->second;
29088626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
29188626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
29288626eedSJames Wright 
293fb455ff0SLeila Ghaffari   T_inf *= Kelvin;
2942518f336SLeila Ghaffari   T_wall *= Kelvin;
29588626eedSJames Wright   P0 *= Pascal;
296fb455ff0SLeila Ghaffari   U_inf *= meter / second;
29788626eedSJames Wright   delta0 *= meter;
29888626eedSJames Wright 
299d2714514SJames Wright   PetscReal *mesh_ynodes  = NULL;
300d2714514SJames Wright   PetscInt   mesh_nynodes = 0;
301d2714514SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
3022b730f8bSJeremy L Thompson     PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
303d2714514SJames Wright   }
3042b730f8bSJeremy L Thompson   PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
30588626eedSJames Wright 
306841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
3072b730f8bSJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
30888626eedSJames Wright 
309ba6664aeSJames Wright   blasius_ctx->weakT         = weakT;
310fb455ff0SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
311fb455ff0SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
3122518f336SLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
313841e4c73SJed Brown   blasius_ctx->delta0        = delta0;
314841e4c73SJed Brown   blasius_ctx->P0            = P0;
3152518f336SLeila Ghaffari   blasius_ctx->n_cheb        = N;
31630e9fa81SJames Wright   newtonian_ig_ctx->P0       = P0;
317841e4c73SJed Brown   blasius_ctx->implicit      = user->phys->implicit;
318841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
319ba6664aeSJames Wright 
320f1122ed0SJames Wright   {
3212518f336SLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3222b730f8bSJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
323f1122ed0SJames Wright     blasius_ctx->x_inflow = domain_min[0];
3242518f336SLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
325f1122ed0SJames Wright   }
3263ce1835eSJames Wright   if (!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
327f1122ed0SJames Wright 
3282b730f8bSJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
32988626eedSJames Wright 
330841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
3312b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx);
3322b730f8bSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc);
33388626eedSJames Wright 
334b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
335841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
336ba6664aeSJames Wright   if (use_stg) {
3372b730f8bSJeremy L Thompson     PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes));
33865dd5cafSJames Wright   } else {
33965dd5cafSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
34065dd5cafSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
341fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
342fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
3432b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context);
3442b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context);
345ba6664aeSJames Wright   }
3462b730f8bSJeremy L Thompson   PetscCall(PetscFree(mesh_ynodes));
34788626eedSJames Wright   PetscFunctionReturn(0);
34888626eedSJames Wright }
349