xref: /honee/problems/blasius.c (revision ae2b091fac884a554e48acc4b4c187524c2a2818)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 /// @file
5 /// Utility functions for setting up Blasius Boundary Layer
6 
7 #include "../qfunctions/blasius.h"
8 
9 #include <ceed.h>
10 #include <petscdm.h>
11 #include <petscdt.h>
12 
13 #include <navierstokes.h>
14 #include "stg_shur14.h"
15 
16 PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) {
17   const BlasiusContext blasius = (BlasiusContext)ctx;
18   const PetscScalar   *Tf, *Th;  // Chebyshev coefficients
19   PetscScalar         *r, f[4], h[4];
20   PetscInt             N       = blasius->n_cheb;
21   State                S_infty = blasius->S_infty;
22   CeedScalar           U_infty = Norm3(S_infty.Y.velocity);
23 
24   PetscFunctionBeginUser;
25   PetscScalar Ma = Mach(&blasius->newtonian_ctx, S_infty.Y.temperature, U_infty), Pr = Prandtl(&blasius->newtonian_ctx),
26               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
27 
28   PetscCall(VecGetArrayRead(X, &Tf));
29   Th = Tf + N;
30   PetscCall(VecGetArray(R, &r));
31 
32   // Left boundary conditions f = f' = 0
33   ChebyshevEval(N, Tf, -1., blasius->eta_max, f);
34   r[0] = f[0];
35   r[1] = f[1];
36 
37   // f - right end boundary condition
38   ChebyshevEval(N, Tf, 1., blasius->eta_max, f);
39   r[2] = f[1] - 1.;
40 
41   for (int i = 0; i < N - 3; i++) {
42     ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f);
43     ChebyshevEval(N - 1, Th, blasius->X[i], blasius->eta_max, h);
44     // mu and rho generally depend on h.
45     // We naively assume constant mu.
46     // For an ideal gas at constant pressure, density is inversely proportional to enthalpy.
47     // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here.
48     const PetscScalar mu_tilde[2]     = {1, 0};
49     const PetscScalar rho_tilde[2]    = {1 / h[0], -h[1] / PetscSqr(h[0])};
50     const PetscScalar mu_rho_tilde[2] = {
51         mu_tilde[0] * rho_tilde[0],
52         mu_tilde[1] * rho_tilde[0] + mu_tilde[0] * rho_tilde[1],
53     };
54     r[3 + i]     = 2 * (mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0];
55     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]);
56   }
57 
58   // h - left end boundary condition
59   ChebyshevEval(N - 1, Th, -1., blasius->eta_max, h);
60   r[N] = h[0] - blasius->T_wall / S_infty.Y.temperature;
61 
62   // h - right end boundary condition
63   ChebyshevEval(N - 1, Th, 1., blasius->eta_max, h);
64   r[N + 1] = h[0] - 1.;
65 
66   // Restore vectors
67   PetscCall(VecRestoreArrayRead(X, &Tf));
68   PetscCall(VecRestoreArray(R, &r));
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) {
73   SNES                snes;
74   Vec                 sol, res;
75   PetscReal          *w;
76   PetscInt            N = blasius->n_cheb;
77   SNESConvergedReason reason;
78   const PetscScalar  *cheb_coefs;
79 
80   PetscFunctionBeginUser;
81   // Allocate memory
82   PetscCall(PetscMalloc2(N - 3, &blasius->X, N - 3, &w));
83   PetscCall(PetscDTGaussQuadrature(N - 3, -1., 1., blasius->X, w));
84 
85   // Snes solve
86   PetscCall(SNESCreate(PETSC_COMM_SELF, &snes));
87   PetscCall(VecCreate(PETSC_COMM_SELF, &sol));
88   PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2 * N - 1));
89   PetscCall(VecSetFromOptions(sol));
90   // Constant relative enthalpy 1 as initial guess
91   PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES));
92   PetscCall(VecDuplicate(sol, &res));
93   PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius));
94   PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_"));
95   PetscCall(SNESSetFromOptions(snes));
96   PetscCall(SNESSolve(snes, NULL, sol));
97   PetscCall(SNESGetConvergedReason(snes, &reason));
98   PetscCheck(reason >= 0, PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, "The Chebyshev solve failed.\n");
99 
100   // Assign Chebyshev coefficients
101   PetscCall(VecGetArrayRead(sol, &cheb_coefs));
102   for (int i = 0; i < N; i++) blasius->Tf_cheb[i] = cheb_coefs[i];
103   for (int i = 0; i < N - 1; i++) blasius->Th_cheb[i] = cheb_coefs[i + N];
104 
105   // Destroy objects
106   PetscCall(PetscFree2(blasius->X, w));
107   PetscCall(VecDestroy(&sol));
108   PetscCall(VecDestroy(&res));
109   PetscCall(SNESDestroy(&snes));
110   PetscFunctionReturn(PETSC_SUCCESS);
111 }
112 
113 static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
114   int            ndims, dims[2] = {0.};
115   FILE          *fp;
116   const PetscInt char_array_len = 512;
117   char           line[char_array_len];
118   PetscReal     *node_locs;
119 
120   PetscFunctionBeginUser;
121   PetscCall(PetscFOpen(comm, path, "r", &fp));
122   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
123   {
124     char **array;
125 
126     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
127     for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
128     PetscCall(PetscStrToArrayDestroy(ndims, array));
129   }
130   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
131   *nynodes = dims[0];
132   PetscCall(PetscMalloc1(*nynodes, &node_locs));
133 
134   for (PetscInt i = 0; i < dims[0]; i++) {
135     char **array;
136 
137     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
138     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
139     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
140                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%d instead of %d)", i, path, ndims, dims[1]);
141 
142     node_locs[i] = (PetscReal)atof(array[0]);
143     PetscCall(PetscStrToArrayDestroy(ndims, array));
144   }
145   PetscCall(PetscFClose(comm, fp));
146   *pynodes = node_locs;
147   PetscFunctionReturn(PETSC_SUCCESS);
148 }
149 
150 /* \brief Modify the domain and mesh for blasius
151  *
152  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
153  * linearly in logspace to the top surface.
154  *
155  * The top surface is also angled downwards, so that it may be used as an outflow.
156  * It's angle is controlled by `top_angle` (in units of degrees).
157  *
158  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
159  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
160  */
161 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
162                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
163   PetscInt     narr, ncoords;
164   PetscReal    domain_min[3], domain_max[3], domain_size[3];
165   PetscScalar *arr_coords;
166   Vec          vec_coords;
167 
168   PetscFunctionBeginUser;
169   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
170   // Get domain boundary information
171   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
172   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
173 
174   // Get coords array from DM
175   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
176   PetscCall(VecGetLocalSize(vec_coords, &narr));
177   PetscCall(VecGetArray(vec_coords, &arr_coords));
178 
179   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
180   ncoords                   = narr / dim;
181 
182   // Get mesh information
183   PetscInt nmax = 3, faces[3];
184   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
185   // Get element size of the box mesh, for indexing each node
186   const PetscReal dybox = domain_size[1] / faces[1];
187 
188   if (!*node_locs) {
189     // Calculate the first element height
190     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
191 
192     // Calculate log of sizing outside BL
193     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
194 
195     *num_node_locs = faces[1] + 1;
196     PetscReal *temp_node_locs;
197     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
198 
199     for (PetscInt i = 0; i < ncoords; i++) {
200       PetscInt y_box_index = round(coords[i][1] / dybox);
201       if (y_box_index <= N) {
202         coords[i][1] =
203             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
204       } else {
205         PetscInt j   = y_box_index - N;
206         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
207       }
208       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
209     }
210 
211     *node_locs = temp_node_locs;
212   } else {
213     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
214                "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given",
215                faces[1] + 1, *num_node_locs);
216     if (*num_node_locs > faces[1] + 1) {
217       PetscCall(PetscPrintf(comm,
218                             "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") "
219                             "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n",
220                             *num_node_locs, faces[1] + 1));
221     }
222     PetscScalar max_y = (*node_locs)[faces[1]];
223 
224     for (PetscInt i = 0; i < ncoords; i++) {
225       // Determine which y-node we're at
226       PetscInt y_box_index = round(coords[i][1] / dybox);
227       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
228     }
229   }
230 
231   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
232   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
233   PetscFunctionReturn(PETSC_SUCCESS);
234 }
235 
236 PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc) {
237   User                     user    = *(User *)ctx;
238   MPI_Comm                 comm    = user->comm;
239   Ceed                     ceed    = user->ceed;
240   PetscBool                use_stg = PETSC_FALSE;
241   BlasiusContext           blasius_ctx;
242   NewtonianIdealGasContext newtonian_ig_ctx;
243   CeedQFunctionContext     blasius_context;
244 
245   PetscFunctionBeginUser;
246   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
247   PetscCall(PetscCalloc1(1, &blasius_ctx));
248 
249   // ------------------------------------------------------
250   //               SET UP Blasius
251   // ------------------------------------------------------
252   problem->ics.qfunction     = ICsBlasius;
253   problem->ics.qfunction_loc = ICsBlasius_loc;
254 
255   CeedScalar U_inf                                = 40;           // m/s
256   CeedScalar T_inf                                = 288.;         // K
257   CeedScalar T_wall                               = 288.;         // K
258   CeedScalar delta0                               = 4.2e-3;       // m
259   CeedScalar P_inf                                = 1.01e5;       // Pa
260   PetscInt   N                                    = 20;           // Number of Chebyshev terms
261   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
262   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
263   PetscReal  mesh_growth                          = 1.08;         // [-]
264   PetscInt   mesh_Ndelta                          = 45;           // [-]
265   PetscReal  mesh_top_angle                       = 5;            // degrees
266   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
267   PetscBool  P0_set;
268 
269   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
270   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
271   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
272   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
273   PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set));  // For maintaining behavior of -P0 flag (which is deprecated)
274   PetscCall(
275       PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0",
276                              "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure"));
277   PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL));
278   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
279   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
280   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
281   PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N,
282              BLASIUS_MAX_N_CHEBYSHEV);
283   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
284     PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
285     PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height,
286                                  &mesh_refine_height, NULL));
287     PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
288     PetscCall(
289         PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
290     PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
291                                  "Path to file with y node locations. "
292                                  "If empty, will use the algorithmic mesh warping.",
293                                  NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
294   }
295   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
296   PetscOptionsEnd();
297 
298   PetscScalar meter  = user->units->meter;
299   PetscScalar second = user->units->second;
300   PetscScalar Kelvin = user->units->Kelvin;
301   PetscScalar Pascal = user->units->Pascal;
302 
303   T_inf *= Kelvin;
304   T_wall *= Kelvin;
305   P_inf *= Pascal;
306   U_inf *= meter / second;
307   delta0 *= meter;
308 
309   if (user->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
310     PetscReal *mesh_ynodes  = NULL;
311     PetscInt   mesh_nynodes = 0;
312     if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
313     PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
314     PetscCall(PetscFree(mesh_ynodes));
315   }
316 
317   // Some properties depend on parameters from NewtonianIdealGas
318   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
319 
320   StatePrimitive Y_inf = {
321       .pressure = P_inf, .velocity = {U_inf, 0, 0},
322            .temperature = T_inf
323   };
324   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
325 
326   blasius_ctx->weakT    = weakT;
327   blasius_ctx->T_wall   = T_wall;
328   blasius_ctx->delta0   = delta0;
329   blasius_ctx->S_infty  = S_infty;
330   blasius_ctx->n_cheb   = N;
331   blasius_ctx->implicit = user->phys->implicit;
332   if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf;  // For maintaining behavior of -P0 flag (which is deprecated)
333   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
334 
335   {
336     PetscReal domain_min[3], domain_max[3];
337     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
338     blasius_ctx->x_inflow = domain_min[0];
339     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
340   }
341   PetscBool diff_filter_mms = PETSC_FALSE;
342   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
343   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
344 
345   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
346 
347   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &blasius_context));
348   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
349   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc));
350 
351   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
352   problem->ics.qfunction_context = blasius_context;
353   if (use_stg) {
354     PetscCall(SetupStg(comm, dm, problem, user, weakT, S_infty.Y.temperature, S_infty.Y.pressure));
355   } else if (diff_filter_mms) {
356     PetscCall(DifferentialFilterMmsICSetup(problem));
357   } else {
358     PetscCheck((user->phys->state_var == STATEVAR_CONSERVATIVE) || (user->app_ctx->test_type == TESTTYPE_DIFF_FILTER), user->comm,
359                PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow");
360     problem->apply_inflow.qfunction              = Blasius_Inflow;
361     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
362     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
363     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
364     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context));
365     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context));
366   }
367   PetscFunctionReturn(PETSC_SUCCESS);
368 }
369