xref: /honee/problems/blasius.c (revision eef2387ddc6232a246fdd621cba3d3d81f24a721)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Utility functions for setting up Blasius Boundary Layer
10 
11 #include "../navierstokes.h"
12 #include "../qfunctions/blasius.h"
13 #include "stg_shur14.h"
14 
15 /* \brief Modify the domain and mesh for blasius
16  *
17  * Modifies mesh such that `N` elements are within `refine_height` with a
18  * geometric growth ratio of `growth`. Excess elements are then distributed
19  * linearly in logspace to the top surface.
20  *
21  * The top surface is also angled downwards, so that it may be used as an
22  * outflow. It's angle is controlled by `top_angle` (in units of degrees).
23  */
24 PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N,
25                           PetscReal refine_height, PetscReal top_angle) {
26 
27   PetscInt ierr, narr, ncoords;
28   PetscReal domain_min[3], domain_max[3], domain_size[3];
29   PetscScalar *arr_coords;
30   Vec vec_coords;
31   PetscFunctionBeginUser;
32 
33   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
34 
35   // Get domain boundary information
36   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
37   for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
38 
39   // Get coords array from DM
40   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
41   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
42   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
43 
44   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
45   ncoords = narr/dim;
46 
47   // Get mesh information
48   PetscInt nmax = 3, faces[3];
49   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
50                                  NULL); CHKERRQ(ierr);
51 
52   // Calculate the first element height
53   PetscReal dybox = domain_size[1]/faces[1];
54   PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
55 
56   // Calculate log of sizing outside BL
57   PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
58 
59   for(PetscInt i=0; i<ncoords; i++) {
60     PetscInt y_box_index = round(coords[i][1]/dybox);
61     if(y_box_index <= N) {
62       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
63                      dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1);
64     } else {
65       PetscInt j = y_box_index - N;
66       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
67                      exp(log(refine_height) + logdy*j);
68     }
69   }
70 
71   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
72   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
73 
74   PetscFunctionReturn(0);
75 }
76 
77 PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) {
78 
79   PetscInt ierr;
80   User           user    = *(User *)ctx;
81   MPI_Comm       comm    = PETSC_COMM_WORLD;
82   PetscBool      use_stg = PETSC_FALSE;
83   BlasiusContext blasius_ctx;
84   NewtonianIdealGasContext newtonian_ig_ctx;
85   CeedQFunctionContext blasius_context;
86 
87   PetscFunctionBeginUser;
88   ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr);
89   ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr);
90 
91   // ------------------------------------------------------
92   //               SET UP Blasius
93   // ------------------------------------------------------
94   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
95   problem->ics.qfunction               = ICsBlasius;
96   problem->ics.qfunction_loc           = ICsBlasius_loc;
97   problem->apply_inflow.qfunction      = Blasius_Inflow;
98   problem->apply_inflow.qfunction_loc  = Blasius_Inflow_loc;
99   problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian;
100   problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc;
101   problem->apply_outflow.qfunction     = Blasius_Outflow;
102   problem->apply_outflow.qfunction_loc = Blasius_Outflow_loc;
103   problem->apply_outflow_jacobian.qfunction = Blasius_Outflow_Jacobian;
104   problem->apply_outflow_jacobian.qfunction_loc = Blasius_Outflow_Jacobian_loc;
105 
106   // CeedScalar mu = .04; // Pa s, dynamic viscosity
107   CeedScalar Uinf          = 40;   // m/s
108   CeedScalar delta0        = 4.2e-4;    // m
109   PetscReal  refine_height = 5.9e-4;    // m
110   PetscReal  growth        = 1.08; // [-]
111   PetscInt   Ndelta        = 45;   // [-]
112   PetscReal  top_angle     = 5;    // degrees
113   CeedScalar theta0        = 288.; // K
114   CeedScalar P0            = 1.01e5; // Pa
115   PetscBool  weakT         = PETSC_FALSE; // weak density or temperature
116 
117   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
118   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
119                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
120   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
121                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
122   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
123                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
124   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
125                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
126   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
127                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
128   ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge",
129                                 NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr);
130   ierr = PetscOptionsScalar("-refine_height",
131                             "Height of boundary layer mesh refinement",
132                             NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr);
133   ierr = PetscOptionsScalar("-growth",
134                             "Geometric growth rate of boundary layer mesh",
135                             NULL, growth, &growth, NULL); CHKERRQ(ierr);
136   ierr = PetscOptionsScalar("-top_angle",
137                             "Geometric top_angle rate of boundary layer mesh",
138                             NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr);
139   ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition",
140                           NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr);
141   PetscOptionsEnd();
142 
143   PetscScalar meter           = user->units->meter;
144   PetscScalar second          = user->units->second;
145   PetscScalar Kelvin          = user->units->Kelvin;
146   PetscScalar Pascal          = user->units->Pascal;
147 
148   theta0 *= Kelvin;
149   P0     *= Pascal;
150   Uinf   *= meter / second;
151   delta0 *= meter;
152 
153   ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle);
154   CHKERRQ(ierr);
155 
156   // Some properties depend on parameters from NewtonianIdealGas
157   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
158                               CEED_MEM_HOST, &newtonian_ig_ctx);
159 
160   blasius_ctx->weakT     = weakT;
161   blasius_ctx->Uinf      = Uinf;
162   blasius_ctx->delta0    = delta0;
163   blasius_ctx->theta0    = theta0;
164   blasius_ctx->P0        = P0;
165   blasius_ctx->implicit  = user->phys->implicit;
166   blasius_ctx->newtonian_ctx = *newtonian_ig_ctx;
167 
168   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
169                                   &newtonian_ig_ctx);
170 
171   CeedQFunctionContextCreate(user->ceed, &blasius_context);
172   CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST,
173                               CEED_USE_POINTER,
174                               sizeof(*blasius_ctx), blasius_ctx);
175   CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST,
176                                      FreeContextPetsc);
177 
178   problem->ics.qfunction_context = blasius_context;
179   CeedQFunctionContextReferenceCopy(blasius_context,
180                                     &problem->apply_inflow.qfunction_context);
181   CeedQFunctionContextReferenceCopy(blasius_context,
182                                     &problem->apply_inflow_jacobian.qfunction_context);
183   CeedQFunctionContextReferenceCopy(blasius_context,
184                                     &problem->apply_outflow.qfunction_context);
185   CeedQFunctionContextReferenceCopy(blasius_context,
186                                     &problem->apply_outflow_jacobian.qfunction_context);
187   if (use_stg) {
188     ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0); CHKERRQ(ierr);
189   }
190   PetscFunctionReturn(0);
191 }
192