xref: /libCEED/examples/fluids/problems/blasius.c (revision 5aed82e4fa97acf4ba24a7f10a35f5303a6798e0)
1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson #include <petscdt.h>
1649aac155SJeremy L Thompson 
172b730f8bSJeremy L Thompson #include "../navierstokes.h"
18ba6664aeSJames Wright #include "stg_shur14.h"
1988626eedSJames Wright 
202518f336SLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
212518f336SLeila Ghaffari   const BlasiusContext blasius = (BlasiusContext)ctx;
222518f336SLeila Ghaffari   const PetscScalar   *Tf, *Th;  // Chebyshev coefficients
232518f336SLeila Ghaffari   PetscScalar         *r, f[4], h[4];
242518f336SLeila Ghaffari   PetscInt             N = blasius->n_cheb;
25f17d818dSJames Wright 
26f17d818dSJames Wright   PetscFunctionBeginUser;
272b730f8bSJeremy L Thompson   PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx),
282518f336SLeila Ghaffari               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
29f17d818dSJames Wright 
302518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
312518f336SLeila Ghaffari   Th = Tf + N;
322518f336SLeila Ghaffari   PetscCall(VecGetArray(R, &r));
332518f336SLeila Ghaffari 
342518f336SLeila Ghaffari   // Left boundary conditions f = f' = 0
352518f336SLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
362518f336SLeila Ghaffari   r[0] = f[0];
372518f336SLeila Ghaffari   r[1] = f[1];
382518f336SLeila Ghaffari 
392518f336SLeila Ghaffari   // f - right end boundary condition
402518f336SLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
412518f336SLeila Ghaffari   r[2] = f[1] - 1.;
422518f336SLeila Ghaffari 
432518f336SLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
442518f336SLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
452518f336SLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
46ea61e9acSJeremy L Thompson     // mu and rho generally depend on h.
47ea61e9acSJeremy L Thompson     // We naively assume constant mu.
4807d14e58SLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
4907d14e58SLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
5007d14e58SLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
5107d14e58SLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
5207d14e58SLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
5307d14e58SLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
5407d14e58SLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
5507d14e58SLeila Ghaffari     };
5607d14e58SLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
572b730f8bSJeremy 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]);
582518f336SLeila Ghaffari   }
592518f336SLeila Ghaffari 
602518f336SLeila Ghaffari   // h - left end boundary condition
612518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
62fb455ff0SLeila Ghaffari   r[N] = h[0] - blasius->T_wall / blasius->T_inf;
632518f336SLeila Ghaffari 
642518f336SLeila Ghaffari   // h - right end boundary condition
652518f336SLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
662518f336SLeila Ghaffari   r[N + 1] = h[0] - 1.;
672518f336SLeila Ghaffari 
682518f336SLeila Ghaffari   // Restore vectors
692518f336SLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
702518f336SLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
71ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
722518f336SLeila Ghaffari }
732518f336SLeila Ghaffari 
742518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
752518f336SLeila Ghaffari   SNES                snes;
762518f336SLeila Ghaffari   Vec                 sol, res;
772518f336SLeila Ghaffari   PetscReal          *w;
782518f336SLeila Ghaffari   PetscInt            N = blasius->n_cheb;
7907d14e58SLeila Ghaffari   SNESConvergedReason reason;
802518f336SLeila Ghaffari   const PetscScalar  *cheb_coefs;
8107d14e58SLeila Ghaffari 
82f17d818dSJames Wright   PetscFunctionBeginUser;
8307d14e58SLeila Ghaffari   // Allocate memory
842518f336SLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
852518f336SLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
8607d14e58SLeila Ghaffari 
8707d14e58SLeila Ghaffari   // Snes solve
882518f336SLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
892518f336SLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
902518f336SLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
912518f336SLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
9207d14e58SLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
9307d14e58SLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
942518f336SLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
952518f336SLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
9607d14e58SLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
972518f336SLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
982518f336SLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
9907d14e58SLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
1000e654f56SJames Wright   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
10107d14e58SLeila Ghaffari 
10207d14e58SLeila Ghaffari   // Assign Chebyshev coefficients
1032518f336SLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
1042518f336SLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
1052518f336SLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
10607d14e58SLeila Ghaffari 
10707d14e58SLeila Ghaffari   // Destroy objects
1082518f336SLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
1092518f336SLeila Ghaffari   PetscCall(VecDestroy(&sol));
1102518f336SLeila Ghaffari   PetscCall(VecDestroy(&res));
1112518f336SLeila Ghaffari   PetscCall(SNESDestroy(&snes));
112ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1132518f336SLeila Ghaffari }
1142518f336SLeila Ghaffari 
1152b730f8bSJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
116457e73b2SJames Wright   int            ndims, dims[2];
117d2714514SJames Wright   FILE          *fp;
118d2714514SJames Wright   const PetscInt char_array_len = 512;
119d2714514SJames Wright   char           line[char_array_len];
120d2714514SJames Wright   char         **array;
121d2714514SJames Wright   PetscReal     *node_locs;
122d2714514SJames Wright 
123f17d818dSJames Wright   PetscFunctionBeginUser;
1242b730f8bSJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1252b730f8bSJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1262b730f8bSJeremy L Thompson   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
127d2714514SJames Wright 
128d2714514SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
129d2714514SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
130d2714514SJames Wright   *nynodes = dims[0];
1312b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
132d2714514SJames Wright 
133d2714514SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
1342b730f8bSJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1352b730f8bSJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1360e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
137457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]);
138d2714514SJames Wright 
139d2714514SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
140d2714514SJames Wright   }
1412b730f8bSJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
142d2714514SJames Wright   *pynodes = node_locs;
143ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
144d2714514SJames Wright }
145d2714514SJames Wright 
14688626eedSJames Wright /* \brief Modify the domain and mesh for blasius
14788626eedSJames Wright  *
148ea61e9acSJeremy L Thompson  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
149ba6664aeSJames Wright  * linearly in logspace to the top surface.
15088626eedSJames Wright  *
151ea61e9acSJeremy L Thompson  * The top surface is also angled downwards, so that it may be used as an outflow.
152ea61e9acSJeremy L Thompson  * It's angle is controlled by `top_angle` (in units of degrees).
153d2714514SJames Wright  *
154ea61e9acSJeremy L Thompson  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
155ea61e9acSJeremy L Thompson  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
15688626eedSJames Wright  */
1572b730f8bSJeremy L Thompson static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
158798d42b8SJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
1592b730f8bSJeremy L Thompson   PetscInt     narr, ncoords;
16088626eedSJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
16188626eedSJames Wright   PetscScalar *arr_coords;
16288626eedSJames Wright   Vec          vec_coords;
163f17d818dSJames Wright 
16488626eedSJames Wright   PetscFunctionBeginUser;
16588626eedSJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
16688626eedSJames Wright   // Get domain boundary information
1672b730f8bSJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
168ba6664aeSJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
16988626eedSJames Wright 
17088626eedSJames Wright   // Get coords array from DM
1712b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1722b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1732b730f8bSJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
17488626eedSJames Wright 
17588626eedSJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
17688626eedSJames Wright   ncoords                   = narr / dim;
17788626eedSJames Wright 
17888626eedSJames Wright   // Get mesh information
17988626eedSJames Wright   PetscInt nmax = 3, faces[3];
1802b730f8bSJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
181d2714514SJames Wright   // Get element size of the box mesh, for indexing each node
182d2714514SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
18388626eedSJames Wright 
184798d42b8SJames Wright   if (!*node_locs) {
18588626eedSJames Wright     // Calculate the first element height
18688626eedSJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
18788626eedSJames Wright 
18888626eedSJames Wright     // Calculate log of sizing outside BL
18988626eedSJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
19088626eedSJames Wright 
191798d42b8SJames Wright     *num_node_locs = faces[1] + 1;
192798d42b8SJames Wright     PetscReal *temp_node_locs;
1932b730f8bSJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
194798d42b8SJames Wright 
195ba6664aeSJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
19688626eedSJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
19788626eedSJames Wright       if (y_box_index <= N) {
1982b730f8bSJeremy L Thompson         coords[i][1] =
1992b730f8bSJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
20088626eedSJames Wright       } else {
20188626eedSJames Wright         PetscInt j   = y_box_index - N;
2022b730f8bSJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
203d2714514SJames Wright       }
2042b730f8bSJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
205d2714514SJames Wright     }
206798d42b8SJames Wright 
207798d42b8SJames Wright     *node_locs = temp_node_locs;
208d2714514SJames Wright   } else {
2090e654f56SJames Wright     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
210457e73b2SJames Wright                "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given",
211457e73b2SJames Wright                faces[1] + 1, *num_node_locs);
212798d42b8SJames Wright     if (*num_node_locs > faces[1] + 1) {
2132b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
214457e73b2SJames Wright                             "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") "
215457e73b2SJames Wright                             "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n",
2162b730f8bSJeremy L Thompson                             *num_node_locs, faces[1] + 1));
217d2714514SJames Wright     }
218798d42b8SJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
219d2714514SJames Wright 
220d2714514SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
221d2714514SJames Wright       // Determine which y-node we're at
222d2714514SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2232b730f8bSJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
22488626eedSJames Wright     }
22588626eedSJames Wright   }
22688626eedSJames Wright 
2272b730f8bSJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2282b730f8bSJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
229ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
23088626eedSJames Wright }
23188626eedSJames Wright 
23246de7363SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
23388626eedSJames Wright   User                     user    = *(User *)ctx;
2340a37e1beSJames Wright   MPI_Comm                 comm    = user->comm;
235a424bcd0SJames Wright   Ceed                     ceed    = user->ceed;
236ba6664aeSJames Wright   PetscBool                use_stg = PETSC_FALSE;
237841e4c73SJed Brown   BlasiusContext           blasius_ctx;
238841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
239841e4c73SJed Brown   CeedQFunctionContext     blasius_context;
240841e4c73SJed Brown 
24188626eedSJames Wright   PetscFunctionBeginUser;
24246de7363SJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
2432b730f8bSJeremy L Thompson   PetscCall(PetscCalloc1(1, &blasius_ctx));
24488626eedSJames Wright 
24588626eedSJames Wright   // ------------------------------------------------------
24688626eedSJames Wright   //               SET UP Blasius
24788626eedSJames Wright   // ------------------------------------------------------
24891e5af17SJed Brown   problem->ics.qfunction     = ICsBlasius;
24991e5af17SJed Brown   problem->ics.qfunction_loc = ICsBlasius_loc;
25088626eedSJames Wright 
251fb455ff0SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
252fb455ff0SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
25307d14e58SLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
2542518f336SLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
25588626eedSJames Wright   CeedScalar P0                                   = 1.01e5;       // Pa
256457e73b2SJames Wright   PetscInt   N                                    = 20;           // Number of Chebyshev terms
257871db79fSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2587b8d3891SJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2597b8d3891SJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2607b8d3891SJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2617b8d3891SJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
262d2714514SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
26388626eedSJames Wright 
2647b37518fSJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2652b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2662b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2672b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
2682b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2692b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2702b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
2712b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
2722b730f8bSJeremy 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,
27307d14e58SLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
2742526956eSJames Wright   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
2752b730f8bSJeremy L Thompson     PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
2762526956eSJames Wright     PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height,
2772526956eSJames Wright                                  &mesh_refine_height, 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));
2852526956eSJames Wright   }
2862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
28788626eedSJames Wright   PetscOptionsEnd();
28888626eedSJames Wright 
28988626eedSJames Wright   PetscScalar meter  = user->units->meter;
29088626eedSJames Wright   PetscScalar second = user->units->second;
29188626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
29288626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
29388626eedSJames Wright 
294fb455ff0SLeila Ghaffari   T_inf *= Kelvin;
2952518f336SLeila Ghaffari   T_wall *= Kelvin;
29688626eedSJames Wright   P0 *= Pascal;
297fb455ff0SLeila Ghaffari   U_inf *= meter / second;
29888626eedSJames Wright   delta0 *= meter;
29988626eedSJames Wright 
3002526956eSJames Wright   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
301d2714514SJames Wright     PetscReal *mesh_ynodes  = NULL;
302d2714514SJames Wright     PetscInt   mesh_nynodes = 0;
3039309e21cSJames Wright     if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
3042b730f8bSJeremy L Thompson     PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
305f8839eb4SJames Wright     PetscCall(PetscFree(mesh_ynodes));
3069309e21cSJames Wright   }
30788626eedSJames Wright 
308841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
309a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
31088626eedSJames Wright 
311ba6664aeSJames Wright   blasius_ctx->weakT         = weakT;
312fb455ff0SLeila Ghaffari   blasius_ctx->U_inf         = U_inf;
313fb455ff0SLeila Ghaffari   blasius_ctx->T_inf         = T_inf;
3142518f336SLeila Ghaffari   blasius_ctx->T_wall        = T_wall;
315841e4c73SJed Brown   blasius_ctx->delta0        = delta0;
316841e4c73SJed Brown   blasius_ctx->P0            = P0;
3172518f336SLeila Ghaffari   blasius_ctx->n_cheb        = N;
31830e9fa81SJames Wright   newtonian_ig_ctx->P0       = P0;
319841e4c73SJed Brown   blasius_ctx->implicit      = user->phys->implicit;
320841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
321ba6664aeSJames Wright 
322f1122ed0SJames Wright   {
3232518f336SLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3242b730f8bSJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
325f1122ed0SJames Wright     blasius_ctx->x_inflow = domain_min[0];
3262518f336SLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
327f1122ed0SJames Wright   }
3282029a677SJames Wright   PetscBool diff_filter_mms = PETSC_FALSE;
3292029a677SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
3302029a677SJames Wright   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
331f1122ed0SJames Wright 
332a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
33388626eedSJames Wright 
334a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &blasius_context));
335a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
336a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc));
33788626eedSJames Wright 
338a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
339841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
340ba6664aeSJames Wright   if (use_stg) {
341cbef7084SJames Wright     PetscCall(SetupStg(comm, dm, problem, user, weakT, T_inf, P0));
3422029a677SJames Wright   } else if (diff_filter_mms) {
343cbef7084SJames Wright     PetscCall(DifferentialFilterMmsICSetup(problem));
34465dd5cafSJames Wright   } else {
34565dd5cafSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
34665dd5cafSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
347fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
348fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
349a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context));
350a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context));
351ba6664aeSJames Wright   }
352ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
35388626eedSJames Wright }
354