// Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. // // SPDX-License-Identifier: BSD-2-Clause // // This file is part of CEED: http://github.com/ceed /// @file /// Utility functions for setting up Blasius Boundary Layer #include "../navierstokes.h" #include "../qfunctions/blasius.h" #include "stg_shur14.h" static PetscErrorCode GetYNodeLocs(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal **pynodes, PetscInt *nynodes) { PetscErrorCode ierr; PetscInt ndims, dims[2]; FILE *fp; const PetscInt char_array_len = 512; char line[char_array_len]; char **array; PetscReal *node_locs; PetscFunctionBeginUser; ierr = PetscFOpen(comm, path, "r", &fp); CHKERRQ(ierr); ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); for (PetscInt i=0; i faces[1] +1) { ierr = PetscPrintf(comm, "WARNING: y_node_locs_path has more locations (%d) " "than the mesh has nodes (%d). This maybe unintended.", *num_node_locs, faces[1]+1); CHKERRQ(ierr); } PetscScalar max_y = (*node_locs)[faces[1]]; for (PetscInt i=0; iics.qfunction = ICsBlasius; problem->ics.qfunction_loc = ICsBlasius_loc; CeedScalar Uinf = 40; // m/s CeedScalar delta0 = 4.2e-4; // m CeedScalar theta0 = 288.; // K CeedScalar P0 = 1.01e5; // Pa PetscBool weakT = PETSC_FALSE; // weak density or temperature PetscReal mesh_refine_height = 5.9e-4; // m PetscReal mesh_growth = 1.08; // [-] PetscInt mesh_Ndelta = 45; // [-] PetscReal mesh_top_angle = 5; // degrees char mesh_ynodes_path[PETSC_MAX_PATH_LEN] = ""; PetscOptionsBegin(comm, NULL, "Options for CHANNEL problem", NULL); ierr = PetscOptionsBool("-weakT", "Change from rho weak to T weak at inflow", NULL, weakT, &weakT, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-Uinf", "Velocity at boundary layer edge", NULL, Uinf, &Uinf, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-delta0", "Boundary layer height at inflow", NULL, delta0, &delta0, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-theta0", "Wall temperature", NULL, theta0, &theta0, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-P0", "Pressure at outflow", NULL, P0, &P0, NULL); CHKERRQ(ierr); ierr = PetscOptionsBoundedInt("-platemesh_Ndelta", "Velocity at boundary layer edge", NULL, mesh_Ndelta, &mesh_Ndelta, NULL, 1); CHKERRQ(ierr); ierr = PetscOptionsScalar("-platemesh_refine_height", "Height of boundary layer mesh refinement", NULL, mesh_refine_height, &mesh_refine_height, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-platemesh_growth", "Geometric growth rate of boundary layer mesh", NULL, mesh_growth, &mesh_growth, NULL); CHKERRQ(ierr); ierr = PetscOptionsScalar("-platemesh_top_angle", "Geometric top_angle rate of boundary layer mesh", NULL, mesh_top_angle, &mesh_top_angle, NULL); CHKERRQ(ierr); ierr = PetscOptionsString("-platemesh_y_node_locs_path", "Path to file with y node locations. " "If empty, will use the algorithmic mesh warping.", NULL, mesh_ynodes_path, mesh_ynodes_path, sizeof(mesh_ynodes_path), NULL); CHKERRQ(ierr); ierr = PetscOptionsBool("-stg_use", "Use STG inflow boundary condition", NULL, use_stg, &use_stg, NULL); CHKERRQ(ierr); PetscOptionsEnd(); PetscScalar meter = user->units->meter; PetscScalar second = user->units->second; PetscScalar Kelvin = user->units->Kelvin; PetscScalar Pascal = user->units->Pascal; theta0 *= Kelvin; P0 *= Pascal; Uinf *= meter / second; delta0 *= meter; PetscReal *mesh_ynodes = NULL; PetscInt mesh_nynodes = 0; if (strcmp(mesh_ynodes_path, "")) { ierr = GetYNodeLocs(comm, mesh_ynodes_path, &mesh_ynodes, &mesh_nynodes); CHKERRQ(ierr); } ierr = ModifyMesh(comm, dm, problem->dim, mesh_growth, mesh_Ndelta, mesh_refine_height, mesh_top_angle, &mesh_ynodes, &mesh_nynodes); CHKERRQ(ierr); // Some properties depend on parameters from NewtonianIdealGas CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); blasius_ctx->weakT = weakT; blasius_ctx->Uinf = Uinf; blasius_ctx->delta0 = delta0; blasius_ctx->theta0 = theta0; blasius_ctx->P0 = P0; newtonian_ig_ctx->P0 = P0; blasius_ctx->implicit = user->phys->implicit; blasius_ctx->newtonian_ctx = *newtonian_ig_ctx; { PetscReal domain_min[3]; ierr = DMGetBoundingBox(dm, domain_min, NULL); CHKERRQ(ierr); blasius_ctx->x_inflow = domain_min[0]; } CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); CeedQFunctionContextCreate(user->ceed, &blasius_context); CeedQFunctionContextSetData(blasius_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*blasius_ctx), blasius_ctx); CeedQFunctionContextSetDataDestroy(blasius_context, CEED_MEM_HOST, FreeContextPetsc); CeedQFunctionContextDestroy(&problem->ics.qfunction_context); problem->ics.qfunction_context = blasius_context; if (use_stg) { ierr = SetupSTG(comm, dm, problem, user, weakT, theta0, P0, mesh_ynodes, mesh_nynodes); CHKERRQ(ierr); } else { CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context); problem->apply_inflow.qfunction = Blasius_Inflow; problem->apply_inflow.qfunction_loc = Blasius_Inflow_loc; problem->apply_inflow_jacobian.qfunction = Blasius_Inflow_Jacobian; problem->apply_inflow_jacobian.qfunction_loc = Blasius_Inflow_Jacobian_loc; CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow.qfunction_context); CeedQFunctionContextReferenceCopy(blasius_context, &problem->apply_inflow_jacobian.qfunction_context); } ierr = PetscFree(mesh_ynodes); CHKERRQ(ierr); PetscFunctionReturn(0); }