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) { 1643cec1b8dSMohammed Amin PetscInt nmodes = 0, nprofs; 165493642f1SJames Wright 16606f41313SJames Wright PetscFunctionBeginUser; 16742454adaSJames Wright PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs)); 1683cec1b8dSMohammed Amin const PetscBool need_rand = (!(*stg_ctx)->mean_only) || (*stg_ctx)->use_fluctuating_IC; 1693cec1b8dSMohammed Amin if (need_rand) { 1703cec1b8dSMohammed Amin PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes)); 1715d28dccaSJames Wright PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP, 1723cec1b8dSMohammed Amin "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", 1733cec1b8dSMohammed Amin stg_rand_path, nmodes, STG_NMODES_MAX); 1743cec1b8dSMohammed Amin } 175493642f1SJames Wright 176493642f1SJames Wright { 17742454adaSJames Wright StgShur14Context temp_ctx; 1783cec1b8dSMohammed Amin PetscCall(PetscNew(&temp_ctx)); 17968d8e747SJames Wright *temp_ctx = **stg_ctx; 18068d8e747SJames Wright temp_ctx->nmodes = nmodes; 18168d8e747SJames Wright temp_ctx->nprofs = nprofs; 1823cec1b8dSMohammed Amin // nmode = 0 if random numbers are not read from file, therefore offsets will be correctly handled 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)); 2013cec1b8dSMohammed Amin if (need_rand) { 20242454adaSJames Wright PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx)); 20368d8e747SJames Wright CeedScalar *kappa = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa]; 20468d8e747SJames Wright CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist]; 20568d8e747SJames Wright CeedScalar *lt = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt]; 206493642f1SJames Wright CeedScalar le, le_max = 0; 207493642f1SJames Wright 20868d8e747SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) { 209c77f3192SJames Wright le = PetscMin(2 * wall_dist[i], 3 * lt[i]); 210493642f1SJames Wright if (le_max < le) le_max = le; 211493642f1SJames Wright } 212493642f1SJames Wright CeedScalar kmin = M_PI / le_max; 213493642f1SJames Wright 21414bd2a07SJames Wright CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) kappa[i] = kmin * pow((*stg_ctx)->alpha, i); 21568d8e747SJames Wright } 216d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 217493642f1SJames Wright } 218493642f1SJames Wright 219d6cac220SJames Wright static PetscErrorCode STGWeakInflowBCSetup_CreateIFunctionQF(BCDefinition bc_def, CeedQFunction *qf) { 220d6cac220SJames Wright HoneeBCStruct honee_bc; 221d6cac220SJames Wright 222d6cac220SJames Wright PetscFunctionBeginUser; 223d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 224d6cac220SJames Wright PetscCall(HoneeBCCreateIFunctionQF(bc_def, StgShur14Inflow, StgShur14Inflow_loc, honee_bc->qfctx, qf)); 225d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 226d6cac220SJames Wright } 227d6cac220SJames Wright 228d6cac220SJames Wright static PetscErrorCode STGWeakInflowBCSetup_CreateIJacobianQF(BCDefinition bc_def, CeedQFunction *qf) { 229d6cac220SJames Wright HoneeBCStruct honee_bc; 230d6cac220SJames Wright 231d6cac220SJames Wright PetscFunctionBeginUser; 232d6cac220SJames Wright PetscCall(BCDefinitionGetContext(bc_def, &honee_bc)); 233d6cac220SJames Wright PetscCall(HoneeBCCreateIJacobianQF(bc_def, StgShur14Inflow_Jacobian, StgShur14Inflow_Jacobian_loc, honee_bc->qfctx, qf)); 234d6cac220SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 235d6cac220SJames Wright } 236d6cac220SJames Wright 2370c373b74SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, Honee honee, const bool prescribe_T, const CeedScalar theta0, 2389ef62cddSJames Wright const CeedScalar P0) { 2390c373b74SJames Wright Ceed ceed = honee->ceed; 240493642f1SJames Wright char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat"; 241493642f1SJames Wright char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat"; 242b01c649fSJames Wright PetscBool mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE; 2430c373b74SJames Wright CeedScalar u0 = 0.0, alpha = 1.01, stg_dx = -1, stg_h_scale_factor = 1 / honee->app_ctx->degree; 244e07531f7SJames Wright CeedQFunctionContext stg_qfctx; 245493642f1SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 246493642f1SJames Wright 24706f41313SJames Wright PetscFunctionBeginUser; 248493642f1SJames Wright PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL); 2492b916ea7SJeremy L Thompson PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL)); 2503cec1b8dSMohammed Amin PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGRand.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL)); 2512b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL)); 2522b916ea7SJeremy L Thompson PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL)); 2532b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL)); 2542b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL)); 2552b916ea7SJeremy L Thompson PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL, 2562b916ea7SJeremy L Thompson use_fluctuating_IC, &use_fluctuating_IC, NULL)); 25784b557acSJames Wright PetscCall(PetscOptionsReal("-stg_dx", "Element length in x direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx)); 25821ba7ba4SJames Wright if (given_stg_dx && use_stgstrong) PetscCall(PetscPrintf(comm, "WARNING: -stg_dx is ignored for -stg_strong\n")); 25984b557acSJames Wright PetscCall(PetscOptionsReal("-stg_h_scale_factor", "Scale element size for cutoff frequency calculation", NULL, stg_h_scale_factor, 26084b557acSJames Wright &stg_h_scale_factor, NULL)); 26184b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dyScale", NULL, "libCEED 0.12.0", "Use -stg_h_scale_factor to scale all the element dimensions")); 26284b557acSJames Wright PetscCall(PetscOptionsDeprecated("-stg_dz", NULL, "libCEED 0.12.0", NULL)); 263493642f1SJames Wright PetscOptionsEnd(); 264493642f1SJames Wright 2652d898fa6SJames Wright PetscCall(PetscNew(&global_stg_ctx)); 266f5dc303cSJames Wright *global_stg_ctx = (struct STGShur14Context_){ 267f5dc303cSJames Wright .alpha = alpha, 268f5dc303cSJames Wright .u0 = u0, 269f5dc303cSJames Wright .is_implicit = honee->phys->implicit, 270f5dc303cSJames Wright .prescribe_T = prescribe_T, 271f5dc303cSJames Wright .mean_only = mean_only, 272f5dc303cSJames Wright .use_fluctuating_IC = use_fluctuating_IC, 273f5dc303cSJames Wright .theta0 = theta0, 274f5dc303cSJames Wright .P0 = P0, 275f5dc303cSJames Wright .h_scale_factor = stg_h_scale_factor, 276f5dc303cSJames Wright }; 277493642f1SJames Wright 278f9d6418bSJames Wright if (!use_stgstrong) { // Calculate dx assuming constant spacing 279493642f1SJames Wright PetscReal domain_min[3], domain_max[3], domain_size[3]; 28034b5deb1SJames Wright PetscCall(DMGetBoundingBox(dm, domain_min, domain_max)); 281493642f1SJames Wright for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i]; 282493642f1SJames Wright 283493642f1SJames Wright PetscInt nmax = 3, faces[3]; 2842b916ea7SJeremy L Thompson PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL)); 285b01c649fSJames Wright global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0]; 28684b557acSJames Wright PetscCheck((global_stg_ctx->dx > 0) && PetscIsNormalReal((PetscReal)global_stg_ctx->dx), comm, PETSC_ERR_LIB, 28784b557acSJames Wright "STG dx must be positive normal number, got %g", global_stg_ctx->dx); 288493642f1SJames Wright } 289493642f1SJames Wright 290e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 291cde3d787SJames Wright global_stg_ctx->newt_ctx = *newtonian_ig_ctx; 292e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 293493642f1SJames Wright 29442454adaSJames Wright PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx)); 295493642f1SJames Wright 2960c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &stg_qfctx)); 297e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx)); 298e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 299e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_qfctx, "solution time", offsetof(struct STGShur14Context_, time), 1, 300b4c37c5cSJames Wright "Physical time of the solution")); 301493642f1SJames Wright 302e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 303f5dc303cSJames Wright problem->ics = (HoneeQFSpec){.qf_func_ptr = ICsStg, .qf_loc = ICsStg_loc, .qfctx = stg_qfctx}; 30443bff553SJames Wright 305e6098bcdSJames Wright if (use_stgstrong) { 3068085925cSJames Wright // Use default boundary integral QF (BoundaryIntegral) in newtonian.h 3072eb7bf1fSJames Wright problem->use_strong_bc_ceed = PETSC_TRUE; 30858ce1233SJames Wright problem->set_bc_from_ics = PETSC_FALSE; 309e6098bcdSJames Wright } else { 310d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 311d6cac220SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 312d6cac220SJames Wright const char *name; 313d6cac220SJames Wright 314d6cac220SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 315d6cac220SJames Wright if (!strcmp(name, "inflow")) { 316d6cac220SJames Wright HoneeBCStruct honee_bc; 317d6cac220SJames Wright 318d6cac220SJames Wright PetscCall(PetscNew(&honee_bc)); 319d6cac220SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_qfctx, &honee_bc->qfctx)); 320d6cac220SJames Wright honee_bc->honee = honee; 3211abc2837SJames Wright honee_bc->num_comps_jac_data = honee->phys->implicit ? 11 : 0; 322*26d401f3SJames Wright PetscCall(BCDefinitionSetContext(bc_def, (PetscCtxDestroyFn *)HoneeBCDestroy, honee_bc)); 323d6cac220SJames Wright 324d6cac220SJames Wright PetscCall(BCDefinitionSetIFunction(bc_def, STGWeakInflowBCSetup_CreateIFunctionQF, HoneeBCAddIFunctionOp)); 325d6cac220SJames Wright PetscCall(BCDefinitionSetIJacobian(bc_def, STGWeakInflowBCSetup_CreateIJacobianQF, HoneeBCAddIJacobianOp)); 326d6cac220SJames Wright } 327d6cac220SJames Wright } 32858ce1233SJames Wright problem->set_bc_from_ics = PETSC_TRUE; 329e6098bcdSJames Wright } 330d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 331493642f1SJames Wright } 332e6098bcdSJames Wright 3339ef62cddSJames Wright // @brief Set STG strongly enforce components using DMAddBoundary 334d6cac220SJames Wright PetscErrorCode SetupStrongStg(DM dm, ProblemData problem, Physics phys) { 335e6098bcdSJames Wright DMLabel label; 336d374fb47SJames Wright PetscInt comps[5], num_comps = 4; 33706f41313SJames Wright 33806f41313SJames Wright PetscFunctionBeginUser; 3393636f6a4SJames Wright switch (phys->state_var) { 3403636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 341d374fb47SJames Wright // {0,1,2,3} for rho, rho*u, rho*v, rho*w 342d374fb47SJames Wright for (int i = 0; i < 4; i++) comps[i] = i; 3433636f6a4SJames Wright break; 3443636f6a4SJames Wright 3453636f6a4SJames Wright case STATEVAR_PRIMITIVE: 3463636f6a4SJames Wright // {1,2,3,4} for u, v, w, T 3473636f6a4SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3483636f6a4SJames Wright break; 3499b103f75SJames Wright 3509b103f75SJames Wright case STATEVAR_ENTROPY: 3519b103f75SJames Wright // {1,2,3,4} 3529b103f75SJames Wright for (int i = 0; i < 4; i++) comps[i] = i + 1; 3539b103f75SJames Wright break; 354d374fb47SJames Wright } 355d374fb47SJames Wright 35634b5deb1SJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label)); 357d6cac220SJames Wright for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 358d6cac220SJames Wright BCDefinition bc_def = problem->bc_defs[b]; 359d6cac220SJames Wright const char *name; 360d6cac220SJames Wright 361d6cac220SJames Wright PetscCall(BCDefinitionGetInfo(bc_def, &name, NULL, NULL)); 362d6cac220SJames Wright if (!strcmp(name, "inflow")) PetscCall(BCDefinitionSetEssential(bc_def, num_comps, comps)); 363e6098bcdSJames Wright } 364d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 365e6098bcdSJames Wright } 3666d0190e2SJames Wright 367991aef52SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size, 3686f188493SJames Wright CeedQFunction *qf_strongbc) { 3696d0190e2SJames Wright PetscFunctionBeginUser; 37006f41313SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc)); 3716f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 372b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 373b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE)); 374b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 375b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE)); 3766d0190e2SJames Wright 377e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfctx)); 378d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3796d0190e2SJames Wright } 3809eeef72bSJames Wright 381991aef52SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size, 38268d8e747SJames Wright CeedQFunction *qf_strongbc) { 3839eeef72bSJames Wright PetscFunctionBeginUser; 38442454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc)); 3856f188493SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE)); 386b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE)); 387b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE)); 3889eeef72bSJames Wright 389e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfctx)); 390d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3919eeef72bSJames Wright } 392