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