xref: /libCEED/examples/fluids/problems/blasius.c (revision 871db79f2f3e2ed5b6659a8b21aea7cbbaad5cb5)
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"
1388626eedSJames Wright 
1488626eedSJames Wright #ifndef blasius_context_struct
1588626eedSJames Wright #define blasius_context_struct
1688626eedSJames Wright typedef struct BlasiusContext_ *BlasiusContext;
1788626eedSJames Wright struct BlasiusContext_ {
1888626eedSJames Wright   bool       implicit;  // !< Using implicit timesteping or not
1988626eedSJames Wright   CeedScalar delta0;    // !< Boundary layer height at inflow
2088626eedSJames Wright   CeedScalar Uinf;      // !< Velocity at boundary layer edge
2188626eedSJames Wright   CeedScalar P0;        // !< Pressure at outflow
2288626eedSJames Wright   CeedScalar theta0;    // !< Temperature at inflow
23*871db79fSKenneth E. Jansen   CeedInt weakT;        // !< flag to weakly set Temperature at inflow if not set weak rho instead
2488626eedSJames Wright   struct NewtonianIdealGasContext_ newtonian_ctx;
2588626eedSJames Wright };
2688626eedSJames Wright #endif
2788626eedSJames Wright 
2888626eedSJames Wright #ifndef M_PI
2988626eedSJames Wright #define M_PI    3.14159265358979323846
3088626eedSJames Wright #endif
3188626eedSJames Wright 
3288626eedSJames Wright /* \brief Modify the domain and mesh for blasius
3388626eedSJames Wright  *
3488626eedSJames Wright  * Modifies mesh such that `N` elements are within 1.2*`delta0` with a geometric
3588626eedSJames Wright  * growth ratio of `growth`. Excess elements are then geometrically distributed
3688626eedSJames Wright  * to the top surface.
3788626eedSJames Wright  *
3888626eedSJames Wright  * The top surface is also angled downwards, so that it may be used as an
3988626eedSJames Wright  * outflow. It's angle is controlled by top_angle (in units of degrees).
4088626eedSJames Wright  */
4188626eedSJames Wright PetscErrorCode modifyMesh(DM dm, PetscInt dim, PetscReal growth, PetscInt N,
4288626eedSJames Wright                           PetscReal refine_height, PetscReal top_angle) {
4388626eedSJames Wright 
4488626eedSJames Wright   PetscInt ierr, narr, ncoords;
4588626eedSJames Wright   PetscReal domain_min[3], domain_max[3], domain_size[3];
4688626eedSJames Wright   PetscScalar *arr_coords;
4788626eedSJames Wright   Vec vec_coords;
4888626eedSJames Wright   PetscFunctionBeginUser;
4988626eedSJames Wright 
5088626eedSJames Wright   PetscReal angle_coeff = tan(top_angle*(M_PI/180));
5188626eedSJames Wright 
5288626eedSJames Wright   // Get domain boundary information
5388626eedSJames Wright   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
5488626eedSJames Wright   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
5588626eedSJames Wright 
5688626eedSJames Wright   // Get coords array from DM
5788626eedSJames Wright   ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr);
5888626eedSJames Wright   ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr);
5988626eedSJames Wright   ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr);
6088626eedSJames Wright 
6188626eedSJames Wright   PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords;
6288626eedSJames Wright   ncoords = narr/dim;
6388626eedSJames Wright 
6488626eedSJames Wright   // Get mesh information
6588626eedSJames Wright   PetscInt nmax = 3, faces[3];
6688626eedSJames Wright   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
6788626eedSJames Wright                                  NULL); CHKERRQ(ierr);
6888626eedSJames Wright 
6988626eedSJames Wright   // Calculate the first element height
7088626eedSJames Wright   PetscReal dybox = domain_size[1]/faces[1];
7188626eedSJames Wright   PetscReal dy1   = refine_height*(growth-1)/(pow(growth, N)-1);
7288626eedSJames Wright 
7388626eedSJames Wright   // Calculate log of sizing outside BL
7488626eedSJames Wright   PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N);
7588626eedSJames Wright 
7688626eedSJames Wright   for(int i=0; i<ncoords; i++) {
7788626eedSJames Wright     PetscInt y_box_index = round(coords[i][1]/dybox);
7888626eedSJames Wright     if(y_box_index <= N) {
7988626eedSJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
8088626eedSJames Wright                      dy1*(pow(growth, coords[i][1]/dybox)-1)/(growth-1);
8188626eedSJames Wright     } else {
8288626eedSJames Wright       PetscInt j = y_box_index - N;
8388626eedSJames Wright       coords[i][1] = (1 - (coords[i][0]/domain_max[0])*angle_coeff) *
8488626eedSJames Wright                      exp(log(refine_height) + logdy*j);
8588626eedSJames Wright     }
8688626eedSJames Wright   }
8788626eedSJames Wright 
8888626eedSJames Wright   ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr);
8988626eedSJames Wright   ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr);
9088626eedSJames Wright 
9188626eedSJames Wright   PetscFunctionReturn(0);
9288626eedSJames Wright }
9388626eedSJames Wright 
9488626eedSJames Wright PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *setup_ctx,
9588626eedSJames Wright                           void *ctx) {
9688626eedSJames Wright 
9788626eedSJames Wright   PetscInt ierr;
9888626eedSJames Wright   ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr);
9988626eedSJames Wright   User              user = *(User *)ctx;
10088626eedSJames Wright   MPI_Comm          comm = PETSC_COMM_WORLD;
10188626eedSJames Wright   PetscFunctionBeginUser;
10288626eedSJames Wright   ierr = PetscCalloc1(1, &user->phys->blasius_ctx); CHKERRQ(ierr);
10388626eedSJames Wright 
10488626eedSJames Wright   // ------------------------------------------------------
10588626eedSJames Wright   //               SET UP Blasius
10688626eedSJames Wright   // ------------------------------------------------------
10788626eedSJames Wright   problem->ics                     = ICsBlasius;
10888626eedSJames Wright   problem->ics_loc                 = ICsBlasius_loc;
10988626eedSJames Wright   problem->apply_inflow            = Blasius_Inflow;
11088626eedSJames Wright   problem->apply_inflow_loc        = Blasius_Inflow_loc;
11188626eedSJames Wright   problem->apply_outflow           = Blasius_Outflow;
11288626eedSJames Wright   problem->apply_outflow_loc       = Blasius_Outflow_loc;
11388626eedSJames Wright   problem->setup_ctx               = SetupContext_BLASIUS;
11488626eedSJames Wright 
11588626eedSJames Wright   // CeedScalar mu = .04; // Pa s, dynamic viscosity
11688626eedSJames Wright   CeedScalar mu            = 1.8e-5;   // Pa s, dynamic viscosity
11788626eedSJames Wright   CeedScalar Uinf          = 40;   // m/s
11888626eedSJames Wright   CeedScalar delta0        = 4.2e-4;    // m
11988626eedSJames Wright   PetscReal  refine_height = 5.9e-4;    // m
12088626eedSJames Wright   PetscReal  growth        = 1.08; // [-]
12188626eedSJames Wright   PetscInt   Ndelta        = 45;   // [-]
12288626eedSJames Wright   PetscReal  top_angle     = 5;    // degrees
12388626eedSJames Wright   CeedScalar theta0        = 288.; // K
12488626eedSJames Wright   CeedScalar P0            = 1.01e5; // Pa
125*871db79fSKenneth E. Jansen   PetscBool  weakT         = PETSC_FALSE; // weak density or temperature
12688626eedSJames Wright 
12788626eedSJames Wright   PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL);
128*871db79fSKenneth E. Jansen   ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow",
129*871db79fSKenneth E. Jansen                           NULL, weakT, &weakT, NULL); CHKERRQ(ierr);
13088626eedSJames Wright   ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge",
13188626eedSJames Wright                             NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr);
13288626eedSJames Wright   ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow",
13388626eedSJames Wright                             NULL, delta0, &delta0, NULL); CHKERRQ(ierr);
13488626eedSJames Wright   ierr = PetscOptionsScalar("-theta0", "Wall temperature",
13588626eedSJames Wright                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
13688626eedSJames Wright   ierr = PetscOptionsScalar("-P0", "Pressure at outflow",
13788626eedSJames Wright                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
13888626eedSJames Wright   ierr = PetscOptionsBoundedInt("-Ndelta", "Velocity at boundary layer edge",
13988626eedSJames Wright                                 NULL, Ndelta, &Ndelta, NULL, 1); CHKERRQ(ierr);
14088626eedSJames Wright   ierr = PetscOptionsScalar("-refine_height",
14188626eedSJames Wright                             "Height of boundary layer mesh refinement",
14288626eedSJames Wright                             NULL, refine_height, &refine_height, NULL); CHKERRQ(ierr);
14388626eedSJames Wright   ierr = PetscOptionsScalar("-growth",
14488626eedSJames Wright                             "Geometric growth rate of boundary layer mesh",
14588626eedSJames Wright                             NULL, growth, &growth, NULL); CHKERRQ(ierr);
14688626eedSJames Wright   ierr = PetscOptionsScalar("-top_angle",
14788626eedSJames Wright                             "Geometric top_angle rate of boundary layer mesh",
14888626eedSJames Wright                             NULL, top_angle, &top_angle, NULL); CHKERRQ(ierr);
14988626eedSJames Wright   PetscOptionsEnd();
15088626eedSJames Wright 
15188626eedSJames Wright   PetscScalar meter           = user->units->meter;
15288626eedSJames Wright   PetscScalar second          = user->units->second;
15388626eedSJames Wright   PetscScalar Kelvin          = user->units->Kelvin;
15488626eedSJames Wright   PetscScalar Pascal          = user->units->Pascal;
15588626eedSJames Wright 
15688626eedSJames Wright   mu     *= Pascal * second;
15788626eedSJames Wright   theta0 *= Kelvin;
15888626eedSJames Wright   P0     *= Pascal;
15988626eedSJames Wright   Uinf   *= meter / second;
16088626eedSJames Wright   delta0 *= meter;
16188626eedSJames Wright 
16288626eedSJames Wright   ierr = modifyMesh(dm, problem->dim, growth, Ndelta, refine_height, top_angle);
16388626eedSJames Wright   CHKERRQ(ierr);
16488626eedSJames Wright 
165*871db79fSKenneth E. Jansen   user->phys->blasius_ctx->weakT     = !!weakT;
16688626eedSJames Wright   user->phys->blasius_ctx->Uinf      = Uinf;
16788626eedSJames Wright   user->phys->blasius_ctx->delta0    = delta0;
16888626eedSJames Wright   user->phys->blasius_ctx->theta0    = theta0;
16988626eedSJames Wright   user->phys->blasius_ctx->P0        = P0;
17088626eedSJames Wright   user->phys->blasius_ctx->implicit  = user->phys->implicit;
17188626eedSJames Wright 
17288626eedSJames Wright   user->phys->newtonian_ig_ctx->mu = mu;
17388626eedSJames Wright   user->phys->blasius_ctx->newtonian_ctx = *user->phys->newtonian_ig_ctx;
17488626eedSJames Wright   PetscFunctionReturn(0);
17588626eedSJames Wright }
17688626eedSJames Wright 
17788626eedSJames Wright PetscErrorCode SetupContext_BLASIUS(Ceed ceed, CeedData ceed_data,
17888626eedSJames Wright                                     AppCtx app_ctx, SetupContext setup_ctx, Physics phys) {
17988626eedSJames Wright   PetscFunctionBeginUser;
18088626eedSJames Wright   PetscInt ierr;
18188626eedSJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->setup_context);
18288626eedSJames Wright   CeedQFunctionContextSetData(ceed_data->setup_context, CEED_MEM_HOST,
18388626eedSJames Wright                               CEED_USE_POINTER,
18488626eedSJames Wright                               sizeof(*setup_ctx), setup_ctx);
18588626eedSJames Wright   ierr = SetupContext_NEWTONIAN_IG(ceed, ceed_data, app_ctx, setup_ctx, phys);
18688626eedSJames Wright   CHKERRQ(ierr);
18788626eedSJames Wright 
18888626eedSJames Wright   CeedQFunctionContextCreate(ceed, &ceed_data->blasius_context);
18988626eedSJames Wright   CeedQFunctionContextSetData(ceed_data->blasius_context, CEED_MEM_HOST,
19088626eedSJames Wright                               CEED_USE_POINTER, sizeof(*phys->blasius_ctx), phys->blasius_ctx);
19188626eedSJames Wright   phys->has_neumann = PETSC_TRUE;
19288626eedSJames Wright   if (ceed_data->qf_ics)
19388626eedSJames Wright     CeedQFunctionSetContext(ceed_data->qf_ics, ceed_data->blasius_context);
19488626eedSJames Wright   if (ceed_data->qf_apply_inflow)
19588626eedSJames Wright     CeedQFunctionSetContext(ceed_data->qf_apply_inflow, ceed_data->blasius_context);
19688626eedSJames Wright   if (ceed_data->qf_apply_outflow)
19788626eedSJames Wright     CeedQFunctionSetContext(ceed_data->qf_apply_outflow,
19888626eedSJames Wright                             ceed_data->blasius_context);
19988626eedSJames Wright   PetscFunctionReturn(0);
20088626eedSJames Wright }
20188626eedSJames Wright 
202