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