xref: /honee/problems/stg_shur14.c (revision 84b557ac940557fa3d4fece3311510aed39a37d2)
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 
19493642f1SJames 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);
102493642f1SJames Wright   }
1032b916ea7SJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
10434b5deb1SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
10534b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
106d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
107493642f1SJames Wright }
108493642f1SJames Wright 
109493642f1SJames Wright /*
110493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
111493642f1SJames Wright  *
11204e40bb6SJeremy 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.
11304e40bb6SJeremy L Thompson  * Assumes there are 7 columns in the file.
114493642f1SJames Wright  *
115493642f1SJames Wright  * @param[in]     comm    MPI_Comm for the program
116493642f1SJames Wright  * @param[in]     path    Path to the STGRand.dat file
11704e40bb6SJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
118493642f1SJames Wright  */
11942454adaSJames Wright static PetscErrorCode ReadStgRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], StgShur14Context stg_ctx) {
120defe8520SJames Wright   PetscInt       dims[2];
121defe8520SJames Wright   int            ndims;
122493642f1SJames Wright   FILE          *fp;
123493642f1SJames Wright   const PetscInt char_array_len = 512;
124493642f1SJames Wright   char           line[char_array_len];
125493642f1SJames Wright   char         **array;
126493642f1SJames Wright 
127493642f1SJames Wright   PetscFunctionBeginUser;
12842454adaSJames Wright   PetscCall(PhastaDatFileOpen(comm, path, char_array_len, dims, &fp));
129493642f1SJames Wright 
130493642f1SJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
1312b916ea7SJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
1322b916ea7SJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
133493642f1SJames Wright 
134493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
13534b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
13634b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1375d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
138defe8520SJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
139493642f1SJames Wright 
140493642f1SJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
141493642f1SJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
142493642f1SJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
143493642f1SJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
144493642f1SJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
145493642f1SJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
146493642f1SJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
147493642f1SJames Wright   }
14834b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
149d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
150493642f1SJames Wright }
151493642f1SJames Wright 
152493642f1SJames Wright /*
153493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
154493642f1SJames Wright  *
155493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
15668d8e747SJames Wright  * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance.
157493642f1SJames Wright  *
158493642f1SJames Wright  * @param[in]     comm            MPI_Comm for the program
159493642f1SJames Wright  * @param[in]     dm              DM for the program
160493642f1SJames Wright  * @param[in]     stg_inflow_path Path to STGInflow.dat file
161493642f1SJames Wright  * @param[in]     stg_rand_path   Path to STGRand.dat file
16268d8e747SJames Wright  * @param[in,out] stg_ctx         Pointer to STGShur14Context where the data will be loaded into
163493642f1SJames Wright  */
16442454adaSJames 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],
16542454adaSJames Wright                                  StgShur14Context *stg_ctx) {
166493642f1SJames Wright   PetscInt nmodes, nprofs;
167493642f1SJames Wright 
16806f41313SJames Wright   PetscFunctionBeginUser;
16942454adaSJames Wright   PetscCall(PhastaDatFileGetNRows(comm, stg_rand_path, &nmodes));
17042454adaSJames Wright   PetscCall(PhastaDatFileGetNRows(comm, stg_inflow_path, &nprofs));
1715d28dccaSJames Wright   PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP,
172defe8520SJames Wright              "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path,
173defe8520SJames Wright              nmodes, STG_NMODES_MAX);
174493642f1SJames Wright 
175493642f1SJames Wright   {
17642454adaSJames Wright     StgShur14Context temp_ctx;
17768d8e747SJames Wright     PetscCall(PetscCalloc1(1, &temp_ctx));
17868d8e747SJames Wright     *temp_ctx                   = **stg_ctx;
17968d8e747SJames Wright     temp_ctx->nmodes            = nmodes;
18068d8e747SJames Wright     temp_ctx->nprofs            = nprofs;
18168d8e747SJames Wright     temp_ctx->offsets.sigma     = 0;
18268d8e747SJames Wright     temp_ctx->offsets.d         = nmodes * 3;
18368d8e747SJames Wright     temp_ctx->offsets.phi       = temp_ctx->offsets.d + nmodes * 3;
18468d8e747SJames Wright     temp_ctx->offsets.kappa     = temp_ctx->offsets.phi + nmodes;
18568d8e747SJames Wright     temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes;
18668d8e747SJames Wright     temp_ctx->offsets.ubar      = temp_ctx->offsets.wall_dist + nprofs;
18768d8e747SJames Wright     temp_ctx->offsets.cij       = temp_ctx->offsets.ubar + nprofs * 3;
18868d8e747SJames Wright     temp_ctx->offsets.eps       = temp_ctx->offsets.cij + nprofs * 6;
18968d8e747SJames Wright     temp_ctx->offsets.lt        = temp_ctx->offsets.eps + nprofs;
1909ef62cddSJames Wright     PetscInt total_num_scalars  = temp_ctx->offsets.lt + nprofs;
19168d8e747SJames Wright     temp_ctx->total_bytes       = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]);
19268d8e747SJames Wright     PetscCall(PetscFree(*stg_ctx));
19368d8e747SJames Wright     PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx));
19468d8e747SJames Wright     **stg_ctx = *temp_ctx;
19568d8e747SJames Wright     PetscCall(PetscFree(temp_ctx));
196493642f1SJames Wright   }
197493642f1SJames Wright 
19842454adaSJames Wright   PetscCall(ReadStgInflow(comm, stg_inflow_path, *stg_ctx));
19942454adaSJames Wright   PetscCall(ReadStgRand(comm, stg_rand_path, *stg_ctx));
200493642f1SJames Wright 
20168d8e747SJames Wright   {  // -- Calculate kappa
20268d8e747SJames Wright     CeedScalar *kappa     = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa];
20368d8e747SJames Wright     CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist];
20468d8e747SJames Wright     CeedScalar *lt        = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt];
205493642f1SJames Wright     CeedScalar  le, le_max = 0;
206493642f1SJames Wright 
20768d8e747SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) {
208c77f3192SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
209493642f1SJames Wright       if (le_max < le) le_max = le;
210493642f1SJames Wright     }
211493642f1SJames Wright     CeedScalar kmin = M_PI / le_max;
212493642f1SJames Wright 
21368d8e747SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); }
21468d8e747SJames Wright   }
215d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
216493642f1SJames Wright }
217493642f1SJames Wright 
218991aef52SJames Wright PetscErrorCode SetupStg(const MPI_Comm comm, const DM dm, ProblemData problem, User user, const bool prescribe_T, const CeedScalar theta0,
2199ef62cddSJames Wright                         const CeedScalar P0) {
220b4c37c5cSJames Wright   Ceed                     ceed                                = user->ceed;
221493642f1SJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
222493642f1SJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
223b01c649fSJames Wright   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE;
224*84b557acSJames Wright   CeedScalar               u0 = 0.0, alpha = 1.01, stg_dx = -1, stg_h_scale_factor = 1 / user->app_ctx->degree;
225493642f1SJames Wright   CeedQFunctionContext     stg_context;
226493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
227493642f1SJames Wright 
22806f41313SJames Wright   PetscFunctionBeginUser;
229493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
2302b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
2312b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
2322b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
2332b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
2342b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
2352b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
2362b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
2372b916ea7SJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
238*84b557acSJames Wright   PetscCall(PetscOptionsReal("-stg_dx", "Element length in x direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx));
239*84b557acSJames Wright   PetscCall(PetscOptionsReal("-stg_h_scale_factor", "Scale element size for cutoff frequency calculation", NULL, stg_h_scale_factor,
240*84b557acSJames Wright                              &stg_h_scale_factor, NULL));
241*84b557acSJames Wright   PetscCall(PetscOptionsDeprecated("-stg_dyScale", NULL, "libCEED 0.12.0", "Use -stg_h_scale_factor to scale all the element dimensions"));
242*84b557acSJames Wright   PetscCall(PetscOptionsDeprecated("-stg_dz", NULL, "libCEED 0.12.0", NULL));
243493642f1SJames Wright   PetscOptionsEnd();
244493642f1SJames Wright 
24534b5deb1SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
2468085925cSJames Wright   global_stg_ctx->alpha              = alpha;
2478085925cSJames Wright   global_stg_ctx->u0                 = u0;
2488085925cSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
2498085925cSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
2508085925cSJames Wright   global_stg_ctx->mean_only          = mean_only;
251d4e0f297SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
2528085925cSJames Wright   global_stg_ctx->theta0             = theta0;
2538085925cSJames Wright   global_stg_ctx->P0                 = P0;
254*84b557acSJames Wright   global_stg_ctx->h_scale_factor     = stg_h_scale_factor;
255493642f1SJames Wright 
25606f41313SJames Wright   {  // Calculate dx assuming constant spacing
257493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
25834b5deb1SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
259493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
260493642f1SJames Wright 
261493642f1SJames Wright     PetscInt nmax = 3, faces[3];
2622b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
263b01c649fSJames Wright     global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0];
264*84b557acSJames Wright     PetscCheck((global_stg_ctx->dx > 0) && PetscIsNormalReal((PetscReal)global_stg_ctx->dx), comm, PETSC_ERR_LIB,
265*84b557acSJames Wright                "STG dx must be positive normal number, got %g", global_stg_ctx->dx);
266493642f1SJames Wright   }
267493642f1SJames Wright 
268b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
2698085925cSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
270b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
271493642f1SJames Wright 
27242454adaSJames Wright   PetscCall(GetStgContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx));
273493642f1SJames Wright 
274b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context));
275b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx));
276b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc));
277b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1,
278b4c37c5cSJames Wright                                                          "Physical time of the solution"));
279493642f1SJames Wright 
280b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
28142454adaSJames Wright   problem->ics.qfunction         = ICsStg;
28242454adaSJames Wright   problem->ics.qfunction_loc     = ICsStg_loc;
28343bff553SJames Wright   problem->ics.qfunction_context = stg_context;
28443bff553SJames Wright 
285e6098bcdSJames Wright   if (use_stgstrong) {
2868085925cSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
2872eb7bf1fSJames Wright     problem->use_strong_bc_ceed = PETSC_TRUE;
288d7244455SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
289e6098bcdSJames Wright   } else {
29042454adaSJames Wright     problem->apply_inflow.qfunction              = StgShur14Inflow;
29142454adaSJames Wright     problem->apply_inflow.qfunction_loc          = StgShur14Inflow_loc;
29242454adaSJames Wright     problem->apply_inflow_jacobian.qfunction     = StgShur14Inflow_Jacobian;
29342454adaSJames Wright     problem->apply_inflow_jacobian.qfunction_loc = StgShur14Inflow_Jacobian_loc;
294b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context));
295b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context));
296d7244455SJames Wright     problem->bc_from_ics = PETSC_TRUE;
297e6098bcdSJames Wright   }
298d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
299493642f1SJames Wright }
300e6098bcdSJames Wright 
3019ef62cddSJames Wright // @brief Set STG strongly enforce components using DMAddBoundary
302991aef52SJames Wright PetscErrorCode SetupStrongStg(DM dm, SimpleBC bc, ProblemData problem, Physics phys) {
303e6098bcdSJames Wright   DMLabel  label;
304d374fb47SJames Wright   PetscInt comps[5], num_comps = 4;
30506f41313SJames Wright 
30606f41313SJames Wright   PetscFunctionBeginUser;
3073636f6a4SJames Wright   switch (phys->state_var) {
3083636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
309d374fb47SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
310d374fb47SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
3113636f6a4SJames Wright       break;
3123636f6a4SJames Wright 
3133636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
3143636f6a4SJames Wright       // {1,2,3,4} for u, v, w, T
3153636f6a4SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
3163636f6a4SJames Wright       break;
317d374fb47SJames Wright   }
318d374fb47SJames Wright 
31934b5deb1SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
320e6098bcdSJames Wright   if (bc->num_inflow > 0) {
3219ef62cddSJames Wright     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL));
322e6098bcdSJames Wright   }
323d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
324e6098bcdSJames Wright }
3256d0190e2SJames Wright 
326991aef52SJames Wright PetscErrorCode SetupStrongStg_QF(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size, CeedInt dXdx_size,
3276f188493SJames Wright                                  CeedQFunction *qf_strongbc) {
3286d0190e2SJames Wright   PetscFunctionBeginUser;
32906f41313SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14InflowStrongQF, StgShur14InflowStrongQF_loc, qf_strongbc));
3306f188493SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE));
331b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE));
332b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE));
333b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE));
334b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE));
3356d0190e2SJames Wright 
336b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context));
337d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3386d0190e2SJames Wright }
3399eeef72bSJames Wright 
340991aef52SJames Wright PetscErrorCode SetupStrongStg_PreProcessing(Ceed ceed, ProblemData problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt dXdx_size,
34168d8e747SJames Wright                                             CeedQFunction *qf_strongbc) {
3429eeef72bSJames Wright   PetscFunctionBeginUser;
34342454adaSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, StgShur14Preprocess, StgShur14Preprocess_loc, qf_strongbc));
3446f188493SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "dXdx", dXdx_size, CEED_EVAL_NONE));
345b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE));
346b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE));
3479eeef72bSJames Wright 
348b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context));
349d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3509eeef72bSJames Wright }
351