1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 /// @file 9 /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 10 /// presented in Shur et al. 2014 11 12 #include "stg_shur14.h" 13 14 #include <ceed.h> 15 #include <math.h> 16 #include <petscdm.h> 17 #include <stdlib.h> 18 19 #include "../navierstokes.h" 20 #include "../qfunctions/stg_shur14.h" 21 22 StgShur14Context global_stg_ctx; 23 24 /* 25 * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 26 * 27 * This assumes the input matrices are in order [11,22,33,12,13,23]. 28 * This format is also used for the output. 29 * 30 * @param[in] comm MPI_Comm 31 * @param[in] nprofs Number of matrices in Rij 32 * @param[in] Rij Array of the symmetric matrices [6,nprofs] 33 * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 34 */ 35 PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 36 PetscFunctionBeginUser; 37 for (PetscInt i = 0; i < nprofs; i++) { 38 Cij[0][i] = sqrt(Rij[0][i]); 39 Cij[3][i] = Rij[3][i] / Cij[0][i]; 40 Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2)); 41 Cij[4][i] = Rij[4][i] / Cij[0][i]; 42 Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i]; 43 Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 44 45 PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP, 46 "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1); 47 } 48 PetscFunctionReturn(PETSC_SUCCESS); 49 } 50 51 /* 52 * @brief Read the STGInflow file and load the contents into stg_ctx 53 * 54 * 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. 55 * Assumes there are 14 columns in the file. 56 * 57 * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 58 * 59 * @param[in] comm MPI_Comm for the program 60 * @param[in] path Path to the STGInflow.dat file 61 * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 62 */ 63 static PetscErrorCode ReadStgInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 64 PetscInt dims[2]; 65 int ndims; 66 FILE *fp; 67 const PetscInt char_array_len = 512; 68 char line[char_array_len]; 69 char **array; 70 71 PetscFunctionBeginUser; 72 PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 73 74 CeedScalar rij[6][stg_ctx->nprofs]; 75 CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 76 CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 77 CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 78 CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 79 80 for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 81 PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 82 PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 83 PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 84 "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 85 86 wall_dist[i] = (CeedScalar)atof(array[0]); 87 ubar[0][i] = (CeedScalar)atof(array[1]); 88 ubar[1][i] = (CeedScalar)atof(array[2]); 89 ubar[2][i] = (CeedScalar)atof(array[3]); 90 rij[0][i] = (CeedScalar)atof(array[4]); 91 rij[1][i] = (CeedScalar)atof(array[5]); 92 rij[2][i] = (CeedScalar)atof(array[6]); 93 rij[3][i] = (CeedScalar)atof(array[7]); 94 rij[4][i] = (CeedScalar)atof(array[8]); 95 rij[5][i] = (CeedScalar)atof(array[9]); 96 lt[i] = (CeedScalar)atof(array[12]); 97 eps[i] = (CeedScalar)atof(array[13]); 98 99 PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path); 100 PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path); 101 PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path); 102 } 103 CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 104 PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 105 PetscCall(PetscFClose(comm, fp)); 106 PetscFunctionReturn(PETSC_SUCCESS); 107 } 108 109 /* 110 * @brief Read the STGRand file and load the contents into stg_ctx 111 * 112 * 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. 113 * Assumes there are 7 columns in the file. 114 * 115 * @param[in] comm MPI_Comm for the program 116 * @param[in] path Path to the STGRand.dat file 117 * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 118 */ 119 static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 120 PetscInt dims[2]; 121 int ndims; 122 FILE *fp; 123 const PetscInt char_array_len = 512; 124 char line[char_array_len]; 125 char **array; 126 127 PetscFunctionBeginUser; 128 PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 129 130 CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 131 CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 132 CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 133 134 for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 135 PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 136 PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 137 PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 138 "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 139 140 d[0][i] = (CeedScalar)atof(array[0]); 141 d[1][i] = (CeedScalar)atof(array[1]); 142 d[2][i] = (CeedScalar)atof(array[2]); 143 phi[i] = (CeedScalar)atof(array[3]); 144 sigma[0][i] = (CeedScalar)atof(array[4]); 145 sigma[1][i] = (CeedScalar)atof(array[5]); 146 sigma[2][i] = (CeedScalar)atof(array[6]); 147 } 148 PetscCall(PetscFClose(comm, fp)); 149 PetscFunctionReturn(PETSC_SUCCESS); 150 } 151 152 /* 153 * @brief Read STG data from input paths and put in STGShur14Context 154 * 155 * Reads data from input paths and puts them into a STGShur14Context object. 156 * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance. 157 * 158 * @param[in] comm MPI_Comm for the program 159 * @param[in] dm DM for the program 160 * @param[in] stg_inflow_path Path to STGInflow.dat file 161 * @param[in] stg_rand_path Path to STGRand.dat file 162 * @param[in,out] stg_ctx Pointer to STGShur14Context where the data will be loaded into 163 */ 164 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], 165 StgShur14Context *stg_ctx) { 166 PetscInt nmodes, nprofs; 167 168 PetscFunctionBeginUser; 169 PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 170 PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 171 PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 172 "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path, 173 nmodes, STG_NMODES_MAX); 174 175 { 176 StgShur14Context temp_ctx; 177 PetscCall(PetscCalloc1(1, &temp_ctx)); 178 *temp_ctx = **stg_ctx; 179 temp_ctx->nmodes = nmodes; 180 temp_ctx->nprofs = nprofs; 181 temp_ctx->offsets.sigma = 0; 182 temp_ctx->offsets.d = nmodes * 3; 183 temp_ctx->offsets.phi = temp_ctx->offsets.d + nmodes * 3; 184 temp_ctx->offsets.kappa = temp_ctx->offsets.phi + nmodes; 185 temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes; 186 temp_ctx->offsets.ubar = temp_ctx->offsets.wall_dist + nprofs; 187 temp_ctx->offsets.cij = temp_ctx->offsets.ubar + nprofs * 3; 188 temp_ctx->offsets.eps = temp_ctx->offsets.cij + nprofs * 6; 189 temp_ctx->offsets.lt = temp_ctx->offsets.eps + nprofs; 190 PetscInt total_num_scalars = temp_ctx->offsets.lt + nprofs; 191 temp_ctx->total_bytes = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]); 192 PetscCall(PetscFree(*stg_ctx)); 193 PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx)); 194 **stg_ctx = *temp_ctx; 195 PetscCall(PetscFree(temp_ctx)); 196 } 197 198 PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx)); 199 PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 200 201 { // -- Calculate kappa 202 CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 203 CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 204 CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 205 CeedScalar le, le_max = 0; 206 207 CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 208 le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 209 if (le_max < le) le_max = le; 210 } 211 CeedScalar kmin = M_PI / le_max; 212 213 CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); } 214 } 215 PetscFunctionReturn(PETSC_SUCCESS); 216 } 217 218 PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0, 219 const CeedScalar P0) { 220 Ceed ceed = user->ceed; 221 char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 222 char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 223 PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 224 CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = 1.0e-3; 225 CeedQFunctionContext stg_context; 226 NewtonianIdealGasContext newtonian_ig_ctx; 227 228 PetscFunctionBeginUser; 229 PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 230 PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 231 PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 232 PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 233 PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 234 PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 235 PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 236 PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 237 use_fluctuating_IC, &use_fluctuating_IC, NULL)); 238 PetscCall(PetscOptionsReal("-stg_dx", "Element size in streamwise direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 239 PetscOptionsEnd(); 240 241 PetscCall(PetscCalloc1(1, &global_stg_ctx)); 242 global_stg_ctx->alpha = alpha; 243 global_stg_ctx->u0 = u0; 244 global_stg_ctx->is_implicit = user->phys->implicit; 245 global_stg_ctx->prescribe_T = prescribe_T; 246 global_stg_ctx->mean_only = mean_only; 247 global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 248 global_stg_ctx->theta0 = theta0; 249 global_stg_ctx->P0 = P0; 250 251 { // Calculate dx assuming constant spacing 252 PetscReal domain_min[3], domain_max[3], domain_size[3]; 253 PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 254 for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 255 256 PetscInt nmax = 3, faces[3]; 257 PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 258 global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 259 } 260 261 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 262 global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 263 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 264 265 PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 266 267 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context)); 268 PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 269 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc)); 270 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, 271 "Physical time of the solution")); 272 273 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 274 problem->ics.qfunction = ICsStg; 275 problem->ics.qfunction_loc = ICsStg_loc; 276 problem->ics.qfunction_context = stg_context; 277 278 if (use_stgstrong) { 279 // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 280 problem->use_strong_bc_ceed = PETSC_TRUE; 281 problem->bc_from_ics = PETSC_FALSE; 282 } else { 283 problem->apply_inflow.qfunction = StgShur14Inflow; 284 problem->apply_inflow.qfunction_loc = StgShur14Inflow_loc; 285 problem->apply_inflow_jacobian.qfunction = StgShur14Inflow_Jacobian; 286 problem->apply_inflow_jacobian.qfunction_loc = StgShur14Inflow_Jacobian_loc; 287 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context)); 288 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context)); 289 problem->bc_from_ics = PETSC_TRUE; 290 } 291 PetscFunctionReturn(PETSC_SUCCESS); 292 } 293 294 // @brief Set STG strongly enforce components using DMAddBoundary 295 PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) { 296 DMLabel label; 297 PetscInt comps[5], num_comps = 4; 298 299 PetscFunctionBeginUser; 300 switch (phys->state_var) { 301 case STATEVAR_CONSERVATIVE: 302 // {0,1,2,3} for rho, rho*u, rho*v, rho*w 303 for (int i = 0; i < 4; i++) comps[i] = i; 304 break; 305 306 case STATEVAR_PRIMITIVE: 307 // {1,2,3,4} for u, v, w, T 308 for (int i = 0; i < 4; i++) comps[i] = i + 1; 309 break; 310 } 311 312 PetscCall(DMGetLabel(dm, "Face Sets", &label)); 313 if (bc->num_inflow > 0) { 314 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL)); 315 } 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318 319 PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 320 CeedQFunction *qf_strongbc) { 321 PetscFunctionBeginUser; 322 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 323 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 324 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 325 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 326 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 327 PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 328 329 PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 330 PetscFunctionReturn(PETSC_SUCCESS); 331 } 332 333 PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 334 CeedQFunction *qf_strongbc) { 335 PetscFunctionBeginUser; 336 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 337 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 338 PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 339 PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 340 341 PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344