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