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 14e419654dSJeremy L Thompson #include <ceed.h> 15493642f1SJames Wright #include <math.h> 16e419654dSJeremy 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 45*5d28dccaSJames Wright PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP, 46*5d28dccaSJames Wright "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1); 47493642f1SJames Wright } 48493642f1SJames Wright PetscFunctionReturn(0); 49493642f1SJames Wright } 50493642f1SJames Wright 51493642f1SJames Wright /* 52493642f1SJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 53493642f1SJames Wright * 5404e40bb6SJeremy 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. 5504e40bb6SJeremy L Thompson * Assumes there are 14 columns in the file. 56493642f1SJames Wright * 5704e40bb6SJeremy L Thompson * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 58493642f1SJames Wright * 59493642f1SJames Wright * @param[in] comm MPI_Comm for the program 60493642f1SJames Wright * @param[in] path Path to the STGInflow.dat file 6104e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 62493642f1SJames Wright */ 632b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) { 64493642f1SJames Wright PetscInt ndims, dims[2]; 65493642f1SJames Wright FILE *fp; 66493642f1SJames Wright const PetscInt char_array_len = 512; 67493642f1SJames Wright char line[char_array_len]; 68493642f1SJames Wright char **array; 69493642f1SJames Wright 70493642f1SJames Wright PetscFunctionBeginUser; 71493642f1SJames Wright 72676080b4SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 73493642f1SJames Wright 74493642f1SJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 75c77f3192SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 76493642f1SJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 77493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 782b916ea7SJeremy L Thompson CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 79493642f1SJames Wright 80493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 8134b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 8234b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 83*5d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 84*5d28dccaSJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, ndims, 85*5d28dccaSJames Wright dims[1]); 86493642f1SJames Wright 87c77f3192SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 88493642f1SJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 89493642f1SJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 90493642f1SJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 91493642f1SJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 92493642f1SJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 93493642f1SJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 94493642f1SJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 95493642f1SJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 96493642f1SJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 97493642f1SJames Wright lt[i] = (CeedScalar)atof(array[12]); 98493642f1SJames Wright eps[i] = (CeedScalar)atof(array[13]); 99493642f1SJames Wright 100*5d28dccaSJames Wright PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path); 101*5d28dccaSJames Wright PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path); 102*5d28dccaSJames Wright PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path); 103493642f1SJames Wright } 1042b916ea7SJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 10534b5deb1SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 10634b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 1072b916ea7SJeremy L Thompson 108493642f1SJames Wright PetscFunctionReturn(0); 109493642f1SJames Wright } 110493642f1SJames Wright 111493642f1SJames Wright /* 112493642f1SJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 113493642f1SJames Wright * 11404e40bb6SJeremy 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. 11504e40bb6SJeremy L Thompson * Assumes there are 7 columns in the file. 116493642f1SJames Wright * 117493642f1SJames Wright * @param[in] comm MPI_Comm for the program 118493642f1SJames Wright * @param[in] path Path to the STGRand.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 ReadSTGRand(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; 129676080b4SJames Wright PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp)); 130493642f1SJames Wright 131493642f1SJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1322b916ea7SJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1332b916ea7SJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 134493642f1SJames Wright 135493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 13634b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 13734b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 138*5d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 139*5d28dccaSJames Wright "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path, ndims, 140*5d28dccaSJames Wright dims[1]); 141493642f1SJames Wright 142493642f1SJames Wright d[0][i] = (CeedScalar)atof(array[0]); 143493642f1SJames Wright d[1][i] = (CeedScalar)atof(array[1]); 144493642f1SJames Wright d[2][i] = (CeedScalar)atof(array[2]); 145493642f1SJames Wright phi[i] = (CeedScalar)atof(array[3]); 146493642f1SJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 147493642f1SJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 148493642f1SJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 149493642f1SJames Wright } 15034b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 1512b916ea7SJeremy L Thompson 152493642f1SJames Wright PetscFunctionReturn(0); 153493642f1SJames Wright } 154493642f1SJames Wright 155493642f1SJames Wright /* 156493642f1SJames Wright * @brief Read STG data from input paths and put in STGShur14Context 157493642f1SJames Wright * 158493642f1SJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 15904e40bb6SJeremy L Thompson * Data stored initially in `*pstg_ctx` will be copied over to the new STGShur14Context instance. 160493642f1SJames Wright * 161493642f1SJames Wright * @param[in] comm MPI_Comm for the program 162493642f1SJames Wright * @param[in] dm DM for the program 163493642f1SJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 164493642f1SJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 16504e40bb6SJeremy L Thompson * @param[in,out] pstg_ctx Pointer to STGShur14Context where the data will be loaded into 166493642f1SJames Wright */ 1672b916ea7SJeremy 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], 1682b916ea7SJeremy L Thompson STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) { 169493642f1SJames Wright PetscInt nmodes, nprofs; 170493642f1SJames Wright STGShur14Context stg_ctx; 171493642f1SJames Wright PetscFunctionBeginUser; 172493642f1SJames Wright 173493642f1SJames Wright // Get options 174676080b4SJames Wright PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes)); 175676080b4SJames Wright PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs)); 176*5d28dccaSJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 177*5d28dccaSJames Wright "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT "). Change size of STG_NMODES_MAX and recompile", 1782b916ea7SJeremy L Thompson stg_rand_path, nmodes, STG_NMODES_MAX); 179493642f1SJames Wright 180493642f1SJames Wright { 181493642f1SJames Wright STGShur14Context s; 18234b5deb1SJames Wright PetscCall(PetscCalloc1(1, &s)); 183493642f1SJames Wright *s = **pstg_ctx; 184493642f1SJames Wright s->nmodes = nmodes; 185493642f1SJames Wright s->nprofs = nprofs; 186493642f1SJames Wright s->offsets.sigma = 0; 187493642f1SJames Wright s->offsets.d = nmodes * 3; 188493642f1SJames Wright s->offsets.phi = s->offsets.d + nmodes * 3; 189493642f1SJames Wright s->offsets.kappa = s->offsets.phi + nmodes; 190c77f3192SJames Wright s->offsets.wall_dist = s->offsets.kappa + nmodes; 191c77f3192SJames Wright s->offsets.ubar = s->offsets.wall_dist + nprofs; 192493642f1SJames Wright s->offsets.cij = s->offsets.ubar + nprofs * 3; 193493642f1SJames Wright s->offsets.eps = s->offsets.cij + nprofs * 6; 194493642f1SJames Wright s->offsets.lt = s->offsets.eps + nprofs; 195e6098bcdSJames Wright s->offsets.ynodes = s->offsets.lt + nprofs; 196e6098bcdSJames Wright PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes; 197493642f1SJames Wright s->total_bytes = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]); 19834b5deb1SJames Wright PetscCall(PetscMalloc(s->total_bytes, &stg_ctx)); 199493642f1SJames Wright *stg_ctx = *s; 20034b5deb1SJames Wright PetscCall(PetscFree(s)); 201493642f1SJames Wright } 202493642f1SJames Wright 20334b5deb1SJames Wright PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx)); 20434b5deb1SJames Wright PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx)); 205493642f1SJames Wright 206e6098bcdSJames Wright if (stg_ctx->nynodes > 0) { 207e6098bcdSJames Wright CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes]; 208e6098bcdSJames Wright for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i]; 209e6098bcdSJames Wright } 210e6098bcdSJames Wright 211493642f1SJames Wright // -- Calculate kappa 212493642f1SJames Wright { 213493642f1SJames Wright CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa]; 214c77f3192SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 215493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 216493642f1SJames Wright CeedScalar le, le_max = 0; 217493642f1SJames Wright 2182b916ea7SJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 219c77f3192SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 220493642f1SJames Wright if (le_max < le) le_max = le; 221493642f1SJames Wright } 222493642f1SJames Wright CeedScalar kmin = M_PI / le_max; 223493642f1SJames Wright 2242b916ea7SJeremy L Thompson CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); } 225493642f1SJames Wright } // end calculate kappa 226493642f1SJames Wright 22734b5deb1SJames Wright PetscCall(PetscFree(*pstg_ctx)); 228493642f1SJames Wright *pstg_ctx = stg_ctx; 229493642f1SJames Wright PetscFunctionReturn(0); 230493642f1SJames Wright } 231493642f1SJames Wright 2322b916ea7SJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0, 2332b916ea7SJeremy L Thompson const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) { 234493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 235493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 2362b916ea7SJeremy L Thompson PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE; 2372b916ea7SJeremy L Thompson CeedScalar u0 = 0.0, alpha = 1.01; 238493642f1SJames Wright CeedQFunctionContext stg_context; 239493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 240493642f1SJames Wright PetscFunctionBeginUser; 241493642f1SJames Wright 242493642f1SJames Wright // Get options 243493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2442b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2452b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2462b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2472b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2482b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2492b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2502b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2512b916ea7SJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 252493642f1SJames Wright PetscOptionsEnd(); 253493642f1SJames Wright 25434b5deb1SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 2558085925cSJames Wright global_stg_ctx->alpha = alpha; 2568085925cSJames Wright global_stg_ctx->u0 = u0; 2578085925cSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 2588085925cSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 2598085925cSJames Wright global_stg_ctx->mean_only = mean_only; 260d4e0f297SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 2618085925cSJames Wright global_stg_ctx->theta0 = theta0; 2628085925cSJames Wright global_stg_ctx->P0 = P0; 2638085925cSJames Wright global_stg_ctx->nynodes = nynodes; 264493642f1SJames Wright 265493642f1SJames Wright { 266493642f1SJames Wright // Calculate dx assuming constant spacing 267493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 26834b5deb1SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 269493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 270493642f1SJames Wright 271493642f1SJames Wright PetscInt nmax = 3, faces[3]; 2722b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 2738085925cSJames Wright global_stg_ctx->dx = domain_size[0] / faces[0]; 2748085925cSJames Wright global_stg_ctx->dz = domain_size[2] / faces[2]; 275493642f1SJames Wright } 276493642f1SJames Wright 2772b916ea7SJeremy L Thompson CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 2788085925cSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 2792b916ea7SJeremy L Thompson CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 280493642f1SJames Wright 2812b916ea7SJeremy L Thompson PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes)); 282493642f1SJames Wright 283493642f1SJames Wright CeedQFunctionContextCreate(user->ceed, &stg_context); 2842b916ea7SJeremy L Thompson CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx); 2852b916ea7SJeremy L Thompson CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc); 2862b916ea7SJeremy L Thompson CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution"); 287493642f1SJames Wright 28843bff553SJames Wright CeedQFunctionContextDestroy(&problem->ics.qfunction_context); 28943bff553SJames Wright problem->ics.qfunction = ICsSTG; 29043bff553SJames Wright problem->ics.qfunction_loc = ICsSTG_loc; 29143bff553SJames Wright problem->ics.qfunction_context = stg_context; 29243bff553SJames Wright 293e6098bcdSJames Wright if (use_stgstrong) { 2948085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 2956d0190e2SJames Wright problem->use_dirichlet_ceed = PETSC_TRUE; 296d7244455SJames Wright problem->bc_from_ics = PETSC_FALSE; 297e6098bcdSJames Wright } else { 298493642f1SJames Wright problem->apply_inflow.qfunction = STGShur14_Inflow; 299493642f1SJames Wright problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc; 300a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction = STGShur14_Inflow_Jacobian; 301a6e8f989SJames Wright problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc; 3022b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context); 3032b916ea7SJeremy L Thompson CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context); 304d7244455SJames Wright problem->bc_from_ics = PETSC_TRUE; 305e6098bcdSJames Wright } 306493642f1SJames Wright 307493642f1SJames Wright PetscFunctionReturn(0); 308493642f1SJames Wright } 309e6098bcdSJames Wright 3102b916ea7SJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) { 311e6098bcdSJames Wright const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]); 312e6098bcdSJames Wright // ^^assuming min(dy) is first element off the wall 313e6098bcdSJames Wright PetscInt idx = -1; // Index 314e6098bcdSJames Wright 315e6098bcdSJames Wright for (PetscInt i = 0; i < nynodes; i++) { 316e6098bcdSJames Wright if (y < ynodes[i] + half_mindy) { 3172b916ea7SJeremy L Thompson idx = i; 3182b916ea7SJeremy L Thompson break; 319e6098bcdSJames Wright } 320e6098bcdSJames Wright } 321e6098bcdSJames Wright if (idx == 0) return ynodes[1] - ynodes[0]; 322e6098bcdSJames Wright else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1]; 323e6098bcdSJames Wright else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]); 324e6098bcdSJames Wright } 325e6098bcdSJames Wright 326e6098bcdSJames Wright // Function passed to DMAddBoundary 3276d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method 3282b916ea7SJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) { 329e6098bcdSJames Wright PetscFunctionBeginUser; 3302b916ea7SJeremy L Thompson 331e6098bcdSJames Wright const STGShur14Context stg_ctx = (STGShur14Context)ctx; 332e6098bcdSJames Wright PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt; 333e6098bcdSJames Wright const bool mean_only = stg_ctx->mean_only; 334e6098bcdSJames Wright const PetscScalar dx = stg_ctx->dx; 335e6098bcdSJames Wright const PetscScalar dz = stg_ctx->dz; 336e6098bcdSJames Wright const PetscScalar mu = stg_ctx->newtonian_ctx.mu; 337e6098bcdSJames Wright const PetscScalar theta0 = stg_ctx->theta0; 338e6098bcdSJames Wright const PetscScalar P0 = stg_ctx->P0; 339e6098bcdSJames Wright const PetscScalar cv = stg_ctx->newtonian_ctx.cv; 340e6098bcdSJames Wright const PetscScalar cp = stg_ctx->newtonian_ctx.cp; 341e6098bcdSJames Wright const PetscScalar Rd = cp - cv; 342e6098bcdSJames Wright 343e6098bcdSJames Wright const CeedScalar rho = P0 / (Rd * theta0); 344e6098bcdSJames Wright InterpolateProfile(x[1], ubar, cij, &eps, <, stg_ctx); 345e6098bcdSJames Wright if (!mean_only) { 346e6098bcdSJames Wright const PetscInt nynodes = stg_ctx->nynodes; 347e6098bcdSJames Wright const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes]; 348e6098bcdSJames Wright const PetscScalar h[3] = {dx, FindDy(ynodes, nynodes, x[1]), dz}; 349e6098bcdSJames Wright CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx); 350e6098bcdSJames Wright STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx); 351e6098bcdSJames Wright } else { 352e6098bcdSJames Wright for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j]; 353e6098bcdSJames Wright } 354e6098bcdSJames Wright 355e6098bcdSJames Wright bcval[0] = rho; 356e6098bcdSJames Wright bcval[1] = rho * u[0]; 357e6098bcdSJames Wright bcval[2] = rho * u[1]; 358e6098bcdSJames Wright bcval[3] = rho * u[2]; 359e6098bcdSJames Wright PetscFunctionReturn(0); 360e6098bcdSJames Wright } 361e6098bcdSJames Wright 3622b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) { 363e6098bcdSJames Wright DMLabel label; 364e6098bcdSJames Wright PetscFunctionBeginUser; 365e6098bcdSJames Wright 366d374fb47SJames Wright PetscInt comps[5], num_comps = 4; 3673636f6a4SJames Wright switch (phys->state_var) { 3683636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 369d374fb47SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 370d374fb47SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 3713636f6a4SJames Wright break; 3723636f6a4SJames Wright 3733636f6a4SJames Wright case STATEVAR_PRIMITIVE: 3743636f6a4SJames Wright // {1,2,3,4} for u, v, w, T 3753636f6a4SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3763636f6a4SJames Wright break; 377d374fb47SJames Wright } 378d374fb47SJames Wright 37934b5deb1SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 380e6098bcdSJames Wright // Set wall BCs 381e6098bcdSJames Wright if (bc->num_inflow > 0) { 3822b916ea7SJeremy L Thompson PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc, 38334b5deb1SJames Wright NULL, global_stg_ctx, NULL)); 384e6098bcdSJames Wright } 385e6098bcdSJames Wright 386e6098bcdSJames Wright PetscFunctionReturn(0); 387e6098bcdSJames Wright } 3886d0190e2SJames Wright 3892b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, 3909eeef72bSJames Wright CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) { 3916d0190e2SJames Wright CeedQFunction qf_strongbc; 3926d0190e2SJames Wright PetscFunctionBeginUser; 3932b916ea7SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc); 3942b916ea7SJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 3956d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 3966d0190e2SJames Wright CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE); 3979eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 3986d0190e2SJames Wright CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE); 3996d0190e2SJames Wright 4006d0190e2SJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 4016d0190e2SJames Wright *pqf_strongbc = qf_strongbc; 4026d0190e2SJames Wright PetscFunctionReturn(0); 4036d0190e2SJames Wright } 4049eeef72bSJames Wright 4052b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur, 4069eeef72bSJames Wright CeedQFunction *pqf_strongbc) { 4079eeef72bSJames Wright CeedQFunction qf_strongbc; 4089eeef72bSJames Wright PetscFunctionBeginUser; 4092b916ea7SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc); 4102b916ea7SJeremy L Thompson CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE); 4119eeef72bSJames Wright CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE); 4129eeef72bSJames Wright CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE); 4139eeef72bSJames Wright 4149eeef72bSJames Wright CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context); 4159eeef72bSJames Wright *pqf_strongbc = qf_strongbc; 4169eeef72bSJames Wright PetscFunctionReturn(0); 4179eeef72bSJames Wright } 418