xref: /libCEED/examples/fluids/problems/blasius.c (revision 990fdeb6bb8fc9af2da4472bdc0d1f57da5da0e5)
188626eedSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
288626eedSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
388626eedSJames Wright //
488626eedSJames Wright // SPDX-License-Identifier: BSD-2-Clause
588626eedSJames Wright //
688626eedSJames Wright // This file is part of CEED:  http://github.com/ceed
788626eedSJames Wright 
888626eedSJames Wright /// @file
988626eedSJames Wright /// Utility functions for setting up Blasius Boundary Layer
1088626eedSJames Wright 
1188626eedSJames Wright #include "../navierstokes.h"
1288626eedSJames Wright #include "../qfunctions/blasius.h"
13ba6664aeSJames Wright #include "stg_shur14.h"
1488626eedSJames Wright 
15d2714514SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm,
16d2714514SJames Wright                                    const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes,
17d2714514SJames Wright                                    PetscInt *nynodes) {
18d2714514SJames Wright   PetscErrorCode ierr;
19d2714514SJames Wright   PetscInt ndims, dims[2];
20d2714514SJames Wright   FILE *fp;
21d2714514SJames Wright   const PetscInt char_array_len = 512;
22d2714514SJames Wright   char line[char_array_len];
23d2714514SJames Wright   char **array;
24d2714514SJames Wright   PetscReal *node_locs;
25d2714514SJames Wright   PetscFunctionBeginUser;
26d2714514SJames Wright 
27d2714514SJames Wright   ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr);
28d2714514SJames Wright   ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
29d2714514SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
30d2714514SJames Wright 
31d2714514SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
32d2714514SJames Wright   if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified
33d2714514SJames Wright   *nynodes = dims[0];
34d2714514SJames Wright   ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr);
35d2714514SJames Wright 
36d2714514SJames Wright   for (PetscInt i=0; i<dims[0]; i++) {
37d2714514SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
38d2714514SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
39d2714514SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
40*990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
41*990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT ")", i,
42d2714514SJames Wright                                    path, ndims, dims[1]);
43d2714514SJames Wright 
44d2714514SJames Wright     node_locs[i] = (PetscReal) atof(array[0]);
45d2714514SJames Wright   }
46d2714514SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
47d2714514SJames Wright   *pynodes = node_locs;
48d2714514SJames Wright   PetscFunctionReturn(0);
49d2714514SJames Wright }
50d2714514SJames Wright 
5188626eedSJames Wright /* \brief Modify the domain and mesh for blasius
5288626eedSJames Wright  *
53ba6664aeSJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
54ba6664aeSJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
55ba6664aeSJames Wright  * linearly in logspace to the top surface.
5688626eedSJames Wright  *
5788626eedSJames Wright  * The top surface is also angled downwards, so that it may be used as an
58ba6664aeSJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
59d2714514SJames Wright  *
60d2714514SJames Wright  * If `node_locs` is not NULL, then the nodes will be placed at `node_locs`
61798d42b8SJames Wright  * locations. If it is NULL, then the modified coordinate values will be set in
62798d42b8SJames Wright  * the array, along with `num_node_locs`.
6388626eedSJames Wright  */
64d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim,
65d2714514SJames Wright                                  PetscReal growth, PetscInt N,
66d2714514SJames Wright                                  PetscReal refine_height, PetscReal top_angle,
67798d42b8SJames Wright                                  PetscReal *node_locs[], PetscInt *num_node_locs) {
6888626eedSJames Wright   PetscInt ierr, narr, ncoords;
6988626eedSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
7088626eedSJames Wright   PetscScalar *arr_coords;
7188626eedSJames Wright   Vec vec_coords;
7288626eedSJames Wright   PetscFunctionBeginUser;
7388626eedSJames Wright 
7488626eedSJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
7588626eedSJames Wright 
7688626eedSJames Wright   // Get domain boundary information
7788626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
78ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
7988626eedSJames Wright 
8088626eedSJames Wright   // Get coords array from DM
8188626eedSJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
8288626eedSJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
8388626eedSJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
8488626eedSJames Wright 
8588626eedSJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
8688626eedSJames Wright   ncoords = narr/dim;
8788626eedSJames Wright 
8888626eedSJames Wright   // Get mesh information
8988626eedSJames Wright   PetscInt nmax = 3, faces[3];
9088626eedSJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
9188626eedSJames Wright                                  NULL); CHKERRQ(ierr);
92d2714514SJames Wright   // Get element size of the box mesh, for indexing each node
93d2714514SJames Wright   const PetscReal dybox = domain_size[1]/faces[1];
9488626eedSJames Wright 
95798d42b8SJames Wright   if (!*node_locs) {
9688626eedSJames Wright     // Calculate the first element height
9788626eedSJames Wright     PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
9888626eedSJames Wright 
9988626eedSJames Wright     // Calculate log of sizing outside BL
10088626eedSJames Wright     PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
10188626eedSJames Wright 
102798d42b8SJames Wright     *num_node_locs = faces[1] + 1;
103798d42b8SJames Wright     PetscReal *temp_node_locs;
104798d42b8SJames Wright     ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr);
105798d42b8SJames Wright 
106ba6664aeSJames Wright     for (PetscInt i=0; i<ncoords; i++) {
10788626eedSJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
10888626eedSJames Wright       if (y_box_index <= N) {
109855641eaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
110d2714514SJames Wright                        * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1);
11188626eedSJames Wright       } else {
11288626eedSJames Wright         PetscInt j = y_box_index - N;
113855641eaSJames Wright         coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1])
114d2714514SJames Wright                        * exp(log(refine_height) + logdy*j);
115d2714514SJames Wright       }
116798d42b8SJames Wright       if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2])
117798d42b8SJames Wright         temp_node_locs[y_box_index] = coords[i][1];
118d2714514SJames Wright     }
119798d42b8SJames Wright 
120798d42b8SJames Wright     *node_locs = temp_node_locs;
121d2714514SJames Wright   } else {
122d2714514SJames Wright     // Error checking
123798d42b8SJames Wright     if (*num_node_locs < faces[1] +1)
124d2714514SJames Wright       SETERRQ(comm, -1, "The y_node_locs_path has too few locations; "
125d2714514SJames Wright               "There are %d + 1 nodes, but only %d locations given",
126798d42b8SJames Wright               faces[1]+1, *num_node_locs);
127798d42b8SJames Wright     if (*num_node_locs > faces[1] +1) {
128d2714514SJames Wright       ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) "
129d2714514SJames Wright                          "than the mesh has nodes (%d). This maybe unintended.",
130798d42b8SJames Wright                          *num_node_locs, faces[1]+1); CHKERRQ(ierr);
131d2714514SJames Wright     }
132798d42b8SJames Wright     PetscScalar max_y = (*node_locs)[faces[1]];
133d2714514SJames Wright 
134d2714514SJames Wright     for (PetscInt i=0; i<ncoords; i++) {
135d2714514SJames Wright       // Determine which y-node we're at
136d2714514SJames Wright       PetscInt y_box_index = round(coords[i][1]/dybox);
137855641eaSJames Wright       coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y)
138798d42b8SJames Wright                      * (*node_locs)[y_box_index];
13988626eedSJames Wright     }
14088626eedSJames Wright   }
14188626eedSJames Wright 
14288626eedSJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
14388626eedSJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
14488626eedSJames Wright 
14588626eedSJames Wright   PetscFunctionReturn(0);
14688626eedSJames Wright }
14788626eedSJames Wright 
148a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
14988626eedSJames Wright 
15088626eedSJames Wright   PetscInt ierr;
15188626eedSJames Wright   User      user    = *(User *)ctx;
15288626eedSJames Wright   MPI_Comm  comm    = PETSC_COMM_WORLD;
153ba6664aeSJames Wright   PetscBool use_stg = PETSC_FALSE;
154841e4c73SJed Brown   BlasiusContext blasius_ctx;
155841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
156841e4c73SJed Brown   CeedQFunctionContext blasius_context;
157841e4c73SJed Brown 
15888626eedSJames Wright   PetscFunctionBeginUser;
159a0add3c9SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
160841e4c73SJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
16188626eedSJames Wright 
16288626eedSJames Wright   // ------------------------------------------------------
16388626eedSJames Wright   //               SET UP Blasius
16488626eedSJames Wright   // ------------------------------------------------------
16591e5af17SJed Brown   problem->ics.qfunction                       = ICsBlasius;
16691e5af17SJed Brown   problem->ics.qfunction_loc                   = ICsBlasius_loc;
16788626eedSJames Wright 
16888626eedSJames Wright   CeedScalar Uinf   = 40;          // m/s
16988626eedSJames Wright   CeedScalar delta0 = 4.2e-4;      // m
17088626eedSJames Wright   CeedScalar theta0 = 288.;        // K
17188626eedSJames Wright   CeedScalar P0     = 1.01e5;      // Pa
172871db79fSKenneth E. Jansen   PetscBool  weakT  = PETSC_FALSE; // weak density or temperature
1737b8d3891SJames Wright   PetscReal  mesh_refine_height = 5.9e-4; // m
1747b8d3891SJames Wright   PetscReal  mesh_growth        = 1.08;   // [-]
1757b8d3891SJames Wright   PetscInt   mesh_Ndelta        = 45;     // [-]
1767b8d3891SJames Wright   PetscReal  mesh_top_angle     = 5;      // degrees
177d2714514SJames Wright   char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = "";
17888626eedSJames Wright 
17988626eedSJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
180871db79fSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
181871db79fSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
18288626eedSJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
18388626eedSJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
18488626eedSJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
18588626eedSJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
18688626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
18788626eedSJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
18888626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
18988626eedSJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
1907b8d3891SJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
1917b8d3891SJames Wright                                 "Velocity at boundary layer edge",
1927b8d3891SJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
1937b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
19488626eedSJames Wright                             "Height of boundary layer mesh refinement",
1957b8d3891SJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
1967b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
19788626eedSJames Wright                             "Geometric growth rate of boundary layer mesh",
1987b8d3891SJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
1997b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
20088626eedSJames Wright                             "Geometric top_angle rate of boundary layer mesh",
2017b8d3891SJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
202d2714514SJames Wright   ierr = PetscOptionsString("-platemesh_y_node_locs_path",
203d2714514SJames Wright                             "Path to file with y node locations. "
204d2714514SJames Wright                             "If empty, will use the algorithmic mesh warping.", NULL,
205d2714514SJames Wright                             mesh_ynodes_path, mesh_ynodes_path,
206d2714514SJames Wright                             sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr);
207ba6664aeSJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
208ba6664aeSJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
20988626eedSJames Wright   PetscOptionsEnd();
21088626eedSJames Wright 
21188626eedSJames Wright   PetscScalar meter  = user->units->meter;
21288626eedSJames Wright   PetscScalar second = user->units->second;
21388626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
21488626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
21588626eedSJames Wright 
21688626eedSJames Wright   theta0 *= Kelvin;
21788626eedSJames Wright   P0     *= Pascal;
21888626eedSJames Wright   Uinf   *= meter / second;
21988626eedSJames Wright   delta0 *= meter;
22088626eedSJames Wright 
221d2714514SJames Wright   PetscReal *mesh_ynodes = NULL;
222d2714514SJames Wright   PetscInt  mesh_nynodes = 0;
223d2714514SJames Wright   if (strcmp(mesh_ynodes_path, "")) {
224d2714514SJames Wright     ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes);
22588626eedSJames Wright     CHKERRQ(ierr);
226d2714514SJames Wright   }
227d2714514SJames Wright   ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta,
228798d42b8SJames Wright                     mesh_refine_height, mesh_top_angle, &mesh_ynodes,
229798d42b8SJames Wright                     &mesh_nynodes); CHKERRQ(ierr);
23088626eedSJames Wright 
231841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
232841e4c73SJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
233841e4c73SJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
23488626eedSJames Wright 
235ba6664aeSJames Wright   blasius_ctx->weakT         = weakT;
236841e4c73SJed Brown   blasius_ctx->Uinf          = Uinf;
237841e4c73SJed Brown   blasius_ctx->delta0        = delta0;
238841e4c73SJed Brown   blasius_ctx->theta0        = theta0;
239841e4c73SJed Brown   blasius_ctx->P0            = P0;
24030e9fa81SJames Wright   newtonian_ig_ctx->P0       = P0;
241841e4c73SJed Brown   blasius_ctx->implicit      = user->phys->implicit;
242841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
243ba6664aeSJames Wright 
244f1122ed0SJames Wright   {
245f1122ed0SJames Wright     PetscReal domain_min[3];
246f1122ed0SJames Wright     ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr);
247f1122ed0SJames Wright     blasius_ctx->x_inflow = domain_min[0];
248f1122ed0SJames Wright   }
249f1122ed0SJames Wright 
250841e4c73SJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
251841e4c73SJed Brown                                   &newtonian_ig_ctx);
25288626eedSJames Wright 
253841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
254841e4c73SJed Brown   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
25588626eedSJames Wright                               CEED_USE_POINTER,
256841e4c73SJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
257841e4c73SJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
258841e4c73SJed Brown                                      FreeContextPetsc);
25988626eedSJames Wright 
260b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
261841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
262ba6664aeSJames Wright   if (use_stg) {
2634e139266SJames Wright     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes,
2644e139266SJames Wright                     mesh_nynodes); CHKERRQ(ierr);
26565dd5cafSJames Wright   } else {
26665dd5cafSJames Wright     problem->apply_inflow.qfunction              = Blasius_Inflow;
26765dd5cafSJames Wright     problem->apply_inflow.qfunction_loc          = Blasius_Inflow_loc;
268fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction     = Blasius_Inflow_Jacobian;
269fc02c281SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
27065dd5cafSJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
27165dd5cafSJames Wright                                       &problem->apply_inflow.qfunction_context);
272fc02c281SJames Wright     CeedQFunctionContextReferenceCopy(blasius_context,
273fc02c281SJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
274ba6664aeSJames Wright   }
2754e139266SJames Wright   ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr);
27688626eedSJames Wright   PetscFunctionReturn(0);
27788626eedSJames Wright }
278