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