xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision ee4ca9cbfe2be39196684117442f3ce8d00b6506)
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   }
48*ee4ca9cbSJames 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) {
64ba6664aeSJames Wright   PetscInt       ndims, dims[2];
65ba6664aeSJames Wright   FILE          *fp;
66ba6664aeSJames Wright   const PetscInt char_array_len = 512;
67ba6664aeSJames Wright   char           line[char_array_len];
68ba6664aeSJames Wright   char         **array;
69ba6664aeSJames Wright 
70ba6664aeSJames Wright   PetscFunctionBeginUser;
71ba6664aeSJames Wright 
728b892a05SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
73ba6664aeSJames Wright 
74ba6664aeSJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
75175f00a6SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
76ba6664aeSJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
77ba6664aeSJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
782b730f8bSJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
79ba6664aeSJames Wright 
80ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
8124941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
8224941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
830e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
840a37e1beSJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
850a37e1beSJames Wright                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 
108*ee4ca9cbSJames 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) {
122ba6664aeSJames Wright   PetscInt       ndims, dims[2];
123ba6664aeSJames Wright   FILE          *fp;
124ba6664aeSJames Wright   const PetscInt char_array_len = 512;
125ba6664aeSJames Wright   char           line[char_array_len];
126ba6664aeSJames Wright   char         **array;
127ba6664aeSJames Wright 
128ba6664aeSJames Wright   PetscFunctionBeginUser;
1298b892a05SJames Wright   PetscCall(PHASTADatFileOpen(comm, path, char_array_len, dims, &fp));
130ba6664aeSJames Wright 
131ba6664aeSJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
1322b730f8bSJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
1332b730f8bSJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
134ba6664aeSJames Wright 
135ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
13624941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
13724941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
1380e654f56SJames Wright     PetscCheck(ndims == dims[1], comm, PETSC_ERR_FILE_UNEXPECTED,
1390a37e1beSJames Wright                "Line %" PetscInt_FMT " of %s does not have correct number of columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
1400a37e1beSJames Wright                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 
152*ee4ca9cbSJames 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],
16805675bc2SJames Wright                                  STGShur14Context *stg_ctx, const CeedScalar ynodes[]) {
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,
1760e654f56SJames Wright              "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT "). Change size of STG_NMODES_MAX and recompile",
1772b730f8bSJeremy L Thompson              stg_rand_path, 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;
19405675bc2SJames Wright     temp_ctx->offsets.ynodes    = temp_ctx->offsets.lt + nprofs;
19505675bc2SJames Wright     PetscInt total_num_scalars  = temp_ctx->offsets.ynodes + temp_ctx->nynodes;
19605675bc2SJames Wright     temp_ctx->total_bytes       = sizeof(*temp_ctx) + total_num_scalars * sizeof(temp_ctx->data[0]);
19705675bc2SJames Wright     PetscCall(PetscFree(*stg_ctx));
19805675bc2SJames Wright     PetscCall(PetscMalloc(temp_ctx->total_bytes, stg_ctx));
19905675bc2SJames Wright     **stg_ctx = *temp_ctx;
20005675bc2SJames Wright     PetscCall(PetscFree(temp_ctx));
201ba6664aeSJames Wright   }
202ba6664aeSJames Wright 
20305675bc2SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, *stg_ctx));
20405675bc2SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, *stg_ctx));
205ba6664aeSJames Wright 
20605675bc2SJames Wright   if ((*stg_ctx)->nynodes > 0) {
20705675bc2SJames Wright     CeedScalar *ynodes_ctx = &(*stg_ctx)->data[(*stg_ctx)->offsets.ynodes];
20805675bc2SJames Wright     for (PetscInt i = 0; i < (*stg_ctx)->nynodes; i++) ynodes_ctx[i] = ynodes[i];
2094e139266SJames Wright   }
2104e139266SJames Wright 
21105675bc2SJames Wright   {  // -- Calculate kappa
21205675bc2SJames Wright     CeedScalar *kappa     = &(*stg_ctx)->data[(*stg_ctx)->offsets.kappa];
21305675bc2SJames Wright     CeedScalar *wall_dist = &(*stg_ctx)->data[(*stg_ctx)->offsets.wall_dist];
21405675bc2SJames Wright     CeedScalar *lt        = &(*stg_ctx)->data[(*stg_ctx)->offsets.lt];
215ba6664aeSJames Wright     CeedScalar  le, le_max = 0;
216ba6664aeSJames Wright 
21705675bc2SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nprofs; i++) {
218175f00a6SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
219ba6664aeSJames Wright       if (le_max < le) le_max = le;
220ba6664aeSJames Wright     }
221ba6664aeSJames Wright     CeedScalar kmin = M_PI / le_max;
222ba6664aeSJames Wright 
22305675bc2SJames Wright     CeedPragmaSIMD for (PetscInt i = 0; i < (*stg_ctx)->nmodes; i++) { kappa[i] = kmin * pow((*stg_ctx)->alpha, i); }
22405675bc2SJames Wright   }
225ba6664aeSJames Wright 
226*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
227ba6664aeSJames Wright }
228ba6664aeSJames Wright 
2292b730f8bSJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
2302b730f8bSJeremy L Thompson                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
231ba6664aeSJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
232ba6664aeSJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
2332b730f8bSJeremy L Thompson   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
2342b730f8bSJeremy L Thompson   CeedScalar               u0 = 0.0, alpha = 1.01;
235ba6664aeSJames Wright   CeedQFunctionContext     stg_context;
236ba6664aeSJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
237ba6664aeSJames Wright   PetscFunctionBeginUser;
238ba6664aeSJames Wright 
239ba6664aeSJames Wright   // Get options
240ba6664aeSJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
2412b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
2422b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
2432b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
2442b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
2452b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
2462b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
2472b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
2482b730f8bSJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
249ba6664aeSJames Wright   PetscOptionsEnd();
250ba6664aeSJames Wright 
25124941a69SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
25265dd5cafSJames Wright   global_stg_ctx->alpha              = alpha;
25365dd5cafSJames Wright   global_stg_ctx->u0                 = u0;
25465dd5cafSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
25565dd5cafSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
25665dd5cafSJames Wright   global_stg_ctx->mean_only          = mean_only;
25789060322SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
25865dd5cafSJames Wright   global_stg_ctx->theta0             = theta0;
25965dd5cafSJames Wright   global_stg_ctx->P0                 = P0;
26065dd5cafSJames Wright   global_stg_ctx->nynodes            = nynodes;
261ba6664aeSJames Wright 
262ba6664aeSJames Wright   {
263ba6664aeSJames Wright     // Calculate dx assuming constant spacing
264ba6664aeSJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
26524941a69SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
266ba6664aeSJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
267ba6664aeSJames Wright 
268ba6664aeSJames Wright     PetscInt nmax = 3, faces[3];
2692b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
27065dd5cafSJames Wright     global_stg_ctx->dx = domain_size[0] / faces[0];
27165dd5cafSJames Wright     global_stg_ctx->dz = domain_size[2] / faces[2];
272ba6664aeSJames Wright   }
273ba6664aeSJames Wright 
2742b730f8bSJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
27565dd5cafSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
2762b730f8bSJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
277ba6664aeSJames Wright 
2782b730f8bSJeremy L Thompson   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
279ba6664aeSJames Wright 
280ba6664aeSJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
2812b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
2822b730f8bSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
2832b730f8bSJeremy L Thompson   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
284ba6664aeSJames Wright 
285b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
286b77c53c9SJames Wright   problem->ics.qfunction         = ICsSTG;
287b77c53c9SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
288b77c53c9SJames Wright   problem->ics.qfunction_context = stg_context;
289b77c53c9SJames Wright 
2904e139266SJames Wright   if (use_stgstrong) {
29165dd5cafSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
2923722cd23SJames Wright     problem->use_strong_bc_ceed = PETSC_TRUE;
29336b31e27SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
2944e139266SJames Wright   } else {
295ba6664aeSJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
296ba6664aeSJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
2974dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
2984dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
2992b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
3002b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
30136b31e27SJames Wright     problem->bc_from_ics = PETSC_TRUE;
3024e139266SJames Wright   }
303ba6664aeSJames Wright 
304*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
305ba6664aeSJames Wright }
3064e139266SJames Wright 
3072b730f8bSJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
3084e139266SJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
3094e139266SJames Wright   // ^^assuming min(dy) is first element off the wall
3104e139266SJames Wright   PetscInt idx = -1;  // Index
3114e139266SJames Wright 
3124e139266SJames Wright   for (PetscInt i = 0; i < nynodes; i++) {
3134e139266SJames Wright     if (y < ynodes[i] + half_mindy) {
3142b730f8bSJeremy L Thompson       idx = i;
3152b730f8bSJeremy L Thompson       break;
3164e139266SJames Wright     }
3174e139266SJames Wright   }
3184e139266SJames Wright   if (idx == 0) return ynodes[1] - ynodes[0];
3194e139266SJames Wright   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
3204e139266SJames Wright   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
3214e139266SJames Wright }
3224e139266SJames Wright 
3234e139266SJames Wright // Function passed to DMAddBoundary
324dada6cc0SJames Wright // NOTE: Not used in favor of QFunction-based method
3252b730f8bSJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
3264e139266SJames Wright   PetscFunctionBeginUser;
3272b730f8bSJeremy L Thompson 
3284e139266SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
3294e139266SJames Wright   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
3304e139266SJames Wright   const bool             mean_only = stg_ctx->mean_only;
3314e139266SJames Wright   const PetscScalar      dx        = stg_ctx->dx;
3324e139266SJames Wright   const PetscScalar      dz        = stg_ctx->dz;
3334e139266SJames Wright   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
3344e139266SJames Wright   const PetscScalar      theta0    = stg_ctx->theta0;
3354e139266SJames Wright   const PetscScalar      P0        = stg_ctx->P0;
3364e139266SJames Wright   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
3374e139266SJames Wright   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
3384e139266SJames Wright   const PetscScalar      Rd        = cp - cv;
3394e139266SJames Wright 
3404e139266SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
3414e139266SJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
3424e139266SJames Wright   if (!mean_only) {
3434e139266SJames Wright     const PetscInt     nynodes = stg_ctx->nynodes;
3444e139266SJames Wright     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
3454e139266SJames Wright     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
3464e139266SJames Wright     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
3474e139266SJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
3484e139266SJames Wright   } else {
3494e139266SJames Wright     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
3504e139266SJames Wright   }
3514e139266SJames Wright 
3524e139266SJames Wright   bcval[0] = rho;
3534e139266SJames Wright   bcval[1] = rho * u[0];
3544e139266SJames Wright   bcval[2] = rho * u[1];
3554e139266SJames Wright   bcval[3] = rho * u[2];
356*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3574e139266SJames Wright }
3584e139266SJames Wright 
3592b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
3604e139266SJames Wright   DMLabel label;
3614e139266SJames Wright   PetscFunctionBeginUser;
3624e139266SJames Wright 
363192a7459SJames Wright   PetscInt comps[5], num_comps = 4;
36497baf651SJames Wright   switch (phys->state_var) {
36597baf651SJames Wright     case STATEVAR_CONSERVATIVE:
366192a7459SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
367192a7459SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
36897baf651SJames Wright       break;
36997baf651SJames Wright 
37097baf651SJames Wright     case STATEVAR_PRIMITIVE:
37197baf651SJames Wright       // {1,2,3,4} for u, v, w, T
37297baf651SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
37397baf651SJames Wright       break;
374192a7459SJames Wright   }
375192a7459SJames Wright 
37624941a69SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
3774e139266SJames Wright   // Set wall BCs
3784e139266SJames Wright   if (bc->num_inflow > 0) {
3792b730f8bSJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
38024941a69SJames Wright                             NULL, global_stg_ctx, NULL));
3814e139266SJames Wright   }
3824e139266SJames Wright 
383*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3844e139266SJames Wright }
385dada6cc0SJames Wright 
3862b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
38705675bc2SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *qf_strongbc) {
388dada6cc0SJames Wright   PetscFunctionBeginUser;
38905675bc2SJames Wright   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, qf_strongbc);
39005675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
39105675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
39205675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "scale", 1, CEED_EVAL_NONE);
39305675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
39405675bc2SJames Wright   CeedQFunctionAddOutput(*qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
395dada6cc0SJames Wright 
39605675bc2SJames Wright   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
397*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
398dada6cc0SJames Wright }
3995dc40723SJames Wright 
4002b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
40105675bc2SJames Wright                                             CeedQFunction *qf_strongbc) {
4025dc40723SJames Wright   PetscFunctionBeginUser;
40305675bc2SJames Wright   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, qf_strongbc);
40405675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
40505675bc2SJames Wright   CeedQFunctionAddInput(*qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
40605675bc2SJames Wright   CeedQFunctionAddOutput(*qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4075dc40723SJames Wright 
40805675bc2SJames Wright   CeedQFunctionSetContext(*qf_strongbc, problem->ics.qfunction_context);
409*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4105dc40723SJames Wright }
411