xref: /honee/problems/stg_shur14.c (revision d949ddfc92e3b3e8aad90700e3fb60ca7f4585db)
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 
455d28dccaSJames Wright     PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP,
465d28dccaSJames Wright                "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1);
47493642f1SJames Wright   }
48*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
49493642f1SJames Wright }
50493642f1SJames Wright 
51493642f1SJames Wright /*
52493642f1SJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
53493642f1SJames Wright  *
5404e40bb6SJeremy 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.
5504e40bb6SJeremy L Thompson  * Assumes there are 14 columns in the file.
56493642f1SJames Wright  *
5704e40bb6SJeremy L Thompson  * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file.
58493642f1SJames Wright  *
59493642f1SJames Wright  * @param[in]     comm    MPI_Comm for the program
60493642f1SJames Wright  * @param[in]     path    Path to the STGInflow.dat file
6104e40bb6SJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
62493642f1SJames Wright  */
632b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
64493642f1SJames Wright   PetscInt       ndims, dims[2];
65493642f1SJames Wright   FILE          *fp;
66493642f1SJames Wright   const PetscInt char_array_len = 512;
67493642f1SJames Wright   char           line[char_array_len];
68493642f1SJames Wright   char         **array;
69493642f1SJames Wright 
70493642f1SJames Wright   PetscFunctionBeginUser;
71493642f1SJames Wright 
72676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
73493642f1SJames Wright 
74493642f1SJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
75c77f3192SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
76493642f1SJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
77493642f1SJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
782b916ea7SJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
79493642f1SJames Wright 
80493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
8134b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
8234b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
835d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
844b0f6111SJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
854b0f6111SJames Wright                ndims, dims[1]);
86493642f1SJames Wright 
87c77f3192SJames Wright     wall_dist[i] = (CeedScalar)atof(array[0]);
88493642f1SJames Wright     ubar[0][i]   = (CeedScalar)atof(array[1]);
89493642f1SJames Wright     ubar[1][i]   = (CeedScalar)atof(array[2]);
90493642f1SJames Wright     ubar[2][i]   = (CeedScalar)atof(array[3]);
91493642f1SJames Wright     rij[0][i]    = (CeedScalar)atof(array[4]);
92493642f1SJames Wright     rij[1][i]    = (CeedScalar)atof(array[5]);
93493642f1SJames Wright     rij[2][i]    = (CeedScalar)atof(array[6]);
94493642f1SJames Wright     rij[3][i]    = (CeedScalar)atof(array[7]);
95493642f1SJames Wright     rij[4][i]    = (CeedScalar)atof(array[8]);
96493642f1SJames Wright     rij[5][i]    = (CeedScalar)atof(array[9]);
97493642f1SJames Wright     lt[i]        = (CeedScalar)atof(array[12]);
98493642f1SJames Wright     eps[i]       = (CeedScalar)atof(array[13]);
99493642f1SJames Wright 
1005d28dccaSJames Wright     PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path);
1015d28dccaSJames Wright     PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path);
1025d28dccaSJames Wright     PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path);
103493642f1SJames Wright   }
1042b916ea7SJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
10534b5deb1SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
10634b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
1072b916ea7SJeremy L Thompson 
108*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
109493642f1SJames Wright }
110493642f1SJames Wright 
111493642f1SJames Wright /*
112493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
113493642f1SJames Wright  *
11404e40bb6SJeremy 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.
11504e40bb6SJeremy L Thompson  * Assumes there are 7 columns in the file.
116493642f1SJames Wright  *
117493642f1SJames Wright  * @param[in]     comm    MPI_Comm for the program
118493642f1SJames Wright  * @param[in]     path    Path to the STGRand.dat file
11904e40bb6SJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
120493642f1SJames Wright  */
1212b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
122493642f1SJames Wright   PetscInt       ndims, dims[2];
123493642f1SJames Wright   FILE          *fp;
124493642f1SJames Wright   const PetscInt char_array_len = 512;
125493642f1SJames Wright   char           line[char_array_len];
126493642f1SJames Wright   char         **array;
127493642f1SJames Wright 
128493642f1SJames Wright   PetscFunctionBeginUser;
129676080b4SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
130493642f1SJames Wright 
131493642f1SJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
1322b916ea7SJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
1332b916ea7SJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
134493642f1SJames Wright 
135493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
13634b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
13734b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1385d28dccaSJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
1394b0f6111SJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1404b0f6111SJames Wright                ndims, dims[1]);
141493642f1SJames Wright 
142493642f1SJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
143493642f1SJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
144493642f1SJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
145493642f1SJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
146493642f1SJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
147493642f1SJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
148493642f1SJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
149493642f1SJames Wright   }
15034b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
1512b916ea7SJeremy L Thompson 
152*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153493642f1SJames Wright }
154493642f1SJames Wright 
155493642f1SJames Wright /*
156493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
157493642f1SJames Wright  *
158493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
15968d8e747SJames Wright  * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance.
160493642f1SJames Wright  *
161493642f1SJames Wright  * @param[in]     comm            MPI_Comm for the program
162493642f1SJames Wright  * @param[in]     dm              DM for the program
163493642f1SJames Wright  * @param[in]     stg_inflow_path Path to STGInflow.dat file
164493642f1SJames Wright  * @param[in]     stg_rand_path   Path to STGRand.dat file
16568d8e747SJames Wright  * @param[in,out] stg_ctx         Pointer to STGShur14Context where the data will be loaded into
166493642f1SJames Wright  */
1672b916ea7SJeremy 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],
16868d8e747SJames Wright                                  STGShur14Context *stg_ctx, const CeedScalar ynodes[]) {
169493642f1SJames Wright   PetscInt nmodes, nprofs;
170493642f1SJames Wright   PetscFunctionBeginUser;
171493642f1SJames Wright 
172493642f1SJames Wright   // Get options
173676080b4SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes));
174676080b4SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs));
1755d28dccaSJames Wright   PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP,
1765d28dccaSJames Wright              "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT "). Change size of STG_NMODES_MAX and recompile",
1772b916ea7SJeremy L Thompson              stg_rand_path, nmodes, STG_NMODES_MAX);
178493642f1SJames Wright 
179493642f1SJames Wright   {
18068d8e747SJames Wright     STGShur14Context temp_ctx;
18168d8e747SJames Wright     PetscCall(PetscCalloc1(1, &temp_ctx));
18268d8e747SJames Wright     *temp_ctx                   = **stg_ctx;
18368d8e747SJames Wright     temp_ctx->nmodes            = nmodes;
18468d8e747SJames Wright     temp_ctx->nprofs            = nprofs;
18568d8e747SJames Wright     temp_ctx->offsets.sigma     = 0;
18668d8e747SJames Wright     temp_ctx->offsets.d         = nmodes * 3;
18768d8e747SJames Wright     temp_ctx->offsets.phi       = temp_ctx->offsets.d + nmodes * 3;
18868d8e747SJames Wright     temp_ctx->offsets.kappa     = temp_ctx->offsets.phi + nmodes;
18968d8e747SJames Wright     temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes;
19068d8e747SJames Wright     temp_ctx->offsets.ubar      = temp_ctx->offsets.wall_dist + nprofs;
19168d8e747SJames Wright     temp_ctx->offsets.cij       = temp_ctx->offsets.ubar + nprofs * 3;
19268d8e747SJames Wright     temp_ctx->offsets.eps       = temp_ctx->offsets.cij + nprofs * 6;
19368d8e747SJames Wright     temp_ctx->offsets.lt        = temp_ctx->offsets.eps + nprofs;
19468d8e747SJames Wright     temp_ctx->offsets.ynodes    = temp_ctx->offsets.lt + nprofs;
19568d8e747SJames Wright     PetscInt total_num_scalars  = temp_ctx->offsets.ynodes + temp_ctx->nynodes;
19668d8e747SJames Wright     temp_ctx->total_bytes       = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]);
19768d8e747SJames Wright     PetscCall(PetscFree(*stg_ctx));
19868d8e747SJames Wright     PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx));
19968d8e747SJames Wright     **stg_ctx = *temp_ctx;
20068d8e747SJames Wright     PetscCall(PetscFree(temp_ctx));
201493642f1SJames Wright   }
202493642f1SJames Wright 
20368d8e747SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, *stg_ctx));
20468d8e747SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, *stg_ctx));
205493642f1SJames Wright 
20668d8e747SJames Wright   if ((*stg_ctx)->nynodes > 0) {
20768d8e747SJames Wright     CeedScalar *ynodes_ctx = &(*stg_ctx)->data[(*stg_ctx)->offsets.ynodes];
20868d8e747SJames Wright     for (PetscInt i = 0; i < (*stg_ctx)->nynodes; i++) ynodes_ctx[i] = ynodes[i];
209e6098bcdSJames Wright   }
210e6098bcdSJames Wright 
21168d8e747SJames Wright   {  // -- Calculate kappa
21268d8e747SJames Wright     CeedScalar *kappa     = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa];
21368d8e747SJames Wright     CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist];
21468d8e747SJames Wright     CeedScalar *lt        = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt];
215493642f1SJames Wright     CeedScalar  le, le_max = 0;
216493642f1SJames Wright 
21768d8e747SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) {
218c77f3192SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
219493642f1SJames Wright       if (le_max < le) le_max = le;
220493642f1SJames Wright     }
221493642f1SJames Wright     CeedScalar kmin = M_PI / le_max;
222493642f1SJames Wright 
22368d8e747SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); }
22468d8e747SJames Wright   }
225493642f1SJames Wright 
226*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
227493642f1SJames Wright }
228493642f1SJames Wright 
2292b916ea7SJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
2302b916ea7SJeremy L Thompson                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
231493642f1SJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
232493642f1SJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
2332b916ea7SJeremy L Thompson   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
2342b916ea7SJeremy L Thompson   CeedScalar               u0 = 0.0, alpha = 1.01;
235493642f1SJames Wright   CeedQFunctionContext     stg_context;
236493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
237493642f1SJames Wright   PetscFunctionBeginUser;
238493642f1SJames Wright 
239493642f1SJames Wright   // Get options
240493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
2412b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
2422b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
2432b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
2442b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
2452b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
2462b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
2472b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
2482b916ea7SJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
249493642f1SJames Wright   PetscOptionsEnd();
250493642f1SJames Wright 
25134b5deb1SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
2528085925cSJames Wright   global_stg_ctx->alpha              = alpha;
2538085925cSJames Wright   global_stg_ctx->u0                 = u0;
2548085925cSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
2558085925cSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
2568085925cSJames Wright   global_stg_ctx->mean_only          = mean_only;
257d4e0f297SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
2588085925cSJames Wright   global_stg_ctx->theta0             = theta0;
2598085925cSJames Wright   global_stg_ctx->P0                 = P0;
2608085925cSJames Wright   global_stg_ctx->nynodes            = nynodes;
261493642f1SJames Wright 
262493642f1SJames Wright   {
263493642f1SJames Wright     // Calculate dx assuming constant spacing
264493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
26534b5deb1SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
266493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
267493642f1SJames Wright 
268493642f1SJames Wright     PetscInt nmax = 3, faces[3];
2692b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
2708085925cSJames Wright     global_stg_ctx->dx = domain_size[0] / faces[0];
2718085925cSJames Wright     global_stg_ctx->dz = domain_size[2] / faces[2];
272493642f1SJames Wright   }
273493642f1SJames Wright 
2742b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
2758085925cSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
2762b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
277493642f1SJames Wright 
2782b916ea7SJeremy L Thompson   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
279493642f1SJames Wright 
280493642f1SJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
2812b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
2822b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
2832b916ea7SJeremy L Thompson   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
284493642f1SJames Wright 
28543bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
28643bff553SJames Wright   problem->ics.qfunction         = ICsSTG;
28743bff553SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
28843bff553SJames Wright   problem->ics.qfunction_context = stg_context;
28943bff553SJames Wright 
290e6098bcdSJames Wright   if (use_stgstrong) {
2918085925cSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
2922eb7bf1fSJames Wright     problem->use_strong_bc_ceed = PETSC_TRUE;
293d7244455SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
294e6098bcdSJames Wright   } else {
295493642f1SJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
296493642f1SJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
297a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
298a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
2992b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
3002b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
301d7244455SJames Wright     problem->bc_from_ics = PETSC_TRUE;
302e6098bcdSJames Wright   }
303493642f1SJames Wright 
304*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
305493642f1SJames Wright }
306e6098bcdSJames Wright 
3072b916ea7SJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
308e6098bcdSJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
309e6098bcdSJames Wright   // ^^assuming min(dy) is first element off the wall
310e6098bcdSJames Wright   PetscInt idx = -1;  // Index
311e6098bcdSJames Wright 
312e6098bcdSJames Wright   for (PetscInt i = 0; i < nynodes; i++) {
313e6098bcdSJames Wright     if (y < ynodes[i] + half_mindy) {
3142b916ea7SJeremy L Thompson       idx = i;
3152b916ea7SJeremy L Thompson       break;
316e6098bcdSJames Wright     }
317e6098bcdSJames Wright   }
318e6098bcdSJames Wright   if (idx == 0) return ynodes[1] - ynodes[0];
319e6098bcdSJames Wright   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
320e6098bcdSJames Wright   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
321e6098bcdSJames Wright }
322e6098bcdSJames Wright 
323e6098bcdSJames Wright // Function passed to DMAddBoundary
3246d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method
3252b916ea7SJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
326e6098bcdSJames Wright   PetscFunctionBeginUser;
3272b916ea7SJeremy L Thompson 
328e6098bcdSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
329e6098bcdSJames Wright   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
330e6098bcdSJames Wright   const bool             mean_only = stg_ctx->mean_only;
331e6098bcdSJames Wright   const PetscScalar      dx        = stg_ctx->dx;
332e6098bcdSJames Wright   const PetscScalar      dz        = stg_ctx->dz;
333e6098bcdSJames Wright   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
334e6098bcdSJames Wright   const PetscScalar      theta0    = stg_ctx->theta0;
335e6098bcdSJames Wright   const PetscScalar      P0        = stg_ctx->P0;
336e6098bcdSJames Wright   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
337e6098bcdSJames Wright   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
338e6098bcdSJames Wright   const PetscScalar      Rd        = cp - cv;
339e6098bcdSJames Wright 
340e6098bcdSJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
341e6098bcdSJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
342e6098bcdSJames Wright   if (!mean_only) {
343e6098bcdSJames Wright     const PetscInt     nynodes = stg_ctx->nynodes;
344e6098bcdSJames Wright     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
345e6098bcdSJames Wright     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
346e6098bcdSJames Wright     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
347e6098bcdSJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
348e6098bcdSJames Wright   } else {
349e6098bcdSJames Wright     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
350e6098bcdSJames Wright   }
351e6098bcdSJames Wright 
352e6098bcdSJames Wright   bcval[0] = rho;
353e6098bcdSJames Wright   bcval[1] = rho * u[0];
354e6098bcdSJames Wright   bcval[2] = rho * u[1];
355e6098bcdSJames Wright   bcval[3] = rho * u[2];
356*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
357e6098bcdSJames Wright }
358e6098bcdSJames Wright 
3592b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
360e6098bcdSJames Wright   DMLabel label;
361e6098bcdSJames Wright   PetscFunctionBeginUser;
362e6098bcdSJames Wright 
363d374fb47SJames Wright   PetscInt comps[5], num_comps = 4;
3643636f6a4SJames Wright   switch (phys->state_var) {
3653636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
366d374fb47SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
367d374fb47SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
3683636f6a4SJames Wright       break;
3693636f6a4SJames Wright 
3703636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
3713636f6a4SJames Wright       // {1,2,3,4} for u, v, w, T
3723636f6a4SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
3733636f6a4SJames Wright       break;
374d374fb47SJames Wright   }
375d374fb47SJames Wright 
37634b5deb1SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
377e6098bcdSJames Wright   // Set wall BCs
378e6098bcdSJames Wright   if (bc->num_inflow > 0) {
3792b916ea7SJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
38034b5deb1SJames Wright                             NULL, global_stg_ctx, NULL));
381e6098bcdSJames Wright   }
382e6098bcdSJames Wright 
383*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
384e6098bcdSJames Wright }
3856d0190e2SJames Wright 
3862b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
38768d8e747SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *qf_strongbc) {
3886d0190e2SJames Wright   PetscFunctionBeginUser;
38968d8e747SJames Wright   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, qf_strongbc);
39068d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
39168d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
39268d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE);
39368d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
39468d8e747SJames Wright   CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
3956d0190e2SJames Wright 
39668d8e747SJames Wright   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
397*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3986d0190e2SJames Wright }
3999eeef72bSJames Wright 
4002b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
40168d8e747SJames Wright                                             CeedQFunction *qf_strongbc) {
4029eeef72bSJames Wright   PetscFunctionBeginUser;
40368d8e747SJames Wright   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, qf_strongbc);
40468d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
40568d8e747SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
40668d8e747SJames Wright   CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4079eeef72bSJames Wright 
40868d8e747SJames Wright   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
409*d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4109eeef72bSJames Wright }
411