xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision b2e3f8ecbfa285d0f4ffde9b24c57cc13f0319fb)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10 /// presented in Shur et al. 2014
11 
12 #include "stg_shur14.h"
13 
14 #include <ceed.h>
15 #include <math.h>
16 #include <petscdm.h>
17 #include <stdlib.h>
18 
19 #include "../navierstokes.h"
20 #include "../qfunctions/stg_shur14.h"
21 
22 STGShur14Context global_stg_ctx;
23 
24 /*
25  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
26  *
27  * This assumes the input matrices are in order [11,22,33,12,13,23].
28  * This format is also used for the output.
29  *
30  * @param[in]  comm   MPI_Comm
31  * @param[in]  nprofs Number of matrices in Rij
32  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
33  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
34  */
35 PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
36   PetscFunctionBeginUser;
37   for (PetscInt i = 0; i < nprofs; i++) {
38     Cij[0][i] = sqrt(Rij[0][i]);
39     Cij[3][i] = Rij[3][i] / Cij[0][i];
40     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2));
41     Cij[4][i] = Rij[4][i] / Cij[0][i];
42     Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i];
43     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
44 
45     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) {
46       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.",
47               i + 1);
48     }
49   }
50   PetscFunctionReturn(0);
51 }
52 
53 /*
54  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
55  *
56  * This function opens the file specified by `path` using `PetscFOpen` and passes the file pointer in `fp`.
57  * It is not closed in this function, thus `fp` must be closed sometime after this function has been called (using `PetscFClose` for example).
58  *
59  * 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.
60  *
61  * @param[in]  comm           MPI_Comm for the program
62  * @param[in]  path           Path to the file
63  * @param[in]  char_array_len Length of the character array that should contain each line
64  * @param[out] dims           Dimensions of the file, taken from the first line of the file
65  * @param[out] fp File        pointer to the opened file
66  */
67 static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
68                                         FILE **fp) {
69   PetscInt ndims;
70   char     line[char_array_len];
71   char   **array;
72 
73   PetscFunctionBeginUser;
74   PetscCall(PetscFOpen(comm, path, "r", fp));
75   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
76   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
77   if (ndims != 2) {
78     SETERRQ(comm, -1, "Found %" PetscInt_FMT " dimensions instead of 2 on the first line of %s", ndims, path);
79   }
80 
81   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
82   PetscCall(PetscStrToArrayDestroy(ndims, array));
83 
84   PetscFunctionReturn(0);
85 }
86 
87 /*
88  * @brief Get the number of rows for the PHASTA file at path.
89  *
90  * 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.
91  *
92  * @param[in]  comm  MPI_Comm for the program
93  * @param[in]  path  Path to the file
94  * @param[out] nrows Number of rows
95  */
96 static PetscErrorCode GetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
97   const PetscInt char_array_len = 512;
98   PetscInt       dims[2];
99   FILE          *fp;
100 
101   PetscFunctionBeginUser;
102   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
103   *nrows = dims[0];
104   PetscCall(PetscFClose(comm, fp));
105 
106   PetscFunctionReturn(0);
107 }
108 
109 /*
110  * @brief Read the STGInflow file and load the contents into stg_ctx
111  *
112  * 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.
113  * Assumes there are 14 columns in the file.
114  *
115  * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file.
116  *
117  * @param[in]     comm    MPI_Comm for the program
118  * @param[in]     path    Path to the STGInflow.dat file
119  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
120  */
121 static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
122   PetscInt       ndims, dims[2];
123   FILE          *fp;
124   const PetscInt char_array_len = 512;
125   char           line[char_array_len];
126   char         **array;
127 
128   PetscFunctionBeginUser;
129 
130   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
131 
132   CeedScalar  rij[6][stg_ctx->nprofs];
133   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
134   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
135   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
136   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
137 
138   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
139     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
140     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
141     if (ndims < dims[1]) {
142       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
143               ndims, dims[1]);
144     }
145 
146     wall_dist[i] = (CeedScalar)atof(array[0]);
147     ubar[0][i]   = (CeedScalar)atof(array[1]);
148     ubar[1][i]   = (CeedScalar)atof(array[2]);
149     ubar[2][i]   = (CeedScalar)atof(array[3]);
150     rij[0][i]    = (CeedScalar)atof(array[4]);
151     rij[1][i]    = (CeedScalar)atof(array[5]);
152     rij[2][i]    = (CeedScalar)atof(array[6]);
153     rij[3][i]    = (CeedScalar)atof(array[7]);
154     rij[4][i]    = (CeedScalar)atof(array[8]);
155     rij[5][i]    = (CeedScalar)atof(array[9]);
156     lt[i]        = (CeedScalar)atof(array[12]);
157     eps[i]       = (CeedScalar)atof(array[13]);
158 
159     if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path);
160     if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path);
161     if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path);
162   }
163   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
164   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
165   PetscCall(PetscFClose(comm, fp));
166 
167   PetscFunctionReturn(0);
168 }
169 
170 /*
171  * @brief Read the STGRand file and load the contents into stg_ctx
172  *
173  * 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.
174  * Assumes there are 7 columns in the file.
175  *
176  * @param[in]     comm    MPI_Comm for the program
177  * @param[in]     path    Path to the STGRand.dat file
178  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
179  */
180 static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
181   PetscInt       ndims, dims[2];
182   FILE          *fp;
183   const PetscInt char_array_len = 512;
184   char           line[char_array_len];
185   char         **array;
186 
187   PetscFunctionBeginUser;
188   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
189 
190   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
191   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
192   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
193 
194   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
195     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
196     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
197     if (ndims < dims[1]) {
198       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
199               ndims, dims[1]);
200     }
201 
202     d[0][i]     = (CeedScalar)atof(array[0]);
203     d[1][i]     = (CeedScalar)atof(array[1]);
204     d[2][i]     = (CeedScalar)atof(array[2]);
205     phi[i]      = (CeedScalar)atof(array[3]);
206     sigma[0][i] = (CeedScalar)atof(array[4]);
207     sigma[1][i] = (CeedScalar)atof(array[5]);
208     sigma[2][i] = (CeedScalar)atof(array[6]);
209   }
210   PetscCall(PetscFClose(comm, fp));
211 
212   PetscFunctionReturn(0);
213 }
214 
215 /*
216  * @brief Read STG data from input paths and put in STGShur14Context
217  *
218  * Reads data from input paths and puts them into a STGShur14Context object.
219  * Data stored initially in `*pstg_ctx` will be copied over to the new STGShur14Context instance.
220  *
221  * @param[in]     comm            MPI_Comm for the program
222  * @param[in]     dm              DM for the program
223  * @param[in]     stg_inflow_path Path to STGInflow.dat file
224  * @param[in]     stg_rand_path   Path to STGRand.dat file
225  * @param[in,out] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
226  */
227 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],
228                                  STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) {
229   PetscInt         nmodes, nprofs;
230   STGShur14Context stg_ctx;
231   PetscFunctionBeginUser;
232 
233   // Get options
234   PetscCall(GetNRows(comm, stg_rand_path, &nmodes));
235   PetscCall(GetNRows(comm, stg_inflow_path, &nprofs));
236   if (nmodes > STG_NMODES_MAX) {
237     SETERRQ(comm, 1,
238             "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT
239             "). "
240             "Change size of STG_NMODES_MAX and recompile",
241             stg_rand_path, nmodes, STG_NMODES_MAX);
242   }
243 
244   {
245     STGShur14Context s;
246     PetscCall(PetscCalloc1(1, &s));
247     *s                         = **pstg_ctx;
248     s->nmodes                  = nmodes;
249     s->nprofs                  = nprofs;
250     s->offsets.sigma           = 0;
251     s->offsets.d               = nmodes * 3;
252     s->offsets.phi             = s->offsets.d + nmodes * 3;
253     s->offsets.kappa           = s->offsets.phi + nmodes;
254     s->offsets.wall_dist       = s->offsets.kappa + nmodes;
255     s->offsets.ubar            = s->offsets.wall_dist + nprofs;
256     s->offsets.cij             = s->offsets.ubar + nprofs * 3;
257     s->offsets.eps             = s->offsets.cij + nprofs * 6;
258     s->offsets.lt              = s->offsets.eps + nprofs;
259     s->offsets.ynodes          = s->offsets.lt + nprofs;
260     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
261     s->total_bytes             = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]);
262     PetscCall(PetscMalloc(s->total_bytes, &stg_ctx));
263     *stg_ctx = *s;
264     PetscCall(PetscFree(s));
265   }
266 
267   PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx));
268   PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx));
269 
270   if (stg_ctx->nynodes > 0) {
271     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
272     for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
273   }
274 
275   // -- Calculate kappa
276   {
277     CeedScalar *kappa     = &stg_ctx->data[stg_ctx->offsets.kappa];
278     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
279     CeedScalar *lt        = &stg_ctx->data[stg_ctx->offsets.lt];
280     CeedScalar  le, le_max = 0;
281 
282     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
283       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
284       if (le_max < le) le_max = le;
285     }
286     CeedScalar kmin = M_PI / le_max;
287 
288     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); }
289   }  // end calculate kappa
290 
291   PetscCall(PetscFree(*pstg_ctx));
292   *pstg_ctx = stg_ctx;
293   PetscFunctionReturn(0);
294 }
295 
296 PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
297                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
298   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
299   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
300   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
301   CeedScalar               u0 = 0.0, alpha = 1.01;
302   CeedQFunctionContext     stg_context;
303   NewtonianIdealGasContext newtonian_ig_ctx;
304   PetscFunctionBeginUser;
305 
306   // Get options
307   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
308   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
309   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
310   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
311   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
312   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
313   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
314   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
315                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
316   PetscOptionsEnd();
317 
318   PetscCall(PetscCalloc1(1, &global_stg_ctx));
319   global_stg_ctx->alpha              = alpha;
320   global_stg_ctx->u0                 = u0;
321   global_stg_ctx->is_implicit        = user->phys->implicit;
322   global_stg_ctx->prescribe_T        = prescribe_T;
323   global_stg_ctx->mean_only          = mean_only;
324   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
325   global_stg_ctx->theta0             = theta0;
326   global_stg_ctx->P0                 = P0;
327   global_stg_ctx->nynodes            = nynodes;
328 
329   {
330     // Calculate dx assuming constant spacing
331     PetscReal domain_min[3], domain_max[3], domain_size[3];
332     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
333     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
334 
335     PetscInt nmax = 3, faces[3];
336     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
337     global_stg_ctx->dx = domain_size[0] / faces[0];
338     global_stg_ctx->dz = domain_size[2] / faces[2];
339   }
340 
341   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
342   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
343   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
344 
345   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
346 
347   CeedQFunctionContextCreate(user->ceed, &stg_context);
348   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
349   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
350   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
351 
352   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
353   problem->ics.qfunction         = ICsSTG;
354   problem->ics.qfunction_loc     = ICsSTG_loc;
355   problem->ics.qfunction_context = stg_context;
356 
357   if (use_stgstrong) {
358     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
359     problem->use_dirichlet_ceed = PETSC_TRUE;
360     problem->bc_from_ics        = PETSC_FALSE;
361   } else {
362     problem->apply_inflow.qfunction              = STGShur14_Inflow;
363     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
364     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
365     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
366     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
367     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
368     problem->bc_from_ics = PETSC_TRUE;
369   }
370 
371   PetscFunctionReturn(0);
372 }
373 
374 static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
375   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
376   // ^^assuming min(dy) is first element off the wall
377   PetscInt idx = -1;  // Index
378 
379   for (PetscInt i = 0; i < nynodes; i++) {
380     if (y < ynodes[i] + half_mindy) {
381       idx = i;
382       break;
383     }
384   }
385   if (idx == 0) return ynodes[1] - ynodes[0];
386   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
387   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
388 }
389 
390 // Function passed to DMAddBoundary
391 // NOTE: Not used in favor of QFunction-based method
392 PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
393   PetscFunctionBeginUser;
394 
395   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
396   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
397   const bool             mean_only = stg_ctx->mean_only;
398   const PetscScalar      dx        = stg_ctx->dx;
399   const PetscScalar      dz        = stg_ctx->dz;
400   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
401   const PetscScalar      theta0    = stg_ctx->theta0;
402   const PetscScalar      P0        = stg_ctx->P0;
403   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
404   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
405   const PetscScalar      Rd        = cp - cv;
406 
407   const CeedScalar rho = P0 / (Rd * theta0);
408   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
409   if (!mean_only) {
410     const PetscInt     nynodes = stg_ctx->nynodes;
411     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
412     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
413     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
414     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
415   } else {
416     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
417   }
418 
419   bcval[0] = rho;
420   bcval[1] = rho * u[0];
421   bcval[2] = rho * u[1];
422   bcval[3] = rho * u[2];
423   PetscFunctionReturn(0);
424 }
425 
426 PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
427   DMLabel label;
428   PetscFunctionBeginUser;
429 
430   PetscInt comps[5], num_comps = 4;
431   switch (phys->state_var) {
432     case STATEVAR_CONSERVATIVE:
433       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
434       for (int i = 0; i < 4; i++) comps[i] = i;
435       break;
436 
437     case STATEVAR_PRIMITIVE:
438       // {1,2,3,4} for u, v, w, T
439       for (int i = 0; i < 4; i++) comps[i] = i + 1;
440       break;
441   }
442 
443   PetscCall(DMGetLabel(dm, "Face Sets", &label));
444   // Set wall BCs
445   if (bc->num_inflow > 0) {
446     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
447                             NULL, global_stg_ctx, NULL));
448   }
449 
450   PetscFunctionReturn(0);
451 }
452 
453 PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
454                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
455   CeedQFunction qf_strongbc;
456   PetscFunctionBeginUser;
457   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
458   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
459   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
460   CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE);
461   CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
462   CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
463 
464   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
465   *pqf_strongbc = qf_strongbc;
466   PetscFunctionReturn(0);
467 }
468 
469 PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
470                                             CeedQFunction *pqf_strongbc) {
471   CeedQFunction qf_strongbc;
472   PetscFunctionBeginUser;
473   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc);
474   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
475   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
476   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
477 
478   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
479   *pqf_strongbc = qf_strongbc;
480   PetscFunctionReturn(0);
481 }
482