1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3493642f1SJames Wright 4493642f1SJames Wright /// @file 5493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm 6493642f1SJames Wright /// presented in Shur et al. 2014 7493642f1SJames Wright 82b916ea7SJeremy L Thompson #include "stg_shur14.h" 92b916ea7SJeremy L Thompson 10e419654dSJeremy L Thompson #include <ceed.h> 11493642f1SJames Wright #include <math.h> 12e419654dSJeremy L Thompson #include <petscdm.h> 132b916ea7SJeremy L Thompson #include <stdlib.h> 142b916ea7SJeremy L Thompson 15149fb536SJames Wright #include <navierstokes.h> 16493642f1SJames Wright #include "../qfunctions/stg_shur14.h" 17493642f1SJames Wright 1842454adaSJames Wright StgShur14Context global_stg_ctx; 198085925cSJames Wright 20493642f1SJames Wright /* 21493642f1SJames Wright * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices 22493642f1SJames Wright * 2304e40bb6SJeremy L Thompson * This assumes the input matrices are in order [11,22,33,12,13,23]. 2404e40bb6SJeremy L Thompson * This format is also used for the output. 25493642f1SJames Wright * 26493642f1SJames Wright * @param[in] comm MPI_Comm 27493642f1SJames Wright * @param[in] nprofs Number of matrices in Rij 28493642f1SJames Wright * @param[in] Rij Array of the symmetric matrices [6,nprofs] 29493642f1SJames Wright * @param[out] Cij Array of the Cholesky Decomposition matrices, [6,nprofs] 30493642f1SJames Wright */ 312b916ea7SJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) { 32493642f1SJames Wright PetscFunctionBeginUser; 33493642f1SJames Wright for (PetscInt i = 0; i < nprofs; i++) { 34493642f1SJames Wright Cij[0][i] = sqrt(Rij[0][i]); 35493642f1SJames Wright Cij[3][i] = Rij[3][i] / Cij[0][i]; 3674960ff5SJames Wright Cij[1][i] = sqrt(Rij[1][i] - Square(Cij[3][i])); 37493642f1SJames Wright Cij[4][i] = Rij[4][i] / Cij[0][i]; 38493642f1SJames Wright Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i]; 3974960ff5SJames Wright Cij[2][i] = sqrt(Rij[2][i] - Square(Cij[4][i]) - Square(Cij[5][i])); 40493642f1SJames Wright 415d28dccaSJames Wright PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP, 425d28dccaSJames Wright "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1); 43493642f1SJames Wright } 44d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 45493642f1SJames Wright } 46493642f1SJames Wright 47493642f1SJames Wright /* 48493642f1SJames Wright * @brief Read the STGInflow file and load the contents into stg_ctx 49493642f1SJames Wright * 5004e40bb6SJeremy 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. 5104e40bb6SJeremy L Thompson * Assumes there are 14 columns in the file. 52493642f1SJames Wright * 5304e40bb6SJeremy L Thompson * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file. 54493642f1SJames Wright * 55493642f1SJames Wright * @param[in] comm MPI_Comm for the program 56493642f1SJames Wright * @param[in] path Path to the STGInflow.dat file 5704e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 58493642f1SJames Wright */ 5942454adaSJames Wright static PetscErrorCode ReadStgInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 60defe8520SJames Wright PetscInt dims[2]; 61defe8520SJames Wright int ndims; 62493642f1SJames Wright FILE *fp; 63493642f1SJames Wright const PetscInt char_array_len = 512; 64493642f1SJames Wright char line[char_array_len]; 65493642f1SJames Wright char **array; 66493642f1SJames Wright 67493642f1SJames Wright PetscFunctionBeginUser; 6842454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 69493642f1SJames Wright 70493642f1SJames Wright CeedScalar rij[6][stg_ctx->nprofs]; 71c77f3192SJames Wright CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist]; 72493642f1SJames Wright CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps]; 73493642f1SJames Wright CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt]; 742b916ea7SJeremy L Thompson CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar]; 75493642f1SJames Wright 76493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nprofs; i++) { 7734b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 7834b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 795d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 80defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 81493642f1SJames Wright 82c77f3192SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 83493642f1SJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 84493642f1SJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 85493642f1SJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 86493642f1SJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 87493642f1SJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 88493642f1SJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 89493642f1SJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 90493642f1SJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 91493642f1SJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 92493642f1SJames Wright lt[i] = (CeedScalar)atof(array[12]); 93493642f1SJames Wright eps[i] = (CeedScalar)atof(array[13]); 94493642f1SJames Wright 955d28dccaSJames Wright PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path); 965d28dccaSJames Wright PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path); 975d28dccaSJames Wright PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path); 98039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 99493642f1SJames Wright } 1002b916ea7SJeremy L Thompson CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij]; 10134b5deb1SJames Wright PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij)); 10234b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 103d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 104493642f1SJames Wright } 105493642f1SJames Wright 106493642f1SJames Wright /* 107493642f1SJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 108493642f1SJames Wright * 10904e40bb6SJeremy 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. 11004e40bb6SJeremy L Thompson * Assumes there are 7 columns in the file. 111493642f1SJames Wright * 112493642f1SJames Wright * @param[in] comm MPI_Comm for the program 113493642f1SJames Wright * @param[in] path Path to the STGRand.dat file 11404e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 115493642f1SJames Wright */ 11642454adaSJames Wright static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 117defe8520SJames Wright PetscInt dims[2]; 118defe8520SJames Wright int ndims; 119493642f1SJames Wright FILE *fp; 120493642f1SJames Wright const PetscInt char_array_len = 512; 121493642f1SJames Wright char line[char_array_len]; 122493642f1SJames Wright char **array; 123493642f1SJames Wright 124493642f1SJames Wright PetscFunctionBeginUser; 12542454adaSJames Wright PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp)); 126493642f1SJames Wright 127493642f1SJames Wright CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi]; 1282b916ea7SJeremy L Thompson CeedScalar(*d)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d]; 1292b916ea7SJeremy L Thompson CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma]; 130493642f1SJames Wright 131493642f1SJames Wright for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { 13234b5deb1SJames Wright PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line)); 13334b5deb1SJames Wright PetscCall(PetscStrToArray(line, ' ', &ndims, &array)); 1345d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 135defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 136493642f1SJames Wright 137493642f1SJames Wright d[0][i] = (CeedScalar)atof(array[0]); 138493642f1SJames Wright d[1][i] = (CeedScalar)atof(array[1]); 139493642f1SJames Wright d[2][i] = (CeedScalar)atof(array[2]); 140493642f1SJames Wright phi[i] = (CeedScalar)atof(array[3]); 141493642f1SJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 142493642f1SJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 143493642f1SJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 144039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 145493642f1SJames Wright } 14634b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 147d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 148493642f1SJames Wright } 149493642f1SJames Wright 150493642f1SJames Wright /* 151493642f1SJames Wright * @brief Read STG data from input paths and put in STGShur14Context 152493642f1SJames Wright * 153493642f1SJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 15468d8e747SJames Wright * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance. 155493642f1SJames Wright * 156493642f1SJames Wright * @param[in] comm MPI_Comm for the program 157493642f1SJames Wright * @param[in] dm DM for the program 158493642f1SJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 159493642f1SJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 16068d8e747SJames Wright * @param[in,out] stg_ctx Pointer to STGShur14Context where the data will be loaded into 161493642f1SJames Wright */ 16242454adaSJames 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], 16342454adaSJames Wright StgShur14Context *stg_ctx) { 164493642f1SJames Wright PetscInt nmodes, nprofs; 165493642f1SJames Wright 16606f41313SJames Wright PetscFunctionBeginUser; 16742454adaSJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 16842454adaSJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 1695d28dccaSJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 170defe8520SJames Wright "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path, 171defe8520SJames Wright nmodes, STG_NMODES_MAX); 172493642f1SJames Wright 173493642f1SJames Wright { 17442454adaSJames Wright StgShur14Context temp_ctx; 17568d8e747SJames Wright PetscCall(PetscCalloc1(1, &temp_ctx)); 17668d8e747SJames Wright *temp_ctx = **stg_ctx; 17768d8e747SJames Wright temp_ctx->nmodes = nmodes; 17868d8e747SJames Wright temp_ctx->nprofs = nprofs; 17968d8e747SJames Wright temp_ctx->offsets.sigma = 0; 18068d8e747SJames Wright temp_ctx->offsets.d = nmodes * 3; 18168d8e747SJames Wright temp_ctx->offsets.phi = temp_ctx->offsets.d + nmodes * 3; 18268d8e747SJames Wright temp_ctx->offsets.kappa = temp_ctx->offsets.phi + nmodes; 18368d8e747SJames Wright temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes; 18468d8e747SJames Wright temp_ctx->offsets.ubar = temp_ctx->offsets.wall_dist + nprofs; 18568d8e747SJames Wright temp_ctx->offsets.cij = temp_ctx->offsets.ubar + nprofs * 3; 18668d8e747SJames Wright temp_ctx->offsets.eps = temp_ctx->offsets.cij + nprofs * 6; 18768d8e747SJames Wright temp_ctx->offsets.lt = temp_ctx->offsets.eps + nprofs; 1889ef62cddSJames Wright PetscInt total_num_scalars = temp_ctx->offsets.lt + nprofs; 18968d8e747SJames Wright temp_ctx->total_bytes = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]); 19068d8e747SJames Wright PetscCall(PetscFree(*stg_ctx)); 19168d8e747SJames Wright PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx)); 19268d8e747SJames Wright **stg_ctx = *temp_ctx; 19368d8e747SJames Wright PetscCall(PetscFree(temp_ctx)); 194493642f1SJames Wright } 195493642f1SJames Wright 19642454adaSJames Wright PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx)); 19742454adaSJames Wright PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 198493642f1SJames Wright 19968d8e747SJames Wright { // -- Calculate kappa 20068d8e747SJames Wright CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 20168d8e747SJames Wright CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 20268d8e747SJames Wright CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 203493642f1SJames Wright CeedScalar le, le_max = 0; 204493642f1SJames Wright 20568d8e747SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 206c77f3192SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 207493642f1SJames Wright if (le_max < le) le_max = le; 208493642f1SJames Wright } 209493642f1SJames Wright CeedScalar kmin = M_PI / le_max; 210493642f1SJames Wright 21168d8e747SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); } 21268d8e747SJames Wright } 213d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 214493642f1SJames Wright } 215493642f1SJames Wright 2160c373b74SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, Honee honee, const bool prescribe_T, const CeedScalar theta0, 2179ef62cddSJames Wright const CeedScalar P0) { 2180c373b74SJames Wright Ceed ceed = honee->ceed; 219493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 220493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 221b01c649fSJames Wright PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 2220c373b74SJames Wright CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = -1, stg_h_scale_factor = 1 / honee->app_ctx->degree; 223e07531f7SJames Wright CeedQFunctionContext stg_qfctx; 224493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 225493642f1SJames Wright 22606f41313SJames Wright PetscFunctionBeginUser; 227493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2282b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2292b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2302b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2312b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2322b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2332b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2342b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2352b916ea7SJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 23684b557acSJames Wright PetscCall(PetscOptionsReal("-stg_dx", "Element length in x direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 23721ba7ba4SJames Wright if (given_stg_dx && use_stgstrong) PetscCall(PetscPrintf(comm, "WARNING: -stg_dx is ignored for -stg_strong\n")); 23884b557acSJames Wright PetscCall(PetscOptionsReal("-stg_h_scale_factor", "Scale element size for cutoff frequency calculation", NULL, stg_h_scale_factor, 23984b557acSJames Wright &stg_h_scale_factor, NULL)); 24084b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dyScale", NULL, "libCEED 0.12.0", "Use -stg_h_scale_factor to scale all the element dimensions")); 24184b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dz", NULL, "libCEED 0.12.0", NULL)); 242493642f1SJames Wright PetscOptionsEnd(); 243493642f1SJames Wright 24434b5deb1SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 2458085925cSJames Wright global_stg_ctx->alpha = alpha; 2468085925cSJames Wright global_stg_ctx->u0 = u0; 2470c373b74SJames Wright global_stg_ctx->is_implicit = honee->phys->implicit; 2488085925cSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 2498085925cSJames Wright global_stg_ctx->mean_only = mean_only; 250d4e0f297SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 2518085925cSJames Wright global_stg_ctx->theta0 = theta0; 2528085925cSJames Wright global_stg_ctx->P0 = P0; 25384b557acSJames Wright global_stg_ctx->h_scale_factor = stg_h_scale_factor; 254493642f1SJames Wright 255*f9d6418bSJames Wright if (!use_stgstrong) { // Calculate dx assuming constant spacing 256493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 25734b5deb1SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 258493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 259493642f1SJames Wright 260493642f1SJames Wright PetscInt nmax = 3, faces[3]; 2612b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 262b01c649fSJames Wright global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 26384b557acSJames Wright PetscCheck((global_stg_ctx->dx > 0) && PetscIsNormalReal((PetscReal)global_stg_ctx->dx), comm, PETSC_ERR_LIB, 26484b557acSJames Wright "STG dx must be positive normal number, got %g", global_stg_ctx->dx); 265493642f1SJames Wright } 266493642f1SJames Wright 267e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 2688085925cSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 269e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 270493642f1SJames Wright 27142454adaSJames Wright PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 272493642f1SJames Wright 2730c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &stg_qfctx)); 274e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 275e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 276e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_qfctx, "solution time", offsetof(struct STGShur14Context_, time), 1, 277b4c37c5cSJames Wright "Physical time of the solution")); 278493642f1SJames Wright 279e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 280e07531f7SJames Wright problem->ics.qf_func_ptr = ICsStg; 281e07531f7SJames Wright problem->ics.qf_loc = ICsStg_loc; 282e07531f7SJames Wright problem->ics.qfctx = stg_qfctx; 28343bff553SJames Wright 284e6098bcdSJames Wright if (use_stgstrong) { 2858085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 2862eb7bf1fSJames Wright problem->use_strong_bc_ceed = PETSC_TRUE; 28758ce1233SJames Wright problem->set_bc_from_ics = PETSC_FALSE; 288e6098bcdSJames Wright } else { 289e07531f7SJames Wright problem->apply_inflow.qf_func_ptr = StgShur14Inflow; 290e07531f7SJames Wright problem->apply_inflow.qf_loc = StgShur14Inflow_loc; 291e07531f7SJames Wright problem->apply_inflow_jacobian.qf_func_ptr = StgShur14Inflow_Jacobian; 292e07531f7SJames Wright problem->apply_inflow_jacobian.qf_loc = StgShur14Inflow_Jacobian_loc; 293e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_qfctx, &problem->apply_inflow.qfctx)); 294e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_qfctx, &problem->apply_inflow_jacobian.qfctx)); 29558ce1233SJames Wright problem->set_bc_from_ics = PETSC_TRUE; 296e6098bcdSJames Wright } 297d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 298493642f1SJames Wright } 299e6098bcdSJames Wright 3009ef62cddSJames Wright // @brief Set STG strongly enforce components using DMAddBoundary 301991aef52SJames Wright PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData problem, Physics phys) { 302e6098bcdSJames Wright DMLabel label; 303d374fb47SJames Wright PetscInt comps[5], num_comps = 4; 30406f41313SJames Wright 30506f41313SJames Wright PetscFunctionBeginUser; 3063636f6a4SJames Wright switch (phys->state_var) { 3073636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 308d374fb47SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 309d374fb47SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 3103636f6a4SJames Wright break; 3113636f6a4SJames Wright 3123636f6a4SJames Wright case STATEVAR_PRIMITIVE: 3133636f6a4SJames Wright // {1,2,3,4} for u, v, w, T 3143636f6a4SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3153636f6a4SJames Wright break; 3169b103f75SJames Wright 3179b103f75SJames Wright case STATEVAR_ENTROPY: 3189b103f75SJames Wright // {1,2,3,4} 3199b103f75SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3209b103f75SJames Wright break; 321d374fb47SJames Wright } 322d374fb47SJames Wright 32334b5deb1SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 324e6098bcdSJames Wright if (bc->num_inflow > 0) { 3259ef62cddSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL)); 326e6098bcdSJames Wright } 327d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 328e6098bcdSJames Wright } 3296d0190e2SJames Wright 330991aef52SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 3316f188493SJames Wright CeedQFunction *qf_strongbc) { 3326d0190e2SJames Wright PetscFunctionBeginUser; 33306f41313SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 3346f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 335b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 336b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 337b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 338b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 3396d0190e2SJames Wright 340e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfctx)); 341d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3426d0190e2SJames Wright } 3439eeef72bSJames Wright 344991aef52SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 34568d8e747SJames Wright CeedQFunction *qf_strongbc) { 3469eeef72bSJames Wright PetscFunctionBeginUser; 34742454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 3486f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 349b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 350b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 3519eeef72bSJames Wright 352e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfctx)); 353d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3549eeef72bSJames Wright } 355