1493642f1SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3493642f1SJames Wright // 4493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause 5493642f1SJames Wright // 6493642f1SJames Wright // This file is part of CEED: http://github.com/ceed 7493642f1SJames Wright 8493642f1SJames Wright /// @file 9493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10493642f1SJames Wright /// presented in Shur et al. 2014 11493642f1SJames Wright 122b916ea7SJeremy L Thompson #include "stg_shur14.h" 132b916ea7SJeremy L Thompson 14*e419654dSJeremy L Thompson #include <ceed.h> 15493642f1SJames Wright #include <math.h> 16*e419654dSJeremy L Thompson #include <petscdm.h> 172b916ea7SJeremy L Thompson #include <stdlib.h> 182b916ea7SJeremy L Thompson 19493642f1SJames Wright #include "../navierstokes.h" 20493642f1SJames Wright #include "../qfunctions/stg_shur14.h" 21493642f1SJames Wright 228085925cSJames Wright STGShur14Context global_stg_ctx; 238085925cSJames Wright 24493642f1SJames Wright /* 25493642f1SJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 26493642f1SJames Wright * 2704e40bb6SJeremy L Thompson * This assumes the input matrices are in order [11,22,33,12,13,23]. 2804e40bb6SJeremy L Thompson * This format is also used for the output. 29493642f1SJames Wright * 30493642f1SJames Wright * @param[in] comm MPI_Comm 31493642f1SJames Wright * @param[in] nprofs Number of matrices in Rij 32493642f1SJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 33493642f1SJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 34493642f1SJames Wright */ 352b916ea7SJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 36493642f1SJames Wright PetscFunctionBeginUser; 37493642f1SJames Wright for (PetscInt i = 0; i < nprofs; i++) { 38493642f1SJames Wright Cij[0][i] = sqrt(Rij[0][i]); 39493642f1SJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 40493642f1SJames Wright Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2)); 41493642f1SJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 42493642f1SJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i]; 43493642f1SJames Wright Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 44493642f1SJames Wright 452b916ea7SJeremy L Thompson if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) { 462b916ea7SJeremy L Thompson SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", 472b916ea7SJeremy L Thompson i + 1); 482b916ea7SJeremy L Thompson } 49493642f1SJames Wright } 50493642f1SJames Wright PetscFunctionReturn(0); 51493642f1SJames Wright } 52493642f1SJames Wright 53493642f1SJames Wright /* 54493642f1SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 55493642f1SJames Wright * 5604e40bb6SJeremy L Thompson * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`. 5704e40bb6SJeremy L Thompson * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example). 58493642f1SJames Wright * 5904e40bb6SJeremy 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. 60493642f1SJames Wright * 61493642f1SJames Wright * @param[in] comm MPI_Comm for the program 62493642f1SJames Wright * @param[in] path Path to the file 63493642f1SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 64493642f1SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 65493642f1SJames Wright * @param[out] fp File pointer to the opened file 66493642f1SJames Wright */ 672b916ea7SJeremy L Thompson static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2], 682b916ea7SJeremy L Thompson FILE **fp) { 69493642f1SJames Wright PetscInt ndims; 70493642f1SJames Wright char line[char_array_len]; 71493642f1SJames Wright char **array; 72493642f1SJames Wright 73493642f1SJames Wright PetscFunctionBeginUser; 7434b5deb1SJames Wright PetscCall(PetscFOpen(comm, path, "r", fp)); 7534b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line)); 7634b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 772b916ea7SJeremy L Thompson if (ndims != 2) { 782b916ea7SJeremy L Thompson SETERRQ(comm, -1, "Found %" PetscInt_FMT " dimensions instead of 2 on the first line of %s", ndims, path); 792b916ea7SJeremy L Thompson } 80493642f1SJames Wright 81493642f1SJames Wright for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]); 8234b5deb1SJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 832b916ea7SJeremy L Thompson 84493642f1SJames Wright PetscFunctionReturn(0); 85493642f1SJames Wright } 86493642f1SJames Wright 87493642f1SJames Wright /* 8804e40bb6SJeremy L Thompson * @brief Get the number of rows for the PHASTA file at path. 89493642f1SJames Wright * 9004e40bb6SJeremy 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. 91493642f1SJames Wright * 92493642f1SJames Wright * @param[in] comm MPI_Comm for the program 93493642f1SJames Wright * @param[in] path Path to the file 94493642f1SJames Wright * @param[out] nrows Number of rows 95493642f1SJames Wright */ 962b916ea7SJeremy L Thompson static PetscErrorCode GetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 97493642f1SJames Wright const PetscInt char_array_len = 512; 98493642f1SJames Wright PetscInt dims[2]; 99493642f1SJames Wright FILE *fp; 100493642f1SJames Wright 101493642f1SJames Wright PetscFunctionBeginUser; 10234b5deb1SJames Wright PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp)); 103493642f1SJames Wright *nrows = dims[0]; 10434b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 1052b916ea7SJeremy L Thompson 106493642f1SJames Wright PetscFunctionReturn(0); 107493642f1SJames Wright } 108493642f1SJames Wright 109493642f1SJames Wright /* 110493642f1SJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 111493642f1SJames Wright * 11204e40bb6SJeremy 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. 11304e40bb6SJeremy L Thompson * Assumes there are 14 columns in the file. 114493642f1SJames Wright * 11504e40bb6SJeremy L Thompson * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 116493642f1SJames Wright * 117493642f1SJames Wright * @param[in] comm MPI_Comm for the program 118493642f1SJames Wright * @param[in] path Path to the STGInflow.dat file 11904e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 120493642f1SJames Wright */ 1212b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 122493642f1SJames Wright PetscInt ndims, dims[2]; 123493642f1SJames Wright FILE *fp; 124493642f1SJames Wright const PetscInt char_array_len = 512; 125493642f1SJames Wright char line[char_array_len]; 126493642f1SJames Wright char **array; 127493642f1SJames Wright 128493642f1SJames Wright PetscFunctionBeginUser; 129493642f1SJames Wright 13034b5deb1SJames Wright PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp)); 131493642f1SJames Wright 132493642f1SJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 133c77f3192SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 134493642f1SJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 135493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 1362b916ea7SJeremy L Thompson CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 137493642f1SJames Wright 138493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 13934b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 14034b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1412b916ea7SJeremy L Thompson if (ndims < dims[1]) { 1422b916ea7SJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 1432b916ea7SJeremy L Thompson ndims, dims[1]); 1442b916ea7SJeremy L Thompson } 145493642f1SJames Wright 146c77f3192SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 147493642f1SJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 148493642f1SJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 149493642f1SJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 150493642f1SJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 151493642f1SJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 152493642f1SJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 153493642f1SJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 154493642f1SJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 155493642f1SJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 156493642f1SJames Wright lt[i] = (CeedScalar)atof(array[12]); 157493642f1SJames Wright eps[i] = (CeedScalar)atof(array[13]); 158493642f1SJames Wright 1592b916ea7SJeremy L Thompson if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path); 1602b916ea7SJeremy L Thompson if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path); 1612b916ea7SJeremy L Thompson if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path); 162493642f1SJames Wright } 1632b916ea7SJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 16434b5deb1SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 16534b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 1662b916ea7SJeremy L Thompson 167493642f1SJames Wright PetscFunctionReturn(0); 168493642f1SJames Wright } 169493642f1SJames Wright 170493642f1SJames Wright /* 171493642f1SJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 172493642f1SJames Wright * 17304e40bb6SJeremy 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. 17404e40bb6SJeremy L Thompson * Assumes there are 7 columns in the file. 175493642f1SJames Wright * 176493642f1SJames Wright * @param[in] comm MPI_Comm for the program 177493642f1SJames Wright * @param[in] path Path to the STGRand.dat file 17804e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 179493642f1SJames Wright */ 1802b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 181493642f1SJames Wright PetscInt ndims, dims[2]; 182493642f1SJames Wright FILE *fp; 183493642f1SJames Wright const PetscInt char_array_len = 512; 184493642f1SJames Wright char line[char_array_len]; 185493642f1SJames Wright char **array; 186493642f1SJames Wright 187493642f1SJames Wright PetscFunctionBeginUser; 18834b5deb1SJames Wright PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp)); 189493642f1SJames Wright 190493642f1SJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1912b916ea7SJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1922b916ea7SJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 193493642f1SJames Wright 194493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 19534b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 19634b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1972b916ea7SJeremy L Thompson if (ndims < dims[1]) { 1982b916ea7SJeremy L Thompson SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, 1992b916ea7SJeremy L Thompson ndims, dims[1]); 2002b916ea7SJeremy L Thompson } 201493642f1SJames Wright 202493642f1SJames Wright d[0][i] = (CeedScalar)atof(array[0]); 203493642f1SJames Wright d[1][i] = (CeedScalar)atof(array[1]); 204493642f1SJames Wright d[2][i] = (CeedScalar)atof(array[2]); 205493642f1SJames Wright phi[i] = (CeedScalar)atof(array[3]); 206493642f1SJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 207493642f1SJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 208493642f1SJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 209493642f1SJames Wright } 21034b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 2112b916ea7SJeremy L Thompson 212493642f1SJames Wright PetscFunctionReturn(0); 213493642f1SJames Wright } 214493642f1SJames Wright 215493642f1SJames Wright /* 216493642f1SJames Wright * @brief Read STG data from input paths and put in STGShur14Context 217493642f1SJames Wright * 218493642f1SJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 21904e40bb6SJeremy L Thompson * Data stored initially in `*pstg_ctx` will be copied over to the new STGShur14Context instance. 220493642f1SJames Wright * 221493642f1SJames Wright * @param[in] comm MPI_Comm for the program 222493642f1SJames Wright * @param[in] dm DM for the program 223493642f1SJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 224493642f1SJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 22504e40bb6SJeremy L Thompson * @param[in,out] pstg_ctx Pointer to STGShur14Context where the data will be loaded into 226493642f1SJames Wright */ 2272b916ea7SJeremy 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], 2282b916ea7SJeremy L Thompson STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) { 229493642f1SJames Wright PetscInt nmodes, nprofs; 230493642f1SJames Wright STGShur14Context stg_ctx; 231493642f1SJames Wright PetscFunctionBeginUser; 232493642f1SJames Wright 233493642f1SJames Wright // Get options 23434b5deb1SJames Wright PetscCall(GetNRows(comm, stg_rand_path, &nmodes)); 23534b5deb1SJames Wright PetscCall(GetNRows(comm, stg_inflow_path, &nprofs)); 236*e419654dSJeremy L Thompson if (nmodes > STG_NMODES_MAX) { 2372b916ea7SJeremy L Thompson SETERRQ(comm, 1, 2382b916ea7SJeremy L Thompson "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT 2392b916ea7SJeremy L Thompson "). " 2402b916ea7SJeremy L Thompson "Change size of STG_NMODES_MAX and recompile", 2412b916ea7SJeremy L Thompson stg_rand_path, nmodes, STG_NMODES_MAX); 242*e419654dSJeremy L Thompson } 243493642f1SJames Wright 244493642f1SJames Wright { 245493642f1SJames Wright STGShur14Context s; 24634b5deb1SJames Wright PetscCall(PetscCalloc1(1, &s)); 247493642f1SJames Wright *s = **pstg_ctx; 248493642f1SJames Wright s->nmodes = nmodes; 249493642f1SJames Wright s->nprofs = nprofs; 250493642f1SJames Wright s->offsets.sigma = 0; 251493642f1SJames Wright s->offsets.d = nmodes * 3; 252493642f1SJames Wright s->offsets.phi = s->offsets.d + nmodes * 3; 253493642f1SJames Wright s->offsets.kappa = s->offsets.phi + nmodes; 254c77f3192SJames Wright s->offsets.wall_dist = s->offsets.kappa + nmodes; 255c77f3192SJames Wright s->offsets.ubar = s->offsets.wall_dist + nprofs; 256493642f1SJames Wright s->offsets.cij = s->offsets.ubar + nprofs * 3; 257493642f1SJames Wright s->offsets.eps = s->offsets.cij + nprofs * 6; 258493642f1SJames Wright s->offsets.lt = s->offsets.eps + nprofs; 259e6098bcdSJames Wright s->offsets.ynodes = s->offsets.lt + nprofs; 260e6098bcdSJames Wright PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes; 261493642f1SJames Wright s->total_bytes = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]); 26234b5deb1SJames Wright PetscCall(PetscMalloc(s->total_bytes, &stg_ctx)); 263493642f1SJames Wright *stg_ctx = *s; 26434b5deb1SJames Wright PetscCall(PetscFree(s)); 265493642f1SJames Wright } 266493642f1SJames Wright 26734b5deb1SJames Wright PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx)); 26834b5deb1SJames Wright PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx)); 269493642f1SJames Wright 270e6098bcdSJames Wright if (stg_ctx->nynodes > 0) { 271e6098bcdSJames Wright CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes]; 272e6098bcdSJames Wright for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i]; 273e6098bcdSJames Wright } 274e6098bcdSJames Wright 275493642f1SJames Wright // -- Calculate kappa 276493642f1SJames Wright { 277493642f1SJames Wright CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 278c77f3192SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 279493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 280493642f1SJames Wright CeedScalar le, le_max = 0; 281493642f1SJames Wright 2822b916ea7SJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 283c77f3192SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 284493642f1SJames Wright if (le_max < le) le_max = le; 285493642f1SJames Wright } 286493642f1SJames Wright CeedScalar kmin = M_PI / le_max; 287493642f1SJames Wright 2882b916ea7SJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); } 289493642f1SJames Wright } // end calculate kappa 290493642f1SJames Wright 29134b5deb1SJames Wright PetscCall(PetscFree(*pstg_ctx)); 292493642f1SJames Wright *pstg_ctx = stg_ctx; 293493642f1SJames Wright PetscFunctionReturn(0); 294493642f1SJames Wright } 295493642f1SJames Wright 2962b916ea7SJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0, 2972b916ea7SJeremy L Thompson const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) { 298493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 299493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 3002b916ea7SJeremy L Thompson PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE; 3012b916ea7SJeremy L Thompson CeedScalar u0 = 0.0, alpha = 1.01; 302493642f1SJames Wright CeedQFunctionContext stg_context; 303493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 304493642f1SJames Wright PetscFunctionBeginUser; 305493642f1SJames Wright 306493642f1SJames Wright // Get options 307493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 3082b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 3092b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 3102b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 3112b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 3122b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 3132b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 3142b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 3152b916ea7SJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 316493642f1SJames Wright PetscOptionsEnd(); 317493642f1SJames Wright 31834b5deb1SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 3198085925cSJames Wright global_stg_ctx->alpha = alpha; 3208085925cSJames Wright global_stg_ctx->u0 = u0; 3218085925cSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 3228085925cSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 3238085925cSJames Wright global_stg_ctx->mean_only = mean_only; 324d4e0f297SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 3258085925cSJames Wright global_stg_ctx->theta0 = theta0; 3268085925cSJames Wright global_stg_ctx->P0 = P0; 3278085925cSJames Wright global_stg_ctx->nynodes = nynodes; 328493642f1SJames Wright 329493642f1SJames Wright { 330493642f1SJames Wright // Calculate dx assuming constant spacing 331493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 33234b5deb1SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 333493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 334493642f1SJames Wright 335493642f1SJames Wright PetscInt nmax = 3, faces[3]; 3362b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 3378085925cSJames Wright global_stg_ctx->dx = domain_size[0] / faces[0]; 3388085925cSJames Wright global_stg_ctx->dz = domain_size[2] / faces[2]; 339493642f1SJames Wright } 340493642f1SJames Wright 3412b916ea7SJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 3428085925cSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 3432b916ea7SJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 344493642f1SJames Wright 3452b916ea7SJeremy L Thompson PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes)); 346493642f1SJames Wright 347493642f1SJames Wright CeedQFunctionContextCreate(user->ceed, &stg_context); 3482b916ea7SJeremy L Thompson CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx); 3492b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc); 3502b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution"); 351493642f1SJames Wright 35243bff553SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 35343bff553SJames Wright problem->ics.qfunction = ICsSTG; 35443bff553SJames Wright problem->ics.qfunction_loc = ICsSTG_loc; 35543bff553SJames Wright problem->ics.qfunction_context = stg_context; 35643bff553SJames Wright 357e6098bcdSJames Wright if (use_stgstrong) { 3588085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 3596d0190e2SJames Wright problem->use_dirichlet_ceed = PETSC_TRUE; 360d7244455SJames Wright problem->bc_from_ics = PETSC_FALSE; 361e6098bcdSJames Wright } else { 362493642f1SJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow; 363493642f1SJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc; 364a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction = STGShur14_Inflow_Jacobian; 365a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc; 3662b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context); 3672b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context); 368d7244455SJames Wright problem->bc_from_ics = PETSC_TRUE; 369e6098bcdSJames Wright } 370493642f1SJames Wright 371493642f1SJames Wright PetscFunctionReturn(0); 372493642f1SJames Wright } 373e6098bcdSJames Wright 3742b916ea7SJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) { 375e6098bcdSJames Wright const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]); 376e6098bcdSJames Wright // ^^assuming min(dy) is first element off the wall 377e6098bcdSJames Wright PetscInt idx = -1; // Index 378e6098bcdSJames Wright 379e6098bcdSJames Wright for (PetscInt i = 0; i < nynodes; i++) { 380e6098bcdSJames Wright if (y < ynodes[i] + half_mindy) { 3812b916ea7SJeremy L Thompson idx = i; 3822b916ea7SJeremy L Thompson break; 383e6098bcdSJames Wright } 384e6098bcdSJames Wright } 385e6098bcdSJames Wright if (idx == 0) return ynodes[1] - ynodes[0]; 386e6098bcdSJames Wright else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1]; 387e6098bcdSJames Wright else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]); 388e6098bcdSJames Wright } 389e6098bcdSJames Wright 390e6098bcdSJames Wright // Function passed to DMAddBoundary 3916d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method 3922b916ea7SJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) { 393e6098bcdSJames Wright PetscFunctionBeginUser; 3942b916ea7SJeremy L Thompson 395e6098bcdSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 396e6098bcdSJames Wright PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt; 397e6098bcdSJames Wright const bool mean_only = stg_ctx->mean_only; 398e6098bcdSJames Wright const PetscScalar dx = stg_ctx->dx; 399e6098bcdSJames Wright const PetscScalar dz = stg_ctx->dz; 400e6098bcdSJames Wright const PetscScalar mu = stg_ctx->newtonian_ctx.mu; 401e6098bcdSJames Wright const PetscScalar theta0 = stg_ctx->theta0; 402e6098bcdSJames Wright const PetscScalar P0 = stg_ctx->P0; 403e6098bcdSJames Wright const PetscScalar cv = stg_ctx->newtonian_ctx.cv; 404e6098bcdSJames Wright const PetscScalar cp = stg_ctx->newtonian_ctx.cp; 405e6098bcdSJames Wright const PetscScalar Rd = cp - cv; 406e6098bcdSJames Wright 407e6098bcdSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 408e6098bcdSJames Wright InterpolateProfile(x[1], ubar, cij, &eps, <, stg_ctx); 409e6098bcdSJames Wright if (!mean_only) { 410e6098bcdSJames Wright const PetscInt nynodes = stg_ctx->nynodes; 411e6098bcdSJames Wright const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes]; 412e6098bcdSJames Wright const PetscScalar h[3] = {dx, FindDy(ynodes, nynodes, x[1]), dz}; 413e6098bcdSJames Wright CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx); 414e6098bcdSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 415e6098bcdSJames Wright } else { 416e6098bcdSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 417e6098bcdSJames Wright } 418e6098bcdSJames Wright 419e6098bcdSJames Wright bcval[0] = rho; 420e6098bcdSJames Wright bcval[1] = rho * u[0]; 421e6098bcdSJames Wright bcval[2] = rho * u[1]; 422e6098bcdSJames Wright bcval[3] = rho * u[2]; 423e6098bcdSJames Wright PetscFunctionReturn(0); 424e6098bcdSJames Wright } 425e6098bcdSJames Wright 4262b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) { 427e6098bcdSJames Wright DMLabel label; 428e6098bcdSJames Wright PetscFunctionBeginUser; 429e6098bcdSJames Wright 430d374fb47SJames Wright PetscInt comps[5], num_comps = 4; 4313636f6a4SJames Wright switch (phys->state_var) { 4323636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 433d374fb47SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 434d374fb47SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 4353636f6a4SJames Wright break; 4363636f6a4SJames Wright 4373636f6a4SJames Wright case STATEVAR_PRIMITIVE: 4383636f6a4SJames Wright // {1,2,3,4} for u, v, w, T 4393636f6a4SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 4403636f6a4SJames Wright break; 441d374fb47SJames Wright } 442d374fb47SJames Wright 44334b5deb1SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 444e6098bcdSJames Wright // Set wall BCs 445e6098bcdSJames Wright if (bc->num_inflow > 0) { 4462b916ea7SJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc, 44734b5deb1SJames Wright NULL, global_stg_ctx, NULL)); 448e6098bcdSJames Wright } 449e6098bcdSJames Wright 450e6098bcdSJames Wright PetscFunctionReturn(0); 451e6098bcdSJames Wright } 4526d0190e2SJames Wright 4532b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, 4549eeef72bSJames Wright CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) { 4556d0190e2SJames Wright CeedQFunction qf_strongbc; 4566d0190e2SJames Wright PetscFunctionBeginUser; 4572b916ea7SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc); 4582b916ea7SJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 4596d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 4606d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE); 4619eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 4626d0190e2SJames Wright CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE); 4636d0190e2SJames Wright 4646d0190e2SJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 4656d0190e2SJames Wright *pqf_strongbc = qf_strongbc; 4666d0190e2SJames Wright PetscFunctionReturn(0); 4676d0190e2SJames Wright } 4689eeef72bSJames Wright 4692b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur, 4709eeef72bSJames Wright CeedQFunction *pqf_strongbc) { 4719eeef72bSJames Wright CeedQFunction qf_strongbc; 4729eeef72bSJames Wright PetscFunctionBeginUser; 4732b916ea7SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc); 4742b916ea7SJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 4759eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 4769eeef72bSJames Wright CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 4779eeef72bSJames Wright 4789eeef72bSJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 4799eeef72bSJames Wright *pqf_strongbc = qf_strongbc; 4809eeef72bSJames Wright PetscFunctionReturn(0); 4819eeef72bSJames Wright } 482