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); 102*38690fecSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 103ba6664aeSJames Wright } 1042b730f8bSJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 10524941a69SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 10624941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 107ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 108ba6664aeSJames Wright } 109ba6664aeSJames Wright 110ba6664aeSJames Wright /* 111ba6664aeSJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 112ba6664aeSJames Wright * 113ea61e9acSJeremy 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. 114ea61e9acSJeremy L Thompson * Assumes there are 7 columns in the file. 115ba6664aeSJames Wright * 116ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 117ba6664aeSJames Wright * @param[in] path Path to the STGRand.dat file 118ea61e9acSJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 119ba6664aeSJames Wright */ 120cbef7084SJames Wright static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 121457e73b2SJames Wright PetscInt dims[2]; 122457e73b2SJames Wright int ndims; 123ba6664aeSJames Wright FILE *fp; 124ba6664aeSJames Wright const PetscInt char_array_len = 512; 125ba6664aeSJames Wright char line[char_array_len]; 126ba6664aeSJames Wright char **array; 127ba6664aeSJames Wright 128ba6664aeSJames Wright PetscFunctionBeginUser; 129cbef7084SJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 130ba6664aeSJames Wright 131ba6664aeSJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1322b730f8bSJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1332b730f8bSJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 134ba6664aeSJames Wright 135ba6664aeSJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 13624941a69SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 13724941a69SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1380e654f56SJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 139457e73b2SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 140ba6664aeSJames Wright 141ba6664aeSJames Wright d[0][i] = (CeedScalar)atof(array[0]); 142ba6664aeSJames Wright d[1][i] = (CeedScalar)atof(array[1]); 143ba6664aeSJames Wright d[2][i] = (CeedScalar)atof(array[2]); 144ba6664aeSJames Wright phi[i] = (CeedScalar)atof(array[3]); 145ba6664aeSJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 146ba6664aeSJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 147ba6664aeSJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 148*38690fecSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 149ba6664aeSJames Wright } 15024941a69SJames Wright PetscCall(PetscFClose(comm, fp)); 151ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 152ba6664aeSJames Wright } 153ba6664aeSJames Wright 154ba6664aeSJames Wright /* 155ba6664aeSJames Wright * @brief Read STG data from input paths and put in STGShur14Context 156ba6664aeSJames Wright * 157ba6664aeSJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 15805675bc2SJames Wright * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance. 159ba6664aeSJames Wright * 160ba6664aeSJames Wright * @param[in] comm MPI_Comm for the program 161ba6664aeSJames Wright * @param[in] dm DM for the program 162ba6664aeSJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 163ba6664aeSJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 16405675bc2SJames Wright * @param[in,out] stg_ctx Pointer to STGShur14Context where the data will be loaded into 165ba6664aeSJames Wright */ 166cbef7084SJames 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], 167cbef7084SJames Wright StgShur14Context *stg_ctx) { 168ba6664aeSJames Wright PetscInt nmodes, nprofs; 169ba6664aeSJames Wright 170f17d818dSJames Wright PetscFunctionBeginUser; 171cbef7084SJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 172cbef7084SJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 1730e654f56SJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 174457e73b2SJames Wright "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path, 175457e73b2SJames Wright nmodes, STG_NMODES_MAX); 176ba6664aeSJames Wright 177ba6664aeSJames Wright { 178cbef7084SJames Wright StgShur14Context temp_ctx; 17905675bc2SJames Wright PetscCall(PetscCalloc1(1, &temp_ctx)); 18005675bc2SJames Wright *temp_ctx = **stg_ctx; 18105675bc2SJames Wright temp_ctx->nmodes = nmodes; 18205675bc2SJames Wright temp_ctx->nprofs = nprofs; 18305675bc2SJames Wright temp_ctx->offsets.sigma = 0; 18405675bc2SJames Wright temp_ctx->offsets.d = nmodes * 3; 18505675bc2SJames Wright temp_ctx->offsets.phi = temp_ctx->offsets.d + nmodes * 3; 18605675bc2SJames Wright temp_ctx->offsets.kappa = temp_ctx->offsets.phi + nmodes; 18705675bc2SJames Wright temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes; 18805675bc2SJames Wright temp_ctx->offsets.ubar = temp_ctx->offsets.wall_dist + nprofs; 18905675bc2SJames Wright temp_ctx->offsets.cij = temp_ctx->offsets.ubar + nprofs * 3; 19005675bc2SJames Wright temp_ctx->offsets.eps = temp_ctx->offsets.cij + nprofs * 6; 19105675bc2SJames Wright temp_ctx->offsets.lt = temp_ctx->offsets.eps + nprofs; 192f8839eb4SJames Wright PetscInt total_num_scalars = temp_ctx->offsets.lt + nprofs; 19305675bc2SJames Wright temp_ctx->total_bytes = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]); 19405675bc2SJames Wright PetscCall(PetscFree(*stg_ctx)); 19505675bc2SJames Wright PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx)); 19605675bc2SJames Wright **stg_ctx = *temp_ctx; 19705675bc2SJames Wright PetscCall(PetscFree(temp_ctx)); 198ba6664aeSJames Wright } 199ba6664aeSJames Wright 200cbef7084SJames Wright PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx)); 201cbef7084SJames Wright PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 202ba6664aeSJames Wright 20305675bc2SJames Wright { // -- Calculate kappa 20405675bc2SJames Wright CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 20505675bc2SJames Wright CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 20605675bc2SJames Wright CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 207ba6664aeSJames Wright CeedScalar le, le_max = 0; 208ba6664aeSJames Wright 20905675bc2SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 210175f00a6SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 211ba6664aeSJames Wright if (le_max < le) le_max = le; 212ba6664aeSJames Wright } 213ba6664aeSJames Wright CeedScalar kmin = M_PI / le_max; 214ba6664aeSJames Wright 21505675bc2SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); } 21605675bc2SJames Wright } 217ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 218ba6664aeSJames Wright } 219ba6664aeSJames Wright 220731c13d7SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, User user, const bool prescribe_T, const CeedScalar theta0, 221f8839eb4SJames Wright const CeedScalar P0) { 222a424bcd0SJames Wright Ceed ceed = user->ceed; 223ba6664aeSJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 224ba6664aeSJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 22573557d35SJames Wright PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 226831dbe9eSJames Wright CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = -1, stg_h_scale_factor = 1 / user->app_ctx->degree; 227ba6664aeSJames Wright CeedQFunctionContext stg_context; 228ba6664aeSJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 229ba6664aeSJames Wright 230f17d818dSJames Wright PetscFunctionBeginUser; 231ba6664aeSJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2322b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2332b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2342b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2352b730f8bSJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2362b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2372b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2382b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2392b730f8bSJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 240831dbe9eSJames Wright PetscCall(PetscOptionsReal("-stg_dx", "Element length in x direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 241831dbe9eSJames Wright PetscCall(PetscOptionsReal("-stg_h_scale_factor", "Scale element size for cutoff frequency calculation", NULL, stg_h_scale_factor, 242831dbe9eSJames Wright &stg_h_scale_factor, NULL)); 243831dbe9eSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dyScale", NULL, "libCEED 0.12.0", "Use -stg_h_scale_factor to scale all the element dimensions")); 244831dbe9eSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dz", NULL, "libCEED 0.12.0", NULL)); 245ba6664aeSJames Wright PetscOptionsEnd(); 246ba6664aeSJames Wright 24724941a69SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 24865dd5cafSJames Wright global_stg_ctx->alpha = alpha; 24965dd5cafSJames Wright global_stg_ctx->u0 = u0; 25065dd5cafSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 25165dd5cafSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 25265dd5cafSJames Wright global_stg_ctx->mean_only = mean_only; 25389060322SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 25465dd5cafSJames Wright global_stg_ctx->theta0 = theta0; 25565dd5cafSJames Wright global_stg_ctx->P0 = P0; 256831dbe9eSJames Wright global_stg_ctx->h_scale_factor = stg_h_scale_factor; 257ba6664aeSJames Wright 258f17d818dSJames Wright { // Calculate dx assuming constant spacing 259ba6664aeSJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 26024941a69SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 261ba6664aeSJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 262ba6664aeSJames Wright 263ba6664aeSJames Wright PetscInt nmax = 3, faces[3]; 2642b730f8bSJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 26573557d35SJames Wright global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 266831dbe9eSJames Wright PetscCheck((global_stg_ctx->dx > 0) && PetscIsNormalReal((PetscReal)global_stg_ctx->dx), comm, PETSC_ERR_LIB, 267831dbe9eSJames Wright "STG dx must be positive normal number, got %g", global_stg_ctx->dx); 268ba6664aeSJames Wright } 269ba6664aeSJames Wright 270a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 27165dd5cafSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 272a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 273ba6664aeSJames Wright 274cbef7084SJames Wright PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 275ba6664aeSJames Wright 276a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context)); 277a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 278a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc)); 279a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, 280a424bcd0SJames Wright "Physical time of the solution")); 281ba6664aeSJames Wright 282a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 283cbef7084SJames Wright problem->ics.qfunction = ICsStg; 284cbef7084SJames Wright problem->ics.qfunction_loc = ICsStg_loc; 285b77c53c9SJames Wright problem->ics.qfunction_context = stg_context; 286b77c53c9SJames Wright 2874e139266SJames Wright if (use_stgstrong) { 28865dd5cafSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 2893722cd23SJames Wright problem->use_strong_bc_ceed = PETSC_TRUE; 290de327db4SJames Wright problem->set_bc_from_ics = PETSC_FALSE; 2914e139266SJames Wright } else { 292cbef7084SJames Wright problem->apply_inflow.qfunction = StgShur14Inflow; 293cbef7084SJames Wright problem->apply_inflow.qfunction_loc = StgShur14Inflow_loc; 294cbef7084SJames Wright problem->apply_inflow_jacobian.qfunction = StgShur14Inflow_Jacobian; 295cbef7084SJames Wright problem->apply_inflow_jacobian.qfunction_loc = StgShur14Inflow_Jacobian_loc; 296a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context)); 297a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context)); 298de327db4SJames Wright problem->set_bc_from_ics = PETSC_TRUE; 2994e139266SJames Wright } 300ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 301ba6664aeSJames Wright } 3024e139266SJames Wright 303f8839eb4SJames Wright // @brief Set STG strongly enforce components using DMAddBoundary 304731c13d7SJames Wright PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData problem, Physics phys) { 3054e139266SJames Wright DMLabel label; 306192a7459SJames Wright PetscInt comps[5], num_comps = 4; 307f17d818dSJames Wright 308f17d818dSJames Wright PetscFunctionBeginUser; 30997baf651SJames Wright switch (phys->state_var) { 31097baf651SJames Wright case STATEVAR_CONSERVATIVE: 311192a7459SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 312192a7459SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 31397baf651SJames Wright break; 31497baf651SJames Wright 31597baf651SJames Wright case STATEVAR_PRIMITIVE: 31697baf651SJames Wright // {1,2,3,4} for u, v, w, T 31797baf651SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 31897baf651SJames Wright break; 319a2d72b6fSJames Wright 320a2d72b6fSJames Wright case STATEVAR_ENTROPY: 321a2d72b6fSJames Wright // {1,2,3,4} 322a2d72b6fSJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 323a2d72b6fSJames Wright break; 324192a7459SJames Wright } 325192a7459SJames Wright 32624941a69SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3274e139266SJames Wright if (bc->num_inflow > 0) { 328f8839eb4SJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL)); 3294e139266SJames Wright } 330ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3314e139266SJames Wright } 332dada6cc0SJames Wright 333731c13d7SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 33466170c20SJames Wright CeedQFunction *qf_strongbc) { 335dada6cc0SJames Wright PetscFunctionBeginUser; 336f17d818dSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 33766170c20SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 338a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 339a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 340a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 341a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 342dada6cc0SJames Wright 343a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 344ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 345dada6cc0SJames Wright } 3465dc40723SJames Wright 347731c13d7SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 34805675bc2SJames Wright CeedQFunction *qf_strongbc) { 3495dc40723SJames Wright PetscFunctionBeginUser; 350cbef7084SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 35166170c20SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 352a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 353a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 3545dc40723SJames Wright 355a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 356ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3575dc40723SJames Wright } 358