xref: /honee/problems/blasius.c (revision ef2c71fd650449953c22bdf4b1c61c9c6ed44cff)
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 static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
16                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
17                                    PetscInt *nynodes) {
18   PetscErrorCode ierr;
19   PetscInt ndims, dims[2];
20   FILE *fp;
21   const PetscInt char_array_len = 512;
22   char line[char_array_len];
23   char **array;
24   PetscReal *node_locs;
25   PetscFunctionBeginUser;
26 
27   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
28   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
29   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
30 
31   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
32   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
33   *nynodes = dims[0];
34   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
35 
36   for (PetscInt i=0; i<dims[0]; i++) {
37     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
38     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
39     if (ndims < dims[1]) SETERRQ(comm, -1,
40                                    "Line %d of %s does not contain enough columns (%d instead of %d)", i,
41                                    path, ndims, dims[1]);
42 
43     node_locs[i] = (PetscReal) atof(array[0]);
44   }
45   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
46   *pynodes = node_locs;
47   PetscFunctionReturn(0);
48 }
49 
50 /* \brief Modify the domain and mesh for blasius
51  *
52  * Modifies mesh such that `N` elements are within `refine_height` with a
53  * geometric growth ratio of `growth`. Excess elements are then distributed
54  * linearly in logspace to the top surface.
55  *
56  * The top surface is also angled downwards, so that it may be used as an
57  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
58  *
59  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
60  * locations.
61  */
62 static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
63                                  PetscReal growth, PetscInt N,
64                                  PetscReal refine_height, PetscReal top_angle,
65                                  PetscReal node_locs[], PetscInt num_node_locs) {
66 
67   PetscInt ierr, narr, ncoords;
68   PetscReal domain_min[3], domain_max[3], domain_size[3];
69   PetscScalar *arr_coords;
70   Vec vec_coords;
71   PetscFunctionBeginUser;
72 
73   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
74 
75   // Get domain boundary information
76   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
77   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
78 
79   // Get coords array from DM
80   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
81   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
82   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
83 
84   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
85   ncoords = narr/dim;
86 
87   // Get mesh information
88   PetscInt nmax = 3, faces[3];
89   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
90                                  NULL); CHKERRQ(ierr);
91   // Get element size of the box mesh, for indexing each node
92   const PetscReal dybox = domain_size[1]/faces[1];
93 
94   if (!node_locs) {
95     // Calculate the first element height
96     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
97 
98     // Calculate log of sizing outside BL
99     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
100 
101     for (PetscInt i=0; i<ncoords; i++) {
102       PetscInt y_box_index = round(coords[i][1]/dybox);
103       if (y_box_index <= N) {
104         coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff)
105                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
106       } else {
107         PetscInt j = y_box_index - N;
108         coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff)
109                        * exp(log(refine_height) + logdy*j);
110       }
111     }
112   } else {
113     // Error checking
114     if (num_node_locs < faces[1] +1)
115       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
116               "There are %d + 1 nodes, but only %d locations given",
117               faces[1]+1, num_node_locs);
118     if (num_node_locs > faces[1] +1) {
119       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
120                          "than the mesh has nodes (%d). This maybe unintended.",
121                          num_node_locs, faces[1]+1); CHKERRQ(ierr);
122     }
123 
124     for (PetscInt i=0; i<ncoords; i++) {
125       // Determine which y-node we're at
126       PetscInt y_box_index = round(coords[i][1]/dybox);
127       coords[i][1] = (1 - ((coords[i][0] - domain_min[0])/domain_size[0])*angle_coeff)
128                      * node_locs[y_box_index];
129     }
130   }
131 
132   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
133   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
134 
135   PetscFunctionReturn(0);
136 }
137 
138 PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
139 
140   PetscInt ierr;
141   User      user    = *(User *)ctx;
142   MPI_Comm  comm    = PETSC_COMM_WORLD;
143   PetscBool use_stg = PETSC_FALSE;
144   BlasiusContext blasius_ctx;
145   NewtonianIdealGasContext newtonian_ig_ctx;
146   CeedQFunctionContext blasius_context;
147 
148   PetscFunctionBeginUser;
149   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
150   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
151 
152   // ------------------------------------------------------
153   //               SET UP Blasius
154   // ------------------------------------------------------
155   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
156   problem->ics.qfunction               = ICsBlasius;
157   problem->ics.qfunction_loc           = ICsBlasius_loc;
158   problem->apply_outflow.qfunction     = Blasius_Outflow;
159   problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc;
160   problem->apply_inflow.qfunction      = Blasius_Inflow;
161   problem->apply_inflow.qfunction_loc  = Blasius_Inflow_loc;
162 
163   CeedScalar Uinf   = 40;          // m/s
164   CeedScalar delta0 = 4.2e-4;      // m
165   CeedScalar theta0 = 288.;        // K
166   CeedScalar P0     = 1.01e5;      // Pa
167   PetscBool  weakT  = PETSC_FALSE; // weak density or temperature
168   PetscReal  mesh_refine_height = 5.9e-4; // m
169   PetscReal  mesh_growth        = 1.08;   // [-]
170   PetscInt   mesh_Ndelta        = 45;     // [-]
171   PetscReal  mesh_top_angle     = 5;      // degrees
172   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
173 
174   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
175   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
176                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
177   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
178                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
179   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
180                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
181   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
182                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
183   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
184                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
185   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
186                                 "Velocity at boundary layer edge",
187                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
188   ierr = PetscOptionsScalar("-platemesh_refine_height",
189                             "Height of boundary layer mesh refinement",
190                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
191   ierr = PetscOptionsScalar("-platemesh_growth",
192                             "Geometric growth rate of boundary layer mesh",
193                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
194   ierr = PetscOptionsScalar("-platemesh_top_angle",
195                             "Geometric top_angle rate of boundary layer mesh",
196                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
197   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
198                             "Path to file with y node locations. "
199                             "If empty, will use the algorithmic mesh warping.", NULL,
200                             mesh_ynodes_path, mesh_ynodes_path,
201                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
202   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
203                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
204   PetscOptionsEnd();
205 
206   PetscScalar meter  = user->units->meter;
207   PetscScalar second = user->units->second;
208   PetscScalar Kelvin = user->units->Kelvin;
209   PetscScalar Pascal = user->units->Pascal;
210 
211   theta0 *= Kelvin;
212   P0     *= Pascal;
213   Uinf   *= meter / second;
214   delta0 *= meter;
215 
216   PetscReal *mesh_ynodes = NULL;
217   PetscInt  mesh_nynodes = 0;
218   if (strcmp(mesh_ynodes_path, "")) {
219     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
220     CHKERRQ(ierr);
221   }
222   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
223                     mesh_refine_height, mesh_top_angle, mesh_ynodes,
224                     mesh_nynodes); CHKERRQ(ierr);
225 
226   // Some properties depend on parameters from NewtonianIdealGas
227   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
228                               CEED_MEM_HOST, &newtonian_ig_ctx);
229 
230   blasius_ctx->weakT     = weakT;
231   blasius_ctx->Uinf      = Uinf;
232   blasius_ctx->delta0    = delta0;
233   blasius_ctx->theta0    = theta0;
234   blasius_ctx->P0        = P0;
235   blasius_ctx->implicit  = user->phys->implicit;
236   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
237 
238   {
239     PetscReal domain_min[3];
240     ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr);
241     blasius_ctx->x_inflow = domain_min[0];
242   }
243 
244   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
245                                   &newtonian_ig_ctx);
246 
247   CeedQFunctionContextCreate(user->ceed, &blasius_context);
248   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
249                               CEED_USE_POINTER,
250                               sizeof(*blasius_ctx), blasius_ctx);
251   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
252                                      FreeContextPetsc);
253 
254   problem->ics.qfunction_context = blasius_context;
255   CeedQFunctionContextReferenceCopy(blasius_context,
256                                     &problem->apply_inflow.qfunction_context);
257   CeedQFunctionContextReferenceCopy(blasius_context,
258                                     &problem->apply_outflow.qfunction_context);
259   if (use_stg) {
260     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes,
261                     mesh_nynodes); CHKERRQ(ierr);
262   }
263   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
264   PetscFunctionReturn(0);
265 }
266