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; 20fb455ff0SLeila 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 ChebyshevEval(N-1, Th, blasius->X[i], blasius->eta_max, h); 4007d14e58SLeila Ghaffari // mu and rho generally depend on h. We naively assume constant mu. 4107d14e58SLeila Ghaffari // For an ideal gas at constant pressure, density is inversely proportional to enthalpy. 4207d14e58SLeila Ghaffari // The *_tilde values are *relative* to their freestream values, and we proved first derivatives here. 4307d14e58SLeila Ghaffari const PetscScalar mu_tilde[2] = {1, 0}; 4407d14e58SLeila Ghaffari const PetscScalar rho_tilde[2] = {1 / h[0], -h[1] / PetscSqr(h[0])}; 4507d14e58SLeila Ghaffari const PetscScalar mu_rho_tilde[2] = { 4607d14e58SLeila Ghaffari mu_tilde[0] *rho_tilde[0], 4707d14e58SLeila Ghaffari mu_tilde[1] *rho_tilde[0] + mu_tilde[0] *rho_tilde[1], 4807d14e58SLeila Ghaffari }; 4907d14e58SLeila Ghaffari r[3+i] = 2*(mu_rho_tilde[0] * f[3] + mu_rho_tilde[1] * f[2]) + f[2] * f[0]; 5007d14e58SLeila Ghaffari r[N+2+i] = (mu_rho_tilde[0] * h[2] + mu_rho_tilde[1] * h[1]) + Pr * f[0] * h[1] 5107d14e58SLeila Ghaffari + Pr * (gamma - 1) * mu_rho_tilde[0] * PetscSqr(Ma * f[2]); 522518f336SLeila Ghaffari } 532518f336SLeila Ghaffari 542518f336SLeila Ghaffari // h - left end boundary condition 552518f336SLeila Ghaffari ChebyshevEval(N-1, Th, -1., blasius->eta_max, h); 56fb455ff0SLeila Ghaffari r[N] = h[0] - blasius->T_wall / blasius->T_inf; 572518f336SLeila Ghaffari 582518f336SLeila Ghaffari // h - right end boundary condition 592518f336SLeila Ghaffari ChebyshevEval(N-1, Th, 1., blasius->eta_max, h); 602518f336SLeila Ghaffari r[N+1] = h[0] - 1.; 612518f336SLeila Ghaffari 622518f336SLeila Ghaffari // Restore vectors 632518f336SLeila Ghaffari PetscCall(VecRestoreArrayRead(X, &Tf)); 642518f336SLeila Ghaffari PetscCall(VecRestoreArray(R, &r)); 652518f336SLeila Ghaffari PetscFunctionReturn(0); 662518f336SLeila Ghaffari } 672518f336SLeila Ghaffari 682518f336SLeila Ghaffari PetscErrorCode ComputeChebyshevCoefficients(BlasiusContext blasius) { 692518f336SLeila Ghaffari SNES snes; 702518f336SLeila Ghaffari Vec sol, res; 712518f336SLeila Ghaffari PetscReal *w; 722518f336SLeila Ghaffari PetscInt N = blasius->n_cheb; 7307d14e58SLeila Ghaffari SNESConvergedReason reason; 742518f336SLeila Ghaffari const PetscScalar *cheb_coefs; 752518f336SLeila Ghaffari PetscFunctionBegin; 7607d14e58SLeila Ghaffari 7707d14e58SLeila Ghaffari // Allocate memory 782518f336SLeila Ghaffari PetscCall(PetscMalloc2(N-3, &blasius->X, N-3, &w)); 792518f336SLeila Ghaffari PetscCall(PetscDTGaussQuadrature(N-3, -1., 1., blasius->X, w)); 8007d14e58SLeila Ghaffari 8107d14e58SLeila Ghaffari // Snes solve 822518f336SLeila Ghaffari PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 832518f336SLeila Ghaffari PetscCall(VecCreate(PETSC_COMM_SELF, &sol)); 842518f336SLeila Ghaffari PetscCall(VecSetSizes(sol, PETSC_DECIDE, 2*N-1)); 852518f336SLeila Ghaffari PetscCall(VecSetFromOptions(sol)); 8607d14e58SLeila Ghaffari // Constant relative enthalpy 1 as initial guess 8707d14e58SLeila Ghaffari PetscCall(VecSetValue(sol, N, 1., INSERT_VALUES)); 882518f336SLeila Ghaffari PetscCall(VecDuplicate(sol, &res)); 892518f336SLeila Ghaffari PetscCall(SNESSetFunction(snes, res, CompressibleBlasiusResidual, blasius)); 9007d14e58SLeila Ghaffari PetscCall(SNESSetOptionsPrefix(snes, "chebyshev_")); 912518f336SLeila Ghaffari PetscCall(SNESSetFromOptions(snes)); 922518f336SLeila Ghaffari PetscCall(SNESSolve(snes, NULL, sol)); 9307d14e58SLeila Ghaffari PetscCall(SNESGetConvergedReason(snes, &reason)); 9407d14e58SLeila Ghaffari if (reason < 0) 9507d14e58SLeila Ghaffari SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_CONV_FAILED, 9607d14e58SLeila Ghaffari "The Chebyshev solve failed.\n"); 9707d14e58SLeila Ghaffari 9807d14e58SLeila Ghaffari // Assign Chebyshev coefficients 992518f336SLeila Ghaffari PetscCall(VecGetArrayRead(sol, &cheb_coefs)); 1002518f336SLeila Ghaffari for (int i=0; i<N; i++) blasius->Tf_cheb[i] = cheb_coefs[i]; 1012518f336SLeila Ghaffari for (int i=0; i<N-1; i++) blasius->Th_cheb[i] = cheb_coefs[i+N]; 10207d14e58SLeila Ghaffari 10307d14e58SLeila Ghaffari // Destroy objects 1042518f336SLeila Ghaffari PetscCall(PetscFree2(blasius->X, w)); 1052518f336SLeila Ghaffari PetscCall(VecDestroy(&sol)); 1062518f336SLeila Ghaffari PetscCall(VecDestroy(&res)); 1072518f336SLeila Ghaffari PetscCall(SNESDestroy(&snes)); 1082518f336SLeila Ghaffari PetscFunctionReturn(0); 1092518f336SLeila Ghaffari } 1102518f336SLeila Ghaffari 111d2714514SJames Wright static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, 112d2714514SJames Wright const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, 113d2714514SJames Wright PetscInt *nynodes) { 114d2714514SJames Wright PetscErrorCode ierr; 115d2714514SJames Wright PetscInt ndims, dims[2]; 116d2714514SJames Wright FILE *fp; 117d2714514SJames Wright const PetscInt char_array_len = 512; 118d2714514SJames Wright char line[char_array_len]; 119d2714514SJames Wright char **array; 120d2714514SJames Wright PetscReal *node_locs; 121d2714514SJames Wright PetscFunctionBeginUser; 122d2714514SJames Wright 123d2714514SJames Wright ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); 124d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 125d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 126d2714514SJames Wright 127d2714514SJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 128d2714514SJames Wright if (ndims<2) dims[1] = 1; // Assume 1 column of data is not otherwise specified 129d2714514SJames Wright *nynodes = dims[0]; 130d2714514SJames Wright ierr = PetscMalloc1(*nynodes, &node_locs); CHKERRQ(ierr); 131d2714514SJames Wright 132d2714514SJames Wright for (PetscInt i=0; i<dims[0]; i++) { 133d2714514SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 134d2714514SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 135d2714514SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 136990fdeb6SJeremy L Thompson "Line %" PetscInt_FMT" of %s does not contain enough columns (%" 137990fdeb6SJeremy L Thompson PetscInt_FMT" instead of %" PetscInt_FMT ")", i, 138d2714514SJames Wright path, ndims, dims[1]); 139d2714514SJames Wright 140d2714514SJames Wright node_locs[i] = (PetscReal) atof(array[0]); 141d2714514SJames Wright } 142d2714514SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 143d2714514SJames Wright *pynodes = node_locs; 144d2714514SJames Wright PetscFunctionReturn(0); 145d2714514SJames Wright } 146d2714514SJames Wright 14788626eedSJames Wright /* \brief Modify the domain and mesh for blasius 14888626eedSJames Wright * 149ba6664aeSJames Wright * Modifies mesh such that `N` elements are within `refine_height` with a 150ba6664aeSJames Wright * geometric growth ratio of `growth`. Excess elements are then distributed 151ba6664aeSJames Wright * linearly in logspace to the top surface. 15288626eedSJames Wright * 15388626eedSJames Wright * The top surface is also angled downwards, so that it may be used as an 154ba6664aeSJames Wright * outflow. It's angle is controlled by `top_angle` (in units of degrees). 155d2714514SJames Wright * 156d2714514SJames Wright * If `node_locs` is not NULL, then the nodes will be placed at `node_locs` 157798d42b8SJames Wright * locations. If it is NULL, then the modified coordinate values will be set in 158798d42b8SJames Wright * the array, along with `num_node_locs`. 15988626eedSJames Wright */ 160d2714514SJames Wright static PetscErrorCode ModifyMesh(MPI_Comm comm, DM dm, PetscInt dim, 161d2714514SJames Wright PetscReal growth, PetscInt N, 162d2714514SJames Wright PetscReal refine_height, PetscReal top_angle, 163798d42b8SJames Wright PetscReal *node_locs[], PetscInt *num_node_locs) { 16488626eedSJames Wright PetscInt ierr, narr, ncoords; 16588626eedSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 16688626eedSJames Wright PetscScalar *arr_coords; 16788626eedSJames Wright Vec vec_coords; 16888626eedSJames Wright PetscFunctionBeginUser; 16988626eedSJames Wright 17088626eedSJames Wright PetscReal angle_coeff = tan(top_angle*(M_PI/180)); 17188626eedSJames Wright 17288626eedSJames Wright // Get domain boundary information 17388626eedSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 174ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 17588626eedSJames Wright 17688626eedSJames Wright // Get coords array from DM 17788626eedSJames Wright ierr = DMGetCoordinatesLocal(dm, &vec_coords); CHKERRQ(ierr); 17888626eedSJames Wright ierr = VecGetLocalSize(vec_coords, &narr); CHKERRQ(ierr); 17988626eedSJames Wright ierr = VecGetArray(vec_coords, &arr_coords); CHKERRQ(ierr); 18088626eedSJames Wright 18188626eedSJames Wright PetscScalar (*coords)[dim] = (PetscScalar(*)[dim]) arr_coords; 18288626eedSJames Wright ncoords = narr/dim; 18388626eedSJames Wright 18488626eedSJames Wright // Get mesh information 18588626eedSJames Wright PetscInt nmax = 3, faces[3]; 18688626eedSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 18788626eedSJames Wright NULL); CHKERRQ(ierr); 188d2714514SJames Wright // Get element size of the box mesh, for indexing each node 189d2714514SJames Wright const PetscReal dybox = domain_size[1]/faces[1]; 19088626eedSJames Wright 191798d42b8SJames Wright if (!*node_locs) { 19288626eedSJames Wright // Calculate the first element height 19388626eedSJames Wright PetscReal dy1 = refine_height*(growth-1)/(pow(growth, N)-1); 19488626eedSJames Wright 19588626eedSJames Wright // Calculate log of sizing outside BL 19688626eedSJames Wright PetscReal logdy = (log(domain_max[1]) - log(refine_height)) / (faces[1] - N); 19788626eedSJames Wright 198798d42b8SJames Wright *num_node_locs = faces[1] + 1; 199798d42b8SJames Wright PetscReal *temp_node_locs; 200798d42b8SJames Wright ierr = PetscMalloc1(*num_node_locs, &temp_node_locs); CHKERRQ(ierr); 201798d42b8SJames Wright 202ba6664aeSJames Wright for (PetscInt i=0; i<ncoords; i++) { 20388626eedSJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 20488626eedSJames Wright if (y_box_index <= N) { 205855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 206d2714514SJames Wright * dy1 * (pow(growth, coords[i][1]/dybox)-1)/(growth-1); 20788626eedSJames Wright } else { 20888626eedSJames Wright PetscInt j = y_box_index - N; 209855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/domain_max[1]) 210d2714514SJames Wright * exp(log(refine_height) + logdy*j); 211d2714514SJames Wright } 212798d42b8SJames Wright if (coords[i][0] == domain_min[0] && coords[i][2] == domain_min[2]) 213798d42b8SJames Wright temp_node_locs[y_box_index] = coords[i][1]; 214d2714514SJames Wright } 215798d42b8SJames Wright 216798d42b8SJames Wright *node_locs = temp_node_locs; 217d2714514SJames Wright } else { 218d2714514SJames Wright // Error checking 219798d42b8SJames Wright if (*num_node_locs < faces[1] +1) 220d2714514SJames Wright SETERRQ(comm, -1, "The y_node_locs_path has too few locations; " 221d2714514SJames Wright "There are %d + 1 nodes, but only %d locations given", 222798d42b8SJames Wright faces[1]+1, *num_node_locs); 223798d42b8SJames Wright if (*num_node_locs > faces[1] +1) { 224d2714514SJames Wright ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " 225c10dcd27SJames Wright "than the mesh has nodes (%d). This maybe unintended.\n", 226798d42b8SJames Wright *num_node_locs, faces[1]+1); CHKERRQ(ierr); 227d2714514SJames Wright } 228798d42b8SJames Wright PetscScalar max_y = (*node_locs)[faces[1]]; 229d2714514SJames Wright 230d2714514SJames Wright for (PetscInt i=0; i<ncoords; i++) { 231d2714514SJames Wright // Determine which y-node we're at 232d2714514SJames Wright PetscInt y_box_index = round(coords[i][1]/dybox); 233855641eaSJames Wright coords[i][1] = (1 - (coords[i][0] - domain_min[0])*angle_coeff/max_y) 234798d42b8SJames Wright * (*node_locs)[y_box_index]; 23588626eedSJames Wright } 23688626eedSJames Wright } 23788626eedSJames Wright 23888626eedSJames Wright ierr = VecRestoreArray(vec_coords, &arr_coords); CHKERRQ(ierr); 23988626eedSJames Wright ierr = DMSetCoordinatesLocal(dm, vec_coords); CHKERRQ(ierr); 24088626eedSJames Wright 24188626eedSJames Wright PetscFunctionReturn(0); 24288626eedSJames Wright } 24388626eedSJames Wright 244a0add3c9SJed Brown PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx) { 24588626eedSJames Wright 24688626eedSJames Wright PetscInt ierr; 24788626eedSJames Wright User user = *(User *)ctx; 24888626eedSJames Wright MPI_Comm comm = PETSC_COMM_WORLD; 249ba6664aeSJames Wright PetscBool use_stg = PETSC_FALSE; 250841e4c73SJed Brown BlasiusContext blasius_ctx; 251841e4c73SJed Brown NewtonianIdealGasContext newtonian_ig_ctx; 252841e4c73SJed Brown CeedQFunctionContext blasius_context; 253841e4c73SJed Brown 25488626eedSJames Wright PetscFunctionBeginUser; 255a0add3c9SJed Brown ierr = NS_NEWTONIAN_IG(problem, dm, ctx); CHKERRQ(ierr); 256841e4c73SJed Brown ierr = PetscCalloc1(1, &blasius_ctx); CHKERRQ(ierr); 25788626eedSJames Wright 25888626eedSJames Wright // ------------------------------------------------------ 25988626eedSJames Wright // SET UP Blasius 26088626eedSJames Wright // ------------------------------------------------------ 26191e5af17SJed Brown problem->ics.qfunction = ICsBlasius; 26291e5af17SJed Brown problem->ics.qfunction_loc = ICsBlasius_loc; 26388626eedSJames Wright 264fb455ff0SLeila Ghaffari CeedScalar U_inf = 40; // m/s 265fb455ff0SLeila Ghaffari CeedScalar T_inf = 288.; // K 26607d14e58SLeila Ghaffari CeedScalar T_wall = 288.; // K 2672518f336SLeila Ghaffari CeedScalar delta0 = 4.2e-3; // m 26888626eedSJames Wright CeedScalar P0 = 1.01e5; // Pa 2692518f336SLeila Ghaffari CeedInt N = 20; // Number of Chebyshev terms 270871db79fSKenneth E. Jansen PetscBool weakT = PETSC_FALSE; // weak density or temperature 2717b8d3891SJames Wright PetscReal mesh_refine_height = 5.9e-4; // m 2727b8d3891SJames Wright PetscReal mesh_growth = 1.08; // [-] 2737b8d3891SJames Wright PetscInt mesh_Ndelta = 45; // [-] 2747b8d3891SJames Wright PetscReal mesh_top_angle = 5; // degrees 275d2714514SJames Wright char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; 27688626eedSJames Wright 2777b37518fSJames Wright PetscOptionsBegin(comm, NULL, "Options for BLASIUS problem", NULL); 278871db79fSKenneth E. Jansen ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", 279871db79fSKenneth E. Jansen NULL, weakT, &weakT, NULL); CHKERRQ(ierr); 280fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-velocity_infinity", 281fb455ff0SLeila Ghaffari "Velocity at boundary layer edge", 282fb455ff0SLeila Ghaffari NULL, U_inf, &U_inf, NULL); CHKERRQ(ierr); 283fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-temperature_infinity", 284fb455ff0SLeila Ghaffari "Temperature at boundary layer edge", 285fb455ff0SLeila Ghaffari NULL, T_inf, &T_inf, NULL); CHKERRQ(ierr); 286fb455ff0SLeila Ghaffari ierr = PetscOptionsScalar("-temperature_wall", "Temperature at wall", 2872518f336SLeila Ghaffari NULL, T_wall, &T_wall, NULL); CHKERRQ(ierr); 28888626eedSJames Wright ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", 28988626eedSJames Wright NULL, delta0, &delta0, NULL); CHKERRQ(ierr); 29088626eedSJames Wright ierr = PetscOptionsScalar("-P0", "Pressure at outflow", 29188626eedSJames Wright NULL, P0, &P0, NULL); CHKERRQ(ierr); 29207d14e58SLeila Ghaffari ierr = PetscOptionsInt("-n_chebyshev", "Number of Chebyshev terms", 2932518f336SLeila Ghaffari NULL, N, &N, NULL); CHKERRQ(ierr); 29407d14e58SLeila Ghaffari PetscCheck(3 <= N && N <= BLASIUS_MAX_N_CHEBYSHEV, 29507d14e58SLeila Ghaffari comm, PETSC_ERR_ARG_OUTOFRANGE, 29607d14e58SLeila Ghaffari "-n_chebyshev %" PetscInt_FMT " must be in range [3, %d]", N, 29707d14e58SLeila Ghaffari BLASIUS_MAX_N_CHEBYSHEV); 2987b8d3891SJames Wright ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", 2997b8d3891SJames Wright "Velocity at boundary layer edge", 3007b8d3891SJames Wright NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); 3017b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_refine_height", 30288626eedSJames Wright "Height of boundary layer mesh refinement", 3037b8d3891SJames Wright NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); 3047b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_growth", 30588626eedSJames Wright "Geometric growth rate of boundary layer mesh", 3067b8d3891SJames Wright NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); 3077b8d3891SJames Wright ierr = PetscOptionsScalar("-platemesh_top_angle", 30888626eedSJames Wright "Geometric top_angle rate of boundary layer mesh", 3097b8d3891SJames Wright NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); 310d2714514SJames Wright ierr = PetscOptionsString("-platemesh_y_node_locs_path", 311d2714514SJames Wright "Path to file with y node locations. " 312d2714514SJames Wright "If empty, will use the algorithmic mesh warping.", NULL, 313d2714514SJames Wright mesh_ynodes_path, mesh_ynodes_path, 314d2714514SJames Wright sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); 315ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", 316ba6664aeSJames Wright NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); 31788626eedSJames Wright PetscOptionsEnd(); 31888626eedSJames Wright 31988626eedSJames Wright PetscScalar meter = user->units->meter; 32088626eedSJames Wright PetscScalar second = user->units->second; 32188626eedSJames Wright PetscScalar Kelvin = user->units->Kelvin; 32288626eedSJames Wright PetscScalar Pascal = user->units->Pascal; 32388626eedSJames Wright 324fb455ff0SLeila Ghaffari T_inf *= Kelvin; 3252518f336SLeila Ghaffari T_wall *= Kelvin; 32688626eedSJames Wright P0 *= Pascal; 327fb455ff0SLeila Ghaffari U_inf *= meter / second; 32888626eedSJames Wright delta0 *= meter; 32988626eedSJames Wright 330d2714514SJames Wright PetscReal *mesh_ynodes = NULL; 331d2714514SJames Wright PetscInt mesh_nynodes = 0; 332d2714514SJames Wright if (strcmp(mesh_ynodes_path, "")) { 333d2714514SJames Wright ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); 33488626eedSJames Wright CHKERRQ(ierr); 335d2714514SJames Wright } 336d2714514SJames Wright ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, 337798d42b8SJames Wright mesh_refine_height, mesh_top_angle, &mesh_ynodes, 338798d42b8SJames Wright &mesh_nynodes); CHKERRQ(ierr); 33988626eedSJames Wright 340841e4c73SJed Brown // Some properties depend on parameters from NewtonianIdealGas 341841e4c73SJed Brown CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 342841e4c73SJed Brown CEED_MEM_HOST, &newtonian_ig_ctx); 34388626eedSJames Wright 344ba6664aeSJames Wright blasius_ctx->weakT = weakT; 345fb455ff0SLeila Ghaffari blasius_ctx->U_inf = U_inf; 346fb455ff0SLeila Ghaffari blasius_ctx->T_inf = T_inf; 3472518f336SLeila Ghaffari blasius_ctx->T_wall = T_wall; 348841e4c73SJed Brown blasius_ctx->delta0 = delta0; 349841e4c73SJed Brown blasius_ctx->P0 = P0; 3502518f336SLeila Ghaffari blasius_ctx->n_cheb = N; 35130e9fa81SJames Wright newtonian_ig_ctx->P0 = P0; 352841e4c73SJed Brown blasius_ctx->implicit = user->phys->implicit; 353841e4c73SJed Brown blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; 354ba6664aeSJames Wright 355f1122ed0SJames Wright { 3562518f336SLeila Ghaffari PetscReal domain_min[3], domain_max[3]; 3572518f336SLeila Ghaffari ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 358f1122ed0SJames Wright blasius_ctx->x_inflow = domain_min[0]; 3592518f336SLeila Ghaffari blasius_ctx->eta_max = 5 * domain_max[1] / blasius_ctx->delta0; 360f1122ed0SJames Wright } 361*3ce1835eSJames Wright if(!use_stg) PetscCall(ComputeChebyshevCoefficients(blasius_ctx)); 362f1122ed0SJames Wright 363841e4c73SJed Brown CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 364841e4c73SJed Brown &newtonian_ig_ctx); 36588626eedSJames Wright 366841e4c73SJed Brown CeedQFunctionContextCreate(user->ceed, &blasius_context); 3672518f336SLeila Ghaffari CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, 368841e4c73SJed Brown sizeof(*blasius_ctx), blasius_ctx); 369841e4c73SJed Brown CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, 370841e4c73SJed Brown FreeContextPetsc); 37188626eedSJames Wright 372b77c53c9SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 373841e4c73SJed Brown problem->ics.qfunction_context = blasius_context; 374ba6664aeSJames Wright if (use_stg) { 375fb455ff0SLeila Ghaffari ierr = SetupSTG(comm, dm, problem, user, weakT, T_inf, P0, mesh_ynodes, 3764e139266SJames Wright mesh_nynodes); CHKERRQ(ierr); 37765dd5cafSJames Wright } else { 37865dd5cafSJames Wright problem->apply_inflow.qfunction = Blasius_Inflow; 37965dd5cafSJames Wright problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; 380fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; 381fc02c281SJames Wright problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; 38265dd5cafSJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 38365dd5cafSJames Wright &problem->apply_inflow.qfunction_context); 384fc02c281SJames Wright CeedQFunctionContextReferenceCopy(blasius_context, 385fc02c281SJames Wright &problem->apply_inflow_jacobian.qfunction_context); 386ba6664aeSJames Wright } 3874e139266SJames Wright ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); 38888626eedSJames Wright PetscFunctionReturn(0); 38988626eedSJames Wright } 390