1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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 19149fb536SJames Wright #include <navierstokes.h> 20493642f1SJames Wright #include "../qfunctions/stg_shur14.h" 21493642f1SJames Wright 2242454adaSJames 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 455d28dccaSJames Wright PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP, 465d28dccaSJames Wright "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1); 47493642f1SJames Wright } 48d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 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 */ 6342454adaSJames Wright static PetscErrorCode ReadStgInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 64defe8520SJames Wright PetscInt dims[2]; 65defe8520SJames Wright int ndims; 66493642f1SJames Wright FILE *fp; 67493642f1SJames Wright const PetscInt char_array_len = 512; 68493642f1SJames Wright char line[char_array_len]; 69493642f1SJames Wright char **array; 70493642f1SJames Wright 71493642f1SJames Wright PetscFunctionBeginUser; 7242454adaSJames 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)); 835d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 84defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 85493642f1SJames Wright 86c77f3192SJames Wright wall_dist[i] = (CeedScalar)atof(array[0]); 87493642f1SJames Wright ubar[0][i] = (CeedScalar)atof(array[1]); 88493642f1SJames Wright ubar[1][i] = (CeedScalar)atof(array[2]); 89493642f1SJames Wright ubar[2][i] = (CeedScalar)atof(array[3]); 90493642f1SJames Wright rij[0][i] = (CeedScalar)atof(array[4]); 91493642f1SJames Wright rij[1][i] = (CeedScalar)atof(array[5]); 92493642f1SJames Wright rij[2][i] = (CeedScalar)atof(array[6]); 93493642f1SJames Wright rij[3][i] = (CeedScalar)atof(array[7]); 94493642f1SJames Wright rij[4][i] = (CeedScalar)atof(array[8]); 95493642f1SJames Wright rij[5][i] = (CeedScalar)atof(array[9]); 96493642f1SJames Wright lt[i] = (CeedScalar)atof(array[12]); 97493642f1SJames Wright eps[i] = (CeedScalar)atof(array[13]); 98493642f1SJames Wright 995d28dccaSJames Wright PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path); 1005d28dccaSJames Wright PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path); 1015d28dccaSJames Wright PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path); 102039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 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)); 107d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 108493642f1SJames Wright } 109493642f1SJames Wright 110493642f1SJames Wright /* 111493642f1SJames Wright * @brief Read the STGRand file and load the contents into stg_ctx 112493642f1SJames Wright * 11304e40bb6SJeremy 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. 11404e40bb6SJeremy L Thompson * Assumes there are 7 columns in the file. 115493642f1SJames Wright * 116493642f1SJames Wright * @param[in] comm MPI_Comm for the program 117493642f1SJames Wright * @param[in] path Path to the STGRand.dat file 11804e40bb6SJeremy L Thompson * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into 119493642f1SJames Wright */ 12042454adaSJames Wright static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) { 121defe8520SJames Wright PetscInt dims[2]; 122defe8520SJames Wright int ndims; 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; 12942454adaSJames 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)); 1385d28dccaSJames Wright PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED, 139defe8520SJames Wright "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]); 140493642f1SJames Wright 141493642f1SJames Wright d[0][i] = (CeedScalar)atof(array[0]); 142493642f1SJames Wright d[1][i] = (CeedScalar)atof(array[1]); 143493642f1SJames Wright d[2][i] = (CeedScalar)atof(array[2]); 144493642f1SJames Wright phi[i] = (CeedScalar)atof(array[3]); 145493642f1SJames Wright sigma[0][i] = (CeedScalar)atof(array[4]); 146493642f1SJames Wright sigma[1][i] = (CeedScalar)atof(array[5]); 147493642f1SJames Wright sigma[2][i] = (CeedScalar)atof(array[6]); 148039caf0dSJames Wright PetscCall(PetscStrToArrayDestroy(ndims, array)); 149493642f1SJames Wright } 15034b5deb1SJames Wright PetscCall(PetscFClose(comm, fp)); 151d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 152493642f1SJames Wright } 153493642f1SJames Wright 154493642f1SJames Wright /* 155493642f1SJames Wright * @brief Read STG data from input paths and put in STGShur14Context 156493642f1SJames Wright * 157493642f1SJames Wright * Reads data from input paths and puts them into a STGShur14Context object. 15868d8e747SJames Wright * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance. 159493642f1SJames Wright * 160493642f1SJames Wright * @param[in] comm MPI_Comm for the program 161493642f1SJames Wright * @param[in] dm DM for the program 162493642f1SJames Wright * @param[in] stg_inflow_path Path to STGInflow.dat file 163493642f1SJames Wright * @param[in] stg_rand_path Path to STGRand.dat file 16468d8e747SJames Wright * @param[in,out] stg_ctx Pointer to STGShur14Context where the data will be loaded into 165493642f1SJames Wright */ 16642454adaSJames 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], 16742454adaSJames Wright StgShur14Context *stg_ctx) { 168493642f1SJames Wright PetscInt nmodes, nprofs; 169493642f1SJames Wright 17006f41313SJames Wright PetscFunctionBeginUser; 17142454adaSJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 17242454adaSJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 1735d28dccaSJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 174defe8520SJames Wright "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path, 175defe8520SJames Wright nmodes, STG_NMODES_MAX); 176493642f1SJames Wright 177493642f1SJames Wright { 17842454adaSJames Wright StgShur14Context temp_ctx; 17968d8e747SJames Wright PetscCall(PetscCalloc1(1, &temp_ctx)); 18068d8e747SJames Wright *temp_ctx = **stg_ctx; 18168d8e747SJames Wright temp_ctx->nmodes = nmodes; 18268d8e747SJames Wright temp_ctx->nprofs = nprofs; 18368d8e747SJames Wright temp_ctx->offsets.sigma = 0; 18468d8e747SJames Wright temp_ctx->offsets.d = nmodes * 3; 18568d8e747SJames Wright temp_ctx->offsets.phi = temp_ctx->offsets.d + nmodes * 3; 18668d8e747SJames Wright temp_ctx->offsets.kappa = temp_ctx->offsets.phi + nmodes; 18768d8e747SJames Wright temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes; 18868d8e747SJames Wright temp_ctx->offsets.ubar = temp_ctx->offsets.wall_dist + nprofs; 18968d8e747SJames Wright temp_ctx->offsets.cij = temp_ctx->offsets.ubar + nprofs * 3; 19068d8e747SJames Wright temp_ctx->offsets.eps = temp_ctx->offsets.cij + nprofs * 6; 19168d8e747SJames Wright temp_ctx->offsets.lt = temp_ctx->offsets.eps + nprofs; 1929ef62cddSJames Wright PetscInt total_num_scalars = temp_ctx->offsets.lt + nprofs; 19368d8e747SJames Wright temp_ctx->total_bytes = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]); 19468d8e747SJames Wright PetscCall(PetscFree(*stg_ctx)); 19568d8e747SJames Wright PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx)); 19668d8e747SJames Wright **stg_ctx = *temp_ctx; 19768d8e747SJames Wright PetscCall(PetscFree(temp_ctx)); 198493642f1SJames Wright } 199493642f1SJames Wright 20042454adaSJames Wright PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx)); 20142454adaSJames Wright PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 202493642f1SJames Wright 20368d8e747SJames Wright { // -- Calculate kappa 20468d8e747SJames Wright CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 20568d8e747SJames Wright CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 20668d8e747SJames Wright CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 207493642f1SJames Wright CeedScalar le, le_max = 0; 208493642f1SJames Wright 20968d8e747SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 210c77f3192SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 211493642f1SJames Wright if (le_max < le) le_max = le; 212493642f1SJames Wright } 213493642f1SJames Wright CeedScalar kmin = M_PI / le_max; 214493642f1SJames Wright 21568d8e747SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); } 21668d8e747SJames Wright } 217d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 218493642f1SJames Wright } 219493642f1SJames Wright 220991aef52SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, User user, const bool prescribe_T, const CeedScalar theta0, 2219ef62cddSJames Wright const CeedScalar P0) { 222b4c37c5cSJames Wright Ceed ceed = user->ceed; 223493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 224493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 225b01c649fSJames Wright PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 22684b557acSJames Wright CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = -1, stg_h_scale_factor = 1 / user->app_ctx->degree; 227493642f1SJames Wright CeedQFunctionContext stg_context; 228493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 229493642f1SJames Wright 23006f41313SJames Wright PetscFunctionBeginUser; 231493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2322b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2332b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2342b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2352b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2362b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2372b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2382b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2392b916ea7SJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 24084b557acSJames Wright PetscCall(PetscOptionsReal("-stg_dx", "Element length in x direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 241*21ba7ba4SJames Wright if (given_stg_dx && use_stgstrong) PetscCall(PetscPrintf(comm, "WARNING: -stg_dx is ignored for -stg_strong\n")); 24284b557acSJames Wright PetscCall(PetscOptionsReal("-stg_h_scale_factor", "Scale element size for cutoff frequency calculation", NULL, stg_h_scale_factor, 24384b557acSJames Wright &stg_h_scale_factor, NULL)); 24484b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dyScale", NULL, "libCEED 0.12.0", "Use -stg_h_scale_factor to scale all the element dimensions")); 24584b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dz", NULL, "libCEED 0.12.0", NULL)); 246493642f1SJames Wright PetscOptionsEnd(); 247493642f1SJames Wright 24834b5deb1SJames Wright PetscCall(PetscCalloc1(1, &global_stg_ctx)); 2498085925cSJames Wright global_stg_ctx->alpha = alpha; 2508085925cSJames Wright global_stg_ctx->u0 = u0; 2518085925cSJames Wright global_stg_ctx->is_implicit = user->phys->implicit; 2528085925cSJames Wright global_stg_ctx->prescribe_T = prescribe_T; 2538085925cSJames Wright global_stg_ctx->mean_only = mean_only; 254d4e0f297SJames Wright global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC; 2558085925cSJames Wright global_stg_ctx->theta0 = theta0; 2568085925cSJames Wright global_stg_ctx->P0 = P0; 25784b557acSJames Wright global_stg_ctx->h_scale_factor = stg_h_scale_factor; 258493642f1SJames Wright 25906f41313SJames Wright { // Calculate dx assuming constant spacing 260493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 26134b5deb1SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 262493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 263493642f1SJames Wright 264493642f1SJames Wright PetscInt nmax = 3, faces[3]; 2652b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 266b01c649fSJames Wright global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 26784b557acSJames Wright PetscCheck((global_stg_ctx->dx > 0) && PetscIsNormalReal((PetscReal)global_stg_ctx->dx), comm, PETSC_ERR_LIB, 26884b557acSJames Wright "STG dx must be positive normal number, got %g", global_stg_ctx->dx); 269493642f1SJames Wright } 270493642f1SJames Wright 271b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 2728085925cSJames Wright global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx; 273b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 274493642f1SJames Wright 27542454adaSJames Wright PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 276493642f1SJames Wright 277b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context)); 278b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 279b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc)); 280b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, 281b4c37c5cSJames Wright "Physical time of the solution")); 282493642f1SJames Wright 283b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 28442454adaSJames Wright problem->ics.qfunction = ICsStg; 28542454adaSJames Wright problem->ics.qfunction_loc = ICsStg_loc; 28643bff553SJames Wright problem->ics.qfunction_context = stg_context; 28743bff553SJames Wright 288e6098bcdSJames Wright if (use_stgstrong) { 2898085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 2902eb7bf1fSJames Wright problem->use_strong_bc_ceed = PETSC_TRUE; 29158ce1233SJames Wright problem->set_bc_from_ics = PETSC_FALSE; 292e6098bcdSJames Wright } else { 29342454adaSJames Wright problem->apply_inflow.qfunction = StgShur14Inflow; 29442454adaSJames Wright problem->apply_inflow.qfunction_loc = StgShur14Inflow_loc; 29542454adaSJames Wright problem->apply_inflow_jacobian.qfunction = StgShur14Inflow_Jacobian; 29642454adaSJames Wright problem->apply_inflow_jacobian.qfunction_loc = StgShur14Inflow_Jacobian_loc; 297b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context)); 298b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context)); 29958ce1233SJames Wright problem->set_bc_from_ics = PETSC_TRUE; 300e6098bcdSJames Wright } 301d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 302493642f1SJames Wright } 303e6098bcdSJames Wright 3049ef62cddSJames Wright // @brief Set STG strongly enforce components using DMAddBoundary 305991aef52SJames Wright PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData problem, Physics phys) { 306e6098bcdSJames Wright DMLabel label; 307d374fb47SJames Wright PetscInt comps[5], num_comps = 4; 30806f41313SJames Wright 30906f41313SJames Wright PetscFunctionBeginUser; 3103636f6a4SJames Wright switch (phys->state_var) { 3113636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 312d374fb47SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 313d374fb47SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 3143636f6a4SJames Wright break; 3153636f6a4SJames Wright 3163636f6a4SJames Wright case STATEVAR_PRIMITIVE: 3173636f6a4SJames Wright // {1,2,3,4} for u, v, w, T 3183636f6a4SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3193636f6a4SJames Wright break; 3209b103f75SJames Wright 3219b103f75SJames Wright case STATEVAR_ENTROPY: 3229b103f75SJames Wright // {1,2,3,4} 3239b103f75SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3249b103f75SJames Wright break; 325d374fb47SJames Wright } 326d374fb47SJames Wright 32734b5deb1SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 328e6098bcdSJames Wright if (bc->num_inflow > 0) { 3299ef62cddSJames Wright PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL)); 330e6098bcdSJames Wright } 331d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 332e6098bcdSJames Wright } 3336d0190e2SJames Wright 334991aef52SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 3356f188493SJames Wright CeedQFunction *qf_strongbc) { 3366d0190e2SJames Wright PetscFunctionBeginUser; 33706f41313SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 3386f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 339b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 340b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 341b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 342b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 3436d0190e2SJames Wright 344b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 345d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3466d0190e2SJames Wright } 3479eeef72bSJames Wright 348991aef52SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 34968d8e747SJames Wright CeedQFunction *qf_strongbc) { 3509eeef72bSJames Wright PetscFunctionBeginUser; 35142454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 3526f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 353b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 354b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 3559eeef72bSJames Wright 356b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context)); 357d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3589eeef72bSJames Wright } 359