xref: /honee/problems/blasius.c (revision b77f7f1cdab1a1e467caa338141434688aea171c)
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, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
162                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
163   PetscInt     narr, ncoords, dim;
164   PetscReal    domain_min[3], domain_max[3], domain_size[3];
165   PetscScalar *arr_coords;
166   Vec          vec_coords;
167 
168   PetscFunctionBeginUser;
169   PetscCall(DMGetDimension(dm, &dim));
170   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
171   // Get domain boundary information
172   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
173   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
174 
175   // Get coords array from DM
176   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
177   PetscCall(VecGetLocalSize(vec_coords, &narr));
178   PetscCall(VecGetArray(vec_coords, &arr_coords));
179 
180   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
181   ncoords                   = narr / dim;
182 
183   // Get mesh information
184   PetscInt nmax = 3, faces[3];
185   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
186   // Get element size of the box mesh, for indexing each node
187   const PetscReal dybox = domain_size[1] / faces[1];
188 
189   if (!*node_locs) {
190     // Calculate the first element height
191     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
192 
193     // Calculate log of sizing outside BL
194     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
195 
196     *num_node_locs = faces[1] + 1;
197     PetscReal *temp_node_locs;
198     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
199 
200     for (PetscInt i = 0; i < ncoords; i++) {
201       PetscInt y_box_index = round(coords[i][1] / dybox);
202       if (y_box_index <= N) {
203         coords[i][1] =
204             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
205       } else {
206         PetscInt j   = y_box_index - N;
207         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
208       }
209       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
210     }
211 
212     *node_locs = temp_node_locs;
213   } else {
214     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
215                "The y_node_locs_path has too few locations; There are %" PetscInt_FMT " + 1 nodes, but only %" PetscInt_FMT " locations given",
216                faces[1] + 1, *num_node_locs);
217     if (*num_node_locs > faces[1] + 1) {
218       PetscCall(PetscPrintf(comm,
219                             "WARNING: y_node_locs_path has more locations (%" PetscInt_FMT ") "
220                             "than the mesh has nodes (%" PetscInt_FMT "). This maybe unintended.\n",
221                             *num_node_locs, faces[1] + 1));
222     }
223     PetscScalar max_y = (*node_locs)[faces[1]];
224 
225     for (PetscInt i = 0; i < ncoords; i++) {
226       // Determine which y-node we're at
227       PetscInt y_box_index = round(coords[i][1] / dybox);
228       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
229     }
230   }
231 
232   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
233   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
234   PetscFunctionReturn(PETSC_SUCCESS);
235 }
236 
237 static PetscErrorCode BlasiusInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) {
238   HoneeBCStruct honee_bc;
239 
240   PetscFunctionBeginUser;
241   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
242   PetscCall(HoneeBCCreateIFunctionQF(bc_def, Blasius_Inflow, Blasius_Inflow_loc, honee_bc->qfctx, qf));
243   PetscFunctionReturn(PETSC_SUCCESS);
244 }
245 
246 static PetscErrorCode BlasiusInflowBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) {
247   HoneeBCStruct honee_bc;
248 
249   PetscFunctionBeginUser;
250   PetscCall(BCDefinitionGetContext(bc_def, &honee_bc));
251   PetscCall(HoneeBCCreateIJacobianQF(bc_def, Blasius_Inflow_Jacobian, Blasius_Inflow_Jacobian_loc, honee_bc->qfctx, qf));
252   PetscFunctionReturn(PETSC_SUCCESS);
253 }
254 
255 PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx) {
256   Honee                    honee   = *(Honee *)ctx;
257   MPI_Comm                 comm    = honee->comm;
258   Ceed                     ceed    = honee->ceed;
259   PetscBool                use_stg = PETSC_FALSE;
260   BlasiusContext           blasius_ctx;
261   NewtonianIdealGasContext newtonian_ig_ctx;
262   CeedQFunctionContext     blasius_qfctx;
263 
264   PetscFunctionBeginUser;
265   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx));
266   PetscCall(PetscCalloc1(1, &blasius_ctx));
267 
268   // ------------------------------------------------------
269   //               SET UP Blasius
270   // ------------------------------------------------------
271   problem->ics.qf_func_ptr = ICsBlasius;
272   problem->ics.qf_loc      = ICsBlasius_loc;
273 
274   CeedScalar U_inf                                = 40;           // m/s
275   CeedScalar T_inf                                = 288.;         // K
276   CeedScalar T_wall                               = 288.;         // K
277   CeedScalar delta0                               = 4.2e-3;       // m
278   CeedScalar P_inf                                = 1.01e5;       // Pa
279   PetscInt   N                                    = 20;           // Number of Chebyshev terms
280   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
281   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
282   PetscReal  mesh_growth                          = 1.08;         // [-]
283   PetscInt   mesh_Ndelta                          = 45;           // [-]
284   PetscReal  mesh_top_angle                       = 5;            // degrees
285   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
286   PetscBool  P0_set;
287 
288   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
289   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
290   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
291   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
292   PetscCall(PetscOptionsHasName(NULL, NULL, "-P0", &P0_set));  // For maintaining behavior of -P0 flag (which is deprecated)
293   PetscCall(
294       PetscOptionsDeprecated("-P0", "-pressure_infinity", "libCEED 0.12.0",
295                              "Use -pressure_infinity to set pressure at boundary layer edge and -idl_pressure to set the IDL reference pressure"));
296   PetscCall(PetscOptionsScalar("-pressure_infinity", "Pressure at boundary layer edge", NULL, P_inf, &P_inf, NULL));
297   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
298   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
299   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
300   PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N,
301              BLASIUS_MAX_N_CHEBYSHEV);
302   if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
303     PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
304     PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height,
305                                  &mesh_refine_height, NULL));
306     PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
307     PetscCall(
308         PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
309     PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
310                                  "Path to file with y node locations. "
311                                  "If empty, will use the algorithmic mesh warping.",
312                                  NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
313   }
314   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
315   PetscOptionsEnd();
316 
317   PetscScalar meter  = honee->units->meter;
318   PetscScalar second = honee->units->second;
319   PetscScalar Kelvin = honee->units->Kelvin;
320   PetscScalar Pascal = honee->units->Pascal;
321 
322   T_inf *= Kelvin;
323   T_wall *= Kelvin;
324   P_inf *= Pascal;
325   U_inf *= meter / second;
326   delta0 *= meter;
327 
328   if (honee->app_ctx->mesh_transform_type == MESH_TRANSFORM_PLATEMESH) {
329     PetscReal *mesh_ynodes  = NULL;
330     PetscInt   mesh_nynodes = 0;
331     if (strcmp(mesh_ynodes_path, "")) PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
332     PetscCall(ModifyMesh(comm, dm, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
333     PetscCall(PetscFree(mesh_ynodes));
334   }
335 
336   // Some properties depend on parameters from NewtonianIdealGas
337   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
338 
339   StatePrimitive Y_inf = {
340       .pressure = P_inf, .velocity = {U_inf, 0, 0},
341            .temperature = T_inf
342   };
343   State S_infty = StateFromPrimitive(newtonian_ig_ctx, Y_inf);
344 
345   blasius_ctx->weakT    = weakT;
346   blasius_ctx->T_wall   = T_wall;
347   blasius_ctx->delta0   = delta0;
348   blasius_ctx->S_infty  = S_infty;
349   blasius_ctx->n_cheb   = N;
350   blasius_ctx->implicit = honee->phys->implicit;
351   if (P0_set) newtonian_ig_ctx->idl_pressure = P_inf;  // For maintaining behavior of -P0 flag (which is deprecated)
352   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
353 
354   {
355     PetscReal domain_min[3], domain_max[3];
356     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
357     blasius_ctx->x_inflow = domain_min[0];
358     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
359   }
360   PetscBool diff_filter_mms = PETSC_FALSE;
361   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
362   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
363 
364   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
365 
366   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &blasius_qfctx));
367   PetscCallCeed(ceed, CeedQFunctionContextSetData(blasius_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx));
368   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(blasius_qfctx, CEED_MEM_HOST, FreeContextPetsc));
369 
370   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx));
371   problem->ics.qfctx = blasius_qfctx;
372   if (use_stg) {
373     PetscCall(SetupStg(comm, dm, problem, honee, weakT, S_infty.Y.temperature, S_infty.Y.pressure));
374   } else if (diff_filter_mms) {
375     PetscCall(DifferentialFilterMmsICSetup(problem));
376   } else {
377     PetscCheck((honee->phys->state_var == STATEVAR_CONSERVATIVE) || (honee->app_ctx->test_type == TESTTYPE_DIFF_FILTER), honee->comm,
378                PETSC_ERR_ARG_INCOMP, "Can only use conservative variables with Blasius and weak inflow");
379     for (PetscCount b = 0; b < problem->num_bc_defs; b++) {
380       BCDefinition bc_def = problem->bc_defs[b];
381       const char  *name;
382 
383       PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL));
384       if (!strcmp(name, "inflow")) {
385         HoneeBCStruct honee_bc;
386 
387         PetscCall(PetscNew(&honee_bc));
388         PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(blasius_qfctx, &honee_bc->qfctx));
389         honee_bc->honee              = honee;
390         honee_bc->num_comps_jac_data = 0;
391         PetscCall(BCDefinitionSetContext(bc_def, HoneeBCDestroy, honee_bc));
392 
393         PetscCall(BCDefinitionSetIFunction(bc_def, BlasiusInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp));
394         PetscCall(BCDefinitionSetIJacobian(bc_def, BlasiusInflowBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp));
395       }
396     }
397   }
398   PetscFunctionReturn(PETSC_SUCCESS);
399 }
400