xref: /honee/problems/blasius.c (revision 14bd2a07a64d3eb0c7fc395a41fb4052a9a41171)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3bb8a0c61SJames Wright 
4bb8a0c61SJames Wright /// @file
5bb8a0c61SJames Wright /// Utility functions for setting up Blasius Boundary Layer
6bb8a0c61SJames Wright 
7bb8a0c61SJames Wright #include "../qfunctions/blasius.h"
82b916ea7SJeremy L Thompson 
9e419654dSJeremy L Thompson #include <ceed.h>
10e419654dSJeremy L Thompson #include <petscdm.h>
11e419654dSJeremy L Thompson #include <petscdt.h>
12e419654dSJeremy L Thompson 
138d78d7c8SJames Wright #include <differential_filter.h>
14149fb536SJames Wright #include <navierstokes.h>
15493642f1SJames Wright #include "stg_shur14.h"
16bb8a0c61SJames Wright 
17e0d1a4dfSLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
18e0d1a4dfSLeila Ghaffari   const BlasiusContext  blasius = (BlasiusContext)ctx;
19e0d1a4dfSLeila Ghaffari   const PetscScalar    *Tf, *Th;  // Chebyshev coefficients
20e0d1a4dfSLeila Ghaffari   PetscScalar          *r, f[4], h[4];
21e0d1a4dfSLeila Ghaffari   PetscInt              N       = blasius->n_cheb;
22fcb2c22aSJames Wright   State                 S_infty = blasius->S_infty;
2364667825SJames Wright   CeedScalar            U_infty = Norm3(S_infty.Y.velocity);
24cde3d787SJames Wright   NewtonianIGProperties gas     = blasius->newt_ctx.gas;
2506f41313SJames Wright 
2606f41313SJames Wright   PetscFunctionBeginUser;
27cde3d787SJames Wright   PetscScalar Ma = Mach(gas, S_infty.Y.temperature, U_infty), Pr = Prandtl(gas), gamma = HeatCapacityRatio(gas);
2806f41313SJames Wright 
29e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(X, &Tf));
30e0d1a4dfSLeila Ghaffari   Th = Tf + N;
31e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArray(R, &r));
32e0d1a4dfSLeila Ghaffari 
33e0d1a4dfSLeila Ghaffari   // Left boundary conditions f = f' = 0
34e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
35e0d1a4dfSLeila Ghaffari   r[0] = f[0];
36e0d1a4dfSLeila Ghaffari   r[1] = f[1];
37e0d1a4dfSLeila Ghaffari 
38e0d1a4dfSLeila Ghaffari   // f - right end boundary condition
39e0d1a4dfSLeila Ghaffari   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
40e0d1a4dfSLeila Ghaffari   r[2] = f[1] - 1.;
41e0d1a4dfSLeila Ghaffari 
42e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 3; i++) {
43e0d1a4dfSLeila Ghaffari     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
44e0d1a4dfSLeila Ghaffari     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
4504e40bb6SJeremy L Thompson     // mu and rho generally depend on h.
4604e40bb6SJeremy L Thompson     // We naively assume constant mu.
470d850f2eSLeila Ghaffari     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
480d850f2eSLeila Ghaffari     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
490d850f2eSLeila Ghaffari     const PetscScalar mu_tilde[2]     = {1, 0};
500d850f2eSLeila Ghaffari     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
510d850f2eSLeila Ghaffari     const PetscScalar mu_rho_tilde[2] = {
520d850f2eSLeila Ghaffari         mu_tilde[0] * rho_tilde[0],
530d850f2eSLeila Ghaffari         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
540d850f2eSLeila Ghaffari     };
550d850f2eSLeila Ghaffari     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
562b916ea7SJeremy 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]);
57e0d1a4dfSLeila Ghaffari   }
58e0d1a4dfSLeila Ghaffari 
59e0d1a4dfSLeila Ghaffari   // h - left end boundary condition
60e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
61fcb2c22aSJames Wright   r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature;
62e0d1a4dfSLeila Ghaffari 
63e0d1a4dfSLeila Ghaffari   // h - right end boundary condition
64e0d1a4dfSLeila Ghaffari   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
65e0d1a4dfSLeila Ghaffari   r[N + 1] = h[0] - 1.;
66e0d1a4dfSLeila Ghaffari 
67e0d1a4dfSLeila Ghaffari   // Restore vectors
68e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArrayRead(X, &Tf));
69e0d1a4dfSLeila Ghaffari   PetscCall(VecRestoreArray(R, &r));
70d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
71e0d1a4dfSLeila Ghaffari }
72e0d1a4dfSLeila Ghaffari 
73e0d1a4dfSLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
74e0d1a4dfSLeila Ghaffari   SNES                snes;
75e0d1a4dfSLeila Ghaffari   Vec                 sol, res;
76e0d1a4dfSLeila Ghaffari   PetscReal          *w;
77e0d1a4dfSLeila Ghaffari   PetscInt            N = blasius->n_cheb;
780d850f2eSLeila Ghaffari   SNESConvergedReason reason;
79e0d1a4dfSLeila Ghaffari   const PetscScalar  *cheb_coefs;
800d850f2eSLeila Ghaffari 
8106f41313SJames Wright   PetscFunctionBeginUser;
820d850f2eSLeila Ghaffari   // Allocate memory
83e0d1a4dfSLeila Ghaffari   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
84e0d1a4dfSLeila Ghaffari   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
850d850f2eSLeila Ghaffari 
860d850f2eSLeila Ghaffari   // Snes solve
87e0d1a4dfSLeila Ghaffari   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
88e0d1a4dfSLeila Ghaffari   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
89e0d1a4dfSLeila Ghaffari   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
90e0d1a4dfSLeila Ghaffari   PetscCall(VecSetFromOptions(sol));
910d850f2eSLeila Ghaffari   // Constant relative enthalpy 1 as initial guess
920d850f2eSLeila Ghaffari   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
93e0d1a4dfSLeila Ghaffari   PetscCall(VecDuplicate(sol, &res));
94e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
950d850f2eSLeila Ghaffari   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
96e0d1a4dfSLeila Ghaffari   PetscCall(SNESSetFromOptions(snes));
97e0d1a4dfSLeila Ghaffari   PetscCall(SNESSolve(snes, NULL, sol));
980d850f2eSLeila Ghaffari   PetscCall(SNESGetConvergedReason(snes, &reason));
99*14bd2a07SJames Wright   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.");
1000d850f2eSLeila Ghaffari 
1010d850f2eSLeila Ghaffari   // Assign Chebyshev coefficients
102e0d1a4dfSLeila Ghaffari   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
103e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
104e0d1a4dfSLeila Ghaffari   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
1050d850f2eSLeila Ghaffari 
1060d850f2eSLeila Ghaffari   // Destroy objects
107e0d1a4dfSLeila Ghaffari   PetscCall(PetscFree2(blasius->X, w));
108e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&sol));
109e0d1a4dfSLeila Ghaffari   PetscCall(VecDestroy(&res));
110e0d1a4dfSLeila Ghaffari   PetscCall(SNESDestroy(&snes));
111d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
112e0d1a4dfSLeila Ghaffari }
113e0d1a4dfSLeila Ghaffari 
1142b916ea7SJeremy L Thompson static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
11587d3884fSJeremy L Thompson   int            ndims, dims[2] = {0.};
116f362fe62SJames Wright   FILE          *fp;
117f362fe62SJames Wright   const PetscInt char_array_len = 512;
118f362fe62SJames Wright   char           line[char_array_len];
119f362fe62SJames Wright   PetscReal     *node_locs;
120f362fe62SJames Wright 
12106f41313SJames Wright   PetscFunctionBeginUser;
1222b916ea7SJeremy L Thompson   PetscCall(PetscFOpen(comm, path, "r", &fp));
1232b916ea7SJeremy L Thompson   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
124039caf0dSJames Wright   {
125039caf0dSJames Wright     char **array;
126f362fe62SJames Wright 
127039caf0dSJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
128f362fe62SJames Wright     for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
129039caf0dSJames Wright     PetscCall(PetscStrToArrayDestroy(ndims, array));
130039caf0dSJames Wright   }
131f362fe62SJames Wright   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
132f362fe62SJames Wright   *nynodes = dims[0];
1332b916ea7SJeremy L Thompson   PetscCall(PetscMalloc1(*nynodes, &node_locs));
134f362fe62SJames Wright 
135f362fe62SJames Wright   for (PetscInt i = 0; i < dims[0]; i++) {
136039caf0dSJames Wright     char **array;
137039caf0dSJames Wright 
1382b916ea7SJeremy L Thompson     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
1392b916ea7SJeremy L Thompson     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1405d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
141defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]);
142f362fe62SJames Wright 
143f362fe62SJames Wright     node_locs[i] = (PetscReal)atof(array[0]);
144039caf0dSJames Wright     PetscCall(PetscStrToArrayDestroy(ndims, array));
145f362fe62SJames Wright   }
1462b916ea7SJeremy L Thompson   PetscCall(PetscFClose(comm, fp));
147f362fe62SJames Wright   *pynodes = node_locs;
148d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
149f362fe62SJames Wright }
150f362fe62SJames Wright 
151bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
152bb8a0c61SJames Wright  *
15304e40bb6SJeremy L Thompson  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
154493642f1SJames Wright  * linearly in logspace to the top surface.
155bb8a0c61SJames Wright  *
15604e40bb6SJeremy L Thompson  * The top surface is also angled downwards, so that it may be used as an outflow.
15704e40bb6SJeremy L Thompson  * It's angle is controlled by `top_angle` (in units of degrees).
158f362fe62SJames Wright  *
15904e40bb6SJeremy L Thompson  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
16004e40bb6SJeremy L Thompson  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
161bb8a0c61SJames Wright  */
16266d54740SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
1638e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
16466d54740SJames Wright   PetscInt     narr, ncoords, dim;
165bb8a0c61SJames Wright   PetscReal    domain_min[3], domain_max[3], domain_size[3];
166bb8a0c61SJames Wright   PetscScalar *arr_coords;
167bb8a0c61SJames Wright   Vec          vec_coords;
16806f41313SJames Wright 
169bb8a0c61SJames Wright   PetscFunctionBeginUser;
17066d54740SJames Wright   PetscCall(DMGetDimension(dm, &dim));
171bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
172bb8a0c61SJames Wright   // Get domain boundary information
1732b916ea7SJeremy L Thompson   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
174493642f1SJames Wright   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
175bb8a0c61SJames Wright 
176bb8a0c61SJames Wright   // Get coords array from DM
1772b916ea7SJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
1782b916ea7SJeremy L Thompson   PetscCall(VecGetLocalSize(vec_coords, &narr));
1792b916ea7SJeremy L Thompson   PetscCall(VecGetArray(vec_coords, &arr_coords));
180bb8a0c61SJames Wright 
181bb8a0c61SJames Wright   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
182bb8a0c61SJames Wright   ncoords                   = narr / dim;
183bb8a0c61SJames Wright 
184bb8a0c61SJames Wright   // Get mesh information
185bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
1862b916ea7SJeremy L Thompson   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
187f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
188f362fe62SJames Wright   const PetscReal dybox = domain_size[1] / faces[1];
189bb8a0c61SJames Wright 
1908e63b9cbSJames Wright   if (!*node_locs) {
191bb8a0c61SJames Wright     // Calculate the first element height
192bb8a0c61SJames Wright     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
193bb8a0c61SJames Wright 
194bb8a0c61SJames Wright     // Calculate log of sizing outside BL
195bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
196bb8a0c61SJames Wright 
1978e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1988e63b9cbSJames Wright     PetscReal *temp_node_locs;
1992b916ea7SJeremy L Thompson     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
2008e63b9cbSJames Wright 
201493642f1SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
202bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
203bb8a0c61SJames Wright       if (y_box_index <= N) {
2042b916ea7SJeremy L Thompson         coords[i][1] =
2052b916ea7SJeremy L Thompson             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
206bb8a0c61SJames Wright       } else {
207bb8a0c61SJames Wright         PetscInt j   = y_box_index - N;
2082b916ea7SJeremy L Thompson         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
209f362fe62SJames Wright       }
2102b916ea7SJeremy L Thompson       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
211f362fe62SJames Wright     }
2128e63b9cbSJames Wright 
2138e63b9cbSJames Wright     *node_locs = temp_node_locs;
214f362fe62SJames Wright   } else {
2155d28dccaSJames Wright     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
216defe8520SJames Wright                "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given",
217defe8520SJames Wright                faces[1] + 1, *num_node_locs);
2188e63b9cbSJames Wright     if (*num_node_locs > faces[1] + 1) {
2192b916ea7SJeremy L Thompson       PetscCall(PetscPrintf(comm,
220defe8520SJames Wright                             "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") "
221defe8520SJames Wright                             "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n",
2222b916ea7SJeremy L Thompson                             *num_node_locs, faces[1] + 1));
223f362fe62SJames Wright     }
2248e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
225f362fe62SJames Wright 
226f362fe62SJames Wright     for (PetscInt i = 0; i < ncoords; i++) {
227f362fe62SJames Wright       // Determine which y-node we're at
228f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1] / dybox);
2292b916ea7SJeremy L Thompson       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
230bb8a0c61SJames Wright     }
231bb8a0c61SJames Wright   }
232bb8a0c61SJames Wright 
2332b916ea7SJeremy L Thompson   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
2342b916ea7SJeremy L Thompson   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
235d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
236bb8a0c61SJames Wright }
237bb8a0c61SJames Wright 
238d6cac220SJames Wright static PetscErrorCode BlasiusInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
239d6cac220SJames Wright   HoneeBCStruct honee_bc;
240d6cac220SJames Wright 
241d6cac220SJames Wright   PetscFunctionBeginUser;
242d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
243d6cac220SJames Wright   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Blasius_Inflow, Blasius_Inflow_loc, honee_bc->qfctx, qf));
244d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
245d6cac220SJames Wright }
246d6cac220SJames Wright 
247d6cac220SJames Wright static PetscErrorCode BlasiusInflowBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
248d6cac220SJames Wright   HoneeBCStruct honee_bc;
249d6cac220SJames Wright 
250d6cac220SJames Wright   PetscFunctionBeginUser;
251d6cac220SJames Wright   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
252d6cac220SJames Wright   PetscCall(HoneeBCCreateIJacobianQF(bc_def, Blasius_Inflow_Jacobian, Blasius_Inflow_Jacobian_loc, honee_bc->qfctx, qf));
253d6cac220SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
254d6cac220SJames Wright }
255d6cac220SJames Wright 
256d3c60affSJames Wright PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx) {
2570c373b74SJames Wright   Honee                    honee   = *(Honee *)ctx;
2580c373b74SJames Wright   MPI_Comm                 comm    = honee->comm;
2590c373b74SJames Wright   Ceed                     ceed    = honee->ceed;
260493642f1SJames Wright   PetscBool                use_stg = PETSC_FALSE;
26115a3537eSJed Brown   BlasiusContext           blasius_ctx;
26215a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
263e07531f7SJames Wright   CeedQFunctionContext     blasius_qfctx;
26415a3537eSJed Brown 
265bb8a0c61SJames Wright   PetscFunctionBeginUser;
266d3c60affSJames Wright   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx));
2672d898fa6SJames Wright   PetscCall(PetscNew(&blasius_ctx));
268bb8a0c61SJames Wright 
269bb8a0c61SJames Wright   // ------------------------------------------------------
270bb8a0c61SJames Wright   //               SET UP Blasius
271bb8a0c61SJames Wright   // ------------------------------------------------------
272f5dc303cSJames Wright   problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsBlasius, .qf_loc = ICsBlasius_loc, .qfctx = problem->ics.qfctx};
273bb8a0c61SJames Wright 
274aef1eb53SLeila Ghaffari   CeedScalar U_inf                                = 40;           // m/s
275aef1eb53SLeila Ghaffari   CeedScalar T_inf                                = 288.;         // K
2760d850f2eSLeila Ghaffari   CeedScalar T_wall                               = 288.;         // K
277e0d1a4dfSLeila Ghaffari   CeedScalar delta0                               = 4.2e-3;       // m
278fcb2c22aSJames Wright   CeedScalar P_inf                                = 1.01e5;       // Pa
279defe8520SJames Wright   PetscInt   N                                    = 20;           // Number of Chebyshev terms
2802acc7cbcSKenneth E. Jansen   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
2817470235eSJames Wright   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
2827470235eSJames Wright   PetscReal  mesh_growth                          = 1.08;         // [-]
2837470235eSJames Wright   PetscInt   mesh_Ndelta                          = 45;           // [-]
2847470235eSJames Wright   PetscReal  mesh_top_angle                       = 5;            // degrees
285f362fe62SJames Wright   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
286a213b8aaSJames Wright   PetscBool  P0_set;
287bb8a0c61SJames Wright 
2883fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
2892b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
2902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
2912b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
292a213b8aaSJames Wright   PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set));  // For maintaining behavior of -P0 flag (which is deprecated)
293a213b8aaSJames Wright   PetscCall(
294a213b8aaSJames Wright       PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0",
295a213b8aaSJames Wright                              "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure"));
296a213b8aaSJames Wright   PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL));
2972b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
2982b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
2992b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
3002b916ea7SJeremy 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,
3010d850f2eSLeila Ghaffari              BLASIUS_MAX_N_CHEBYSHEV);
3020c373b74SJames Wright   if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
3032b916ea7SJeremy L Thompson     PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
304f31f4833SJames Wright     PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height,
305f31f4833SJames Wright                                  &mesh_refine_height, NULL));
3062b916ea7SJeremy L Thompson     PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
3072b916ea7SJeremy L Thompson     PetscCall(
3082b916ea7SJeremy L Thompson         PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
3092b916ea7SJeremy L Thompson     PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
310f362fe62SJames Wright                                  "Path to file with y node locations. "
3112b916ea7SJeremy L Thompson                                  "If empty, will use the algorithmic mesh warping.",
3122b916ea7SJeremy L Thompson                                  NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
313f31f4833SJames Wright   }
3142b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
315bb8a0c61SJames Wright   PetscOptionsEnd();
316bb8a0c61SJames Wright 
317c9f37605SMohammed Amin   Units units = honee->units;
318bb8a0c61SJames Wright 
319c9f37605SMohammed Amin   T_inf *= units->Kelvin;
320c9f37605SMohammed Amin   T_wall *= units->Kelvin;
321c9f37605SMohammed Amin   P_inf *= units->Pascal;
322c9f37605SMohammed Amin   U_inf *= units->meter / units->second;
323c9f37605SMohammed Amin   delta0 *= units->meter;
324bb8a0c61SJames Wright 
3250c373b74SJames Wright   if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
326f362fe62SJames Wright     PetscReal *mesh_ynodes  = NULL;
327f362fe62SJames Wright     PetscInt   mesh_nynodes = 0;
328c029f0c5SJames Wright     if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
32966d54740SJames Wright     PetscCall(ModifyMesh(comm, dm, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
3309ef62cddSJames Wright     PetscCall(PetscFree(mesh_ynodes));
331c029f0c5SJames Wright   }
332bb8a0c61SJames Wright 
33315a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
334e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
335bb8a0c61SJames Wright 
336fcb2c22aSJames Wright   StatePrimitive Y_inf = {
337fcb2c22aSJames Wright       .pressure = P_inf, .velocity = {U_inf, 0, 0},
338fcb2c22aSJames Wright            .temperature = T_inf
339fcb2c22aSJames Wright   };
340cde3d787SJames Wright   State S_infty = StateFromPrimitive(newtonian_ig_ctx->gas, Y_inf);
341fcb2c22aSJames Wright 
342493642f1SJames Wright   blasius_ctx->weakT    = weakT;
343e0d1a4dfSLeila Ghaffari   blasius_ctx->T_wall   = T_wall;
34415a3537eSJed Brown   blasius_ctx->delta0   = delta0;
345fcb2c22aSJames Wright   blasius_ctx->S_infty  = S_infty;
346e0d1a4dfSLeila Ghaffari   blasius_ctx->n_cheb   = N;
3470c373b74SJames Wright   blasius_ctx->implicit = honee->phys->implicit;
348a213b8aaSJames Wright   if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf;  // For maintaining behavior of -P0 flag (which is deprecated)
349cde3d787SJames Wright   blasius_ctx->newt_ctx = *newtonian_ig_ctx;
350493642f1SJames Wright 
351ef2c71fdSJames Wright   {
352e0d1a4dfSLeila Ghaffari     PetscReal domain_min[3], domain_max[3];
3532b916ea7SJeremy L Thompson     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
354ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
355e0d1a4dfSLeila Ghaffari     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
356ef2c71fdSJames Wright   }
35725c92e8fSJames Wright   PetscBool diff_filter_mms = PETSC_FALSE;
35825c92e8fSJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
35925c92e8fSJames Wright   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
360ef2c71fdSJames Wright 
361e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
362bb8a0c61SJames Wright 
3630c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &blasius_qfctx));
364e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
365e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_qfctx, CEED_MEM_HOST, FreeContextPetsc));
366bb8a0c61SJames Wright 
367e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
368e07531f7SJames Wright   problem->ics.qfctx = blasius_qfctx;
369493642f1SJames Wright   if (use_stg) {
3700c373b74SJames Wright     PetscCall(SetupStg(comm, dm, problem, honee, weakT, S_infty.Y.temperature, S_infty.Y.pressure));
37125c92e8fSJames Wright   } else if (diff_filter_mms) {
3728d78d7c8SJames Wright     PetscCall(DifferentialFilterMmsICSetup(honee));
3738085925cSJames Wright   } else {
3740c373b74SJames Wright     PetscCheck((honee->phys->state_var == STATEVAR_CONSERVATIVE) || (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER), honee->comm,
375727ec889SJames Wright                PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow");
376d6cac220SJames Wright     for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
377d6cac220SJames Wright       BCDefinition bc_def = problem->bc_defs[b];
378d6cac220SJames Wright       const char  *name;
379d6cac220SJames Wright 
380d6cac220SJames Wright       PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
381d6cac220SJames Wright       if (!strcmp(name, "inflow")) {
382d6cac220SJames Wright         HoneeBCStruct honee_bc;
383d6cac220SJames Wright 
384d6cac220SJames Wright         PetscCall(PetscNew(&honee_bc));
385d6cac220SJames Wright         PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_qfctx, &honee_bc->qfctx));
386d6cac220SJames Wright         honee_bc->honee              = honee;
3871abc2837SJames Wright         honee_bc->num_comps_jac_data = 0;
388d6cac220SJames Wright         PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
389d6cac220SJames Wright 
390d6cac220SJames Wright         PetscCall(BCDefinitionSetIFunction(bc_def, BlasiusInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
391d6cac220SJames Wright         PetscCall(BCDefinitionSetIJacobian(bc_def, BlasiusInflowBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
392d6cac220SJames Wright       }
393d6cac220SJames Wright     }
394493642f1SJames Wright   }
395d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
396bb8a0c61SJames Wright }
397