xref: /honee/problems/blasius.c (revision 3fd71269be8b4f6146b0317f2c92fbc428536ccb)
1bb8a0c61SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2bb8a0c61SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3bb8a0c61SJames Wright //
4bb8a0c61SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5bb8a0c61SJames Wright //
6bb8a0c61SJames Wright // This file is part of CEED:  http://github.com/ceed
7bb8a0c61SJames Wright 
8bb8a0c61SJames Wright /// @file
9bb8a0c61SJames Wright /// Utility functions for setting up Blasius Boundary Layer
10bb8a0c61SJames Wright 
11bb8a0c61SJames Wright #include "../navierstokes.h"
12bb8a0c61SJames Wright #include "../qfunctions/blasius.h"
13493642f1SJames Wright #include "stg_shur14.h"
14bb8a0c61SJames Wright 
15f362fe62SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
16f362fe62SJames Wright                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
17f362fe62SJames Wright                                    PetscInt *nynodes) {
18f362fe62SJames Wright   PetscErrorCode ierr;
19f362fe62SJames Wright   PetscInt ndims, dims[2];
20f362fe62SJames Wright   FILE *fp;
21f362fe62SJames Wright   const PetscInt char_array_len = 512;
22f362fe62SJames Wright   char line[char_array_len];
23f362fe62SJames Wright   char **array;
24f362fe62SJames Wright   PetscReal *node_locs;
25f362fe62SJames Wright   PetscFunctionBeginUser;
26f362fe62SJames Wright 
27f362fe62SJames Wright   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
28f362fe62SJames Wright   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
29f362fe62SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
30f362fe62SJames Wright 
31f362fe62SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
32f362fe62SJames Wright   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
33f362fe62SJames Wright   *nynodes = dims[0];
34f362fe62SJames Wright   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
35f362fe62SJames Wright 
36f362fe62SJames Wright   for (PetscInt i=0; i<dims[0]; i++) {
37f362fe62SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
38f362fe62SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
39f362fe62SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
4045aa3cadSJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
4145aa3cadSJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT ")", i,
42f362fe62SJames Wright                                    path, ndims, dims[1]);
43f362fe62SJames Wright 
44f362fe62SJames Wright     node_locs[i] = (PetscReal) atof(array[0]);
45f362fe62SJames Wright   }
46f362fe62SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
47f362fe62SJames Wright   *pynodes = node_locs;
48f362fe62SJames Wright   PetscFunctionReturn(0);
49f362fe62SJames Wright }
50f362fe62SJames Wright 
51bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
52bb8a0c61SJames Wright  *
53493642f1SJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
54493642f1SJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
55493642f1SJames Wright  * linearly in logspace to the top surface.
56bb8a0c61SJames Wright  *
57bb8a0c61SJames Wright  * The top surface is also angled downwards, so that it may be used as an
58493642f1SJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
59f362fe62SJames Wright  *
60f362fe62SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
618e63b9cbSJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
628e63b9cbSJames Wright  * the array, along with `num_node_locs`.
63bb8a0c61SJames Wright  */
64f362fe62SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
65f362fe62SJames Wright                                  PetscReal growth, PetscInt N,
66f362fe62SJames Wright                                  PetscReal refine_height, PetscReal top_angle,
678e63b9cbSJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
68bb8a0c61SJames Wright   PetscInt ierr, narr, ncoords;
69bb8a0c61SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
70bb8a0c61SJames Wright   PetscScalar *arr_coords;
71bb8a0c61SJames Wright   Vec vec_coords;
72bb8a0c61SJames Wright   PetscFunctionBeginUser;
73bb8a0c61SJames Wright 
74bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
75bb8a0c61SJames Wright 
76bb8a0c61SJames Wright   // Get domain boundary information
77bb8a0c61SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
78493642f1SJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
79bb8a0c61SJames Wright 
80bb8a0c61SJames Wright   // Get coords array from DM
81bb8a0c61SJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
82bb8a0c61SJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
83bb8a0c61SJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
84bb8a0c61SJames Wright 
85bb8a0c61SJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
86bb8a0c61SJames Wright   ncoords = narr/dim;
87bb8a0c61SJames Wright 
88bb8a0c61SJames Wright   // Get mesh information
89bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
90bb8a0c61SJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
91bb8a0c61SJames Wright                                  NULL); CHKERRQ(ierr);
92f362fe62SJames Wright   // Get element size of the box mesh, for indexing each node
93f362fe62SJames Wright   const PetscReal dybox = domain_size[1]/faces[1];
94bb8a0c61SJames Wright 
958e63b9cbSJames Wright   if (!*node_locs) {
96bb8a0c61SJames Wright     // Calculate the first element height
97bb8a0c61SJames Wright     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
98bb8a0c61SJames Wright 
99bb8a0c61SJames Wright     // Calculate log of sizing outside BL
100bb8a0c61SJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
101bb8a0c61SJames Wright 
1028e63b9cbSJames Wright     *num_node_locs = faces[1] + 1;
1038e63b9cbSJames Wright     PetscReal *temp_node_locs;
1048e63b9cbSJames Wright     ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr);
1058e63b9cbSJames Wright 
106493642f1SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
107bb8a0c61SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
108bb8a0c61SJames Wright       if (y_box_index <= N) {
109029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
110f362fe62SJames Wright                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
111bb8a0c61SJames Wright       } else {
112bb8a0c61SJames Wright         PetscInt j = y_box_index - N;
113029c6dfaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
114f362fe62SJames Wright                        * exp(log(refine_height) + logdy*j);
115f362fe62SJames Wright       }
1168e63b9cbSJames Wright       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2])
1178e63b9cbSJames Wright         temp_node_locs[y_box_index] = coords[i][1];
118f362fe62SJames Wright     }
1198e63b9cbSJames Wright 
1208e63b9cbSJames Wright     *node_locs = temp_node_locs;
121f362fe62SJames Wright   } else {
122f362fe62SJames Wright     // Error checking
1238e63b9cbSJames Wright     if (*num_node_locs < faces[1] +1)
124f362fe62SJames Wright       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
125f362fe62SJames Wright               "There are %d + 1 nodes, but only %d locations given",
1268e63b9cbSJames Wright               faces[1]+1, *num_node_locs);
1278e63b9cbSJames Wright     if (*num_node_locs > faces[1] +1) {
128f362fe62SJames Wright       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
129f362fe62SJames Wright                          "than the mesh has nodes (%d). This maybe unintended.",
1308e63b9cbSJames Wright                          *num_node_locs, faces[1]+1); CHKERRQ(ierr);
131f362fe62SJames Wright     }
1328e63b9cbSJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
133f362fe62SJames Wright 
134f362fe62SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
135f362fe62SJames Wright       // Determine which y-node we're at
136f362fe62SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
137029c6dfaSJames Wright       coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y)
1388e63b9cbSJames Wright                      * (*node_locs)[y_box_index];
139bb8a0c61SJames Wright     }
140bb8a0c61SJames Wright   }
141bb8a0c61SJames Wright 
142bb8a0c61SJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
143bb8a0c61SJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
144bb8a0c61SJames Wright 
145bb8a0c61SJames Wright   PetscFunctionReturn(0);
146bb8a0c61SJames Wright }
147bb8a0c61SJames Wright 
148b7f03f12SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
149bb8a0c61SJames Wright 
150bb8a0c61SJames Wright   PetscInt ierr;
151bb8a0c61SJames Wright   User      user    = *(User *)ctx;
152bb8a0c61SJames Wright   MPI_Comm  comm    = PETSC_COMM_WORLD;
153493642f1SJames Wright   PetscBool use_stg = PETSC_FALSE;
15415a3537eSJed Brown   BlasiusContext blasius_ctx;
15515a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
15615a3537eSJed Brown   CeedQFunctionContext blasius_context;
15715a3537eSJed Brown 
158bb8a0c61SJames Wright   PetscFunctionBeginUser;
159b7f03f12SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
16015a3537eSJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
161bb8a0c61SJames Wright 
162bb8a0c61SJames Wright   // ------------------------------------------------------
163bb8a0c61SJames Wright   //               SET UP Blasius
164bb8a0c61SJames Wright   // ------------------------------------------------------
1659785fe93SJed Brown   problem->ics.qfunction                       = ICsBlasius;
1669785fe93SJed Brown   problem->ics.qfunction_loc                   = ICsBlasius_loc;
167bb8a0c61SJames Wright 
168bb8a0c61SJames Wright   CeedScalar Uinf   = 40;          // m/s
169bb8a0c61SJames Wright   CeedScalar delta0 = 4.2e-4;      // m
170bb8a0c61SJames Wright   CeedScalar theta0 = 288.;        // K
171bb8a0c61SJames Wright   CeedScalar P0     = 1.01e5;      // Pa
1722acc7cbcSKenneth E. Jansen   PetscBool  weakT  = PETSC_FALSE; // weak density or temperature
1737470235eSJames Wright   PetscReal  mesh_refine_height = 5.9e-4; // m
1747470235eSJames Wright   PetscReal  mesh_growth        = 1.08;   // [-]
1757470235eSJames Wright   PetscInt   mesh_Ndelta        = 45;     // [-]
1767470235eSJames Wright   PetscReal  mesh_top_angle     = 5;      // degrees
177f362fe62SJames Wright   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
178bb8a0c61SJames Wright 
179*3fd71269SJames Wright   PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL);
1802acc7cbcSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
1812acc7cbcSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
182bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
183bb8a0c61SJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
184bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
185bb8a0c61SJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
186bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
187bb8a0c61SJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
188bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
189bb8a0c61SJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
1907470235eSJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
1917470235eSJames Wright                                 "Velocity at boundary layer edge",
1927470235eSJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
1937470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
194bb8a0c61SJames Wright                             "Height of boundary layer mesh refinement",
1957470235eSJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
1967470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
197bb8a0c61SJames Wright                             "Geometric growth rate of boundary layer mesh",
1987470235eSJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
1997470235eSJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
200bb8a0c61SJames Wright                             "Geometric top_angle rate of boundary layer mesh",
2017470235eSJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
202f362fe62SJames Wright   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
203f362fe62SJames Wright                             "Path to file with y node locations. "
204f362fe62SJames Wright                             "If empty, will use the algorithmic mesh warping.", NULL,
205f362fe62SJames Wright                             mesh_ynodes_path, mesh_ynodes_path,
206f362fe62SJames Wright                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
207493642f1SJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
208493642f1SJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
209bb8a0c61SJames Wright   PetscOptionsEnd();
210bb8a0c61SJames Wright 
211bb8a0c61SJames Wright   PetscScalar meter  = user->units->meter;
212bb8a0c61SJames Wright   PetscScalar second = user->units->second;
213bb8a0c61SJames Wright   PetscScalar Kelvin = user->units->Kelvin;
214bb8a0c61SJames Wright   PetscScalar Pascal = user->units->Pascal;
215bb8a0c61SJames Wright 
216bb8a0c61SJames Wright   theta0 *= Kelvin;
217bb8a0c61SJames Wright   P0     *= Pascal;
218bb8a0c61SJames Wright   Uinf   *= meter / second;
219bb8a0c61SJames Wright   delta0 *= meter;
220bb8a0c61SJames Wright 
221f362fe62SJames Wright   PetscReal *mesh_ynodes = NULL;
222f362fe62SJames Wright   PetscInt  mesh_nynodes = 0;
223f362fe62SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
224f362fe62SJames Wright     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
225bb8a0c61SJames Wright     CHKERRQ(ierr);
226f362fe62SJames Wright   }
227f362fe62SJames Wright   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
2288e63b9cbSJames Wright                     mesh_refine_height, mesh_top_angle, &mesh_ynodes,
2298e63b9cbSJames Wright                     &mesh_nynodes); CHKERRQ(ierr);
230bb8a0c61SJames Wright 
23115a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
23215a3537eSJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
23315a3537eSJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
234bb8a0c61SJames Wright 
235493642f1SJames Wright   blasius_ctx->weakT         = weakT;
23615a3537eSJed Brown   blasius_ctx->Uinf          = Uinf;
23715a3537eSJed Brown   blasius_ctx->delta0        = delta0;
23815a3537eSJed Brown   blasius_ctx->theta0        = theta0;
23915a3537eSJed Brown   blasius_ctx->P0            = P0;
24004b9037bSJames Wright   newtonian_ig_ctx->P0       = P0;
24115a3537eSJed Brown   blasius_ctx->implicit      = user->phys->implicit;
24215a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
243493642f1SJames Wright 
244ef2c71fdSJames Wright   {
245ef2c71fdSJames Wright     PetscReal domain_min[3];
246ef2c71fdSJames Wright     ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr);
247ef2c71fdSJames Wright     blasius_ctx->x_inflow = domain_min[0];
248ef2c71fdSJames Wright   }
249ef2c71fdSJames Wright 
25015a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
25115a3537eSJed Brown                                   &newtonian_ig_ctx);
252bb8a0c61SJames Wright 
25315a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
25415a3537eSJed Brown   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
255bb8a0c61SJames Wright                               CEED_USE_POINTER,
25615a3537eSJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
25715a3537eSJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
25815a3537eSJed Brown                                      FreeContextPetsc);
259bb8a0c61SJames Wright 
26043bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
26115a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
262493642f1SJames Wright   if (use_stg) {
263e6098bcdSJames Wright     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes,
264e6098bcdSJames Wright                     mesh_nynodes); CHKERRQ(ierr);
2658085925cSJames Wright   } else {
2668085925cSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
2678085925cSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
2680aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
2690aa41abfSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
2708085925cSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
2718085925cSJames Wright                                       &problem->apply_inflow.qfunction_context);
2720aa41abfSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
2730aa41abfSJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
274493642f1SJames Wright   }
275e6098bcdSJames Wright   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
276bb8a0c61SJames Wright   PetscFunctionReturn(0);
277bb8a0c61SJames Wright }
278