xref: /honee/problems/stg_shur14.c (revision 676080b4cb0b4c593f7c39e26e1c02044ea4890e)
1493642f1SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3493642f1SJames Wright //
4493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5493642f1SJames Wright //
6493642f1SJames Wright // This file is part of CEED:  http://github.com/ceed
7493642f1SJames Wright 
8493642f1SJames Wright /// @file
9493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10493642f1SJames Wright /// presented in Shur et al. 2014
11493642f1SJames Wright 
122b916ea7SJeremy L Thompson #include "stg_shur14.h"
132b916ea7SJeremy L Thompson 
14e419654dSJeremy L Thompson #include <ceed.h>
15493642f1SJames Wright #include <math.h>
16e419654dSJeremy L Thompson #include <petscdm.h>
172b916ea7SJeremy L Thompson #include <stdlib.h>
182b916ea7SJeremy L Thompson 
19493642f1SJames Wright #include "../navierstokes.h"
20493642f1SJames Wright #include "../qfunctions/stg_shur14.h"
21493642f1SJames Wright 
228085925cSJames Wright STGShur14Context global_stg_ctx;
238085925cSJames Wright 
24493642f1SJames Wright /*
25493642f1SJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
26493642f1SJames Wright  *
2704e40bb6SJeremy L Thompson  * This assumes the input matrices are in order [11,22,33,12,13,23].
2804e40bb6SJeremy L Thompson  * This format is also used for the output.
29493642f1SJames Wright  *
30493642f1SJames Wright  * @param[in]  comm   MPI_Comm
31493642f1SJames Wright  * @param[in]  nprofs Number of matrices in Rij
32493642f1SJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
33493642f1SJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
34493642f1SJames Wright  */
352b916ea7SJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
36493642f1SJames Wright   PetscFunctionBeginUser;
37493642f1SJames Wright   for (PetscInt i = 0; i < nprofs; i++) {
38493642f1SJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
39493642f1SJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
40493642f1SJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2));
41493642f1SJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
42493642f1SJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i];
43493642f1SJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
44493642f1SJames Wright 
452b916ea7SJeremy L Thompson     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) {
462b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.",
472b916ea7SJeremy L Thompson               i + 1);
482b916ea7SJeremy L Thompson     }
49493642f1SJames Wright   }
50493642f1SJames Wright   PetscFunctionReturn(0);
51493642f1SJames Wright }
52493642f1SJames Wright 
53493642f1SJames Wright /*
54493642f1SJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
55493642f1SJames Wright  *
5604e40bb6SJeremy L Thompson  * 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.
5704e40bb6SJeremy L Thompson  * Assumes there are 14 columns in the file.
58493642f1SJames Wright  *
5904e40bb6SJeremy L Thompson  * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file.
60493642f1SJames Wright  *
61493642f1SJames Wright  * @param[in]     comm    MPI_Comm for the program
62493642f1SJames Wright  * @param[in]     path    Path to the STGInflow.dat file
6304e40bb6SJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
64493642f1SJames Wright  */
652b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
66493642f1SJames Wright   PetscInt       ndims, dims[2];
67493642f1SJames Wright   FILE          *fp;
68493642f1SJames Wright   const PetscInt char_array_len = 512;
69493642f1SJames Wright   char           line[char_array_len];
70493642f1SJames Wright   char         **array;
71493642f1SJames Wright 
72493642f1SJames Wright   PetscFunctionBeginUser;
73493642f1SJames Wright 
74*676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
75493642f1SJames Wright 
76493642f1SJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
77c77f3192SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
78493642f1SJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
79493642f1SJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
802b916ea7SJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
81493642f1SJames Wright 
82493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
8334b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
8434b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
852b916ea7SJeremy L Thompson     if (ndims < dims[1]) {
862b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
872b916ea7SJeremy L Thompson               ndims, dims[1]);
882b916ea7SJeremy L Thompson     }
89493642f1SJames Wright 
90c77f3192SJames Wright     wall_dist[i] = (CeedScalar)atof(array[0]);
91493642f1SJames Wright     ubar[0][i]   = (CeedScalar)atof(array[1]);
92493642f1SJames Wright     ubar[1][i]   = (CeedScalar)atof(array[2]);
93493642f1SJames Wright     ubar[2][i]   = (CeedScalar)atof(array[3]);
94493642f1SJames Wright     rij[0][i]    = (CeedScalar)atof(array[4]);
95493642f1SJames Wright     rij[1][i]    = (CeedScalar)atof(array[5]);
96493642f1SJames Wright     rij[2][i]    = (CeedScalar)atof(array[6]);
97493642f1SJames Wright     rij[3][i]    = (CeedScalar)atof(array[7]);
98493642f1SJames Wright     rij[4][i]    = (CeedScalar)atof(array[8]);
99493642f1SJames Wright     rij[5][i]    = (CeedScalar)atof(array[9]);
100493642f1SJames Wright     lt[i]        = (CeedScalar)atof(array[12]);
101493642f1SJames Wright     eps[i]       = (CeedScalar)atof(array[13]);
102493642f1SJames Wright 
1032b916ea7SJeremy L Thompson     if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path);
1042b916ea7SJeremy L Thompson     if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path);
1052b916ea7SJeremy L Thompson     if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path);
106493642f1SJames Wright   }
1072b916ea7SJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
10834b5deb1SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
10934b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
1102b916ea7SJeremy L Thompson 
111493642f1SJames Wright   PetscFunctionReturn(0);
112493642f1SJames Wright }
113493642f1SJames Wright 
114493642f1SJames Wright /*
115493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
116493642f1SJames Wright  *
11704e40bb6SJeremy L Thompson  * 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.
11804e40bb6SJeremy L Thompson  * Assumes there are 7 columns in the file.
119493642f1SJames Wright  *
120493642f1SJames Wright  * @param[in]     comm    MPI_Comm for the program
121493642f1SJames Wright  * @param[in]     path    Path to the STGRand.dat file
12204e40bb6SJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
123493642f1SJames Wright  */
1242b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
125493642f1SJames Wright   PetscInt       ndims, dims[2];
126493642f1SJames Wright   FILE          *fp;
127493642f1SJames Wright   const PetscInt char_array_len = 512;
128493642f1SJames Wright   char           line[char_array_len];
129493642f1SJames Wright   char         **array;
130493642f1SJames Wright 
131493642f1SJames Wright   PetscFunctionBeginUser;
132*676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
133493642f1SJames Wright 
134493642f1SJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
1352b916ea7SJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
1362b916ea7SJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
137493642f1SJames Wright 
138493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
13934b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
14034b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1412b916ea7SJeremy L Thompson     if (ndims < dims[1]) {
1422b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1432b916ea7SJeremy L Thompson               ndims, dims[1]);
1442b916ea7SJeremy L Thompson     }
145493642f1SJames Wright 
146493642f1SJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
147493642f1SJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
148493642f1SJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
149493642f1SJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
150493642f1SJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
151493642f1SJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
152493642f1SJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
153493642f1SJames Wright   }
15434b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
1552b916ea7SJeremy L Thompson 
156493642f1SJames Wright   PetscFunctionReturn(0);
157493642f1SJames Wright }
158493642f1SJames Wright 
159493642f1SJames Wright /*
160493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
161493642f1SJames Wright  *
162493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
16304e40bb6SJeremy L Thompson  * Data stored initially in `*pstg_ctx` will be copied over to the new STGShur14Context instance.
164493642f1SJames Wright  *
165493642f1SJames Wright  * @param[in]     comm            MPI_Comm for the program
166493642f1SJames Wright  * @param[in]     dm              DM for the program
167493642f1SJames Wright  * @param[in]     stg_inflow_path Path to STGInflow.dat file
168493642f1SJames Wright  * @param[in]     stg_rand_path   Path to STGRand.dat file
16904e40bb6SJeremy L Thompson  * @param[in,out] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
170493642f1SJames Wright  */
1712b916ea7SJeremy L Thompson 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],
1722b916ea7SJeremy L Thompson                                  STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) {
173493642f1SJames Wright   PetscInt         nmodes, nprofs;
174493642f1SJames Wright   STGShur14Context stg_ctx;
175493642f1SJames Wright   PetscFunctionBeginUser;
176493642f1SJames Wright 
177493642f1SJames Wright   // Get options
178*676080b4SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes));
179*676080b4SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs));
180e419654dSJeremy L Thompson   if (nmodes > STG_NMODES_MAX) {
1812b916ea7SJeremy L Thompson     SETERRQ(comm, 1,
1822b916ea7SJeremy L Thompson             "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT
1832b916ea7SJeremy L Thompson             "). "
1842b916ea7SJeremy L Thompson             "Change size of STG_NMODES_MAX and recompile",
1852b916ea7SJeremy L Thompson             stg_rand_path, nmodes, STG_NMODES_MAX);
186e419654dSJeremy L Thompson   }
187493642f1SJames Wright 
188493642f1SJames Wright   {
189493642f1SJames Wright     STGShur14Context s;
19034b5deb1SJames Wright     PetscCall(PetscCalloc1(1, &s));
191493642f1SJames Wright     *s                         = **pstg_ctx;
192493642f1SJames Wright     s->nmodes                  = nmodes;
193493642f1SJames Wright     s->nprofs                  = nprofs;
194493642f1SJames Wright     s->offsets.sigma           = 0;
195493642f1SJames Wright     s->offsets.d               = nmodes * 3;
196493642f1SJames Wright     s->offsets.phi             = s->offsets.d + nmodes * 3;
197493642f1SJames Wright     s->offsets.kappa           = s->offsets.phi + nmodes;
198c77f3192SJames Wright     s->offsets.wall_dist       = s->offsets.kappa + nmodes;
199c77f3192SJames Wright     s->offsets.ubar            = s->offsets.wall_dist + nprofs;
200493642f1SJames Wright     s->offsets.cij             = s->offsets.ubar + nprofs * 3;
201493642f1SJames Wright     s->offsets.eps             = s->offsets.cij + nprofs * 6;
202493642f1SJames Wright     s->offsets.lt              = s->offsets.eps + nprofs;
203e6098bcdSJames Wright     s->offsets.ynodes          = s->offsets.lt + nprofs;
204e6098bcdSJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
205493642f1SJames Wright     s->total_bytes             = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]);
20634b5deb1SJames Wright     PetscCall(PetscMalloc(s->total_bytes, &stg_ctx));
207493642f1SJames Wright     *stg_ctx = *s;
20834b5deb1SJames Wright     PetscCall(PetscFree(s));
209493642f1SJames Wright   }
210493642f1SJames Wright 
21134b5deb1SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx));
21234b5deb1SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx));
213493642f1SJames Wright 
214e6098bcdSJames Wright   if (stg_ctx->nynodes > 0) {
215e6098bcdSJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
216e6098bcdSJames Wright     for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
217e6098bcdSJames Wright   }
218e6098bcdSJames Wright 
219493642f1SJames Wright   // -- Calculate kappa
220493642f1SJames Wright   {
221493642f1SJames Wright     CeedScalar *kappa     = &stg_ctx->data[stg_ctx->offsets.kappa];
222c77f3192SJames Wright     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
223493642f1SJames Wright     CeedScalar *lt        = &stg_ctx->data[stg_ctx->offsets.lt];
224493642f1SJames Wright     CeedScalar  le, le_max = 0;
225493642f1SJames Wright 
2262b916ea7SJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
227c77f3192SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
228493642f1SJames Wright       if (le_max < le) le_max = le;
229493642f1SJames Wright     }
230493642f1SJames Wright     CeedScalar kmin = M_PI / le_max;
231493642f1SJames Wright 
2322b916ea7SJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); }
233493642f1SJames Wright   }  // end calculate kappa
234493642f1SJames Wright 
23534b5deb1SJames Wright   PetscCall(PetscFree(*pstg_ctx));
236493642f1SJames Wright   *pstg_ctx = stg_ctx;
237493642f1SJames Wright   PetscFunctionReturn(0);
238493642f1SJames Wright }
239493642f1SJames Wright 
2402b916ea7SJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
2412b916ea7SJeremy L Thompson                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
242493642f1SJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
243493642f1SJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
2442b916ea7SJeremy L Thompson   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
2452b916ea7SJeremy L Thompson   CeedScalar               u0 = 0.0, alpha = 1.01;
246493642f1SJames Wright   CeedQFunctionContext     stg_context;
247493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
248493642f1SJames Wright   PetscFunctionBeginUser;
249493642f1SJames Wright 
250493642f1SJames Wright   // Get options
251493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
2522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
2532b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
2542b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
2552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
2562b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
2572b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
2582b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
2592b916ea7SJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
260493642f1SJames Wright   PetscOptionsEnd();
261493642f1SJames Wright 
26234b5deb1SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
2638085925cSJames Wright   global_stg_ctx->alpha              = alpha;
2648085925cSJames Wright   global_stg_ctx->u0                 = u0;
2658085925cSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
2668085925cSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
2678085925cSJames Wright   global_stg_ctx->mean_only          = mean_only;
268d4e0f297SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
2698085925cSJames Wright   global_stg_ctx->theta0             = theta0;
2708085925cSJames Wright   global_stg_ctx->P0                 = P0;
2718085925cSJames Wright   global_stg_ctx->nynodes            = nynodes;
272493642f1SJames Wright 
273493642f1SJames Wright   {
274493642f1SJames Wright     // Calculate dx assuming constant spacing
275493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
27634b5deb1SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
277493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
278493642f1SJames Wright 
279493642f1SJames Wright     PetscInt nmax = 3, faces[3];
2802b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
2818085925cSJames Wright     global_stg_ctx->dx = domain_size[0] / faces[0];
2828085925cSJames Wright     global_stg_ctx->dz = domain_size[2] / faces[2];
283493642f1SJames Wright   }
284493642f1SJames Wright 
2852b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
2868085925cSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
2872b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
288493642f1SJames Wright 
2892b916ea7SJeremy L Thompson   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
290493642f1SJames Wright 
291493642f1SJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
2922b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
2932b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
2942b916ea7SJeremy L Thompson   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
295493642f1SJames Wright 
29643bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
29743bff553SJames Wright   problem->ics.qfunction         = ICsSTG;
29843bff553SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
29943bff553SJames Wright   problem->ics.qfunction_context = stg_context;
30043bff553SJames Wright 
301e6098bcdSJames Wright   if (use_stgstrong) {
3028085925cSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
3036d0190e2SJames Wright     problem->use_dirichlet_ceed = PETSC_TRUE;
304d7244455SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
305e6098bcdSJames Wright   } else {
306493642f1SJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
307493642f1SJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
308a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
309a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
3102b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
3112b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
312d7244455SJames Wright     problem->bc_from_ics = PETSC_TRUE;
313e6098bcdSJames Wright   }
314493642f1SJames Wright 
315493642f1SJames Wright   PetscFunctionReturn(0);
316493642f1SJames Wright }
317e6098bcdSJames Wright 
3182b916ea7SJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
319e6098bcdSJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
320e6098bcdSJames Wright   // ^^assuming min(dy) is first element off the wall
321e6098bcdSJames Wright   PetscInt idx = -1;  // Index
322e6098bcdSJames Wright 
323e6098bcdSJames Wright   for (PetscInt i = 0; i < nynodes; i++) {
324e6098bcdSJames Wright     if (y < ynodes[i] + half_mindy) {
3252b916ea7SJeremy L Thompson       idx = i;
3262b916ea7SJeremy L Thompson       break;
327e6098bcdSJames Wright     }
328e6098bcdSJames Wright   }
329e6098bcdSJames Wright   if (idx == 0) return ynodes[1] - ynodes[0];
330e6098bcdSJames Wright   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
331e6098bcdSJames Wright   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
332e6098bcdSJames Wright }
333e6098bcdSJames Wright 
334e6098bcdSJames Wright // Function passed to DMAddBoundary
3356d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method
3362b916ea7SJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
337e6098bcdSJames Wright   PetscFunctionBeginUser;
3382b916ea7SJeremy L Thompson 
339e6098bcdSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
340e6098bcdSJames Wright   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
341e6098bcdSJames Wright   const bool             mean_only = stg_ctx->mean_only;
342e6098bcdSJames Wright   const PetscScalar      dx        = stg_ctx->dx;
343e6098bcdSJames Wright   const PetscScalar      dz        = stg_ctx->dz;
344e6098bcdSJames Wright   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
345e6098bcdSJames Wright   const PetscScalar      theta0    = stg_ctx->theta0;
346e6098bcdSJames Wright   const PetscScalar      P0        = stg_ctx->P0;
347e6098bcdSJames Wright   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
348e6098bcdSJames Wright   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
349e6098bcdSJames Wright   const PetscScalar      Rd        = cp - cv;
350e6098bcdSJames Wright 
351e6098bcdSJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
352e6098bcdSJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
353e6098bcdSJames Wright   if (!mean_only) {
354e6098bcdSJames Wright     const PetscInt     nynodes = stg_ctx->nynodes;
355e6098bcdSJames Wright     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
356e6098bcdSJames Wright     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
357e6098bcdSJames Wright     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
358e6098bcdSJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
359e6098bcdSJames Wright   } else {
360e6098bcdSJames Wright     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
361e6098bcdSJames Wright   }
362e6098bcdSJames Wright 
363e6098bcdSJames Wright   bcval[0] = rho;
364e6098bcdSJames Wright   bcval[1] = rho * u[0];
365e6098bcdSJames Wright   bcval[2] = rho * u[1];
366e6098bcdSJames Wright   bcval[3] = rho * u[2];
367e6098bcdSJames Wright   PetscFunctionReturn(0);
368e6098bcdSJames Wright }
369e6098bcdSJames Wright 
3702b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
371e6098bcdSJames Wright   DMLabel label;
372e6098bcdSJames Wright   PetscFunctionBeginUser;
373e6098bcdSJames Wright 
374d374fb47SJames Wright   PetscInt comps[5], num_comps = 4;
3753636f6a4SJames Wright   switch (phys->state_var) {
3763636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
377d374fb47SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
378d374fb47SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
3793636f6a4SJames Wright       break;
3803636f6a4SJames Wright 
3813636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
3823636f6a4SJames Wright       // {1,2,3,4} for u, v, w, T
3833636f6a4SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
3843636f6a4SJames Wright       break;
385d374fb47SJames Wright   }
386d374fb47SJames Wright 
38734b5deb1SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
388e6098bcdSJames Wright   // Set wall BCs
389e6098bcdSJames Wright   if (bc->num_inflow > 0) {
3902b916ea7SJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
39134b5deb1SJames Wright                             NULL, global_stg_ctx, NULL));
392e6098bcdSJames Wright   }
393e6098bcdSJames Wright 
394e6098bcdSJames Wright   PetscFunctionReturn(0);
395e6098bcdSJames Wright }
3966d0190e2SJames Wright 
3972b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
3989eeef72bSJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
3996d0190e2SJames Wright   CeedQFunction qf_strongbc;
4006d0190e2SJames Wright   PetscFunctionBeginUser;
4012b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
4022b916ea7SJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
4036d0190e2SJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
4046d0190e2SJames Wright   CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE);
4059eeef72bSJames Wright   CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4066d0190e2SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
4076d0190e2SJames Wright 
4086d0190e2SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
4096d0190e2SJames Wright   *pqf_strongbc = qf_strongbc;
4106d0190e2SJames Wright   PetscFunctionReturn(0);
4116d0190e2SJames Wright }
4129eeef72bSJames Wright 
4132b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
4149eeef72bSJames Wright                                             CeedQFunction *pqf_strongbc) {
4159eeef72bSJames Wright   CeedQFunction qf_strongbc;
4169eeef72bSJames Wright   PetscFunctionBeginUser;
4172b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc);
4182b916ea7SJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
4199eeef72bSJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
4209eeef72bSJames Wright   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4219eeef72bSJames Wright 
4229eeef72bSJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
4239eeef72bSJames Wright   *pqf_strongbc = qf_strongbc;
4249eeef72bSJames Wright   PetscFunctionReturn(0);
4259eeef72bSJames Wright }
426