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 152518f336SLeila Ghaffari PetscErrorCode CompressibleBlasiusResidual(SNES snes, Vec X, Vec R, void *ctx) { 162518f336SLeila Ghaffari const BlasiusContext blasius = (BlasiusContext)ctx; 172518f336SLeila Ghaffari const PetscScalar *Tf, *Th; // Chebyshev coefficients 182518f336SLeila Ghaffari PetscScalar *r, f[4], h[4]; 192518f336SLeila Ghaffari PetscInt N = blasius->n_cheb; 20*fb455ff0SLeila Ghaffari PetscScalar Ma = Mach(&blasius->newtonian_ctx, blasius->T_inf, blasius->U_inf), 212518f336SLeila Ghaffari Pr = Prandtl(&blasius->newtonian_ctx), 222518f336SLeila Ghaffari gamma = HeatCapacityRatio(&blasius->newtonian_ctx); 232518f336SLeila Ghaffari PetscFunctionBegin; 242518f336SLeila Ghaffari PetscCall(VecGetArrayRead(X, &Tf)); 252518f336SLeila Ghaffari Th = Tf + N; 262518f336SLeila Ghaffari PetscCall(VecGetArray(R, &r)); 272518f336SLeila Ghaffari 282518f336SLeila Ghaffari // Left boundary conditions f = f' = 0 292518f336SLeila Ghaffari ChebyshevEval(N, Tf, -1., blasius->eta_max, f); 302518f336SLeila Ghaffari r[0] = f[0]; 312518f336SLeila Ghaffari r[1] = f[1]; 322518f336SLeila Ghaffari 332518f336SLeila Ghaffari // f - right end boundary condition 342518f336SLeila Ghaffari ChebyshevEval(N, Tf, 1., blasius->eta_max, f); 352518f336SLeila Ghaffari r[2] = f[1] - 1.; 362518f336SLeila Ghaffari 372518f336SLeila Ghaffari for (int i=0; i<N-3; i++) { 382518f336SLeila Ghaffari ChebyshevEval(N, Tf, blasius->X[i], blasius->eta_max, f); 392518f336SLeila Ghaffari r[3+i] = 2*f[3] + f[2] * f[0]; 402518f336SLeila Ghaffari ChebyshevEval(N-1, Th, blasius->X[i], blasius->eta_max, h); 412518f336SLeila Ghaffari r[N+2+i] = h[2] + Pr * f[0] * h[1] + 422518f336SLeila Ghaffari Pr * (gamma - 1) * PetscSqr(Ma * f[2]); 432518f336SLeila Ghaffari } 442518f336SLeila Ghaffari 452518f336SLeila Ghaffari // h - left end boundary condition 462518f336SLeila Ghaffari ChebyshevEval(N-1, Th, -1., blasius->eta_max, h); 47*fb455ff0SLeila Ghaffari r[N] = h[0] - blasius->T_wall / blasius->T_inf; 482518f336SLeila Ghaffari 492518f336SLeila Ghaffari // h - right end boundary condition 502518f336SLeila Ghaffari ChebyshevEval(N-1, Th, 1., blasius->eta_max, h); 512518f336SLeila Ghaffari r[N+1] = h[0] - 1.; 522518f336SLeila Ghaffari 532518f336SLeila Ghaffari // Restore vectors 542518f336SLeila Ghaffari PetscCall(VecRestoreArrayRead(X, &Tf)); 552518f336SLeila Ghaffari PetscCall(VecRestoreArray(R, &r)); 562518f336SLeila Ghaffari PetscFunctionReturn(0); 572518f336SLeila Ghaffari } 582518f336SLeila Ghaffari 592518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 602518f336SLeila Ghaffari SNES snes; 612518f336SLeila Ghaffari Vec sol, res; 622518f336SLeila Ghaffari PetscReal *w; 632518f336SLeila Ghaffari PetscInt N = blasius->n_cheb; 642518f336SLeila Ghaffari const PetscScalar *cheb_coefs; 652518f336SLeila Ghaffari PetscFunctionBegin; 662518f336SLeila Ghaffari PetscCall(PetscMalloc2(N-3, &blasius->X, N-3, &w)); 672518f336SLeila Ghaffari PetscCall(PetscDTGaussQuadrature(N-3, -1., 1., blasius->X, w)); 682518f336SLeila Ghaffari PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 692518f336SLeila Ghaffari PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 702518f336SLeila Ghaffari PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2*N-1)); 712518f336SLeila Ghaffari PetscCall(VecSetFromOptions(sol)); 722518f336SLeila Ghaffari PetscCall(VecDuplicate(sol, &res)); 732518f336SLeila Ghaffari PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 742518f336SLeila Ghaffari PetscCall(SNESSetFromOptions(snes)); 752518f336SLeila Ghaffari PetscCall(SNESSolve(snes, NULL, sol)); 762518f336SLeila Ghaffari PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 772518f336SLeila Ghaffari for (int i=0; i<N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 782518f336SLeila Ghaffari for (int i=0; i<N-1; i++) blasius->Th_cheb[i] = cheb_coefs[i+N]; 792518f336SLeila Ghaffari PetscCall(PetscFree2(blasius->X, w)); 802518f336SLeila Ghaffari PetscCall(VecDestroy(&sol)); 812518f336SLeila Ghaffari PetscCall(VecDestroy(&res)); 822518f336SLeila Ghaffari PetscCall(SNESDestroy(&snes)); 832518f336SLeila Ghaffari PetscFunctionReturn(0); 842518f336SLeila Ghaffari } 852518f336SLeila Ghaffari 86d2714514SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, 87d2714514SJames Wright const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, 88d2714514SJames Wright PetscInt *nynodes) { 89d2714514SJames Wright PetscErrorCode ierr; 90d2714514SJames Wright PetscInt ndims, dims[2]; 91d2714514SJames Wright FILE *fp; 92d2714514SJames Wright const PetscInt char_array_len = 512; 93d2714514SJames Wright char line[char_array_len]; 94d2714514SJames Wright char **array; 95d2714514SJames Wright PetscReal *node_locs; 96d2714514SJames Wright PetscFunctionBeginUser; 97d2714514SJames Wright 98d2714514SJames Wright ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); 99d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 100d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 101d2714514SJames Wright 102d2714514SJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 103d2714514SJames Wright if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 104d2714514SJames Wright *nynodes = dims[0]; 105d2714514SJames Wright ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr); 106d2714514SJames Wright 107d2714514SJames Wright for (PetscInt i=0; i<dims[0]; i++) { 108d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 109d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 110d2714514SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 111990fdeb6SJeremy L Thompson "Line %" PetscInt_FMT" of %s does not contain enough columns (%" 112990fdeb6SJeremy L Thompson PetscInt_FMT" instead of %" PetscInt_FMT ")", i, 113d2714514SJames Wright path, ndims, dims[1]); 114d2714514SJames Wright 115d2714514SJames Wright node_locs[i] = (PetscReal) atof(array[0]); 116d2714514SJames Wright } 117d2714514SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 118d2714514SJames Wright *pynodes = node_locs; 119d2714514SJames Wright PetscFunctionReturn(0); 120d2714514SJames Wright } 121d2714514SJames Wright 12288626eedSJames Wright /* \brief Modify the domain and mesh for blasius 12388626eedSJames Wright * 124ba6664aeSJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 125ba6664aeSJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 126ba6664aeSJames Wright * linearly in logspace to the top surface. 12788626eedSJames Wright * 12888626eedSJames Wright * The top surface is also angled downwards, so that it may be used as an 129ba6664aeSJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 130d2714514SJames Wright * 131d2714514SJames Wright * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 132798d42b8SJames Wright * locations. If it is NULL, then the modified coordinate values will be set in 133798d42b8SJames Wright * the array, along with `num_node_locs`. 13488626eedSJames Wright */ 135d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 136d2714514SJames Wright PetscReal growth, PetscInt N, 137d2714514SJames Wright PetscReal refine_height, PetscReal top_angle, 138798d42b8SJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 13988626eedSJames Wright PetscInt ierr, narr, ncoords; 14088626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 14188626eedSJames Wright PetscScalar *arr_coords; 14288626eedSJames Wright Vec vec_coords; 14388626eedSJames Wright PetscFunctionBeginUser; 14488626eedSJames Wright 14588626eedSJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 14688626eedSJames Wright 14788626eedSJames Wright // Get domain boundary information 14888626eedSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 149ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 15088626eedSJames Wright 15188626eedSJames Wright // Get coords array from DM 15288626eedSJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 15388626eedSJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 15488626eedSJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 15588626eedSJames Wright 15688626eedSJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 15788626eedSJames Wright ncoords = narr/dim; 15888626eedSJames Wright 15988626eedSJames Wright // Get mesh information 16088626eedSJames Wright PetscInt nmax = 3, faces[3]; 16188626eedSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 16288626eedSJames Wright NULL); CHKERRQ(ierr); 163d2714514SJames Wright // Get element size of the box mesh, for indexing each node 164d2714514SJames Wright const PetscReal dybox = domain_size[1]/faces[1]; 16588626eedSJames Wright 166798d42b8SJames Wright if (!*node_locs) { 16788626eedSJames Wright // Calculate the first element height 16888626eedSJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 16988626eedSJames Wright 17088626eedSJames Wright // Calculate log of sizing outside BL 17188626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 17288626eedSJames Wright 173798d42b8SJames Wright *num_node_locs = faces[1] + 1; 174798d42b8SJames Wright PetscReal *temp_node_locs; 175798d42b8SJames Wright ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr); 176798d42b8SJames Wright 177ba6664aeSJames Wright for (PetscInt i=0; i<ncoords; i++) { 17888626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 17988626eedSJames Wright if (y_box_index <= N) { 180855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 181d2714514SJames Wright * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 18288626eedSJames Wright } else { 18388626eedSJames Wright PetscInt j = y_box_index - N; 184855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 185d2714514SJames Wright * exp(log(refine_height) + logdy*j); 186d2714514SJames Wright } 187798d42b8SJames Wright if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) 188798d42b8SJames Wright temp_node_locs[y_box_index] = coords[i][1]; 189d2714514SJames Wright } 190798d42b8SJames Wright 191798d42b8SJames Wright *node_locs = temp_node_locs; 192d2714514SJames Wright } else { 193d2714514SJames Wright // Error checking 194798d42b8SJames Wright if (*num_node_locs < faces[1] +1) 195d2714514SJames Wright SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 196d2714514SJames Wright "There are %d + 1 nodes, but only %d locations given", 197798d42b8SJames Wright faces[1]+1, *num_node_locs); 198798d42b8SJames Wright if (*num_node_locs > faces[1] +1) { 199d2714514SJames Wright ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 200c10dcd27SJames Wright "than the mesh has nodes (%d). This maybe unintended.\n", 201798d42b8SJames Wright *num_node_locs, faces[1]+1); CHKERRQ(ierr); 202d2714514SJames Wright } 203798d42b8SJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 204d2714514SJames Wright 205d2714514SJames Wright for (PetscInt i=0; i<ncoords; i++) { 206d2714514SJames Wright // Determine which y-node we're at 207d2714514SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 208855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y) 209798d42b8SJames Wright * (*node_locs)[y_box_index]; 21088626eedSJames Wright } 21188626eedSJames Wright } 21288626eedSJames Wright 21388626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 21488626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 21588626eedSJames Wright 21688626eedSJames Wright PetscFunctionReturn(0); 21788626eedSJames Wright } 21888626eedSJames Wright 219a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 22088626eedSJames Wright 22188626eedSJames Wright PetscInt ierr; 22288626eedSJames Wright User user = *(User *)ctx; 22388626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 224ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 225841e4c73SJed Brown BlasiusContext blasius_ctx; 226841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 227841e4c73SJed Brown CeedQFunctionContext blasius_context; 228841e4c73SJed Brown 22988626eedSJames Wright PetscFunctionBeginUser; 230a0add3c9SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 231841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 23288626eedSJames Wright 23388626eedSJames Wright // ------------------------------------------------------ 23488626eedSJames Wright // SET UP Blasius 23588626eedSJames Wright // ------------------------------------------------------ 23691e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 23791e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 23888626eedSJames Wright 239*fb455ff0SLeila Ghaffari CeedScalar U_inf = 40; // m/s 240*fb455ff0SLeila Ghaffari CeedScalar T_inf = 288.; // K 2412518f336SLeila Ghaffari CeedScalar T_wall = 400.; // K 2422518f336SLeila Ghaffari CeedScalar delta0 = 4.2e-3; // m 24388626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 2442518f336SLeila Ghaffari CeedInt N = 20; // Number of Chebyshev terms 245871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 2467b8d3891SJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 2477b8d3891SJames Wright PetscReal mesh_growth = 1.08; // [-] 2487b8d3891SJames Wright PetscInt mesh_Ndelta = 45; // [-] 2497b8d3891SJames Wright PetscReal mesh_top_angle = 5; // degrees 250d2714514SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 25188626eedSJames Wright 2527b37518fSJames Wright PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 253871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 254871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 255*fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-velocity_infinity", 256*fb455ff0SLeila Ghaffari "Velocity at boundary layer edge", 257*fb455ff0SLeila Ghaffari NULL, U_inf, &U_inf, NULL); CHKERRQ(ierr); 258*fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-temperature_infinity", 259*fb455ff0SLeila Ghaffari "Temperature at boundary layer edge", 260*fb455ff0SLeila Ghaffari NULL, T_inf, &T_inf, NULL); CHKERRQ(ierr); 261*fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-temperature_wall", "Temperature at wall", 2622518f336SLeila Ghaffari NULL, T_wall, &T_wall, NULL); CHKERRQ(ierr); 26388626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 26488626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 26588626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 26688626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 2672518f336SLeila Ghaffari ierr = PetscOptionsInt("-N_Chebyshev", "Number of Chebyshev terms", 2682518f336SLeila Ghaffari NULL, N, &N, NULL); CHKERRQ(ierr); 2697b8d3891SJames Wright ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 2707b8d3891SJames Wright "Velocity at boundary layer edge", 2717b8d3891SJames Wright NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 2727b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_refine_height", 27388626eedSJames Wright "Height of boundary layer mesh refinement", 2747b8d3891SJames Wright NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 2757b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_growth", 27688626eedSJames Wright "Geometric growth rate of boundary layer mesh", 2777b8d3891SJames Wright NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 2787b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_top_angle", 27988626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 2807b8d3891SJames Wright NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 281d2714514SJames Wright ierr = PetscOptionsString("-platemesh_y_node_locs_path", 282d2714514SJames Wright "Path to file with y node locations. " 283d2714514SJames Wright "If empty, will use the algorithmic mesh warping.", NULL, 284d2714514SJames Wright mesh_ynodes_path, mesh_ynodes_path, 285d2714514SJames Wright sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 286ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 287ba6664aeSJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 28888626eedSJames Wright PetscOptionsEnd(); 28988626eedSJames Wright 29088626eedSJames Wright PetscScalar meter = user->units->meter; 29188626eedSJames Wright PetscScalar second = user->units->second; 29288626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 29388626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 29488626eedSJames Wright 295*fb455ff0SLeila Ghaffari T_inf *= Kelvin; 2962518f336SLeila Ghaffari T_wall *= Kelvin; 29788626eedSJames Wright P0 *= Pascal; 298*fb455ff0SLeila Ghaffari U_inf *= meter / second; 29988626eedSJames Wright delta0 *= meter; 30088626eedSJames Wright 301d2714514SJames Wright PetscReal *mesh_ynodes = NULL; 302d2714514SJames Wright PetscInt mesh_nynodes = 0; 303d2714514SJames Wright if (strcmp(mesh_ynodes_path, "")) { 304d2714514SJames Wright ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 30588626eedSJames Wright CHKERRQ(ierr); 306d2714514SJames Wright } 307d2714514SJames Wright ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 308798d42b8SJames Wright mesh_refine_height, mesh_top_angle, &mesh_ynodes, 309798d42b8SJames Wright &mesh_nynodes); CHKERRQ(ierr); 31088626eedSJames Wright 311841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 312841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 313841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 31488626eedSJames Wright 315ba6664aeSJames Wright blasius_ctx->weakT = weakT; 316*fb455ff0SLeila Ghaffari blasius_ctx->U_inf = U_inf; 317*fb455ff0SLeila Ghaffari blasius_ctx->T_inf = T_inf; 3182518f336SLeila Ghaffari blasius_ctx->T_wall = T_wall; 319841e4c73SJed Brown blasius_ctx->delta0 = delta0; 320841e4c73SJed Brown blasius_ctx->P0 = P0; 3212518f336SLeila Ghaffari blasius_ctx->n_cheb = N; 32230e9fa81SJames Wright newtonian_ig_ctx->P0 = P0; 323841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 324841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 325ba6664aeSJames Wright 326f1122ed0SJames Wright { 3272518f336SLeila Ghaffari PetscReal domain_min[3], domain_max[3]; 3282518f336SLeila Ghaffari ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 329f1122ed0SJames Wright blasius_ctx->x_inflow = domain_min[0]; 3302518f336SLeila Ghaffari blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 331f1122ed0SJames Wright } 3322518f336SLeila Ghaffari PetscCall(PetscMalloc2(blasius_ctx->n_cheb, &blasius_ctx->Tf_cheb, 3332518f336SLeila Ghaffari blasius_ctx->n_cheb-1, &blasius_ctx->Th_cheb)); 3342518f336SLeila Ghaffari PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 335f1122ed0SJames Wright 336841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 337841e4c73SJed Brown &newtonian_ig_ctx); 33888626eedSJames Wright 339841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 3402518f336SLeila Ghaffari CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, 341841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 342841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 343841e4c73SJed Brown FreeContextPetsc); 34488626eedSJames Wright 345b77c53c9SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 346841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 347ba6664aeSJames Wright if (use_stg) { 348*fb455ff0SLeila Ghaffari ierr = SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, 3494e139266SJames Wright mesh_nynodes); CHKERRQ(ierr); 35065dd5cafSJames Wright } else { 35165dd5cafSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 35265dd5cafSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 353fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 354fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 35565dd5cafSJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 35665dd5cafSJames Wright &problem->apply_inflow.qfunction_context); 357fc02c281SJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 358fc02c281SJames Wright &problem->apply_inflow_jacobian.qfunction_context); 359ba6664aeSJames Wright } 3604e139266SJames Wright ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 36188626eedSJames Wright PetscFunctionReturn(0); 36288626eedSJames Wright } 363