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 12ba6664aeSJames Wright #include <stdlib.h> 13ba6664aeSJames Wright #include <math.h> 14ba6664aeSJames Wright #include <petsc.h> 15ba6664aeSJames Wright #include "../navierstokes.h" 16ba6664aeSJames Wright #include "stg_shur14.h" 17ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h" 18ba6664aeSJames Wright 19ba6664aeSJames Wright #ifndef M_PI 20ba6664aeSJames Wright #define M_PI 3.14159265358979323846 21ba6664aeSJames Wright #endif 22ba6664aeSJames Wright 23ba6664aeSJames Wright /* 24ba6664aeSJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 25ba6664aeSJames Wright * 26ba6664aeSJames Wright * This assumes the input matrices are in order [11,22,33,12,13,23]. This 27ba6664aeSJames Wright * format is also used for the output. 28ba6664aeSJames Wright * 29ba6664aeSJames Wright * @param[in] comm MPI_Comm 30ba6664aeSJames Wright * @param[in] nprofs Number of matrices in Rij 31ba6664aeSJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 32ba6664aeSJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 33ba6664aeSJames Wright */ 34ba6664aeSJames Wright PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, 35ba6664aeSJames Wright const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 36ba6664aeSJames Wright 37ba6664aeSJames Wright PetscFunctionBeginUser; 38ba6664aeSJames Wright for (PetscInt i=0; i<nprofs; i++) { 39ba6664aeSJames Wright Cij[0][i] = sqrt(Rij[0][i]); 40ba6664aeSJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 41ba6664aeSJames Wright Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2) ); 42ba6664aeSJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 43ba6664aeSJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i]*Cij[4][i]) / Cij[1][i]; 44ba6664aeSJames Wright Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 45ba6664aeSJames Wright 46ba6664aeSJames Wright if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) 47ba6664aeSJames Wright SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %d. " 48ba6664aeSJames Wright "Either STGInflow has non-SPD matrix or contains nan.", i+1); 49ba6664aeSJames Wright } 50ba6664aeSJames Wright PetscFunctionReturn(0); 51ba6664aeSJames Wright } 52ba6664aeSJames Wright 53ba6664aeSJames Wright 54ba6664aeSJames Wright /* 55ba6664aeSJames Wright * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer 56ba6664aeSJames Wright * 57ba6664aeSJames Wright * This function opens the file specified by `path` using `PetscFOpen` and 58ba6664aeSJames Wright * passes the file pointer in `fp`. It is not closed in this function, thus 59ba6664aeSJames Wright * `fp` must be closed sometime after this function has been called (using 60ba6664aeSJames Wright * `PetscFClose` for example). 61ba6664aeSJames Wright * 62ba6664aeSJames Wright * Assumes that the first line of the file has the number of rows and columns 63ba6664aeSJames Wright * as the only two entries, separated by a single space 64ba6664aeSJames Wright * 65ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 66ba6664aeSJames Wright * @param[in] path Path to the file 67ba6664aeSJames Wright * @param[in] char_array_len Length of the character array that should contain each line 68ba6664aeSJames Wright * @param[out] dims Dimensions of the file, taken from the first line of the file 69ba6664aeSJames Wright * @param[out] fp File pointer to the opened file 70ba6664aeSJames Wright */ 71ba6664aeSJames Wright static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm, 72ba6664aeSJames Wright const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, 73ba6664aeSJames Wright PetscInt dims[2], FILE **fp) { 74ba6664aeSJames Wright PetscErrorCode ierr; 75ba6664aeSJames Wright PetscInt ndims; 76ba6664aeSJames Wright char line[char_array_len]; 77ba6664aeSJames Wright char **array; 78ba6664aeSJames Wright 79ba6664aeSJames Wright PetscFunctionBeginUser; 80ba6664aeSJames Wright ierr = PetscFOpen(comm, path, "r", fp); CHKERRQ(ierr); 81ba6664aeSJames Wright ierr = PetscSynchronizedFGets(comm, *fp, char_array_len, line); CHKERRQ(ierr); 82ba6664aeSJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 83ba6664aeSJames Wright if (ndims != 2) SETERRQ(comm, -1, 84ba6664aeSJames Wright "Found %d dimensions instead of 2 on the first line of %s", 85ba6664aeSJames Wright ndims, path); 86ba6664aeSJames Wright 87ba6664aeSJames Wright for (PetscInt i=0; i<ndims; i++) dims[i] = atoi(array[i]); 88ba6664aeSJames Wright ierr = PetscStrToArrayDestroy(ndims, array); CHKERRQ(ierr); 89ba6664aeSJames Wright PetscFunctionReturn(0); 90ba6664aeSJames Wright } 91ba6664aeSJames Wright 92ba6664aeSJames Wright /* 93ba6664aeSJames Wright * @brief Get the number of rows for the PHASTA file at path 94ba6664aeSJames Wright * 95ba6664aeSJames Wright * Assumes that the first line of the file has the number of rows and columns 96ba6664aeSJames Wright * as the only two entries, separated by a single space 97ba6664aeSJames Wright * 98ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 99ba6664aeSJames Wright * @param[in] path Path to the file 100ba6664aeSJames Wright * @param[out] nrows Number of rows 101ba6664aeSJames Wright */ 102ba6664aeSJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm, 103ba6664aeSJames Wright const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) { 104ba6664aeSJames Wright PetscErrorCode ierr; 105ba6664aeSJames Wright const PetscInt char_array_len = 512; 106ba6664aeSJames Wright PetscInt dims[2]; 107ba6664aeSJames Wright FILE *fp; 108ba6664aeSJames Wright 109ba6664aeSJames Wright PetscFunctionBeginUser; 110ba6664aeSJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 111ba6664aeSJames Wright *nrows = dims[0]; 112ba6664aeSJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 113ba6664aeSJames Wright PetscFunctionReturn(0); 114ba6664aeSJames Wright } 115ba6664aeSJames Wright 116ba6664aeSJames Wright /* 117ba6664aeSJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 118ba6664aeSJames Wright * 119ba6664aeSJames Wright * Assumes that the first line of the file has the number of rows and columns 120ba6664aeSJames Wright * as the only two entries, separated by a single space. 121ba6664aeSJames Wright * Assumes there are 14 columns in the file 122ba6664aeSJames Wright * 123ba6664aeSJames Wright * Function calculates the Cholesky decomposition from the Reynolds stress 124ba6664aeSJames Wright * profile found in the file 125ba6664aeSJames Wright * 126ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 127ba6664aeSJames Wright * @param[in] path Path to the STGInflow.dat file 128ba6664aeSJames Wright * @param[inout] stg_ctx STGShur14Context where the data will be loaded into 129ba6664aeSJames Wright */ 130ba6664aeSJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, 131ba6664aeSJames Wright const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 132ba6664aeSJames Wright PetscErrorCode ierr; 133ba6664aeSJames Wright PetscInt ndims, dims[2]; 134ba6664aeSJames Wright FILE *fp; 135ba6664aeSJames Wright const PetscInt char_array_len=512; 136ba6664aeSJames Wright char line[char_array_len]; 137ba6664aeSJames Wright char **array; 138ba6664aeSJames Wright 139ba6664aeSJames Wright PetscFunctionBeginUser; 140ba6664aeSJames Wright 141ba6664aeSJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 142ba6664aeSJames Wright 143ba6664aeSJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 144ba6664aeSJames Wright CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw]; 145ba6664aeSJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 146ba6664aeSJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 147ba6664aeSJames Wright CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs]) 148ba6664aeSJames Wright &stg_ctx->data[stg_ctx->offsets.ubar]; 149ba6664aeSJames Wright 150ba6664aeSJames Wright for (PetscInt i=0; i<stg_ctx->nprofs; i++) { 151ba6664aeSJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 152ba6664aeSJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 153ba6664aeSJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 154ba6664aeSJames Wright "Line %d of %s does not contain enough columns (%d instead of %d)", i, 155ba6664aeSJames Wright path, ndims, dims[1]); 156ba6664aeSJames Wright 157ba6664aeSJames Wright prof_dw[i] = (CeedScalar) atof(array[0]); 158ba6664aeSJames Wright ubar[0][i] = (CeedScalar) atof(array[1]); 159ba6664aeSJames Wright ubar[1][i] = (CeedScalar) atof(array[2]); 160ba6664aeSJames Wright ubar[2][i] = (CeedScalar) atof(array[3]); 161ba6664aeSJames Wright rij[0][i] = (CeedScalar) atof(array[4]); 162ba6664aeSJames Wright rij[1][i] = (CeedScalar) atof(array[5]); 163ba6664aeSJames Wright rij[2][i] = (CeedScalar) atof(array[6]); 164ba6664aeSJames Wright rij[3][i] = (CeedScalar) atof(array[7]); 165ba6664aeSJames Wright rij[4][i] = (CeedScalar) atof(array[8]); 166ba6664aeSJames Wright rij[5][i] = (CeedScalar) atof(array[9]); 167ba6664aeSJames Wright lt[i] = (CeedScalar) atof(array[12]); 168ba6664aeSJames Wright eps[i] = (CeedScalar) atof(array[13]); 169ba6664aeSJames Wright 170ba6664aeSJames Wright if (prof_dw[i] < 0) SETERRQ(comm, -1, 171ba6664aeSJames Wright "Distance to wall in %s cannot be negative", path); 172ba6664aeSJames Wright if (lt[i] < 0) SETERRQ(comm, -1, 173ba6664aeSJames Wright "Turbulent length scale in %s cannot be negative", path); 174ba6664aeSJames Wright if (eps[i] < 0) SETERRQ(comm, -1, 175ba6664aeSJames Wright "Turbulent dissipation in %s cannot be negative", path); 176ba6664aeSJames Wright 177ba6664aeSJames Wright } 178ba6664aeSJames Wright CeedScalar (*cij)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs]) 179ba6664aeSJames Wright &stg_ctx->data[stg_ctx->offsets.cij]; 180ba6664aeSJames Wright ierr = CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij); CHKERRQ(ierr); 181ba6664aeSJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 182ba6664aeSJames Wright PetscFunctionReturn(0); 183ba6664aeSJames Wright } 184ba6664aeSJames Wright 185ba6664aeSJames Wright /* 186ba6664aeSJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 187ba6664aeSJames Wright * 188ba6664aeSJames Wright * Assumes that the first line of the file has the number of rows and columns 189ba6664aeSJames Wright * as the only two entries, separated by a single space. 190ba6664aeSJames Wright * Assumes there are 7 columns in the file 191ba6664aeSJames Wright * 192ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 193ba6664aeSJames Wright * @param[in] path Path to the STGRand.dat file 194ba6664aeSJames Wright * @param[inout] stg_ctx STGShur14Context where the data will be loaded into 195ba6664aeSJames Wright */ 196ba6664aeSJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm, 197ba6664aeSJames Wright const char path[PETSC_MAX_PATH_LEN], 198ba6664aeSJames Wright STGShur14Context stg_ctx) { 199ba6664aeSJames Wright 200ba6664aeSJames Wright PetscErrorCode ierr; 201ba6664aeSJames Wright PetscInt ndims, dims[2]; 202ba6664aeSJames Wright FILE *fp; 203ba6664aeSJames Wright const PetscInt char_array_len = 512; 204ba6664aeSJames Wright char line[char_array_len]; 205ba6664aeSJames Wright char **array; 206ba6664aeSJames Wright 207ba6664aeSJames Wright PetscFunctionBeginUser; 208ba6664aeSJames Wright ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr); 209ba6664aeSJames Wright 210ba6664aeSJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 211ba6664aeSJames Wright CeedScalar (*d)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes]) 212ba6664aeSJames Wright &stg_ctx->data[stg_ctx->offsets.d]; 213ba6664aeSJames Wright CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes]) 214ba6664aeSJames Wright &stg_ctx->data[stg_ctx->offsets.sigma]; 215ba6664aeSJames Wright 216ba6664aeSJames Wright for (PetscInt i=0; i<stg_ctx->nmodes; i++) { 217ba6664aeSJames Wright ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr); 218ba6664aeSJames Wright ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr); 219ba6664aeSJames Wright if (ndims < dims[1]) SETERRQ(comm, -1, 220ba6664aeSJames Wright "Line %d of %s does not contain enough columns (%d instead of %d)", i, 221ba6664aeSJames Wright path, ndims, dims[1]); 222ba6664aeSJames Wright 223ba6664aeSJames Wright d[0][i] = (CeedScalar) atof(array[0]); 224ba6664aeSJames Wright d[1][i] = (CeedScalar) atof(array[1]); 225ba6664aeSJames Wright d[2][i] = (CeedScalar) atof(array[2]); 226ba6664aeSJames Wright phi[i] = (CeedScalar) atof(array[3]); 227ba6664aeSJames Wright sigma[0][i] = (CeedScalar) atof(array[4]); 228ba6664aeSJames Wright sigma[1][i] = (CeedScalar) atof(array[5]); 229ba6664aeSJames Wright sigma[2][i] = (CeedScalar) atof(array[6]); 230ba6664aeSJames Wright } 231ba6664aeSJames Wright ierr = PetscFClose(comm, fp); CHKERRQ(ierr); 232ba6664aeSJames Wright PetscFunctionReturn(0); 233ba6664aeSJames Wright } 234ba6664aeSJames Wright 235ba6664aeSJames Wright /* 236ba6664aeSJames Wright * @brief Read STG data from input paths and put in STGShur14Context 237ba6664aeSJames Wright * 238ba6664aeSJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 239ba6664aeSJames Wright * Data stored initially in `*pstg_ctx` will be copied over to the new 240ba6664aeSJames Wright * STGShur14Context instance. 241ba6664aeSJames Wright * 242ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 243ba6664aeSJames Wright * @param[in] dm DM for the program 244ba6664aeSJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 245ba6664aeSJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 246ba6664aeSJames Wright * @param[inout] pstg_ctx Pointer to STGShur14Context where the data will be loaded into 247ba6664aeSJames Wright */ 248ba6664aeSJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm, 249ba6664aeSJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN], 250ba6664aeSJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN], 2514e139266SJames Wright STGShur14Context *pstg_ctx, 2524e139266SJames Wright const CeedScalar ynodes[]) { 253ba6664aeSJames Wright PetscErrorCode ierr; 254ba6664aeSJames Wright PetscInt nmodes, nprofs; 255ba6664aeSJames Wright STGShur14Context stg_ctx; 256ba6664aeSJames Wright PetscFunctionBeginUser; 257ba6664aeSJames Wright 258ba6664aeSJames Wright // Get options 259ba6664aeSJames Wright ierr = GetNRows(comm, stg_rand_path, &nmodes); CHKERRQ(ierr); 260ba6664aeSJames Wright ierr = GetNRows(comm, stg_inflow_path, &nprofs); CHKERRQ(ierr); 261ba6664aeSJames Wright if (nmodes > STG_NMODES_MAX) 262ba6664aeSJames Wright SETERRQ(comm, 1, "Number of wavemodes in %s (%d) exceeds STG_NMODES_MAX (%d). " 263ba6664aeSJames Wright "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes, 264ba6664aeSJames Wright STG_NMODES_MAX); 265ba6664aeSJames Wright 266ba6664aeSJames Wright { 267ba6664aeSJames Wright STGShur14Context s; 268ba6664aeSJames Wright ierr = PetscCalloc1(1, &s); CHKERRQ(ierr); 269ba6664aeSJames Wright *s = **pstg_ctx; 270ba6664aeSJames Wright s->nmodes = nmodes; 271ba6664aeSJames Wright s->nprofs = nprofs; 272ba6664aeSJames Wright s->offsets.sigma = 0; 273ba6664aeSJames Wright s->offsets.d = nmodes*3; 274ba6664aeSJames Wright s->offsets.phi = s->offsets.d + nmodes*3; 275ba6664aeSJames Wright s->offsets.kappa = s->offsets.phi + nmodes; 276ba6664aeSJames Wright s->offsets.prof_dw = s->offsets.kappa + nmodes; 277ba6664aeSJames Wright s->offsets.ubar = s->offsets.prof_dw + nprofs; 278ba6664aeSJames Wright s->offsets.cij = s->offsets.ubar + nprofs*3; 279ba6664aeSJames Wright s->offsets.eps = s->offsets.cij + nprofs*6; 280ba6664aeSJames Wright s->offsets.lt = s->offsets.eps + nprofs; 2814e139266SJames Wright s->offsets.ynodes = s->offsets.lt + nprofs; 2824e139266SJames Wright PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes; 283ba6664aeSJames Wright s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]); 284ba6664aeSJames Wright ierr = PetscMalloc(s->total_bytes, &stg_ctx); CHKERRQ(ierr); 285ba6664aeSJames Wright *stg_ctx = *s; 286ba6664aeSJames Wright ierr = PetscFree(s); CHKERRQ(ierr); 287ba6664aeSJames Wright } 288ba6664aeSJames Wright 289ba6664aeSJames Wright ierr = ReadSTGInflow(comm, stg_inflow_path, stg_ctx); CHKERRQ(ierr); 290ba6664aeSJames Wright ierr = ReadSTGRand(comm, stg_rand_path, stg_ctx); CHKERRQ(ierr); 291ba6664aeSJames Wright 2924e139266SJames Wright if (stg_ctx->nynodes > 0) { 2934e139266SJames Wright CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes]; 2944e139266SJames Wright for (PetscInt i=0; i<stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i]; 2954e139266SJames Wright } 2964e139266SJames Wright 297ba6664aeSJames Wright // -- Calculate kappa 298ba6664aeSJames Wright { 299ba6664aeSJames Wright CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 300ba6664aeSJames Wright CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw]; 301ba6664aeSJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 302ba6664aeSJames Wright CeedScalar le, le_max=0; 303ba6664aeSJames Wright 304ba6664aeSJames Wright CeedPragmaSIMD 305ba6664aeSJames Wright for (PetscInt i=0; i<stg_ctx->nprofs; i++) { 306ba6664aeSJames Wright le = PetscMin(2*prof_dw[i], 3*lt[i]); 307ba6664aeSJames Wright if (le_max < le) le_max = le; 308ba6664aeSJames Wright } 309ba6664aeSJames Wright CeedScalar kmin = M_PI/le_max; 310ba6664aeSJames Wright 311ba6664aeSJames Wright CeedPragmaSIMD 312ba6664aeSJames Wright for (PetscInt i=0; i<stg_ctx->nmodes; i++) { 313ba6664aeSJames Wright kappa[i] = kmin*pow(stg_ctx->alpha, i); 314ba6664aeSJames Wright } 315ba6664aeSJames Wright } //end calculate kappa 316ba6664aeSJames Wright 317d95c92dbSJames Wright ierr = PetscFree(*pstg_ctx); CHKERRQ(ierr); 318ba6664aeSJames Wright *pstg_ctx = stg_ctx; 319ba6664aeSJames Wright PetscFunctionReturn(0); 320ba6664aeSJames Wright } 321ba6664aeSJames Wright 322ba6664aeSJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, 323ba6664aeSJames Wright User user, const bool prescribe_T, 3244e139266SJames Wright const CeedScalar theta0, const CeedScalar P0, 3254e139266SJames Wright const CeedScalar ynodes[], const CeedInt nynodes) { 326ba6664aeSJames Wright PetscErrorCode ierr; 327ba6664aeSJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 328ba6664aeSJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 3294e139266SJames Wright PetscBool mean_only = PETSC_FALSE, 3304e139266SJames Wright use_stgstrong = PETSC_FALSE; 3314e139266SJames Wright CeedScalar u0 = 0.0, 3324e139266SJames Wright alpha = 1.01; 333ba6664aeSJames Wright STGShur14Context stg_ctx; 334ba6664aeSJames Wright CeedQFunctionContext stg_context; 335ba6664aeSJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 336ba6664aeSJames Wright PetscFunctionBeginUser; 337ba6664aeSJames Wright 338ba6664aeSJames Wright // Get options 339ba6664aeSJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 340ba6664aeSJames Wright ierr = PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, 341ba6664aeSJames Wright stg_inflow_path, stg_inflow_path, 342ba6664aeSJames Wright sizeof(stg_inflow_path), NULL); CHKERRQ(ierr); 343ba6664aeSJames Wright ierr = PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, 344ba6664aeSJames Wright stg_rand_path,stg_rand_path, 345ba6664aeSJames Wright sizeof(stg_rand_path), NULL); CHKERRQ(ierr); 346ba6664aeSJames Wright ierr = PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, 347ba6664aeSJames Wright alpha, &alpha, NULL); CHKERRQ(ierr); 348ba6664aeSJames Wright ierr = PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", 349ba6664aeSJames Wright NULL, u0, &u0, NULL); CHKERRQ(ierr); 350ba6664aeSJames Wright ierr = PetscOptionsBool("-stg_mean_only", "Only apply mean profile", 351ba6664aeSJames Wright NULL, mean_only, &mean_only, NULL); CHKERRQ(ierr); 3524e139266SJames Wright ierr = PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", 3534e139266SJames Wright NULL, use_stgstrong, &use_stgstrong, NULL); CHKERRQ(ierr); 354ba6664aeSJames Wright PetscOptionsEnd(); 355ba6664aeSJames Wright 356ba6664aeSJames Wright ierr = PetscCalloc1(1, &stg_ctx); CHKERRQ(ierr); 357ba6664aeSJames Wright stg_ctx->alpha = alpha; 358ba6664aeSJames Wright stg_ctx->u0 = u0; 359ba6664aeSJames Wright stg_ctx->is_implicit = user->phys->implicit; 360ba6664aeSJames Wright stg_ctx->prescribe_T = prescribe_T; 361ba6664aeSJames Wright stg_ctx->mean_only = mean_only; 362ba6664aeSJames Wright stg_ctx->theta0 = theta0; 363ba6664aeSJames Wright stg_ctx->P0 = P0; 3644e139266SJames Wright stg_ctx->nynodes = nynodes; 365ba6664aeSJames Wright 366ba6664aeSJames Wright { 367ba6664aeSJames Wright // Calculate dx assuming constant spacing 368ba6664aeSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 369ba6664aeSJames Wright ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr); 370ba6664aeSJames Wright for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 371ba6664aeSJames Wright 372ba6664aeSJames Wright PetscInt nmax = 3, faces[3]; 373ba6664aeSJames Wright ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, 374ba6664aeSJames Wright NULL); CHKERRQ(ierr); 375ba6664aeSJames Wright stg_ctx->dx = domain_size[0]/faces[0]; 3764e139266SJames Wright stg_ctx->dz = domain_size[2]/faces[2]; 377ba6664aeSJames Wright } 378ba6664aeSJames Wright 379ba6664aeSJames Wright CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, 380ba6664aeSJames Wright CEED_MEM_HOST, &newtonian_ig_ctx); 381ba6664aeSJames Wright stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 382ba6664aeSJames Wright CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, 383ba6664aeSJames Wright &newtonian_ig_ctx); 384ba6664aeSJames Wright 3854e139266SJames Wright ierr = GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &stg_ctx, 3864e139266SJames Wright ynodes); CHKERRQ(ierr); 387ba6664aeSJames Wright 388ba6664aeSJames Wright CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context); 389ba6664aeSJames Wright CeedQFunctionContextCreate(user->ceed, &stg_context); 390ba6664aeSJames Wright CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, 391ba6664aeSJames Wright CEED_USE_POINTER, stg_ctx->total_bytes, stg_ctx); 392ba6664aeSJames Wright CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, 393ba6664aeSJames Wright FreeContextPetsc); 394ba6664aeSJames Wright CeedQFunctionContextRegisterDouble(stg_context, "solution time", 395ba6664aeSJames Wright offsetof(struct STGShur14Context_, time), 1, 396ba6664aeSJames Wright "Phyiscal time of the solution"); 397ba6664aeSJames Wright 3984e139266SJames Wright if (use_stgstrong) { 3994e139266SJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow_Strong; 4004e139266SJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_Strong_loc; 401*36b31e27SJames Wright problem->bc_from_ics = PETSC_FALSE; 4024e139266SJames Wright } else { 403ba6664aeSJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow; 404ba6664aeSJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc; 405*36b31e27SJames Wright problem->bc_from_ics = PETSC_TRUE; 4064e139266SJames Wright } 407ba6664aeSJames Wright problem->apply_inflow.qfunction_context = stg_context; 408ba6664aeSJames Wright 409ba6664aeSJames Wright PetscFunctionReturn(0); 410ba6664aeSJames Wright } 4114e139266SJames Wright 4124e139266SJames Wright static inline PetscScalar FindDy(const PetscScalar ynodes[], 4134e139266SJames Wright const PetscInt nynodes, const PetscScalar y) { 4144e139266SJames Wright const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]); 4154e139266SJames Wright // ^^assuming min(dy) is first element off the wall 4164e139266SJames Wright PetscInt idx = -1; // Index 4174e139266SJames Wright 4184e139266SJames Wright for (PetscInt i=0; i<nynodes; i++) { 4194e139266SJames Wright if (y < ynodes[i] + half_mindy) { 4204e139266SJames Wright idx = i; break; 4214e139266SJames Wright } 4224e139266SJames Wright } 4234e139266SJames Wright if (idx == 0) return ynodes[1] - ynodes[0]; 4244e139266SJames Wright else if (idx == nynodes-1) return ynodes[nynodes-2] - ynodes[nynodes-1]; 4254e139266SJames Wright else return 0.5 * (ynodes[idx+1] - ynodes[idx-1]); 4264e139266SJames Wright } 4274e139266SJames Wright 4284e139266SJames Wright // Function passed to DMAddBoundary 4294e139266SJames Wright PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, 4304e139266SJames Wright const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) { 4314e139266SJames Wright PetscFunctionBeginUser; 4324e139266SJames Wright 4334e139266SJames Wright const STGShur14Context stg_ctx = (STGShur14Context) ctx; 4344e139266SJames Wright PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt; 4354e139266SJames Wright const bool mean_only = stg_ctx->mean_only; 4364e139266SJames Wright const PetscScalar dx = stg_ctx->dx; 4374e139266SJames Wright const PetscScalar dz = stg_ctx->dz; 4384e139266SJames Wright const PetscScalar mu = stg_ctx->newtonian_ctx.mu; 4394e139266SJames Wright const PetscScalar theta0 = stg_ctx->theta0; 4404e139266SJames Wright const PetscScalar P0 = stg_ctx->P0; 4414e139266SJames Wright const PetscScalar cv = stg_ctx->newtonian_ctx.cv; 4424e139266SJames Wright const PetscScalar cp = stg_ctx->newtonian_ctx.cp; 4434e139266SJames Wright const PetscScalar Rd = cp - cv; 4444e139266SJames Wright 4454e139266SJames Wright const CeedScalar rho = P0 / (Rd * theta0); 4464e139266SJames Wright InterpolateProfile(x[1], ubar, cij, &eps, <, stg_ctx); 4474e139266SJames Wright if (!mean_only) { 4484e139266SJames Wright const PetscInt nynodes = stg_ctx->nynodes; 4494e139266SJames Wright const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes]; 4504e139266SJames Wright const PetscScalar h[3] = {dx, FindDy(ynodes, nynodes, x[1]), dz}; 4514e139266SJames Wright CalcSpectrum(x[1], eps, lt, h, mu/rho, qn, stg_ctx); 4524e139266SJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 4534e139266SJames Wright } else { 4544e139266SJames Wright for (CeedInt j=0; j<3; j++) u[j] = ubar[j]; 4554e139266SJames Wright } 4564e139266SJames Wright 4574e139266SJames Wright bcval[0] = rho; 4584e139266SJames Wright bcval[1] = rho * u[0]; 4594e139266SJames Wright bcval[2] = rho * u[1]; 4604e139266SJames Wright bcval[3] = rho * u[2]; 4614e139266SJames Wright PetscFunctionReturn(0); 4624e139266SJames Wright } 4634e139266SJames Wright 4644e139266SJames Wright PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, 4654e139266SJames Wright STGShur14Context stg_ctx) { 4664e139266SJames Wright 4674e139266SJames Wright PetscErrorCode ierr; 4684e139266SJames Wright DMLabel label; 4694e139266SJames Wright const PetscInt comps[] = {0, 1, 2, 3}; 4704e139266SJames Wright const PetscInt num_comps = 4; 4714e139266SJames Wright PetscFunctionBeginUser; 4724e139266SJames Wright 4734e139266SJames Wright ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr); 4744e139266SJames Wright // Set wall BCs 4754e139266SJames Wright if (bc->num_inflow > 0) { 4764e139266SJames Wright ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, 4774e139266SJames Wright bc->num_inflow, bc->inflows, 0, num_comps, 4784e139266SJames Wright comps, (void(*)(void))StrongSTGbcFunc, 4794e139266SJames Wright NULL, stg_ctx, NULL); CHKERRQ(ierr); 4804e139266SJames Wright } 4814e139266SJames Wright 4824e139266SJames Wright PetscFunctionReturn(0); 4834e139266SJames Wright } 484