xref: /libCEED/examples/fluids/problems/blasius.c (revision 7b8d38910f090ec36d75f5ddf0ab25cdb6e25045)
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 
1588626eedSJames Wright /* \brief Modify the domain and mesh for blasius
1688626eedSJames Wright  *
17ba6664aeSJames Wright  * Modifies mesh such that `N` elements are within `refine_height` with a
18ba6664aeSJames Wright  * geometric growth ratio of `growth`. Excess elements are then distributed
19ba6664aeSJames Wright  * linearly in logspace to the top surface.
2088626eedSJames Wright  *
2188626eedSJames Wright  * The top surface is also angled downwards, so that it may be used as an
22ba6664aeSJames Wright  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
2388626eedSJames Wright  */
24*7b8d3891SJames Wright static PetscErrorCode ModifyMesh(DM dm, PetscInt dim, PetscReal growth,
25*7b8d3891SJames Wright                                  PetscInt N, PetscReal refine_height,
26*7b8d3891SJames Wright                                  PetscReal top_angle) {
2788626eedSJames Wright 
2888626eedSJames Wright   PetscInt ierr, narr, ncoords;
2988626eedSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
3088626eedSJames Wright   PetscScalar *arr_coords;
3188626eedSJames Wright   Vec vec_coords;
3288626eedSJames Wright   PetscFunctionBeginUser;
3388626eedSJames Wright 
3488626eedSJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
3588626eedSJames Wright 
3688626eedSJames Wright   // Get domain boundary information
3788626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
38ba6664aeSJames Wright   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
3988626eedSJames Wright 
4088626eedSJames Wright   // Get coords array from DM
4188626eedSJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
4288626eedSJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
4388626eedSJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
4488626eedSJames Wright 
4588626eedSJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
4688626eedSJames Wright   ncoords = narr/dim;
4788626eedSJames Wright 
4888626eedSJames Wright   // Get mesh information
4988626eedSJames Wright   PetscInt nmax = 3, faces[3];
5088626eedSJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
5188626eedSJames Wright                                  NULL); CHKERRQ(ierr);
5288626eedSJames Wright 
5388626eedSJames Wright   // Calculate the first element height
5488626eedSJames Wright   PetscReal dybox = domain_size[1]/faces[1];
5588626eedSJames Wright   PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
5688626eedSJames Wright 
5788626eedSJames Wright   // Calculate log of sizing outside BL
5888626eedSJames Wright   PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
5988626eedSJames Wright 
60ba6664aeSJames Wright   for(PetscInt i=0; i<ncoords; i++) {
6188626eedSJames Wright     PetscInt y_box_index = round(coords[i][1]/dybox);
6288626eedSJames Wright     if(y_box_index <= N) {
6388626eedSJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
6488626eedSJames Wright                      dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1);
6588626eedSJames Wright     } else {
6688626eedSJames Wright       PetscInt j = y_box_index - N;
6788626eedSJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
6888626eedSJames Wright                      exp(log(refine_height) + logdy*j);
6988626eedSJames Wright     }
7088626eedSJames Wright   }
7188626eedSJames Wright 
7288626eedSJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
7388626eedSJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
7488626eedSJames Wright 
7588626eedSJames Wright   PetscFunctionReturn(0);
7688626eedSJames Wright }
7788626eedSJames Wright 
78a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
7988626eedSJames Wright 
8088626eedSJames Wright   PetscInt ierr;
8188626eedSJames Wright   User           user    = *(User *)ctx;
8288626eedSJames Wright   MPI_Comm       comm    = PETSC_COMM_WORLD;
83ba6664aeSJames Wright   PetscBool      use_stg = PETSC_FALSE;
84841e4c73SJed Brown   BlasiusContext blasius_ctx;
85841e4c73SJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
86841e4c73SJed Brown   CeedQFunctionContext blasius_context;
87841e4c73SJed Brown 
8888626eedSJames Wright   PetscFunctionBeginUser;
89a0add3c9SJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
90841e4c73SJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
9188626eedSJames Wright 
9288626eedSJames Wright   // ------------------------------------------------------
9388626eedSJames Wright   //               SET UP Blasius
9488626eedSJames Wright   // ------------------------------------------------------
95841e4c73SJed Brown   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
9691e5af17SJed Brown   problem->ics.qfunction               = ICsBlasius;
9791e5af17SJed Brown   problem->ics.qfunction_loc           = ICsBlasius_loc;
9891e5af17SJed Brown   problem->apply_outflow.qfunction     = Blasius_Outflow;
9991e5af17SJed Brown   problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc;
100ba6664aeSJames Wright   problem->apply_inflow.qfunction      = Blasius_Inflow;
101ba6664aeSJames Wright   problem->apply_inflow.qfunction_loc  = Blasius_Inflow_loc;
10288626eedSJames Wright 
10388626eedSJames Wright   CeedScalar Uinf   = 40;          // m/s
10488626eedSJames Wright   CeedScalar delta0 = 4.2e-4;      // m
10588626eedSJames Wright   CeedScalar theta0 = 288.;        // K
10688626eedSJames Wright   CeedScalar P0     = 1.01e5;      // Pa
107871db79fSKenneth E. Jansen   PetscBool  weakT  = PETSC_FALSE; // weak density or temperature
108*7b8d3891SJames Wright   PetscReal  mesh_refine_height = 5.9e-4; // m
109*7b8d3891SJames Wright   PetscReal  mesh_growth        = 1.08;   // [-]
110*7b8d3891SJames Wright   PetscInt   mesh_Ndelta        = 45;     // [-]
111*7b8d3891SJames Wright   PetscReal  mesh_top_angle     = 5;      // degrees
11288626eedSJames Wright 
11388626eedSJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
114871db79fSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
115871db79fSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
11688626eedSJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
11788626eedSJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
11888626eedSJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
11988626eedSJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
12088626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
12188626eedSJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
12288626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
12388626eedSJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
124*7b8d3891SJames Wright   ierr = PetscOptionsBoundedInt("-platemesh_Ndelta",
125*7b8d3891SJames Wright                                 "Velocity at boundary layer edge",
126*7b8d3891SJames Wright                                 NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr);
127*7b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_refine_height",
12888626eedSJames Wright                             "Height of boundary layer mesh refinement",
129*7b8d3891SJames Wright                             NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr);
130*7b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_growth",
13188626eedSJames Wright                             "Geometric growth rate of boundary layer mesh",
132*7b8d3891SJames Wright                             NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr);
133*7b8d3891SJames Wright   ierr = PetscOptionsScalar("-platemesh_top_angle",
13488626eedSJames Wright                             "Geometric top_angle rate of boundary layer mesh",
135*7b8d3891SJames Wright                             NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr);
136ba6664aeSJames Wright   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
137ba6664aeSJames Wright                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
13888626eedSJames Wright   PetscOptionsEnd();
13988626eedSJames Wright 
14088626eedSJames Wright   PetscScalar meter  = user->units->meter;
14188626eedSJames Wright   PetscScalar second = user->units->second;
14288626eedSJames Wright   PetscScalar Kelvin = user->units->Kelvin;
14388626eedSJames Wright   PetscScalar Pascal = user->units->Pascal;
14488626eedSJames Wright 
14588626eedSJames Wright   theta0 *= Kelvin;
14688626eedSJames Wright   P0     *= Pascal;
14788626eedSJames Wright   Uinf   *= meter / second;
14888626eedSJames Wright   delta0 *= meter;
14988626eedSJames Wright 
150*7b8d3891SJames Wright   ierr = ModifyMesh(dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle);
15188626eedSJames Wright   CHKERRQ(ierr);
15288626eedSJames Wright 
153841e4c73SJed Brown   // Some properties depend on parameters from NewtonianIdealGas
154841e4c73SJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
155841e4c73SJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
15688626eedSJames Wright 
157ba6664aeSJames Wright   blasius_ctx->weakT     = weakT;
158841e4c73SJed Brown   blasius_ctx->Uinf      = Uinf;
159841e4c73SJed Brown   blasius_ctx->delta0    = delta0;
160841e4c73SJed Brown   blasius_ctx->theta0    = theta0;
161841e4c73SJed Brown   blasius_ctx->P0        = P0;
162841e4c73SJed Brown   blasius_ctx->implicit  = user->phys->implicit;
163841e4c73SJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
164ba6664aeSJames Wright 
165841e4c73SJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
166841e4c73SJed Brown                                   &newtonian_ig_ctx);
16788626eedSJames Wright 
168841e4c73SJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
169841e4c73SJed Brown   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
17088626eedSJames Wright                               CEED_USE_POINTER,
171841e4c73SJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
172841e4c73SJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
173841e4c73SJed Brown                                      FreeContextPetsc);
17488626eedSJames Wright 
175841e4c73SJed Brown   problem->ics.qfunction_context = blasius_context;
176841e4c73SJed Brown   CeedQFunctionContextReferenceCopy(blasius_context,
177841e4c73SJed Brown                                     &problem->apply_inflow.qfunction_context);
178841e4c73SJed Brown   CeedQFunctionContextReferenceCopy(blasius_context,
179841e4c73SJed Brown                                     &problem->apply_outflow.qfunction_context);
180ba6664aeSJames Wright   if (use_stg) {
181ba6664aeSJames Wright     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0); CHKERRQ(ierr);
182ba6664aeSJames Wright   }
18388626eedSJames Wright   PetscFunctionReturn(0);
18488626eedSJames Wright }
185