15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 122b730f8bSJeremy L Thompson #include "stg_shur14.h" 132b730f8bSJeremy L Thompson 1449aac155SJeremy L Thompson #include <ceed.h> 15ba6664aeSJames Wright #include <math.h> 1649aac155SJeremy L Thompson #include <petscdm.h> 172b730f8bSJeremy L Thompson #include <stdlib.h> 182b730f8bSJeremy L Thompson 19ba6664aeSJames Wright #include "../navierstokes.h" 20ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h" 21ba6664aeSJames Wright 22cbef7084SJames Wright StgShur14Context global_stg_ctx; 2365dd5cafSJames Wright 24ba6664aeSJames Wright /* 25ba6664aeSJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 26ba6664aeSJames Wright * 27ea61e9acSJeremy L Thompson * This assumes the input matrices are in order [11,22,33,12,13,23]. 28ea61e9acSJeremy L Thompson * This format is also used for the output. 29ba6664aeSJames Wright * 30ba6664aeSJames Wright * @param[in] comm MPI_Comm 31ba6664aeSJames Wright * @param[in] nprofs Number of matrices in Rij 32ba6664aeSJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 33ba6664aeSJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 34ba6664aeSJames Wright */ 352b730f8bSJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 36ba6664aeSJames Wright PetscFunctionBeginUser; 37ba6664aeSJames Wright for (PetscInt i = 0; i < nprofs; i++) { 38ba6664aeSJames Wright Cij[0][i] = sqrt(Rij[0][i]); 39ba6664aeSJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 40ba6664aeSJames Wright Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2)); 41ba6664aeSJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 42ba6664aeSJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i]; 43ba6664aeSJames Wright Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2)); 44ba6664aeSJames Wright 450e654f56SJames Wright PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP, 460e654f56SJames Wright "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1); 47ba6664aeSJames Wright } 48ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 49ba6664aeSJames Wright } 50ba6664aeSJames Wright 51ba6664aeSJames Wright /* 52ba6664aeSJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 53ba6664aeSJames Wright * 54ea61e9acSJeremy 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. 55ea61e9acSJeremy L Thompson * Assumes there are 14 columns in the file. 56ba6664aeSJames Wright * 57ea61e9acSJeremy L Thompson * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 58ba6664aeSJames Wright * 59ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 60ba6664aeSJames Wright * @param[in] path Path to the STGInflow.dat file 61ea61e9acSJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 62ba6664aeSJames Wright */ 63cbef7084SJames Wright static PetscErrorCode ReadStgInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 64457e73b2SJames Wright PetscInt dims[2]; 65457e73b2SJames Wright int ndims; 66ba6664aeSJames Wright FILE *fp; 67ba6664aeSJames Wright const PetscInt char_array_len = 512; 68ba6664aeSJames Wright char line[char_array_len]; 69ba6664aeSJames Wright char **array; 70ba6664aeSJames Wright 71ba6664aeSJames Wright PetscFunctionBeginUser; 72cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 73ba6664aeSJames Wright 74ba6664aeSJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 75175f00a6SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 76ba6664aeSJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 77ba6664aeSJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 782b730f8bSJeremy L Thompson CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 79ba6664aeSJames Wright 80ba6664aeSJames Wright for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 8124941a69SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 8224941a69SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 830e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 84457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 85ba6664aeSJames Wright 86175f00a6SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 87ba6664aeSJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 88ba6664aeSJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 89ba6664aeSJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 90ba6664aeSJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 91ba6664aeSJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 92ba6664aeSJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 93ba6664aeSJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 94ba6664aeSJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 95ba6664aeSJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 96ba6664aeSJames Wright lt[i] = (CeedScalar)atof(array[12]); 97ba6664aeSJames Wright eps[i] = (CeedScalar)atof(array[13]); 98ba6664aeSJames Wright 990e654f56SJames Wright PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path); 1000e654f56SJames Wright PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path); 1010e654f56SJames Wright PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path); 102ba6664aeSJames Wright } 1032b730f8bSJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 10424941a69SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 10524941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 106ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 107ba6664aeSJames Wright } 108ba6664aeSJames Wright 109ba6664aeSJames Wright /* 110ba6664aeSJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 111ba6664aeSJames Wright * 112ea61e9acSJeremy 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. 113ea61e9acSJeremy L Thompson * Assumes there are 7 columns in the file. 114ba6664aeSJames Wright * 115ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 116ba6664aeSJames Wright * @param[in] path Path to the STGRand.dat file 117ea61e9acSJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 118ba6664aeSJames Wright */ 119cbef7084SJames Wright static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 120457e73b2SJames Wright PetscInt dims[2]; 121457e73b2SJames Wright int ndims; 122ba6664aeSJames Wright FILE *fp; 123ba6664aeSJames Wright const PetscInt char_array_len = 512; 124ba6664aeSJames Wright char line[char_array_len]; 125ba6664aeSJames Wright char **array; 126ba6664aeSJames Wright 127ba6664aeSJames Wright PetscFunctionBeginUser; 128cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 129ba6664aeSJames Wright 130ba6664aeSJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1312b730f8bSJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1322b730f8bSJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 133ba6664aeSJames Wright 134ba6664aeSJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 13524941a69SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 13624941a69SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1370e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 138457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 139ba6664aeSJames Wright 140ba6664aeSJames Wright d[0][i] = (CeedScalar)atof(array[0]); 141ba6664aeSJames Wright d[1][i] = (CeedScalar)atof(array[1]); 142ba6664aeSJames Wright d[2][i] = (CeedScalar)atof(array[2]); 143ba6664aeSJames Wright phi[i] = (CeedScalar)atof(array[3]); 144ba6664aeSJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 145ba6664aeSJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 146ba6664aeSJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 147ba6664aeSJames Wright } 14824941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 149ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 150ba6664aeSJames Wright } 151ba6664aeSJames Wright 152ba6664aeSJames Wright /* 153ba6664aeSJames Wright * @brief Read STG data from input paths and put in STGShur14Context 154ba6664aeSJames Wright * 155ba6664aeSJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 15605675bc2SJames Wright * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance. 157ba6664aeSJames Wright * 158ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 159ba6664aeSJames Wright * @param[in] dm DM for the program 160ba6664aeSJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 161ba6664aeSJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 16205675bc2SJames Wright * @param[in,out] stg_ctx Pointer to STGShur14Context where the data will be loaded into 163ba6664aeSJames Wright */ 164cbef7084SJames Wright 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], 165cbef7084SJames Wright StgShur14Context *stg_ctx) { 166ba6664aeSJames Wright PetscInt nmodes, nprofs; 167ba6664aeSJames Wright 168f17d818dSJames Wright PetscFunctionBeginUser; 169cbef7084SJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 170cbef7084SJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 1710e654f56SJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 172457e73b2SJames Wright "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path, 173457e73b2SJames Wright nmodes, STG_NMODES_MAX); 174ba6664aeSJames Wright 175ba6664aeSJames Wright { 176cbef7084SJames Wright StgShur14Context temp_ctx; 17705675bc2SJames Wright PetscCall(PetscCalloc1(1, &temp_ctx)); 17805675bc2SJames Wright *temp_ctx = **stg_ctx; 17905675bc2SJames Wright temp_ctx->nmodes = nmodes; 18005675bc2SJames Wright temp_ctx->nprofs = nprofs; 18105675bc2SJames Wright temp_ctx->offsets.sigma = 0; 18205675bc2SJames Wright temp_ctx->offsets.d = nmodes * 3; 18305675bc2SJames Wright temp_ctx->offsets.phi = temp_ctx->offsets.d + nmodes * 3; 18405675bc2SJames Wright temp_ctx->offsets.kappa = temp_ctx->offsets.phi + nmodes; 18505675bc2SJames Wright temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes; 18605675bc2SJames Wright temp_ctx->offsets.ubar = temp_ctx->offsets.wall_dist + nprofs; 18705675bc2SJames Wright temp_ctx->offsets.cij = temp_ctx->offsets.ubar + nprofs * 3; 18805675bc2SJames Wright temp_ctx->offsets.eps = temp_ctx->offsets.cij + nprofs * 6; 18905675bc2SJames Wright temp_ctx->offsets.lt = temp_ctx->offsets.eps + nprofs; 190f8839eb4SJames Wright PetscInt total_num_scalars = temp_ctx->offsets.lt + nprofs; 19105675bc2SJames Wright temp_ctx->total_bytes = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]); 19205675bc2SJames Wright PetscCall(PetscFree(*stg_ctx)); 19305675bc2SJames Wright PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx)); 19405675bc2SJames Wright **stg_ctx = *temp_ctx; 19505675bc2SJames Wright PetscCall(PetscFree(temp_ctx)); 196ba6664aeSJames Wright } 197ba6664aeSJames Wright 198cbef7084SJames Wright PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx)); 199cbef7084SJames Wright PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 200ba6664aeSJames Wright 20105675bc2SJames Wright { // -- Calculate kappa 20205675bc2SJames Wright CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 20305675bc2SJames Wright CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 20405675bc2SJames Wright CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 205ba6664aeSJames Wright CeedScalar le, le_max = 0; 206ba6664aeSJames Wright 20705675bc2SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 208175f00a6SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 209ba6664aeSJames Wright if (le_max < le) le_max = le; 210ba6664aeSJames Wright } 211ba6664aeSJames Wright CeedScalar kmin = M_PI / le_max; 212ba6664aeSJames Wright 21305675bc2SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); } 21405675bc2SJames Wright } 215ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 216ba6664aeSJames Wright } 217ba6664aeSJames Wright 218731c13d7SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, User user, const bool prescribe_T, const CeedScalar theta0, 219f8839eb4SJames Wright const CeedScalar P0) { 220a424bcd0SJames Wright Ceed ceed = user->ceed; 221ba6664aeSJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 222ba6664aeSJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 22373557d35SJames Wright PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 22473557d35SJames Wright CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = 1.0e-3; 225ba6664aeSJames Wright CeedQFunctionContext stg_context; 226ba6664aeSJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 227ba6664aeSJames Wright 228f17d818dSJames Wright PetscFunctionBeginUser; 229ba6664aeSJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2302b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2312b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2322b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2332b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2342b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2352b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2362b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2372b730f8bSJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 23873557d35SJames Wright PetscCall(PetscOptionsReal("-stg_dx", "Element size in streamwise direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 239ba6664aeSJames Wright PetscOptionsEnd(); 240ba6664aeSJames Wright 24124941a69SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 24265dd5cafSJames Wright global_stg_ctx->alpha = alpha; 24365dd5cafSJames Wright global_stg_ctx->u0 = u0; 24465dd5cafSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 24565dd5cafSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 24665dd5cafSJames Wright global_stg_ctx->mean_only = mean_only; 24789060322SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 24865dd5cafSJames Wright global_stg_ctx->theta0 = theta0; 24965dd5cafSJames Wright global_stg_ctx->P0 = P0; 250ba6664aeSJames Wright 251f17d818dSJames Wright { // Calculate dx assuming constant spacing 252ba6664aeSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 25324941a69SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 254ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 255ba6664aeSJames Wright 256ba6664aeSJames Wright PetscInt nmax = 3, faces[3]; 2572b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 25873557d35SJames Wright global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 259ba6664aeSJames Wright } 260ba6664aeSJames Wright 261a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 26265dd5cafSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 263a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 264ba6664aeSJames Wright 265cbef7084SJames Wright PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 266ba6664aeSJames Wright 267a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context)); 268a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 269a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc)); 270a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, 271a424bcd0SJames Wright "Physical time of the solution")); 272ba6664aeSJames Wright 273a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 274cbef7084SJames Wright problem->ics.qfunction = ICsStg; 275cbef7084SJames Wright problem->ics.qfunction_loc = ICsStg_loc; 276b77c53c9SJames Wright problem->ics.qfunction_context = stg_context; 277b77c53c9SJames Wright 2784e139266SJames Wright if (use_stgstrong) { 27965dd5cafSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 2803722cd23SJames Wright problem->use_strong_bc_ceed = PETSC_TRUE; 28136b31e27SJames Wright problem->bc_from_ics = PETSC_FALSE; 2824e139266SJames Wright } else { 283cbef7084SJames Wright problem->apply_inflow.qfunction = StgShur14Inflow; 284cbef7084SJames Wright problem->apply_inflow.qfunction_loc = StgShur14Inflow_loc; 285cbef7084SJames Wright problem->apply_inflow_jacobian.qfunction = StgShur14Inflow_Jacobian; 286cbef7084SJames Wright problem->apply_inflow_jacobian.qfunction_loc = StgShur14Inflow_Jacobian_loc; 287a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context)); 288a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context)); 28936b31e27SJames Wright problem->bc_from_ics = PETSC_TRUE; 2904e139266SJames Wright } 291ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 292ba6664aeSJames Wright } 2934e139266SJames Wright 294f8839eb4SJames Wright // @brief Set STG strongly enforce components using DMAddBoundary 295731c13d7SJames Wright PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData problem, Physics phys) { 2964e139266SJames Wright DMLabel label; 297192a7459SJames Wright PetscInt comps[5], num_comps = 4; 298f17d818dSJames Wright 299f17d818dSJames Wright PetscFunctionBeginUser; 30097baf651SJames Wright switch (phys->state_var) { 30197baf651SJames Wright case STATEVAR_CONSERVATIVE: 302192a7459SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 303192a7459SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 30497baf651SJames Wright break; 30597baf651SJames Wright 30697baf651SJames Wright case STATEVAR_PRIMITIVE: 30797baf651SJames Wright // {1,2,3,4} for u, v, w, T 30897baf651SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 30997baf651SJames Wright break; 310*a2d72b6fSJames Wright 311*a2d72b6fSJames Wright case STATEVAR_ENTROPY: 312*a2d72b6fSJames Wright // {1,2,3,4} 313*a2d72b6fSJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 314*a2d72b6fSJames Wright break; 315192a7459SJames Wright } 316192a7459SJames Wright 31724941a69SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3184e139266SJames Wright if (bc->num_inflow > 0) { 319f8839eb4SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL)); 3204e139266SJames Wright } 321ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3224e139266SJames Wright } 323dada6cc0SJames Wright 324731c13d7SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 32566170c20SJames Wright CeedQFunction *qf_strongbc) { 326dada6cc0SJames Wright PetscFunctionBeginUser; 327f17d818dSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 32866170c20SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 329a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 330a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 331a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 332a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 333dada6cc0SJames Wright 334a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 335ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 336dada6cc0SJames Wright } 3375dc40723SJames Wright 338731c13d7SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 33905675bc2SJames Wright CeedQFunction *qf_strongbc) { 3405dc40723SJames Wright PetscFunctionBeginUser; 341cbef7084SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 34266170c20SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 343a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 344a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 3455dc40723SJames Wright 346a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 347ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3485dc40723SJames Wright } 349