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