1ba6664aeSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2ba6664aeSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ba6664aeSJames Wright // 4ba6664aeSJames Wright // SPDX-License-Identifier: BSD-2-Clause 5ba6664aeSJames Wright // 6ba6664aeSJames Wright // This file is part of CEED: http://github.com/ceed 7ba6664aeSJames Wright 8ba6664aeSJames Wright /// @file 9ba6664aeSJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10ba6664aeSJames Wright /// presented in Shur et al. 2014 11ba6664aeSJames Wright 122b730f8bSJeremy L Thompson #include "stg_shur14.h" 132b730f8bSJeremy L Thompson 1449aac155SJeremy L Thompson #include <ceed.h> 15ba6664aeSJames Wright #include <math.h> 1649aac155SJeremy L Thompson #include <petscdm.h> 172b730f8bSJeremy L Thompson #include <stdlib.h> 182b730f8bSJeremy L Thompson 19ba6664aeSJames Wright #include "../navierstokes.h" 20ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h" 21ba6664aeSJames Wright 2265dd5cafSJames Wright STGShur14Context global_stg_ctx; 2365dd5cafSJames Wright 24ba6664aeSJames Wright /* 25ba6664aeSJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 26ba6664aeSJames Wright * 27ea61e9acSJeremy L Thompson * This assumes the input matrices are in order [11,22,33,12,13,23]. 28ea61e9acSJeremy L Thompson * This format is also used for the output. 29ba6664aeSJames Wright * 30ba6664aeSJames Wright * @param[in] comm MPI_Comm 31ba6664aeSJames Wright * @param[in] nprofs Number of matrices in Rij 32ba6664aeSJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 33ba6664aeSJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 34ba6664aeSJames Wright */ 352b730f8bSJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 36ba6664aeSJames Wright PetscFunctionBeginUser; 37ba6664aeSJames Wright for (PetscInt i = 0; i < nprofs; i++) { 38ba6664aeSJames Wright Cij[0][i] = sqrt(Rij[0][i]); 39ba6664aeSJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 40ba6664aeSJames Wright Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2)); 41ba6664aeSJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 42ba6664aeSJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i]; 43ba6664aeSJames Wright Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 44ba6664aeSJames Wright 452b730f8bSJeremy L Thompson if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) { 462b730f8bSJeremy L Thompson SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", 472b730f8bSJeremy L Thompson i + 1); 482b730f8bSJeremy L Thompson } 49ba6664aeSJames Wright } 50ba6664aeSJames Wright PetscFunctionReturn(0); 51ba6664aeSJames Wright } 52ba6664aeSJames Wright 53ba6664aeSJames Wright /* 54ba6664aeSJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 55ba6664aeSJames Wright * 56ea61e9acSJeremy L Thompson * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 57ea61e9acSJeremy L Thompson * Assumes there are 14 columns in the file. 58ba6664aeSJames Wright * 59ea61e9acSJeremy L Thompson * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 60ba6664aeSJames Wright * 61ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 62ba6664aeSJames Wright * @param[in] path Path to the STGInflow.dat file 63ea61e9acSJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 64ba6664aeSJames Wright */ 652b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 66ba6664aeSJames Wright PetscInt ndims, dims[2]; 67ba6664aeSJames Wright FILE *fp; 68ba6664aeSJames Wright const PetscInt char_array_len = 512; 69ba6664aeSJames Wright char line[char_array_len]; 70ba6664aeSJames Wright char **array; 71ba6664aeSJames Wright 72ba6664aeSJames Wright PetscFunctionBeginUser; 73ba6664aeSJames Wright 74*8b892a05SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 75ba6664aeSJames Wright 76ba6664aeSJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 77175f00a6SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 78ba6664aeSJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 79ba6664aeSJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 802b730f8bSJeremy L Thompson CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 81ba6664aeSJames Wright 82ba6664aeSJames Wright for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 8324941a69SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 8424941a69SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 852b730f8bSJeremy L Thompson if (ndims < dims[1]) { 862b730f8bSJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 872b730f8bSJeremy L Thompson ndims, dims[1]); 882b730f8bSJeremy L Thompson } 89ba6664aeSJames Wright 90175f00a6SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 91ba6664aeSJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 92ba6664aeSJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 93ba6664aeSJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 94ba6664aeSJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 95ba6664aeSJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 96ba6664aeSJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 97ba6664aeSJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 98ba6664aeSJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 99ba6664aeSJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 100ba6664aeSJames Wright lt[i] = (CeedScalar)atof(array[12]); 101ba6664aeSJames Wright eps[i] = (CeedScalar)atof(array[13]); 102ba6664aeSJames Wright 1032b730f8bSJeremy L Thompson if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path); 1042b730f8bSJeremy L Thompson if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path); 1052b730f8bSJeremy L Thompson if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path); 106ba6664aeSJames Wright } 1072b730f8bSJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 10824941a69SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 10924941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 1102b730f8bSJeremy L Thompson 111ba6664aeSJames Wright PetscFunctionReturn(0); 112ba6664aeSJames Wright } 113ba6664aeSJames Wright 114ba6664aeSJames Wright /* 115ba6664aeSJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 116ba6664aeSJames Wright * 117ea61e9acSJeremy L Thompson * Assumes that the first line of the file has the number of rows and columns as the only two entries, separated by a single space. 118ea61e9acSJeremy L Thompson * Assumes there are 7 columns in the file. 119ba6664aeSJames Wright * 120ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 121ba6664aeSJames Wright * @param[in] path Path to the STGRand.dat file 122ea61e9acSJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 123ba6664aeSJames Wright */ 1242b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 125ba6664aeSJames Wright PetscInt ndims, dims[2]; 126ba6664aeSJames Wright FILE *fp; 127ba6664aeSJames Wright const PetscInt char_array_len = 512; 128ba6664aeSJames Wright char line[char_array_len]; 129ba6664aeSJames Wright char **array; 130ba6664aeSJames Wright 131ba6664aeSJames Wright PetscFunctionBeginUser; 132*8b892a05SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 133ba6664aeSJames Wright 134ba6664aeSJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1352b730f8bSJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1362b730f8bSJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 137ba6664aeSJames Wright 138ba6664aeSJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 13924941a69SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 14024941a69SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1412b730f8bSJeremy L Thompson if (ndims < dims[1]) { 1422b730f8bSJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 1432b730f8bSJeremy L Thompson ndims, dims[1]); 1442b730f8bSJeremy L Thompson } 145ba6664aeSJames Wright 146ba6664aeSJames Wright d[0][i] = (CeedScalar)atof(array[0]); 147ba6664aeSJames Wright d[1][i] = (CeedScalar)atof(array[1]); 148ba6664aeSJames Wright d[2][i] = (CeedScalar)atof(array[2]); 149ba6664aeSJames Wright phi[i] = (CeedScalar)atof(array[3]); 150ba6664aeSJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 151ba6664aeSJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 152ba6664aeSJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 153ba6664aeSJames Wright } 15424941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 1552b730f8bSJeremy L Thompson 156ba6664aeSJames Wright PetscFunctionReturn(0); 157ba6664aeSJames Wright } 158ba6664aeSJames Wright 159ba6664aeSJames Wright /* 160ba6664aeSJames Wright * @brief Read STG data from input paths and put in STGShur14Context 161ba6664aeSJames Wright * 162ba6664aeSJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 163ea61e9acSJeremy L Thompson * Data stored initially in `*pstg_ctx` will be copied over to the new STGShur14Context instance. 164ba6664aeSJames Wright * 165ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 166ba6664aeSJames Wright * @param[in] dm DM for the program 167ba6664aeSJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 168ba6664aeSJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 169ea61e9acSJeremy L Thompson * @param[in,out] pstg_ctx Pointer to STGShur14Context where the data will be loaded into 170ba6664aeSJames Wright */ 1712b730f8bSJeremy L Thompson PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm, char stg_inflow_path[PETSC_MAX_PATH_LEN], char stg_rand_path[PETSC_MAX_PATH_LEN], 1722b730f8bSJeremy L Thompson STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) { 173ba6664aeSJames Wright PetscInt nmodes, nprofs; 174ba6664aeSJames Wright STGShur14Context stg_ctx; 175ba6664aeSJames Wright PetscFunctionBeginUser; 176ba6664aeSJames Wright 177ba6664aeSJames Wright // Get options 178*8b892a05SJames Wright PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes)); 179*8b892a05SJames Wright PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs)); 18049aac155SJeremy L Thompson if (nmodes > STG_NMODES_MAX) { 1812b730f8bSJeremy L Thompson SETERRQ(comm, 1, 1822b730f8bSJeremy L Thompson "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT 1832b730f8bSJeremy L Thompson "). " 1842b730f8bSJeremy L Thompson "Change size of STG_NMODES_MAX and recompile", 1852b730f8bSJeremy L Thompson stg_rand_path, nmodes, STG_NMODES_MAX); 18649aac155SJeremy L Thompson } 187ba6664aeSJames Wright 188ba6664aeSJames Wright { 189ba6664aeSJames Wright STGShur14Context s; 19024941a69SJames Wright PetscCall(PetscCalloc1(1, &s)); 191ba6664aeSJames Wright *s = **pstg_ctx; 192ba6664aeSJames Wright s->nmodes = nmodes; 193ba6664aeSJames Wright s->nprofs = nprofs; 194ba6664aeSJames Wright s->offsets.sigma = 0; 195ba6664aeSJames Wright s->offsets.d = nmodes * 3; 196ba6664aeSJames Wright s->offsets.phi = s->offsets.d + nmodes * 3; 197ba6664aeSJames Wright s->offsets.kappa = s->offsets.phi + nmodes; 198175f00a6SJames Wright s->offsets.wall_dist = s->offsets.kappa + nmodes; 199175f00a6SJames Wright s->offsets.ubar = s->offsets.wall_dist + nprofs; 200ba6664aeSJames Wright s->offsets.cij = s->offsets.ubar + nprofs * 3; 201ba6664aeSJames Wright s->offsets.eps = s->offsets.cij + nprofs * 6; 202ba6664aeSJames Wright s->offsets.lt = s->offsets.eps + nprofs; 2034e139266SJames Wright s->offsets.ynodes = s->offsets.lt + nprofs; 2044e139266SJames Wright PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes; 205ba6664aeSJames Wright s->total_bytes = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]); 20624941a69SJames Wright PetscCall(PetscMalloc(s->total_bytes, &stg_ctx)); 207ba6664aeSJames Wright *stg_ctx = *s; 20824941a69SJames Wright PetscCall(PetscFree(s)); 209ba6664aeSJames Wright } 210ba6664aeSJames Wright 21124941a69SJames Wright PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx)); 21224941a69SJames Wright PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx)); 213ba6664aeSJames Wright 2144e139266SJames Wright if (stg_ctx->nynodes > 0) { 2154e139266SJames Wright CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes]; 2164e139266SJames Wright for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i]; 2174e139266SJames Wright } 2184e139266SJames Wright 219ba6664aeSJames Wright // -- Calculate kappa 220ba6664aeSJames Wright { 221ba6664aeSJames Wright CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 222175f00a6SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 223ba6664aeSJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 224ba6664aeSJames Wright CeedScalar le, le_max = 0; 225ba6664aeSJames Wright 2262b730f8bSJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 227175f00a6SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 228ba6664aeSJames Wright if (le_max < le) le_max = le; 229ba6664aeSJames Wright } 230ba6664aeSJames Wright CeedScalar kmin = M_PI / le_max; 231ba6664aeSJames Wright 2322b730f8bSJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); } 233ba6664aeSJames Wright } // end calculate kappa 234ba6664aeSJames Wright 23524941a69SJames Wright PetscCall(PetscFree(*pstg_ctx)); 236ba6664aeSJames Wright *pstg_ctx = stg_ctx; 237ba6664aeSJames Wright PetscFunctionReturn(0); 238ba6664aeSJames Wright } 239ba6664aeSJames Wright 2402b730f8bSJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0, 2412b730f8bSJeremy L Thompson const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) { 242ba6664aeSJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 243ba6664aeSJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 2442b730f8bSJeremy L Thompson PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE; 2452b730f8bSJeremy L Thompson CeedScalar u0 = 0.0, alpha = 1.01; 246ba6664aeSJames Wright CeedQFunctionContext stg_context; 247ba6664aeSJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 248ba6664aeSJames Wright PetscFunctionBeginUser; 249ba6664aeSJames Wright 250ba6664aeSJames Wright // Get options 251ba6664aeSJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2522b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2532b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2542b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2552b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2562b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2572b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2582b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2592b730f8bSJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 260ba6664aeSJames Wright PetscOptionsEnd(); 261ba6664aeSJames Wright 26224941a69SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 26365dd5cafSJames Wright global_stg_ctx->alpha = alpha; 26465dd5cafSJames Wright global_stg_ctx->u0 = u0; 26565dd5cafSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 26665dd5cafSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 26765dd5cafSJames Wright global_stg_ctx->mean_only = mean_only; 26889060322SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 26965dd5cafSJames Wright global_stg_ctx->theta0 = theta0; 27065dd5cafSJames Wright global_stg_ctx->P0 = P0; 27165dd5cafSJames Wright global_stg_ctx->nynodes = nynodes; 272ba6664aeSJames Wright 273ba6664aeSJames Wright { 274ba6664aeSJames Wright // Calculate dx assuming constant spacing 275ba6664aeSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 27624941a69SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 277ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 278ba6664aeSJames Wright 279ba6664aeSJames Wright PetscInt nmax = 3, faces[3]; 2802b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 28165dd5cafSJames Wright global_stg_ctx->dx = domain_size[0] / faces[0]; 28265dd5cafSJames Wright global_stg_ctx->dz = domain_size[2] / faces[2]; 283ba6664aeSJames Wright } 284ba6664aeSJames Wright 2852b730f8bSJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 28665dd5cafSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 2872b730f8bSJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 288ba6664aeSJames Wright 2892b730f8bSJeremy L Thompson PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes)); 290ba6664aeSJames Wright 291ba6664aeSJames Wright CeedQFunctionContextCreate(user->ceed, &stg_context); 2922b730f8bSJeremy L Thompson CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx); 2932b730f8bSJeremy L Thompson CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc); 2942b730f8bSJeremy L Thompson CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution"); 295ba6664aeSJames Wright 296b77c53c9SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 297b77c53c9SJames Wright problem->ics.qfunction = ICsSTG; 298b77c53c9SJames Wright problem->ics.qfunction_loc = ICsSTG_loc; 299b77c53c9SJames Wright problem->ics.qfunction_context = stg_context; 300b77c53c9SJames Wright 3014e139266SJames Wright if (use_stgstrong) { 30265dd5cafSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 303dada6cc0SJames Wright problem->use_dirichlet_ceed = PETSC_TRUE; 30436b31e27SJames Wright problem->bc_from_ics = PETSC_FALSE; 3054e139266SJames Wright } else { 306ba6664aeSJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow; 307ba6664aeSJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc; 3084dbab5e5SJames Wright problem->apply_inflow_jacobian.qfunction = STGShur14_Inflow_Jacobian; 3094dbab5e5SJames Wright problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc; 3102b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context); 3112b730f8bSJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context); 31236b31e27SJames Wright problem->bc_from_ics = PETSC_TRUE; 3134e139266SJames Wright } 314ba6664aeSJames Wright 315ba6664aeSJames Wright PetscFunctionReturn(0); 316ba6664aeSJames Wright } 3174e139266SJames Wright 3182b730f8bSJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) { 3194e139266SJames Wright const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]); 3204e139266SJames Wright // ^^assuming min(dy) is first element off the wall 3214e139266SJames Wright PetscInt idx = -1; // Index 3224e139266SJames Wright 3234e139266SJames Wright for (PetscInt i = 0; i < nynodes; i++) { 3244e139266SJames Wright if (y < ynodes[i] + half_mindy) { 3252b730f8bSJeremy L Thompson idx = i; 3262b730f8bSJeremy L Thompson break; 3274e139266SJames Wright } 3284e139266SJames Wright } 3294e139266SJames Wright if (idx == 0) return ynodes[1] - ynodes[0]; 3304e139266SJames Wright else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1]; 3314e139266SJames Wright else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]); 3324e139266SJames Wright } 3334e139266SJames Wright 3344e139266SJames Wright // Function passed to DMAddBoundary 335dada6cc0SJames Wright // NOTE: Not used in favor of QFunction-based method 3362b730f8bSJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) { 3374e139266SJames Wright PetscFunctionBeginUser; 3382b730f8bSJeremy L Thompson 3394e139266SJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 3404e139266SJames Wright PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt; 3414e139266SJames Wright const bool mean_only = stg_ctx->mean_only; 3424e139266SJames Wright const PetscScalar dx = stg_ctx->dx; 3434e139266SJames Wright const PetscScalar dz = stg_ctx->dz; 3444e139266SJames Wright const PetscScalar mu = stg_ctx->newtonian_ctx.mu; 3454e139266SJames Wright const PetscScalar theta0 = stg_ctx->theta0; 3464e139266SJames Wright const PetscScalar P0 = stg_ctx->P0; 3474e139266SJames Wright const PetscScalar cv = stg_ctx->newtonian_ctx.cv; 3484e139266SJames Wright const PetscScalar cp = stg_ctx->newtonian_ctx.cp; 3494e139266SJames Wright const PetscScalar Rd = cp - cv; 3504e139266SJames Wright 3514e139266SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 3524e139266SJames Wright InterpolateProfile(x[1], ubar, cij, &eps, <, stg_ctx); 3534e139266SJames Wright if (!mean_only) { 3544e139266SJames Wright const PetscInt nynodes = stg_ctx->nynodes; 3554e139266SJames Wright const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes]; 3564e139266SJames Wright const PetscScalar h[3] = {dx, FindDy(ynodes, nynodes, x[1]), dz}; 3574e139266SJames Wright CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx); 3584e139266SJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 3594e139266SJames Wright } else { 3604e139266SJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 3614e139266SJames Wright } 3624e139266SJames Wright 3634e139266SJames Wright bcval[0] = rho; 3644e139266SJames Wright bcval[1] = rho * u[0]; 3654e139266SJames Wright bcval[2] = rho * u[1]; 3664e139266SJames Wright bcval[3] = rho * u[2]; 3674e139266SJames Wright PetscFunctionReturn(0); 3684e139266SJames Wright } 3694e139266SJames Wright 3702b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) { 3714e139266SJames Wright DMLabel label; 3724e139266SJames Wright PetscFunctionBeginUser; 3734e139266SJames Wright 374192a7459SJames Wright PetscInt comps[5], num_comps = 4; 37597baf651SJames Wright switch (phys->state_var) { 37697baf651SJames Wright case STATEVAR_CONSERVATIVE: 377192a7459SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 378192a7459SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 37997baf651SJames Wright break; 38097baf651SJames Wright 38197baf651SJames Wright case STATEVAR_PRIMITIVE: 38297baf651SJames Wright // {1,2,3,4} for u, v, w, T 38397baf651SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 38497baf651SJames Wright break; 385192a7459SJames Wright } 386192a7459SJames Wright 38724941a69SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3884e139266SJames Wright // Set wall BCs 3894e139266SJames Wright if (bc->num_inflow > 0) { 3902b730f8bSJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc, 39124941a69SJames Wright NULL, global_stg_ctx, NULL)); 3924e139266SJames Wright } 3934e139266SJames Wright 3944e139266SJames Wright PetscFunctionReturn(0); 3954e139266SJames Wright } 396dada6cc0SJames Wright 3972b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, 3985dc40723SJames Wright CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) { 399dada6cc0SJames Wright CeedQFunction qf_strongbc; 400dada6cc0SJames Wright PetscFunctionBeginUser; 4012b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc); 4022b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 403dada6cc0SJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 404dada6cc0SJames Wright CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE); 4055dc40723SJames Wright CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 406dada6cc0SJames Wright CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE); 407dada6cc0SJames Wright 408dada6cc0SJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 409dada6cc0SJames Wright *pqf_strongbc = qf_strongbc; 410dada6cc0SJames Wright PetscFunctionReturn(0); 411dada6cc0SJames Wright } 4125dc40723SJames Wright 4132b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur, 4145dc40723SJames Wright CeedQFunction *pqf_strongbc) { 4155dc40723SJames Wright CeedQFunction qf_strongbc; 4165dc40723SJames Wright PetscFunctionBeginUser; 4172b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc); 4182b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 4195dc40723SJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 4205dc40723SJames Wright CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 4215dc40723SJames Wright 4225dc40723SJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 4235dc40723SJames Wright *pqf_strongbc = qf_strongbc; 4245dc40723SJames Wright PetscFunctionReturn(0); 4255dc40723SJames Wright } 426