xref: /honee/problems/blasius.c (revision 0e477d6fa28ddf06be4c7b21d26143f501088c32)
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 "../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   PetscScalar          Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), Pr = Prandtl(&blasius->newtonian_ctx),
26               gamma = HeatCapacityRatio(&blasius->newtonian_ctx);
27   PetscFunctionBegin;
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 / blasius->T_inf;
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(0);
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   PetscFunctionBegin;
80 
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(0);
111 }
112 
113 static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) {
114   PetscInt       ndims, dims[2];
115   FILE          *fp;
116   const PetscInt char_array_len = 512;
117   char           line[char_array_len];
118   char         **array;
119   PetscReal     *node_locs;
120   PetscFunctionBeginUser;
121 
122   PetscCall(PetscFOpen(comm, path, "r", &fp));
123   PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
124   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
125 
126   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
127   if (ndims < 2) dims[1] = 1;  // Assume 1 column of data is not otherwise specified
128   *nynodes = dims[0];
129   PetscCall(PetscMalloc1(*nynodes, &node_locs));
130 
131   for (PetscInt i = 0; i < dims[0]; i++) {
132     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
133     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
134     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
135                "Line %" PetscInt_FMT " of %s does not contain correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
136                ndims, dims[1]);
137 
138     node_locs[i] = (PetscReal)atof(array[0]);
139   }
140   PetscCall(PetscFClose(comm, fp));
141   *pynodes = node_locs;
142   PetscFunctionReturn(0);
143 }
144 
145 /* \brief Modify the domain and mesh for blasius
146  *
147  * Modifies mesh such that `N` elements are within `refine_height` with a geometric growth ratio of `growth`. Excess elements are then distributed
148  * linearly in logspace to the top surface.
149  *
150  * The top surface is also angled downwards, so that it may be used as an outflow.
151  * It's angle is controlled by `top_angle` (in units of degrees).
152  *
153  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` locations.
154  * If it is NULL, then the modified coordinate values will be set in the array, along with `num_node_locs`.
155  */
156 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, PetscReal growth, PetscInt N, PetscReal refine_height, PetscReal top_angle,
157                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
158   PetscInt     narr, ncoords;
159   PetscReal    domain_min[3], domain_max[3], domain_size[3];
160   PetscScalar *arr_coords;
161   Vec          vec_coords;
162   PetscFunctionBeginUser;
163 
164   PetscReal angle_coeff = tan(top_angle * (M_PI / 180));
165 
166   // Get domain boundary information
167   PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
168   for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
169 
170   // Get coords array from DM
171   PetscCall(DMGetCoordinatesLocal(dm, &vec_coords));
172   PetscCall(VecGetLocalSize(vec_coords, &narr));
173   PetscCall(VecGetArray(vec_coords, &arr_coords));
174 
175   PetscScalar(*coords)[dim] = (PetscScalar(*)[dim])arr_coords;
176   ncoords                   = narr / dim;
177 
178   // Get mesh information
179   PetscInt nmax = 3, faces[3];
180   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
181   // Get element size of the box mesh, for indexing each node
182   const PetscReal dybox = domain_size[1] / faces[1];
183 
184   if (!*node_locs) {
185     // Calculate the first element height
186     PetscReal dy1 = refine_height * (growth - 1) / (pow(growth, N) - 1);
187 
188     // Calculate log of sizing outside BL
189     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
190 
191     *num_node_locs = faces[1] + 1;
192     PetscReal *temp_node_locs;
193     PetscCall(PetscMalloc1(*num_node_locs, &temp_node_locs));
194 
195     for (PetscInt i = 0; i < ncoords; i++) {
196       PetscInt y_box_index = round(coords[i][1] / dybox);
197       if (y_box_index <= N) {
198         coords[i][1] =
199             (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * dy1 * (pow(growth, coords[i][1] / dybox) - 1) / (growth - 1);
200       } else {
201         PetscInt j   = y_box_index - N;
202         coords[i][1] = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / domain_max[1]) * exp(log(refine_height) + logdy * j);
203       }
204       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) temp_node_locs[y_box_index] = coords[i][1];
205     }
206 
207     *node_locs = temp_node_locs;
208   } else {
209     PetscCheck(*num_node_locs >= faces[1] + 1, comm, PETSC_ERR_FILE_UNEXPECTED,
210                "The y_node_locs_path has too few locations; There are %d + 1 nodes, but only %d locations given", faces[1] + 1, *num_node_locs);
211     if (*num_node_locs > faces[1] + 1) {
212       PetscCall(PetscPrintf(comm,
213                             "WARNING: y_node_locs_path has more locations (%d) "
214                             "than the mesh has nodes (%d). This maybe unintended.\n",
215                             *num_node_locs, faces[1] + 1));
216     }
217     PetscScalar max_y = (*node_locs)[faces[1]];
218 
219     for (PetscInt i = 0; i < ncoords; i++) {
220       // Determine which y-node we're at
221       PetscInt y_box_index = round(coords[i][1] / dybox);
222       coords[i][1]         = (1 - (coords[i][0] - domain_min[0]) * angle_coeff / max_y) * (*node_locs)[y_box_index];
223     }
224   }
225 
226   PetscCall(VecRestoreArray(vec_coords, &arr_coords));
227   PetscCall(DMSetCoordinatesLocal(dm, vec_coords));
228 
229   PetscFunctionReturn(0);
230 }
231 
232 PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
233   User                     user    = *(User *)ctx;
234   MPI_Comm                 comm    = user->comm;
235   PetscBool                use_stg = PETSC_FALSE;
236   BlasiusContext           blasius_ctx;
237   NewtonianIdealGasContext newtonian_ig_ctx;
238   CeedQFunctionContext     blasius_context;
239 
240   PetscFunctionBeginUser;
241   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
242   PetscCall(PetscCalloc1(1, &blasius_ctx));
243 
244   // ------------------------------------------------------
245   //               SET UP Blasius
246   // ------------------------------------------------------
247   problem->ics.qfunction     = ICsBlasius;
248   problem->ics.qfunction_loc = ICsBlasius_loc;
249 
250   CeedScalar U_inf                                = 40;           // m/s
251   CeedScalar T_inf                                = 288.;         // K
252   CeedScalar T_wall                               = 288.;         // K
253   CeedScalar delta0                               = 4.2e-3;       // m
254   CeedScalar P0                                   = 1.01e5;       // Pa
255   CeedInt    N                                    = 20;           // Number of Chebyshev terms
256   PetscBool  weakT                                = PETSC_FALSE;  // weak density or temperature
257   PetscReal  mesh_refine_height                   = 5.9e-4;       // m
258   PetscReal  mesh_growth                          = 1.08;         // [-]
259   PetscInt   mesh_Ndelta                          = 45;           // [-]
260   PetscReal  mesh_top_angle                       = 5;            // degrees
261   char       mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
262 
263   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
264   PetscCall(PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL));
265   PetscCall(PetscOptionsScalar("-velocity_infinity", "Velocity at boundary layer edge", NULL, U_inf, &U_inf, NULL));
266   PetscCall(PetscOptionsScalar("-temperature_infinity", "Temperature at boundary layer edge", NULL, T_inf, &T_inf, NULL));
267   PetscCall(PetscOptionsScalar("-temperature_wall", "Temperature at wall", NULL, T_wall, &T_wall, NULL));
268   PetscCall(PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL));
269   PetscCall(PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL));
270   PetscCall(PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", NULL, N, &N, NULL));
271   PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, comm, PETSC_ERR_ARG_OUTOFRANGE, "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N,
272              BLASIUS_MAX_N_CHEBYSHEV);
273   PetscCall(PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1));
274   PetscCall(PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height,
275                                NULL));
276   PetscCall(PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL));
277   PetscCall(
278       PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL));
279   PetscCall(PetscOptionsString("-platemesh_y_node_locs_path",
280                                "Path to file with y node locations. "
281                                "If empty, will use the algorithmic mesh warping.",
282                                NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL));
283   PetscCall(PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL));
284   PetscOptionsEnd();
285 
286   PetscScalar meter  = user->units->meter;
287   PetscScalar second = user->units->second;
288   PetscScalar Kelvin = user->units->Kelvin;
289   PetscScalar Pascal = user->units->Pascal;
290 
291   T_inf *= Kelvin;
292   T_wall *= Kelvin;
293   P0 *= Pascal;
294   U_inf *= meter / second;
295   delta0 *= meter;
296 
297   PetscReal *mesh_ynodes  = NULL;
298   PetscInt   mesh_nynodes = 0;
299   if (strcmp(mesh_ynodes_path, "")) {
300     PetscCall(GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes));
301   }
302   PetscCall(ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes));
303 
304   // Some properties depend on parameters from NewtonianIdealGas
305   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
306 
307   blasius_ctx->weakT         = weakT;
308   blasius_ctx->U_inf         = U_inf;
309   blasius_ctx->T_inf         = T_inf;
310   blasius_ctx->T_wall        = T_wall;
311   blasius_ctx->delta0        = delta0;
312   blasius_ctx->P0            = P0;
313   blasius_ctx->n_cheb        = N;
314   newtonian_ig_ctx->P0       = P0;
315   blasius_ctx->implicit      = user->phys->implicit;
316   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
317 
318   {
319     PetscReal domain_min[3], domain_max[3];
320     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
321     blasius_ctx->x_inflow = domain_min[0];
322     blasius_ctx->eta_max  = 5 * domain_max[1] / blasius_ctx->delta0;
323   }
324   PetscBool diff_filter_mms = PETSC_FALSE;
325   PetscCall(PetscOptionsGetBool(NULL, NULL, "-diff_filter_mms", &diff_filter_mms, NULL));
326   if (!use_stg && !diff_filter_mms) PetscCall(ComputeChebyshevCoefficients(blasius_ctx));
327 
328   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
329 
330   CeedQFunctionContextCreate(user->ceed, &blasius_context);
331   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx);
332   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc);
333 
334   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
335   problem->ics.qfunction_context = blasius_context;
336   if (use_stg) {
337     PetscCall(SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, mesh_nynodes));
338   } else if (diff_filter_mms) {
339     PetscCall(DifferentialFilter_MMS_ICSetup(problem));
340   } else {
341     problem->apply_inflow.qfunction              = Blasius_Inflow;
342     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
343     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
344     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
345     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context);
346     CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context);
347   }
348   PetscCall(PetscFree(mesh_ynodes));
349   PetscFunctionReturn(0);
350 }
351