xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision f8839eb4707ae2e91eaaa6df77b0ccbfa1ecfd98)
1ba6664aeSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2ba6664aeSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ba6664aeSJames Wright //
4ba6664aeSJames Wright // SPDX-License-Identifier: BSD-2-Clause
5ba6664aeSJames Wright //
6ba6664aeSJames Wright // This file is part of CEED:  http://github.com/ceed
7ba6664aeSJames Wright 
8ba6664aeSJames Wright /// @file
9ba6664aeSJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10ba6664aeSJames Wright /// presented in Shur et al. 2014
11ba6664aeSJames Wright 
122b730f8bSJeremy L Thompson #include "stg_shur14.h"
132b730f8bSJeremy L Thompson 
1449aac155SJeremy L Thompson #include <ceed.h>
15ba6664aeSJames Wright #include <math.h>
1649aac155SJeremy L Thompson #include <petscdm.h>
172b730f8bSJeremy L Thompson #include <stdlib.h>
182b730f8bSJeremy L Thompson 
19ba6664aeSJames Wright #include "../navierstokes.h"
20ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h"
21ba6664aeSJames Wright 
2265dd5cafSJames Wright STGShur14Context global_stg_ctx;
2365dd5cafSJames Wright 
24ba6664aeSJames Wright /*
25ba6664aeSJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
26ba6664aeSJames Wright  *
27ea61e9acSJeremy L Thompson  * This assumes the input matrices are in order [11,22,33,12,13,23].
28ea61e9acSJeremy L Thompson  * This format is also used for the output.
29ba6664aeSJames Wright  *
30ba6664aeSJames Wright  * @param[in]  comm   MPI_Comm
31ba6664aeSJames Wright  * @param[in]  nprofs Number of matrices in Rij
32ba6664aeSJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
33ba6664aeSJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
34ba6664aeSJames Wright  */
352b730f8bSJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
36ba6664aeSJames Wright   PetscFunctionBeginUser;
37ba6664aeSJames Wright   for (PetscInt i = 0; i < nprofs; i++) {
38ba6664aeSJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
39ba6664aeSJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
40ba6664aeSJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2));
41ba6664aeSJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
42ba6664aeSJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i];
43ba6664aeSJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
44ba6664aeSJames Wright 
450e654f56SJames Wright     PetscCheck(!isnan(Cij[0][i]) && !isnan(Cij[1][i]) && !isnan(Cij[2][i]), comm, PETSC_ERR_FP,
460e654f56SJames Wright                "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i + 1);
47ba6664aeSJames Wright   }
48ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
49ba6664aeSJames Wright }
50ba6664aeSJames Wright 
51ba6664aeSJames Wright /*
52ba6664aeSJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
53ba6664aeSJames Wright  *
54ea61e9acSJeremy 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.
55ea61e9acSJeremy L Thompson  * Assumes there are 14 columns in the file.
56ba6664aeSJames Wright  *
57ea61e9acSJeremy L Thompson  * Function calculates the Cholesky decomposition from the Reynolds stress profile found in the file.
58ba6664aeSJames Wright  *
59ba6664aeSJames Wright  * @param[in]     comm    MPI_Comm for the program
60ba6664aeSJames Wright  * @param[in]     path    Path to the STGInflow.dat file
61ea61e9acSJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
62ba6664aeSJames Wright  */
632b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
64457e73b2SJames Wright   PetscInt       dims[2];
65457e73b2SJames Wright   int            ndims;
66ba6664aeSJames Wright   FILE          *fp;
67ba6664aeSJames Wright   const PetscInt char_array_len = 512;
68ba6664aeSJames Wright   char           line[char_array_len];
69ba6664aeSJames Wright   char         **array;
70ba6664aeSJames Wright 
71ba6664aeSJames Wright   PetscFunctionBeginUser;
72ba6664aeSJames Wright 
738b892a05SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
74ba6664aeSJames Wright 
75ba6664aeSJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
76175f00a6SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
77ba6664aeSJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
78ba6664aeSJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
792b730f8bSJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
80ba6664aeSJames Wright 
81ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
8224941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
8324941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
840e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
85457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
86ba6664aeSJames Wright 
87175f00a6SJames Wright     wall_dist[i] = (CeedScalar)atof(array[0]);
88ba6664aeSJames Wright     ubar[0][i]   = (CeedScalar)atof(array[1]);
89ba6664aeSJames Wright     ubar[1][i]   = (CeedScalar)atof(array[2]);
90ba6664aeSJames Wright     ubar[2][i]   = (CeedScalar)atof(array[3]);
91ba6664aeSJames Wright     rij[0][i]    = (CeedScalar)atof(array[4]);
92ba6664aeSJames Wright     rij[1][i]    = (CeedScalar)atof(array[5]);
93ba6664aeSJames Wright     rij[2][i]    = (CeedScalar)atof(array[6]);
94ba6664aeSJames Wright     rij[3][i]    = (CeedScalar)atof(array[7]);
95ba6664aeSJames Wright     rij[4][i]    = (CeedScalar)atof(array[8]);
96ba6664aeSJames Wright     rij[5][i]    = (CeedScalar)atof(array[9]);
97ba6664aeSJames Wright     lt[i]        = (CeedScalar)atof(array[12]);
98ba6664aeSJames Wright     eps[i]       = (CeedScalar)atof(array[13]);
99ba6664aeSJames Wright 
1000e654f56SJames Wright     PetscCheck(wall_dist[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Distance to wall in %s cannot be negative", path);
1010e654f56SJames Wright     PetscCheck(lt[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent length scale in %s cannot be negative", path);
1020e654f56SJames Wright     PetscCheck(eps[i] >= 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Turbulent dissipation in %s cannot be negative", path);
103ba6664aeSJames Wright   }
1042b730f8bSJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
10524941a69SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
10624941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
1072b730f8bSJeremy L Thompson 
108ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
109ba6664aeSJames Wright }
110ba6664aeSJames Wright 
111ba6664aeSJames Wright /*
112ba6664aeSJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
113ba6664aeSJames Wright  *
114ea61e9acSJeremy 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.
115ea61e9acSJeremy L Thompson  * Assumes there are 7 columns in the file.
116ba6664aeSJames Wright  *
117ba6664aeSJames Wright  * @param[in]     comm    MPI_Comm for the program
118ba6664aeSJames Wright  * @param[in]     path    Path to the STGRand.dat file
119ea61e9acSJeremy L Thompson  * @param[in,out] stg_ctx STGShur14Context where the data will be loaded into
120ba6664aeSJames Wright  */
1212b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
122457e73b2SJames Wright   PetscInt       dims[2];
123457e73b2SJames Wright   int            ndims;
124ba6664aeSJames Wright   FILE          *fp;
125ba6664aeSJames Wright   const PetscInt char_array_len = 512;
126ba6664aeSJames Wright   char           line[char_array_len];
127ba6664aeSJames Wright   char         **array;
128ba6664aeSJames Wright 
129ba6664aeSJames Wright   PetscFunctionBeginUser;
1308b892a05SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
131ba6664aeSJames Wright 
132ba6664aeSJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
1332b730f8bSJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
1342b730f8bSJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
135ba6664aeSJames Wright 
136ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
13724941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
13824941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1390e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
140457e73b2SJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%d instead of %" PetscInt_FMT ")", i, path, ndims, dims[1]);
141ba6664aeSJames Wright 
142ba6664aeSJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
143ba6664aeSJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
144ba6664aeSJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
145ba6664aeSJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
146ba6664aeSJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
147ba6664aeSJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
148ba6664aeSJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
149ba6664aeSJames Wright   }
15024941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
1512b730f8bSJeremy L Thompson 
152ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
153ba6664aeSJames Wright }
154ba6664aeSJames Wright 
155ba6664aeSJames Wright /*
156ba6664aeSJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
157ba6664aeSJames Wright  *
158ba6664aeSJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
15905675bc2SJames Wright  * Data stored initially in `*stg_ctx` will be copied over to the new STGShur14Context instance.
160ba6664aeSJames Wright  *
161ba6664aeSJames Wright  * @param[in]     comm            MPI_Comm for the program
162ba6664aeSJames Wright  * @param[in]     dm              DM for the program
163ba6664aeSJames Wright  * @param[in]     stg_inflow_path Path to STGInflow.dat file
164ba6664aeSJames Wright  * @param[in]     stg_rand_path   Path to STGRand.dat file
16505675bc2SJames Wright  * @param[in,out] stg_ctx         Pointer to STGShur14Context where the data will be loaded into
166ba6664aeSJames Wright  */
1672b730f8bSJeremy 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],
168*f8839eb4SJames Wright                                  STGShur14Context *stg_ctx) {
169ba6664aeSJames Wright   PetscInt nmodes, nprofs;
170ba6664aeSJames Wright   PetscFunctionBeginUser;
171ba6664aeSJames Wright 
172ba6664aeSJames Wright   // Get options
1738b892a05SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_rand_path, &nmodes));
1748b892a05SJames Wright   PetscCall(PHASTADatFileGetNRows(comm, stg_inflow_path, &nprofs));
1750e654f56SJames Wright   PetscCheck(nmodes < STG_NMODES_MAX, comm, PETSC_ERR_SUP,
176457e73b2SJames Wright              "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%d). Change size of STG_NMODES_MAX and recompile", stg_rand_path,
177457e73b2SJames Wright              nmodes, STG_NMODES_MAX);
178ba6664aeSJames Wright 
179ba6664aeSJames Wright   {
18005675bc2SJames Wright     STGShur14Context temp_ctx;
18105675bc2SJames Wright     PetscCall(PetscCalloc1(1, &temp_ctx));
18205675bc2SJames Wright     *temp_ctx                   = **stg_ctx;
18305675bc2SJames Wright     temp_ctx->nmodes            = nmodes;
18405675bc2SJames Wright     temp_ctx->nprofs            = nprofs;
18505675bc2SJames Wright     temp_ctx->offsets.sigma     = 0;
18605675bc2SJames Wright     temp_ctx->offsets.d         = nmodes * 3;
18705675bc2SJames Wright     temp_ctx->offsets.phi       = temp_ctx->offsets.d + nmodes * 3;
18805675bc2SJames Wright     temp_ctx->offsets.kappa     = temp_ctx->offsets.phi + nmodes;
18905675bc2SJames Wright     temp_ctx->offsets.wall_dist = temp_ctx->offsets.kappa + nmodes;
19005675bc2SJames Wright     temp_ctx->offsets.ubar      = temp_ctx->offsets.wall_dist + nprofs;
19105675bc2SJames Wright     temp_ctx->offsets.cij       = temp_ctx->offsets.ubar + nprofs * 3;
19205675bc2SJames Wright     temp_ctx->offsets.eps       = temp_ctx->offsets.cij + nprofs * 6;
19305675bc2SJames Wright     temp_ctx->offsets.lt        = temp_ctx->offsets.eps + nprofs;
194*f8839eb4SJames Wright     PetscInt total_num_scalars  = temp_ctx->offsets.lt + nprofs;
19505675bc2SJames Wright     temp_ctx->total_bytes       = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]);
19605675bc2SJames Wright     PetscCall(PetscFree(*stg_ctx));
19705675bc2SJames Wright     PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx));
19805675bc2SJames Wright     **stg_ctx = *temp_ctx;
19905675bc2SJames Wright     PetscCall(PetscFree(temp_ctx));
200ba6664aeSJames Wright   }
201ba6664aeSJames Wright 
20205675bc2SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, *stg_ctx));
20305675bc2SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, *stg_ctx));
204ba6664aeSJames Wright 
20505675bc2SJames Wright   {  // -- Calculate kappa
20605675bc2SJames Wright     CeedScalar *kappa     = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa];
20705675bc2SJames Wright     CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist];
20805675bc2SJames Wright     CeedScalar *lt        = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt];
209ba6664aeSJames Wright     CeedScalar  le, le_max = 0;
210ba6664aeSJames Wright 
21105675bc2SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) {
212175f00a6SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
213ba6664aeSJames Wright       if (le_max < le) le_max = le;
214ba6664aeSJames Wright     }
215ba6664aeSJames Wright     CeedScalar kmin = M_PI / le_max;
216ba6664aeSJames Wright 
21705675bc2SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); }
21805675bc2SJames Wright   }
219ba6664aeSJames Wright 
220ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
221ba6664aeSJames Wright }
222ba6664aeSJames Wright 
2232b730f8bSJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
224*f8839eb4SJames Wright                         const CeedScalar P0) {
225a424bcd0SJames Wright   Ceed                     ceed                                = user->ceed;
226ba6664aeSJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
227ba6664aeSJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
22873557d35SJames Wright   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE, given_stg_dx = PETSC_FALSE;
22973557d35SJames Wright   CeedScalar               u0 = 0.0, alpha = 1.01, stg_dx = 1.0e-3;
230ba6664aeSJames Wright   CeedQFunctionContext     stg_context;
231ba6664aeSJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
232ba6664aeSJames Wright   PetscFunctionBeginUser;
233ba6664aeSJames Wright 
234ba6664aeSJames Wright   // Get options
235ba6664aeSJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
2362b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
2372b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
2382b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
2392b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
2402b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
2412b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
2422b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
2432b730f8bSJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
24473557d35SJames Wright   PetscCall(PetscOptionsReal("-stg_dx", "Element size in streamwise direction at inflow", NULL, stg_dx, &stg_dx, &given_stg_dx));
245ba6664aeSJames Wright   PetscOptionsEnd();
246ba6664aeSJames Wright 
24724941a69SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
24865dd5cafSJames Wright   global_stg_ctx->alpha              = alpha;
24965dd5cafSJames Wright   global_stg_ctx->u0                 = u0;
25065dd5cafSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
25165dd5cafSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
25265dd5cafSJames Wright   global_stg_ctx->mean_only          = mean_only;
25389060322SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
25465dd5cafSJames Wright   global_stg_ctx->theta0             = theta0;
25565dd5cafSJames Wright   global_stg_ctx->P0                 = P0;
256ba6664aeSJames Wright 
257ba6664aeSJames Wright   {
258ba6664aeSJames Wright     // Calculate dx assuming constant spacing
259ba6664aeSJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
26024941a69SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
261ba6664aeSJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
262ba6664aeSJames Wright 
263ba6664aeSJames Wright     PetscInt nmax = 3, faces[3];
2642b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
26573557d35SJames Wright     global_stg_ctx->dx = given_stg_dx ? stg_dx : domain_size[0] / faces[0];
266ba6664aeSJames Wright   }
267ba6664aeSJames Wright 
268a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
26965dd5cafSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
270a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
271ba6664aeSJames Wright 
272*f8839eb4SJames Wright   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx));
273ba6664aeSJames Wright 
274a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &stg_context));
275a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx));
276a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc));
277a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1,
278a424bcd0SJames Wright                                                          "Physical time of the solution"));
279ba6664aeSJames Wright 
280a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context));
281b77c53c9SJames Wright   problem->ics.qfunction         = ICsSTG;
282b77c53c9SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
283b77c53c9SJames Wright   problem->ics.qfunction_context = stg_context;
284b77c53c9SJames Wright 
2854e139266SJames Wright   if (use_stgstrong) {
28665dd5cafSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
2873722cd23SJames Wright     problem->use_strong_bc_ceed = PETSC_TRUE;
28836b31e27SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
2894e139266SJames Wright   } else {
290ba6664aeSJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
291ba6664aeSJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
2924dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
2934dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
294a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context));
295a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context));
29636b31e27SJames Wright     problem->bc_from_ics = PETSC_TRUE;
2974e139266SJames Wright   }
298ba6664aeSJames Wright 
299ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
300ba6664aeSJames Wright }
3014e139266SJames Wright 
302*f8839eb4SJames Wright // @brief Set STG strongly enforce components using DMAddBoundary
3032b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
3044e139266SJames Wright   DMLabel label;
3054e139266SJames Wright   PetscFunctionBeginUser;
3064e139266SJames Wright 
307192a7459SJames Wright   PetscInt comps[5], num_comps = 4;
30897baf651SJames Wright   switch (phys->state_var) {
30997baf651SJames Wright     case STATEVAR_CONSERVATIVE:
310192a7459SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
311192a7459SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
31297baf651SJames Wright       break;
31397baf651SJames Wright 
31497baf651SJames Wright     case STATEVAR_PRIMITIVE:
31597baf651SJames Wright       // {1,2,3,4} for u, v, w, T
31697baf651SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
31797baf651SJames Wright       break;
318192a7459SJames Wright   }
319192a7459SJames Wright 
32024941a69SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
3214e139266SJames Wright   if (bc->num_inflow > 0) {
322*f8839eb4SJames Wright     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, NULL, NULL, NULL, NULL));
3234e139266SJames Wright   }
3244e139266SJames Wright 
325ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3264e139266SJames Wright }
327dada6cc0SJames Wright 
3282b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
32905675bc2SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *qf_strongbc) {
330dada6cc0SJames Wright   PetscFunctionBeginUser;
331a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, qf_strongbc));
332a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
333a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE));
334a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE));
335a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE));
336a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE));
337dada6cc0SJames Wright 
338a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context));
339ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
340dada6cc0SJames Wright }
3415dc40723SJames Wright 
3422b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
34305675bc2SJames Wright                                             CeedQFunction *qf_strongbc) {
3445dc40723SJames Wright   PetscFunctionBeginUser;
345a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, qf_strongbc));
346a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
347a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE));
348a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE));
3495dc40723SJames Wright 
350a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context));
351ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3525dc40723SJames Wright }
353