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