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 12493642f1SJames Wright #include <stdlib.h> 13493642f1SJames Wright #include <math.h> 14493642f1SJames Wright #include <petsc.h> 15493642f1SJames Wright #include "../navierstokes.h" 16493642f1SJames Wright #include "stg_shur14.h" 17493642f1SJames Wright #include "../qfunctions/stg_shur14.h" 18493642f1SJames Wright 198085925cSJames Wright STGShur14Context global_stg_ctx; 208085925cSJames Wright 21493642f1SJames Wright /* 22493642f1SJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 23493642f1SJames Wright * 24493642f1SJames Wright * This assumes the input matrices are in order [11,22,33,12,13,23]. This 25493642f1SJames Wright * format is also used for the output. 26493642f1SJames Wright * 27493642f1SJames Wright * @param[in] comm MPI_Comm 28493642f1SJames Wright * @param[in] nprofs Number of matrices in Rij 29493642f1SJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 30493642f1SJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 31493642f1SJames Wright */ 32493642f1SJames Wright PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, 33493642f1SJames Wright const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 34493642f1SJames Wright PetscFunctionBeginUser; 35493642f1SJames Wright for (PetscInt i=0; i<nprofs; i++) { 36493642f1SJames Wright Cij[0][i] = sqrt(Rij[0][i]); 37493642f1SJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 38493642f1SJames Wright Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2) ); 39493642f1SJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 40493642f1SJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i]*Cij[4][i]) / Cij[1][i]; 41493642f1SJames Wright Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 42493642f1SJames Wright 43493642f1SJames Wright if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) 4445aa3cadSJeremy L Thompson SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" 4545aa3cadSJeremy L Thompson PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i+1); 46493642f1SJames Wright } 47493642f1SJames Wright PetscFunctionReturn(0); 48493642f1SJames Wright } 49493642f1SJames Wright 50493642f1SJames Wright 51493642f1SJames Wright /* 52493642f1SJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 53493642f1SJames Wright * 54493642f1SJames Wright * This function opens the file specified by `path` using `PetscFOpen` and 55493642f1SJames Wright * passes the file pointer in `fp`. It is not closed in this function, thus 56493642f1SJames Wright * `fp` must be closed sometime after this function has been called (using 57493642f1SJames Wright * `PetscFClose` for example). 58493642f1SJames Wright * 59493642f1SJames Wright * Assumes that the first line of the file has the number of rows and columns 60493642f1SJames Wright * as the only two entries, separated by a single space 61493642f1SJames Wright * 62493642f1SJames Wright * @param[in] comm MPI_Comm for the program 63493642f1SJames Wright * @param[in] path Path to the file 64493642f1SJames Wright * @param[in] char_array_len Length of the character array that should contain each line 65493642f1SJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 66493642f1SJames Wright * @param[out] fp File pointer to the opened file 67493642f1SJames Wright */ 68493642f1SJames Wright static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm, 69493642f1SJames Wright const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, 70493642f1SJames Wright PetscInt dims[2], FILE **fp) { 71493642f1SJames Wright PetscErrorCode ierr; 72493642f1SJames Wright PetscInt ndims; 73493642f1SJames Wright char line[char_array_len]; 74493642f1SJames Wright char **array; 75493642f1SJames Wright 76493642f1SJames Wright PetscFunctionBeginUser; 77493642f1SJames Wright ierr = PetscFOpen(comm, path, "r", fp); CHKERRQ(ierr); 78493642f1SJames Wright ierr = PetscSynchronizedFGets(comm, *fp, char_array_len, line); CHKERRQ(ierr); 79493642f1SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 80493642f1SJames Wright if (ndims != 2) SETERRQ(comm, -1, 8145aa3cadSJeremy L Thompson "Found %" PetscInt_FMT" dimensions instead of 2 on the first line of %s", 82493642f1SJames Wright ndims, path); 83493642f1SJames Wright 84493642f1SJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 85493642f1SJames Wright ierr = PetscStrToArrayDestroy(ndims, array); CHKERRQ(ierr); 86493642f1SJames Wright PetscFunctionReturn(0); 87493642f1SJames Wright } 88493642f1SJames Wright 89493642f1SJames Wright /* 90493642f1SJames Wright * @brief Get the number of rows for the PHASTA file at path 91493642f1SJames Wright * 92493642f1SJames Wright * Assumes that the first line of the file has the number of rows and columns 93493642f1SJames Wright * as the only two entries, separated by a single space 94493642f1SJames Wright * 95493642f1SJames Wright * @param[in] comm MPI_Comm for the program 96493642f1SJames Wright * @param[in] path Path to the file 97493642f1SJames Wright * @param[out] nrows Number of rows 98493642f1SJames Wright */ 99493642f1SJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm, 100493642f1SJames Wright const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 101493642f1SJames Wright PetscErrorCode ierr; 102493642f1SJames Wright const PetscInt char_array_len = 512; 103493642f1SJames Wright PetscInt dims[2]; 104493642f1SJames Wright FILE *fp; 105493642f1SJames Wright 106493642f1SJames Wright PetscFunctionBeginUser; 107493642f1SJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 108493642f1SJames Wright *nrows = dims[0]; 109493642f1SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 110493642f1SJames Wright PetscFunctionReturn(0); 111493642f1SJames Wright } 112493642f1SJames Wright 113493642f1SJames Wright /* 114493642f1SJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 115493642f1SJames Wright * 116493642f1SJames Wright * Assumes that the first line of the file has the number of rows and columns 117493642f1SJames Wright * as the only two entries, separated by a single space. 118493642f1SJames Wright * Assumes there are 14 columns in the file 119493642f1SJames Wright * 120493642f1SJames Wright * Function calculates the Cholesky decomposition from the Reynolds stress 121493642f1SJames Wright * profile found in the file 122493642f1SJames Wright * 123493642f1SJames Wright * @param[in] comm MPI_Comm for the program 124493642f1SJames Wright * @param[in] path Path to the STGInflow.dat file 125493642f1SJames Wright * @param[inout] stg_ctx STGShur14Context where the data will be loaded into 126493642f1SJames Wright */ 127493642f1SJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, 128493642f1SJames Wright const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 129493642f1SJames Wright PetscErrorCode ierr; 130493642f1SJames Wright PetscInt ndims, dims[2]; 131493642f1SJames Wright FILE *fp; 132493642f1SJames Wright const PetscInt char_array_len=512; 133493642f1SJames Wright char line[char_array_len]; 134493642f1SJames Wright char **array; 135493642f1SJames Wright 136493642f1SJames Wright PetscFunctionBeginUser; 137493642f1SJames Wright 138493642f1SJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 139493642f1SJames Wright 140493642f1SJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 141493642f1SJames Wright CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw]; 142493642f1SJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 143493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 144493642f1SJames Wright CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs]) 145493642f1SJames Wright &stg_ctx->data[stg_ctx->offsets.ubar]; 146493642f1SJames Wright 147493642f1SJames Wright for (PetscInt i=0; i<stg_ctx->nprofs; i++) { 148493642f1SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 149493642f1SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 150493642f1SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 15145aa3cadSJeremy L Thompson "Line %" PetscInt_FMT" of %s does not contain enough columns (%" 15245aa3cadSJeremy L Thompson PetscInt_FMT" instead of %" PetscInt_FMT")", i, 153493642f1SJames Wright path, ndims, dims[1]); 154493642f1SJames Wright 155493642f1SJames Wright prof_dw[i] = (CeedScalar) atof(array[0]); 156493642f1SJames Wright ubar[0][i] = (CeedScalar) atof(array[1]); 157493642f1SJames Wright ubar[1][i] = (CeedScalar) atof(array[2]); 158493642f1SJames Wright ubar[2][i] = (CeedScalar) atof(array[3]); 159493642f1SJames Wright rij[0][i] = (CeedScalar) atof(array[4]); 160493642f1SJames Wright rij[1][i] = (CeedScalar) atof(array[5]); 161493642f1SJames Wright rij[2][i] = (CeedScalar) atof(array[6]); 162493642f1SJames Wright rij[3][i] = (CeedScalar) atof(array[7]); 163493642f1SJames Wright rij[4][i] = (CeedScalar) atof(array[8]); 164493642f1SJames Wright rij[5][i] = (CeedScalar) atof(array[9]); 165493642f1SJames Wright lt[i] = (CeedScalar) atof(array[12]); 166493642f1SJames Wright eps[i] = (CeedScalar) atof(array[13]); 167493642f1SJames Wright 168493642f1SJames Wright if (prof_dw[i] < 0) SETERRQ(comm, -1, 169493642f1SJames Wright "Distance to wall in %s cannot be negative", path); 170493642f1SJames Wright if (lt[i] < 0) SETERRQ(comm, -1, 171493642f1SJames Wright "Turbulent length scale in %s cannot be negative", path); 172493642f1SJames Wright if (eps[i] < 0) SETERRQ(comm, -1, 173493642f1SJames Wright "Turbulent dissipation in %s cannot be negative", path); 174493642f1SJames Wright 175493642f1SJames Wright } 176493642f1SJames Wright CeedScalar (*cij)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs]) 177493642f1SJames Wright &stg_ctx->data[stg_ctx->offsets.cij]; 178493642f1SJames Wright ierr = CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij); CHKERRQ(ierr); 179493642f1SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 180493642f1SJames Wright PetscFunctionReturn(0); 181493642f1SJames Wright } 182493642f1SJames Wright 183493642f1SJames Wright /* 184493642f1SJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 185493642f1SJames Wright * 186493642f1SJames Wright * Assumes that the first line of the file has the number of rows and columns 187493642f1SJames Wright * as the only two entries, separated by a single space. 188493642f1SJames Wright * Assumes there are 7 columns in the file 189493642f1SJames Wright * 190493642f1SJames Wright * @param[in] comm MPI_Comm for the program 191493642f1SJames Wright * @param[in] path Path to the STGRand.dat file 192493642f1SJames Wright * @param[inout] stg_ctx STGShur14Context where the data will be loaded into 193493642f1SJames Wright */ 194493642f1SJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm, 195493642f1SJames Wright const char path[PETSC_MAX_PATH_LEN], 196493642f1SJames Wright STGShur14Context stg_ctx) { 197493642f1SJames Wright PetscErrorCode ierr; 198493642f1SJames Wright PetscInt ndims, dims[2]; 199493642f1SJames Wright FILE *fp; 200493642f1SJames Wright const PetscInt char_array_len = 512; 201493642f1SJames Wright char line[char_array_len]; 202493642f1SJames Wright char **array; 203493642f1SJames Wright 204493642f1SJames Wright PetscFunctionBeginUser; 205493642f1SJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 206493642f1SJames Wright 207493642f1SJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 208493642f1SJames Wright CeedScalar (*d)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes]) 209493642f1SJames Wright &stg_ctx->data[stg_ctx->offsets.d]; 210493642f1SJames Wright CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes]) 211493642f1SJames Wright &stg_ctx->data[stg_ctx->offsets.sigma]; 212493642f1SJames Wright 213493642f1SJames Wright for (PetscInt i=0; i<stg_ctx->nmodes; i++) { 214493642f1SJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 215493642f1SJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 216493642f1SJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 21745aa3cadSJeremy L Thompson "Line %" PetscInt_FMT" of %s does not contain enough columns (%" 21845aa3cadSJeremy L Thompson PetscInt_FMT" instead of %" PetscInt_FMT")", i, 219493642f1SJames Wright path, ndims, dims[1]); 220493642f1SJames Wright 221493642f1SJames Wright d[0][i] = (CeedScalar) atof(array[0]); 222493642f1SJames Wright d[1][i] = (CeedScalar) atof(array[1]); 223493642f1SJames Wright d[2][i] = (CeedScalar) atof(array[2]); 224493642f1SJames Wright phi[i] = (CeedScalar) atof(array[3]); 225493642f1SJames Wright sigma[0][i] = (CeedScalar) atof(array[4]); 226493642f1SJames Wright sigma[1][i] = (CeedScalar) atof(array[5]); 227493642f1SJames Wright sigma[2][i] = (CeedScalar) atof(array[6]); 228493642f1SJames Wright } 229493642f1SJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 230493642f1SJames Wright PetscFunctionReturn(0); 231493642f1SJames Wright } 232493642f1SJames Wright 233493642f1SJames Wright /* 234493642f1SJames Wright * @brief Read STG data from input paths and put in STGShur14Context 235493642f1SJames Wright * 236493642f1SJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 237493642f1SJames Wright * Data stored initially in `*pstg_ctx` will be copied over to the new 238493642f1SJames Wright * STGShur14Context instance. 239493642f1SJames Wright * 240493642f1SJames Wright * @param[in] comm MPI_Comm for the program 241493642f1SJames Wright * @param[in] dm DM for the program 242493642f1SJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 243493642f1SJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 244493642f1SJames Wright * @param[inout] pstg_ctx Pointer to STGShur14Context where the data will be loaded into 245493642f1SJames Wright */ 246493642f1SJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm, 247493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN], 248493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN], 249e6098bcdSJames Wright STGShur14Context *pstg_ctx, 250e6098bcdSJames Wright const CeedScalar ynodes[]) { 251493642f1SJames Wright PetscErrorCode ierr; 252493642f1SJames Wright PetscInt nmodes, nprofs; 253493642f1SJames Wright STGShur14Context stg_ctx; 254493642f1SJames Wright PetscFunctionBeginUser; 255493642f1SJames Wright 256493642f1SJames Wright // Get options 257493642f1SJames Wright ierr = GetNRows(comm, stg_rand_path, &nmodes); CHKERRQ(ierr); 258493642f1SJames Wright ierr = GetNRows(comm, stg_inflow_path, &nprofs); CHKERRQ(ierr); 259493642f1SJames Wright if (nmodes > STG_NMODES_MAX) 26045aa3cadSJeremy L Thompson SETERRQ(comm, 1, "Number of wavemodes in %s (%" 26145aa3cadSJeremy L Thompson PetscInt_FMT") exceeds STG_NMODES_MAX (%" PetscInt_FMT"). " 262493642f1SJames Wright "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes, 263493642f1SJames Wright STG_NMODES_MAX); 264493642f1SJames Wright 265493642f1SJames Wright { 266493642f1SJames Wright STGShur14Context s; 267493642f1SJames Wright ierr = PetscCalloc1(1, &s); CHKERRQ(ierr); 268493642f1SJames Wright *s = **pstg_ctx; 269493642f1SJames Wright s->nmodes = nmodes; 270493642f1SJames Wright s->nprofs = nprofs; 271493642f1SJames Wright s->offsets.sigma = 0; 272493642f1SJames Wright s->offsets.d = nmodes*3; 273493642f1SJames Wright s->offsets.phi = s->offsets.d + nmodes*3; 274493642f1SJames Wright s->offsets.kappa = s->offsets.phi + nmodes; 275493642f1SJames Wright s->offsets.prof_dw = s->offsets.kappa + nmodes; 276493642f1SJames Wright s->offsets.ubar = s->offsets.prof_dw + nprofs; 277493642f1SJames Wright s->offsets.cij = s->offsets.ubar + nprofs*3; 278493642f1SJames Wright s->offsets.eps = s->offsets.cij + nprofs*6; 279493642f1SJames Wright s->offsets.lt = s->offsets.eps + nprofs; 280e6098bcdSJames Wright s->offsets.ynodes = s->offsets.lt + nprofs; 281e6098bcdSJames Wright PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes; 282493642f1SJames Wright s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]); 283493642f1SJames Wright ierr = PetscMalloc(s->total_bytes, &stg_ctx); CHKERRQ(ierr); 284493642f1SJames Wright *stg_ctx = *s; 285493642f1SJames Wright ierr = PetscFree(s); CHKERRQ(ierr); 286493642f1SJames Wright } 287493642f1SJames Wright 288493642f1SJames Wright ierr = ReadSTGInflow(comm, stg_inflow_path, stg_ctx); CHKERRQ(ierr); 289493642f1SJames Wright ierr = ReadSTGRand(comm, stg_rand_path, stg_ctx); CHKERRQ(ierr); 290493642f1SJames Wright 291e6098bcdSJames Wright if (stg_ctx->nynodes > 0) { 292e6098bcdSJames Wright CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes]; 293e6098bcdSJames Wright for (PetscInt i=0; i<stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i]; 294e6098bcdSJames Wright } 295e6098bcdSJames Wright 296493642f1SJames Wright // -- Calculate kappa 297493642f1SJames Wright { 298493642f1SJames Wright CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 299493642f1SJames Wright CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw]; 300493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 301493642f1SJames Wright CeedScalar le, le_max=0; 302493642f1SJames Wright 303493642f1SJames Wright CeedPragmaSIMD 304493642f1SJames Wright for (PetscInt i=0; i<stg_ctx->nprofs; i++) { 305493642f1SJames Wright le = PetscMin(2*prof_dw[i], 3*lt[i]); 306493642f1SJames Wright if (le_max < le) le_max = le; 307493642f1SJames Wright } 308493642f1SJames Wright CeedScalar kmin = M_PI/le_max; 309493642f1SJames Wright 310493642f1SJames Wright CeedPragmaSIMD 311493642f1SJames Wright for (PetscInt i=0; i<stg_ctx->nmodes; i++) { 312493642f1SJames Wright kappa[i] = kmin*pow(stg_ctx->alpha, i); 313493642f1SJames Wright } 314493642f1SJames Wright } //end calculate kappa 315493642f1SJames Wright 31676a21e77SJames Wright ierr = PetscFree(*pstg_ctx); CHKERRQ(ierr); 317493642f1SJames Wright *pstg_ctx = stg_ctx; 318493642f1SJames Wright PetscFunctionReturn(0); 319493642f1SJames Wright } 320493642f1SJames Wright 321493642f1SJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, 322493642f1SJames Wright User user, const bool prescribe_T, 323e6098bcdSJames Wright const CeedScalar theta0, const CeedScalar P0, 324e6098bcdSJames Wright const CeedScalar ynodes[], const CeedInt nynodes) { 325493642f1SJames Wright PetscErrorCode ierr; 326493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 327493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 328e6098bcdSJames Wright PetscBool mean_only = PETSC_FALSE, 329e6098bcdSJames Wright use_stgstrong = PETSC_FALSE; 330e6098bcdSJames Wright CeedScalar u0 = 0.0, 331e6098bcdSJames Wright alpha = 1.01; 332493642f1SJames Wright CeedQFunctionContext stg_context; 333493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 334493642f1SJames Wright PetscFunctionBeginUser; 335493642f1SJames Wright 336493642f1SJames Wright // Get options 337493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 338493642f1SJames Wright ierr = PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, 339493642f1SJames Wright stg_inflow_path, stg_inflow_path, 340493642f1SJames Wright sizeof(stg_inflow_path), NULL); CHKERRQ(ierr); 341493642f1SJames Wright ierr = PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, 342493642f1SJames Wright stg_rand_path,stg_rand_path, 343493642f1SJames Wright sizeof(stg_rand_path), NULL); CHKERRQ(ierr); 344493642f1SJames Wright ierr = PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, 345493642f1SJames Wright alpha, &alpha, NULL); CHKERRQ(ierr); 346493642f1SJames Wright ierr = PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", 347493642f1SJames Wright NULL, u0, &u0, NULL); CHKERRQ(ierr); 348493642f1SJames Wright ierr = PetscOptionsBool("-stg_mean_only", "Only apply mean profile", 349493642f1SJames Wright NULL, mean_only, &mean_only, NULL); CHKERRQ(ierr); 350e6098bcdSJames Wright ierr = PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", 351e6098bcdSJames Wright NULL, use_stgstrong, &use_stgstrong, NULL); CHKERRQ(ierr); 352493642f1SJames Wright PetscOptionsEnd(); 353493642f1SJames Wright 3548085925cSJames Wright ierr = PetscCalloc1(1, &global_stg_ctx); CHKERRQ(ierr); 3558085925cSJames Wright global_stg_ctx->alpha = alpha; 3568085925cSJames Wright global_stg_ctx->u0 = u0; 3578085925cSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 3588085925cSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 3598085925cSJames Wright global_stg_ctx->mean_only = mean_only; 3608085925cSJames Wright global_stg_ctx->theta0 = theta0; 3618085925cSJames Wright global_stg_ctx->P0 = P0; 3628085925cSJames Wright global_stg_ctx->nynodes = nynodes; 363493642f1SJames Wright 364493642f1SJames Wright { 365493642f1SJames Wright // Calculate dx assuming constant spacing 366493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 367493642f1SJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 368493642f1SJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 369493642f1SJames Wright 370493642f1SJames Wright PetscInt nmax = 3, faces[3]; 371493642f1SJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 372493642f1SJames Wright NULL); CHKERRQ(ierr); 3738085925cSJames Wright global_stg_ctx->dx = domain_size[0]/faces[0]; 3748085925cSJames Wright global_stg_ctx->dz = domain_size[2]/faces[2]; 375493642f1SJames Wright } 376493642f1SJames Wright 377493642f1SJames Wright CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 378493642f1SJames Wright CEED_MEM_HOST, &newtonian_ig_ctx); 3798085925cSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 380493642f1SJames Wright CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 381493642f1SJames Wright &newtonian_ig_ctx); 382493642f1SJames Wright 3838085925cSJames Wright ierr = GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, 3848085925cSJames Wright &global_stg_ctx, ynodes); CHKERRQ(ierr); 385493642f1SJames Wright 386493642f1SJames Wright CeedQFunctionContextCreate(user->ceed, &stg_context); 387493642f1SJames Wright CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, 3888085925cSJames Wright CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx); 389493642f1SJames Wright CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, 390493642f1SJames Wright FreeContextPetsc); 391493642f1SJames Wright CeedQFunctionContextRegisterDouble(stg_context, "solution time", 392493642f1SJames Wright offsetof(struct STGShur14Context_, time), 1, 3932c498363SJed Brown "Physical time of the solution"); 394493642f1SJames Wright 39543bff553SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 39643bff553SJames Wright problem->ics.qfunction = ICsSTG; 39743bff553SJames Wright problem->ics.qfunction_loc = ICsSTG_loc; 39843bff553SJames Wright problem->ics.qfunction_context = stg_context; 39943bff553SJames Wright 400e6098bcdSJames Wright if (use_stgstrong) { 4018085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 4026d0190e2SJames Wright problem->use_dirichlet_ceed = PETSC_TRUE; 403d7244455SJames Wright problem->bc_from_ics = PETSC_FALSE; 404e6098bcdSJames Wright } else { 405493642f1SJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow; 406493642f1SJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc; 407a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction = STGShur14_Inflow_Jacobian; 408a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc; 409a6e8f989SJames Wright CeedQFunctionContextReferenceCopy(stg_context, 410a6e8f989SJames Wright &problem->apply_inflow.qfunction_context); 411a6e8f989SJames Wright CeedQFunctionContextReferenceCopy(stg_context, 412a6e8f989SJames Wright &problem->apply_inflow_jacobian.qfunction_context); 413d7244455SJames Wright problem->bc_from_ics = PETSC_TRUE; 414e6098bcdSJames Wright } 415493642f1SJames Wright 416493642f1SJames Wright PetscFunctionReturn(0); 417493642f1SJames Wright } 418e6098bcdSJames Wright 419e6098bcdSJames Wright static inline PetscScalar FindDy(const PetscScalar ynodes[], 420e6098bcdSJames Wright const PetscInt nynodes, const PetscScalar y) { 4216dd99bedSLeila Ghaffari 422e6098bcdSJames Wright const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]); 423e6098bcdSJames Wright // ^^assuming min(dy) is first element off the wall 424e6098bcdSJames Wright PetscInt idx = -1; // Index 425e6098bcdSJames Wright 426e6098bcdSJames Wright for (PetscInt i=0; i<nynodes; i++) { 427e6098bcdSJames Wright if (y < ynodes[i] + half_mindy) { 428e6098bcdSJames Wright idx = i; break; 429e6098bcdSJames Wright } 430e6098bcdSJames Wright } 431e6098bcdSJames Wright if (idx == 0) return ynodes[1] - ynodes[0]; 432e6098bcdSJames Wright else if (idx == nynodes-1) return ynodes[nynodes-2] - ynodes[nynodes-1]; 433e6098bcdSJames Wright else return 0.5 * (ynodes[idx+1] - ynodes[idx-1]); 434e6098bcdSJames Wright } 435e6098bcdSJames Wright 436e6098bcdSJames Wright // Function passed to DMAddBoundary 4376d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method 438e6098bcdSJames Wright PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, 439e6098bcdSJames Wright const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) { 440e6098bcdSJames Wright PetscFunctionBeginUser; 441e6098bcdSJames Wright 442e6098bcdSJames Wright const STGShur14Context stg_ctx = (STGShur14Context) ctx; 443e6098bcdSJames Wright PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt; 444e6098bcdSJames Wright const bool mean_only = stg_ctx->mean_only; 445e6098bcdSJames Wright const PetscScalar dx = stg_ctx->dx; 446e6098bcdSJames Wright const PetscScalar dz = stg_ctx->dz; 447e6098bcdSJames Wright const PetscScalar mu = stg_ctx->newtonian_ctx.mu; 448e6098bcdSJames Wright const PetscScalar theta0 = stg_ctx->theta0; 449e6098bcdSJames Wright const PetscScalar P0 = stg_ctx->P0; 450e6098bcdSJames Wright const PetscScalar cv = stg_ctx->newtonian_ctx.cv; 451e6098bcdSJames Wright const PetscScalar cp = stg_ctx->newtonian_ctx.cp; 452e6098bcdSJames Wright const PetscScalar Rd = cp - cv; 453e6098bcdSJames Wright 454e6098bcdSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 455e6098bcdSJames Wright InterpolateProfile(x[1], ubar, cij, &eps, <, stg_ctx); 456e6098bcdSJames Wright if (!mean_only) { 457e6098bcdSJames Wright const PetscInt nynodes = stg_ctx->nynodes; 458e6098bcdSJames Wright const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes]; 459e6098bcdSJames Wright const PetscScalar h[3] = {dx, FindDy(ynodes, nynodes, x[1]), dz}; 460e6098bcdSJames Wright CalcSpectrum(x[1], eps, lt, h, mu/rho, qn, stg_ctx); 461e6098bcdSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 462e6098bcdSJames Wright } else { 463e6098bcdSJames Wright for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 464e6098bcdSJames Wright } 465e6098bcdSJames Wright 466e6098bcdSJames Wright bcval[0] = rho; 467e6098bcdSJames Wright bcval[1] = rho * u[0]; 468e6098bcdSJames Wright bcval[2] = rho * u[1]; 469e6098bcdSJames Wright bcval[3] = rho * u[2]; 470e6098bcdSJames Wright PetscFunctionReturn(0); 471e6098bcdSJames Wright } 472e6098bcdSJames Wright 4738085925cSJames Wright PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem) { 474e6098bcdSJames Wright PetscErrorCode ierr; 475e6098bcdSJames Wright DMLabel label; 476e6098bcdSJames Wright const PetscInt comps[] = {0, 1, 2, 3}; 477e6098bcdSJames Wright const PetscInt num_comps = 4; 478e6098bcdSJames Wright PetscFunctionBeginUser; 479e6098bcdSJames Wright 480e6098bcdSJames Wright ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 481e6098bcdSJames Wright // Set wall BCs 482e6098bcdSJames Wright if (bc->num_inflow > 0) { 483e6098bcdSJames Wright ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, 484e6098bcdSJames Wright bc->num_inflow, bc->inflows, 0, num_comps, 485e6098bcdSJames Wright comps, (void(*)(void))StrongSTGbcFunc, 4868085925cSJames Wright NULL, global_stg_ctx, NULL); CHKERRQ(ierr); 487e6098bcdSJames Wright } 488e6098bcdSJames Wright 489e6098bcdSJames Wright PetscFunctionReturn(0); 490e6098bcdSJames Wright } 4916d0190e2SJames Wright 4926d0190e2SJames Wright PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, 493*9eeef72bSJames Wright CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, 494*9eeef72bSJames Wright CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) { 4956d0190e2SJames Wright 4966d0190e2SJames Wright CeedQFunction qf_strongbc; 4976d0190e2SJames Wright PetscFunctionBeginUser; 4986d0190e2SJames Wright CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, 4996d0190e2SJames Wright STGShur14_Inflow_StrongQF_loc, &qf_strongbc); 5006d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, 5016d0190e2SJames Wright CEED_EVAL_NONE); 5026d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 5036d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE); 504*9eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 5056d0190e2SJames Wright CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE); 5066d0190e2SJames Wright 5076d0190e2SJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 5086d0190e2SJames Wright *pqf_strongbc = qf_strongbc; 5096d0190e2SJames Wright PetscFunctionReturn(0); 5106d0190e2SJames Wright } 511*9eeef72bSJames Wright 512*9eeef72bSJames Wright PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, 513*9eeef72bSJames Wright CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur, 514*9eeef72bSJames Wright CeedQFunction *pqf_strongbc) { 515*9eeef72bSJames Wright 516*9eeef72bSJames Wright CeedQFunction qf_strongbc; 517*9eeef72bSJames Wright PetscFunctionBeginUser; 518*9eeef72bSJames Wright CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, 519*9eeef72bSJames Wright Preprocess_STGShur14_loc, &qf_strongbc); 520*9eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, 521*9eeef72bSJames Wright CEED_EVAL_NONE); 522*9eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 523*9eeef72bSJames Wright CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 524*9eeef72bSJames Wright 525*9eeef72bSJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 526*9eeef72bSJames Wright *pqf_strongbc = qf_strongbc; 527*9eeef72bSJames Wright PetscFunctionReturn(0); 528*9eeef72bSJames Wright } 529