xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision 4548da4e4ef44dc0f2704ad6d48ac0ca4a16bc83)
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(PETSC_SUCCESS);
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       dims[2];
65   int            ndims;
66   FILE          *fp;
67   const PetscInt char_array_len = 512;
68   char           line[char_array_len];
69   char         **array;
70 
71   PetscFunctionBeginUser;
72 
73   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
74 
75   CeedScalar  rij[6][stg_ctx->nprofs];
76   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
77   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
78   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
79   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
80 
81   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
82     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
83     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
84     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
85                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, 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(PETSC_SUCCESS);
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       dims[2];
123   int            ndims;
124   FILE          *fp;
125   const PetscInt char_array_len = 512;
126   char           line[char_array_len];
127   char         **array;
128 
129   PetscFunctionBeginUser;
130   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
131 
132   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
133   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
134   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
135 
136   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
137     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
138     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
139     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
140                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, 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(PETSC_SUCCESS);
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 `*stg_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] stg_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 *stg_ctx, const CeedScalar ynodes[]) {
169   PetscInt nmodes, nprofs;
170   PetscFunctionBeginUser;
171 
172   // Get options
173   PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes));
174   PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs));
175   PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP,
176              "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path,
177              nmodes, STG_NMODES_MAX);
178 
179   {
180     STGShur14Context temp_ctx;
181     PetscCall(PetscCalloc1(1, &temp_ctx));
182     *temp_ctx                   = **stg_ctx;
183     temp_ctx->nmodes            = nmodes;
184     temp_ctx->nprofs            = nprofs;
185     temp_ctx->offsets.sigma     = 0;
186     temp_ctx->offsets.d         = nmodes * 3;
187     temp_ctx->offsets.phi       = temp_ctx->offsets.d + nmodes * 3;
188     temp_ctx->offsets.kappa     = temp_ctx->offsets.phi + nmodes;
189     temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes;
190     temp_ctx->offsets.ubar      = temp_ctx->offsets.wall_dist + nprofs;
191     temp_ctx->offsets.cij       = temp_ctx->offsets.ubar + nprofs * 3;
192     temp_ctx->offsets.eps       = temp_ctx->offsets.cij + nprofs * 6;
193     temp_ctx->offsets.lt        = temp_ctx->offsets.eps + nprofs;
194     temp_ctx->offsets.ynodes    = temp_ctx->offsets.lt + nprofs;
195     PetscInt total_num_scalars  = temp_ctx->offsets.ynodes + temp_ctx->nynodes;
196     temp_ctx->total_bytes       = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]);
197     PetscCall(PetscFree(*stg_ctx));
198     PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx));
199     **stg_ctx = *temp_ctx;
200     PetscCall(PetscFree(temp_ctx));
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     CeedScalar *kappa     = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa];
213     CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist];
214     CeedScalar *lt        = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt];
215     CeedScalar  le, le_max = 0;
216 
217     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) {
218       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
219       if (le_max < le) le_max = le;
220     }
221     CeedScalar kmin = M_PI / le_max;
222 
223     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); }
224   }
225 
226   PetscFunctionReturn(PETSC_SUCCESS);
227 }
228 
229 PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
230                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
231   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
232   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
233   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
234   CeedScalar               u0 = 0.0, alpha = 1.01;
235   CeedQFunctionContext     stg_context;
236   NewtonianIdealGasContext newtonian_ig_ctx;
237   PetscFunctionBeginUser;
238 
239   // Get options
240   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
241   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
242   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
243   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
244   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
245   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
246   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
247   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
248                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
249   PetscOptionsEnd();
250 
251   PetscCall(PetscCalloc1(1, &global_stg_ctx));
252   global_stg_ctx->alpha              = alpha;
253   global_stg_ctx->u0                 = u0;
254   global_stg_ctx->is_implicit        = user->phys->implicit;
255   global_stg_ctx->prescribe_T        = prescribe_T;
256   global_stg_ctx->mean_only          = mean_only;
257   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
258   global_stg_ctx->theta0             = theta0;
259   global_stg_ctx->P0                 = P0;
260   global_stg_ctx->nynodes            = nynodes;
261 
262   {
263     // Calculate dx assuming constant spacing
264     PetscReal domain_min[3], domain_max[3], domain_size[3];
265     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
266     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
267 
268     PetscInt nmax = 3, faces[3];
269     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
270     global_stg_ctx->dx = domain_size[0] / faces[0];
271     global_stg_ctx->dz = domain_size[2] / faces[2];
272   }
273 
274   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
275   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
276   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
277 
278   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
279 
280   CeedQFunctionContextCreate(user->ceed, &stg_context);
281   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
282   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
283   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
284 
285   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
286   problem->ics.qfunction         = ICsSTG;
287   problem->ics.qfunction_loc     = ICsSTG_loc;
288   problem->ics.qfunction_context = stg_context;
289 
290   if (use_stgstrong) {
291     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
292     problem->use_strong_bc_ceed = PETSC_TRUE;
293     problem->bc_from_ics        = PETSC_FALSE;
294   } else {
295     problem->apply_inflow.qfunction              = STGShur14_Inflow;
296     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
297     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
298     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
299     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
300     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
301     problem->bc_from_ics = PETSC_TRUE;
302   }
303 
304   PetscFunctionReturn(PETSC_SUCCESS);
305 }
306 
307 static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
308   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
309   // ^^assuming min(dy) is first element off the wall
310   PetscInt idx = -1;  // Index
311 
312   for (PetscInt i = 0; i < nynodes; i++) {
313     if (y < ynodes[i] + half_mindy) {
314       idx = i;
315       break;
316     }
317   }
318   if (idx == 0) return ynodes[1] - ynodes[0];
319   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
320   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
321 }
322 
323 // Function passed to DMAddBoundary
324 // NOTE: Not used in favor of QFunction-based method
325 PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
326   PetscFunctionBeginUser;
327 
328   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
329   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
330   const bool             mean_only = stg_ctx->mean_only;
331   const PetscScalar      dx        = stg_ctx->dx;
332   const PetscScalar      dz        = stg_ctx->dz;
333   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
334   const PetscScalar      theta0    = stg_ctx->theta0;
335   const PetscScalar      P0        = stg_ctx->P0;
336   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
337   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
338   const PetscScalar      Rd        = cp - cv;
339 
340   const CeedScalar rho = P0 / (Rd * theta0);
341   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
342   if (!mean_only) {
343     const PetscInt     nynodes = stg_ctx->nynodes;
344     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
345     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
346     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
347     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
348   } else {
349     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
350   }
351 
352   bcval[0] = rho;
353   bcval[1] = rho * u[0];
354   bcval[2] = rho * u[1];
355   bcval[3] = rho * u[2];
356   PetscFunctionReturn(PETSC_SUCCESS);
357 }
358 
359 PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
360   DMLabel label;
361   PetscFunctionBeginUser;
362 
363   PetscInt comps[5], num_comps = 4;
364   switch (phys->state_var) {
365     case STATEVAR_CONSERVATIVE:
366       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
367       for (int i = 0; i < 4; i++) comps[i] = i;
368       break;
369 
370     case STATEVAR_PRIMITIVE:
371       // {1,2,3,4} for u, v, w, T
372       for (int i = 0; i < 4; i++) comps[i] = i + 1;
373       break;
374   }
375 
376   PetscCall(DMGetLabel(dm, "Face Sets", &label));
377   // Set wall BCs
378   if (bc->num_inflow > 0) {
379     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
380                             NULL, global_stg_ctx, NULL));
381   }
382 
383   PetscFunctionReturn(PETSC_SUCCESS);
384 }
385 
386 PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
387                                  CeedInt q_data_size_sur, CeedQFunction *qf_strongbc) {
388   PetscFunctionBeginUser;
389   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, qf_strongbc);
390   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
391   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
392   CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE);
393   CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
394   CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
395 
396   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
397   PetscFunctionReturn(PETSC_SUCCESS);
398 }
399 
400 PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
401                                             CeedQFunction *qf_strongbc) {
402   PetscFunctionBeginUser;
403   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, qf_strongbc);
404   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
405   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
406   CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
407 
408   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
409   PetscFunctionReturn(PETSC_SUCCESS);
410 }
411