xref: /honee/problems/blasius.c (revision 15a3537e9964f372164f639b6e497665c63d4d0c)
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"
13bb8a0c61SJames Wright 
14bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
15bb8a0c61SJames Wright  *
16bb8a0c61SJames Wright  * Modifies mesh such that `N` elements are within 1.2*`delta0` with a geometric
17bb8a0c61SJames Wright  * growth ratio of `growth`. Excess elements are then geometrically distributed
18bb8a0c61SJames Wright  * to the top surface.
19bb8a0c61SJames Wright  *
20bb8a0c61SJames Wright  * The top surface is also angled downwards, so that it may be used as an
21bb8a0c61SJames Wright  * outflow. It's angle is controlled by top_angle (in units of degrees).
22bb8a0c61SJames Wright  */
23bb8a0c61SJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N,
24bb8a0c61SJames Wright                           PetscReal refine_height, PetscReal top_angle) {
25bb8a0c61SJames Wright 
26bb8a0c61SJames Wright   PetscInt ierr, narr, ncoords;
27bb8a0c61SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
28bb8a0c61SJames Wright   PetscScalar *arr_coords;
29bb8a0c61SJames Wright   Vec vec_coords;
30bb8a0c61SJames Wright   PetscFunctionBeginUser;
31bb8a0c61SJames Wright 
32bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
33bb8a0c61SJames Wright 
34bb8a0c61SJames Wright   // Get domain boundary information
35bb8a0c61SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
36bb8a0c61SJames Wright   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
37bb8a0c61SJames Wright 
38bb8a0c61SJames Wright   // Get coords array from DM
39bb8a0c61SJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
40bb8a0c61SJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
41bb8a0c61SJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
42bb8a0c61SJames Wright 
43bb8a0c61SJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
44bb8a0c61SJames Wright   ncoords = narr/dim;
45bb8a0c61SJames Wright 
46bb8a0c61SJames Wright   // Get mesh information
47bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
48bb8a0c61SJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
49bb8a0c61SJames Wright                                  NULL); CHKERRQ(ierr);
50bb8a0c61SJames Wright 
51bb8a0c61SJames Wright   // Calculate the first element height
52bb8a0c61SJames Wright   PetscReal dybox = domain_size[1]/faces[1];
53bb8a0c61SJames Wright   PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
54bb8a0c61SJames Wright 
55bb8a0c61SJames Wright   // Calculate log of sizing outside BL
56bb8a0c61SJames Wright   PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
57bb8a0c61SJames Wright 
58bb8a0c61SJames Wright   for(int i=0; i<ncoords; i++) {
59bb8a0c61SJames Wright     PetscInt y_box_index = round(coords[i][1]/dybox);
60bb8a0c61SJames Wright     if(y_box_index <= N) {
61bb8a0c61SJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
62bb8a0c61SJames Wright                      dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1);
63bb8a0c61SJames Wright     } else {
64bb8a0c61SJames Wright       PetscInt j = y_box_index - N;
65bb8a0c61SJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
66bb8a0c61SJames Wright                      exp(log(refine_height) + logdy*j);
67bb8a0c61SJames Wright     }
68bb8a0c61SJames Wright   }
69bb8a0c61SJames Wright 
70bb8a0c61SJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
71bb8a0c61SJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
72bb8a0c61SJames Wright 
73bb8a0c61SJames Wright   PetscFunctionReturn(0);
74bb8a0c61SJames Wright }
75bb8a0c61SJames Wright 
76bb8a0c61SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *setup_ctx,
77bb8a0c61SJames Wright                           void *ctx) {
78bb8a0c61SJames Wright 
79bb8a0c61SJames Wright   PetscInt ierr;
80bb8a0c61SJames Wright   User              user = *(User *)ctx;
81bb8a0c61SJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
82*15a3537eSJed Brown   BlasiusContext    blasius_ctx;
83*15a3537eSJed Brown   NewtonianIdealGasContext newtonian_ig_ctx;
84*15a3537eSJed Brown   CeedQFunctionContext blasius_context;
85*15a3537eSJed Brown 
86bb8a0c61SJames Wright   PetscFunctionBeginUser;
87*15a3537eSJed Brown   ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr);
88*15a3537eSJed Brown   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
89bb8a0c61SJames Wright 
90bb8a0c61SJames Wright   // ------------------------------------------------------
91bb8a0c61SJames Wright   //               SET UP Blasius
92bb8a0c61SJames Wright   // ------------------------------------------------------
93*15a3537eSJed Brown   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
949785fe93SJed Brown   problem->ics.qfunction               = ICsBlasius;
959785fe93SJed Brown   problem->ics.qfunction_loc           = ICsBlasius_loc;
969785fe93SJed Brown   problem->apply_inflow.qfunction      = Blasius_Inflow;
979785fe93SJed Brown   problem->apply_inflow.qfunction_loc  = Blasius_Inflow_loc;
989785fe93SJed Brown   problem->apply_outflow.qfunction     = Blasius_Outflow;
999785fe93SJed Brown   problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc;
100bb8a0c61SJames Wright 
101bb8a0c61SJames Wright   // CeedScalar mu = .04; // Pa s, dynamic viscosity
102bb8a0c61SJames Wright   CeedScalar Uinf          = 40;   // m/s
103bb8a0c61SJames Wright   CeedScalar delta0        = 4.2e-4;    // m
104bb8a0c61SJames Wright   PetscReal  refine_height = 5.9e-4;    // m
105bb8a0c61SJames Wright   PetscReal  growth        = 1.08; // [-]
106bb8a0c61SJames Wright   PetscInt   Ndelta        = 45;   // [-]
107bb8a0c61SJames Wright   PetscReal  top_angle     = 5;    // degrees
108bb8a0c61SJames Wright   CeedScalar theta0        = 288.; // K
109bb8a0c61SJames Wright   CeedScalar P0            = 1.01e5; // Pa
1102acc7cbcSKenneth E. Jansen   PetscBool  weakT         = PETSC_FALSE; // weak density or temperature
111bb8a0c61SJames Wright 
112bb8a0c61SJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
1132acc7cbcSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
1142acc7cbcSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
115bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
116bb8a0c61SJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
117bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
118bb8a0c61SJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
119bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
120bb8a0c61SJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
121bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
122bb8a0c61SJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
123bb8a0c61SJames Wright   ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge",
124bb8a0c61SJames Wright                                 NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr);
125bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-refine_height",
126bb8a0c61SJames Wright                             "Height of boundary layer mesh refinement",
127bb8a0c61SJames Wright                             NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr);
128bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-growth",
129bb8a0c61SJames Wright                             "Geometric growth rate of boundary layer mesh",
130bb8a0c61SJames Wright                             NULL, growth, &growth, NULL); CHKERRQ(ierr);
131bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-top_angle",
132bb8a0c61SJames Wright                             "Geometric top_angle rate of boundary layer mesh",
133bb8a0c61SJames Wright                             NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr);
134bb8a0c61SJames Wright   PetscOptionsEnd();
135bb8a0c61SJames Wright 
136bb8a0c61SJames Wright   PetscScalar meter           = user->units->meter;
137bb8a0c61SJames Wright   PetscScalar second          = user->units->second;
138bb8a0c61SJames Wright   PetscScalar Kelvin          = user->units->Kelvin;
139bb8a0c61SJames Wright   PetscScalar Pascal          = user->units->Pascal;
140bb8a0c61SJames Wright 
141bb8a0c61SJames Wright   theta0 *= Kelvin;
142bb8a0c61SJames Wright   P0     *= Pascal;
143bb8a0c61SJames Wright   Uinf   *= meter / second;
144bb8a0c61SJames Wright   delta0 *= meter;
145bb8a0c61SJames Wright 
146bb8a0c61SJames Wright   ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle);
147bb8a0c61SJames Wright   CHKERRQ(ierr);
148bb8a0c61SJames Wright 
149*15a3537eSJed Brown   // Some properties depend on parameters from NewtonianIdealGas
150*15a3537eSJed Brown   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
151*15a3537eSJed Brown                               CEED_MEM_HOST, &newtonian_ig_ctx);
152bb8a0c61SJames Wright 
153*15a3537eSJed Brown   blasius_ctx->weakT     = !!weakT;
154*15a3537eSJed Brown   blasius_ctx->Uinf      = Uinf;
155*15a3537eSJed Brown   blasius_ctx->delta0    = delta0;
156*15a3537eSJed Brown   blasius_ctx->theta0    = theta0;
157*15a3537eSJed Brown   blasius_ctx->P0        = P0;
158*15a3537eSJed Brown   blasius_ctx->implicit  = user->phys->implicit;
159*15a3537eSJed Brown   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
160*15a3537eSJed Brown   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
161*15a3537eSJed Brown                                   &newtonian_ig_ctx);
162bb8a0c61SJames Wright 
163*15a3537eSJed Brown   CeedQFunctionContextCreate(user->ceed, &blasius_context);
164*15a3537eSJed Brown   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
165bb8a0c61SJames Wright                               CEED_USE_POINTER,
166*15a3537eSJed Brown                               sizeof(*blasius_ctx), blasius_ctx);
167*15a3537eSJed Brown   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
168*15a3537eSJed Brown                                      FreeContextPetsc);
169bb8a0c61SJames Wright 
170*15a3537eSJed Brown   problem->ics.qfunction_context = blasius_context;
171*15a3537eSJed Brown   CeedQFunctionContextReferenceCopy(blasius_context,
172*15a3537eSJed Brown                                     &problem->apply_inflow.qfunction_context);
173*15a3537eSJed Brown   CeedQFunctionContextReferenceCopy(blasius_context,
174*15a3537eSJed Brown                                     &problem->apply_outflow.qfunction_context);
175bb8a0c61SJames Wright   PetscFunctionReturn(0);
176bb8a0c61SJames Wright }
177