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