xref: /honee/problems/blasius.c (revision 9785fe939f22dc29a5e9bc719da8d1f0a67d7249)
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 #ifndef blasius_context_struct
15bb8a0c61SJames Wright #define blasius_context_struct
16bb8a0c61SJames Wright typedef struct BlasiusContext_ *BlasiusContext;
17bb8a0c61SJames Wright struct BlasiusContext_ {
18bb8a0c61SJames Wright   bool       implicit;  // !< Using implicit timesteping or not
19bb8a0c61SJames Wright   CeedScalar delta0;    // !< Boundary layer height at inflow
20bb8a0c61SJames Wright   CeedScalar Uinf;      // !< Velocity at boundary layer edge
21bb8a0c61SJames Wright   CeedScalar P0;        // !< Pressure at outflow
22bb8a0c61SJames Wright   CeedScalar theta0;    // !< Temperature at inflow
232acc7cbcSKenneth E. Jansen   CeedInt weakT;        // !< flag to weakly set Temperature at inflow if not set weak rho instead
24bb8a0c61SJames Wright   struct NewtonianIdealGasContext_ newtonian_ctx;
25bb8a0c61SJames Wright };
26bb8a0c61SJames Wright #endif
27bb8a0c61SJames Wright 
28bb8a0c61SJames Wright #ifndef M_PI
29bb8a0c61SJames Wright #define M_PI    3.14159265358979323846
30bb8a0c61SJames Wright #endif
31bb8a0c61SJames Wright 
32bb8a0c61SJames Wright /* \brief Modify the domain and mesh for blasius
33bb8a0c61SJames Wright  *
34bb8a0c61SJames Wright  * Modifies mesh such that `N` elements are within 1.2*`delta0` with a geometric
35bb8a0c61SJames Wright  * growth ratio of `growth`. Excess elements are then geometrically distributed
36bb8a0c61SJames Wright  * to the top surface.
37bb8a0c61SJames Wright  *
38bb8a0c61SJames Wright  * The top surface is also angled downwards, so that it may be used as an
39bb8a0c61SJames Wright  * outflow. It's angle is controlled by top_angle (in units of degrees).
40bb8a0c61SJames Wright  */
41bb8a0c61SJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N,
42bb8a0c61SJames Wright                           PetscReal refine_height, PetscReal top_angle) {
43bb8a0c61SJames Wright 
44bb8a0c61SJames Wright   PetscInt ierr, narr, ncoords;
45bb8a0c61SJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
46bb8a0c61SJames Wright   PetscScalar *arr_coords;
47bb8a0c61SJames Wright   Vec vec_coords;
48bb8a0c61SJames Wright   PetscFunctionBeginUser;
49bb8a0c61SJames Wright 
50bb8a0c61SJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
51bb8a0c61SJames Wright 
52bb8a0c61SJames Wright   // Get domain boundary information
53bb8a0c61SJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
54bb8a0c61SJames Wright   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
55bb8a0c61SJames Wright 
56bb8a0c61SJames Wright   // Get coords array from DM
57bb8a0c61SJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
58bb8a0c61SJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
59bb8a0c61SJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
60bb8a0c61SJames Wright 
61bb8a0c61SJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
62bb8a0c61SJames Wright   ncoords = narr/dim;
63bb8a0c61SJames Wright 
64bb8a0c61SJames Wright   // Get mesh information
65bb8a0c61SJames Wright   PetscInt nmax = 3, faces[3];
66bb8a0c61SJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
67bb8a0c61SJames Wright                                  NULL); CHKERRQ(ierr);
68bb8a0c61SJames Wright 
69bb8a0c61SJames Wright   // Calculate the first element height
70bb8a0c61SJames Wright   PetscReal dybox = domain_size[1]/faces[1];
71bb8a0c61SJames Wright   PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
72bb8a0c61SJames Wright 
73bb8a0c61SJames Wright   // Calculate log of sizing outside BL
74bb8a0c61SJames Wright   PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
75bb8a0c61SJames Wright 
76bb8a0c61SJames Wright   for(int i=0; i<ncoords; i++) {
77bb8a0c61SJames Wright     PetscInt y_box_index = round(coords[i][1]/dybox);
78bb8a0c61SJames Wright     if(y_box_index <= N) {
79bb8a0c61SJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
80bb8a0c61SJames Wright                      dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1);
81bb8a0c61SJames Wright     } else {
82bb8a0c61SJames Wright       PetscInt j = y_box_index - N;
83bb8a0c61SJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
84bb8a0c61SJames Wright                      exp(log(refine_height) + logdy*j);
85bb8a0c61SJames Wright     }
86bb8a0c61SJames Wright   }
87bb8a0c61SJames Wright 
88bb8a0c61SJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
89bb8a0c61SJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
90bb8a0c61SJames Wright 
91bb8a0c61SJames Wright   PetscFunctionReturn(0);
92bb8a0c61SJames Wright }
93bb8a0c61SJames Wright 
94bb8a0c61SJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *setup_ctx,
95bb8a0c61SJames Wright                           void *ctx) {
96bb8a0c61SJames Wright 
97bb8a0c61SJames Wright   PetscInt ierr;
98bb8a0c61SJames Wright   ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr);
99bb8a0c61SJames Wright   User              user = *(User *)ctx;
100bb8a0c61SJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
101bb8a0c61SJames Wright   PetscFunctionBeginUser;
102bb8a0c61SJames Wright   ierr = PetscCalloc1(1, &user->phys->blasius_ctx); CHKERRQ(ierr);
103bb8a0c61SJames Wright 
104bb8a0c61SJames Wright   // ------------------------------------------------------
105bb8a0c61SJames Wright   //               SET UP Blasius
106bb8a0c61SJames Wright   // ------------------------------------------------------
107*9785fe93SJed Brown   problem->ics.qfunction               = ICsBlasius;
108*9785fe93SJed Brown   problem->ics.qfunction_loc           = ICsBlasius_loc;
109*9785fe93SJed Brown   problem->apply_inflow.qfunction      = Blasius_Inflow;
110*9785fe93SJed Brown   problem->apply_inflow.qfunction_loc  = Blasius_Inflow_loc;
111*9785fe93SJed Brown   problem->apply_outflow.qfunction     = Blasius_Outflow;
112*9785fe93SJed Brown   problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc;
113bb8a0c61SJames Wright   problem->setup_ctx                   = SetupContext_BLASIUS;
114bb8a0c61SJames Wright 
115bb8a0c61SJames Wright   // CeedScalar mu = .04; // Pa s, dynamic viscosity
116bb8a0c61SJames Wright   CeedScalar mu            = 1.8e-5;   // Pa s, dynamic viscosity
117bb8a0c61SJames Wright   CeedScalar Uinf          = 40;   // m/s
118bb8a0c61SJames Wright   CeedScalar delta0        = 4.2e-4;    // m
119bb8a0c61SJames Wright   PetscReal  refine_height = 5.9e-4;    // m
120bb8a0c61SJames Wright   PetscReal  growth        = 1.08; // [-]
121bb8a0c61SJames Wright   PetscInt   Ndelta        = 45;   // [-]
122bb8a0c61SJames Wright   PetscReal  top_angle     = 5;    // degrees
123bb8a0c61SJames Wright   CeedScalar theta0        = 288.; // K
124bb8a0c61SJames Wright   CeedScalar P0            = 1.01e5; // Pa
1252acc7cbcSKenneth E. Jansen   PetscBool  weakT         = PETSC_FALSE; // weak density or temperature
126bb8a0c61SJames Wright 
127bb8a0c61SJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
1282acc7cbcSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
1292acc7cbcSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
130bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
131bb8a0c61SJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
132bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
133bb8a0c61SJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
134bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
135bb8a0c61SJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
136bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
137bb8a0c61SJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
138bb8a0c61SJames Wright   ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge",
139bb8a0c61SJames Wright                                 NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr);
140bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-refine_height",
141bb8a0c61SJames Wright                             "Height of boundary layer mesh refinement",
142bb8a0c61SJames Wright                             NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr);
143bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-growth",
144bb8a0c61SJames Wright                             "Geometric growth rate of boundary layer mesh",
145bb8a0c61SJames Wright                             NULL, growth, &growth, NULL); CHKERRQ(ierr);
146bb8a0c61SJames Wright   ierr = PetscOptionsScalar("-top_angle",
147bb8a0c61SJames Wright                             "Geometric top_angle rate of boundary layer mesh",
148bb8a0c61SJames Wright                             NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr);
149bb8a0c61SJames Wright   PetscOptionsEnd();
150bb8a0c61SJames Wright 
151bb8a0c61SJames Wright   PetscScalar meter           = user->units->meter;
152bb8a0c61SJames Wright   PetscScalar second          = user->units->second;
153bb8a0c61SJames Wright   PetscScalar Kelvin          = user->units->Kelvin;
154bb8a0c61SJames Wright   PetscScalar Pascal          = user->units->Pascal;
155bb8a0c61SJames Wright 
156bb8a0c61SJames Wright   mu     *= Pascal * second;
157bb8a0c61SJames Wright   theta0 *= Kelvin;
158bb8a0c61SJames Wright   P0     *= Pascal;
159bb8a0c61SJames Wright   Uinf   *= meter / second;
160bb8a0c61SJames Wright   delta0 *= meter;
161bb8a0c61SJames Wright 
162bb8a0c61SJames Wright   ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle);
163bb8a0c61SJames Wright   CHKERRQ(ierr);
164bb8a0c61SJames Wright 
1652acc7cbcSKenneth E. Jansen   user->phys->blasius_ctx->weakT     = !!weakT;
166bb8a0c61SJames Wright   user->phys->blasius_ctx->Uinf      = Uinf;
167bb8a0c61SJames Wright   user->phys->blasius_ctx->delta0    = delta0;
168bb8a0c61SJames Wright   user->phys->blasius_ctx->theta0    = theta0;
169bb8a0c61SJames Wright   user->phys->blasius_ctx->P0        = P0;
170bb8a0c61SJames Wright   user->phys->blasius_ctx->implicit  = user->phys->implicit;
171bb8a0c61SJames Wright 
172bb8a0c61SJames Wright   user->phys->newtonian_ig_ctx->mu = mu;
173bb8a0c61SJames Wright   user->phys->blasius_ctx->newtonian_ctx = *user->phys->newtonian_ig_ctx;
174bb8a0c61SJames Wright   PetscFunctionReturn(0);
175bb8a0c61SJames Wright }
176bb8a0c61SJames Wright 
177bb8a0c61SJames Wright PetscErrorCode SetupContext_BLASIUS(Ceed ceed, CeedData ceed_data,
178bb8a0c61SJames Wright                                     AppCtx app_ctx, SetupContext setup_ctx, Physics phys) {
179bb8a0c61SJames Wright   PetscFunctionBeginUser;
180bb8a0c61SJames Wright   PetscInt ierr;
181bb8a0c61SJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
182bb8a0c61SJames Wright   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
183bb8a0c61SJames Wright                               CEED_USE_POINTER,
184bb8a0c61SJames Wright                               sizeof(*setup_ctx), setup_ctx);
185bb8a0c61SJames Wright   ierr = SetupContext_NEWTONIAN_IG(ceed, ceed_data, app_ctx, setup_ctx, phys);
186bb8a0c61SJames Wright   CHKERRQ(ierr);
187bb8a0c61SJames Wright 
188bb8a0c61SJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->blasius_context);
189bb8a0c61SJames Wright   CeedQFunctionContextSetData(ceed_data->blasius_context, CEED_MEM_HOST,
190bb8a0c61SJames Wright                               CEED_USE_POINTER, sizeof(*phys->blasius_ctx), phys->blasius_ctx);
191bb8a0c61SJames Wright   phys->has_neumann = PETSC_TRUE;
192bb8a0c61SJames Wright   if (ceed_data->qf_ics)
193bb8a0c61SJames Wright     CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->blasius_context);
194bb8a0c61SJames Wright   if (ceed_data->qf_apply_inflow)
195bb8a0c61SJames Wright     CeedQFunctionSetContext(ceed_data->qf_apply_inflow, ceed_data->blasius_context);
196bb8a0c61SJames Wright   if (ceed_data->qf_apply_outflow)
197bb8a0c61SJames Wright     CeedQFunctionSetContext(ceed_data->qf_apply_outflow,
198bb8a0c61SJames Wright                             ceed_data->blasius_context);
199bb8a0c61SJames Wright   PetscFunctionReturn(0);
200bb8a0c61SJames Wright }
201bb8a0c61SJames Wright 
202