xref: /honee/problems/stg_shur14.c (revision 2b916ea7fa53b5c2584160b9274b1b14ca18ff4f)
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 
12*2b916ea7SJeremy L Thompson #include "stg_shur14.h"
13*2b916ea7SJeremy L Thompson 
14493642f1SJames Wright #include <math.h>
15493642f1SJames Wright #include <petsc.h>
16*2b916ea7SJeremy L Thompson #include <stdlib.h>
17*2b916ea7SJeremy L Thompson 
18493642f1SJames Wright #include "../navierstokes.h"
19493642f1SJames Wright #include "../qfunctions/stg_shur14.h"
20493642f1SJames Wright 
218085925cSJames Wright STGShur14Context global_stg_ctx;
228085925cSJames Wright 
23493642f1SJames Wright /*
24493642f1SJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
25493642f1SJames Wright  *
26493642f1SJames Wright  * This assumes the input matrices are in order [11,22,33,12,13,23]. This
27493642f1SJames Wright  * format is also used for the output.
28493642f1SJames Wright  *
29493642f1SJames Wright  * @param[in]  comm   MPI_Comm
30493642f1SJames Wright  * @param[in]  nprofs Number of matrices in Rij
31493642f1SJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
32493642f1SJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
33493642f1SJames Wright  */
34*2b916ea7SJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
35493642f1SJames Wright   PetscFunctionBeginUser;
36493642f1SJames Wright   for (PetscInt i = 0; i < nprofs; i++) {
37493642f1SJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
38493642f1SJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
39493642f1SJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2));
40493642f1SJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
41493642f1SJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i];
42493642f1SJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
43493642f1SJames Wright 
44*2b916ea7SJeremy L Thompson     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) {
45*2b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.",
46*2b916ea7SJeremy L Thompson               i + 1);
47*2b916ea7SJeremy L Thompson     }
48493642f1SJames Wright   }
49493642f1SJames Wright   PetscFunctionReturn(0);
50493642f1SJames Wright }
51493642f1SJames Wright 
52493642f1SJames Wright /*
53493642f1SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
54493642f1SJames Wright  *
55493642f1SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and
56493642f1SJames Wright  * passes the file pointer in `fp`. It is not closed in this function, thus
57493642f1SJames Wright  * `fp` must be closed sometime after this function has been called (using
58493642f1SJames Wright  * `PetscFClose` for example).
59493642f1SJames Wright  *
60493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
61493642f1SJames Wright  * as the only two entries, separated by a single space
62493642f1SJames Wright  *
63493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
64493642f1SJames Wright  * @param[in] path Path to the file
65493642f1SJames Wright  * @param[in] char_array_len Length of the character array that should contain each line
66493642f1SJames Wright  * @param[out] dims Dimensions of the file, taken from the first line of the file
67493642f1SJames Wright  * @param[out] fp File pointer to the opened file
68493642f1SJames Wright  */
69*2b916ea7SJeremy L Thompson static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
70*2b916ea7SJeremy L Thompson                                         FILE **fp) {
71493642f1SJames Wright   PetscInt ndims;
72493642f1SJames Wright   char     line[char_array_len];
73493642f1SJames Wright   char   **array;
74493642f1SJames Wright 
75493642f1SJames Wright   PetscFunctionBeginUser;
7634b5deb1SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
7734b5deb1SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
7834b5deb1SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
79*2b916ea7SJeremy L Thompson   if (ndims != 2) {
80*2b916ea7SJeremy L Thompson     SETERRQ(comm, -1, "Found %" PetscInt_FMT " dimensions instead of 2 on the first line of %s", ndims, path);
81*2b916ea7SJeremy L Thompson   }
82493642f1SJames Wright 
83493642f1SJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
8434b5deb1SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
85*2b916ea7SJeremy L Thompson 
86493642f1SJames Wright   PetscFunctionReturn(0);
87493642f1SJames Wright }
88493642f1SJames Wright 
89493642f1SJames Wright /*
90493642f1SJames Wright  * @brief Get the number of rows for the PHASTA file at path
91493642f1SJames Wright  *
92493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
93493642f1SJames Wright  * as the only two entries, separated by a single space
94493642f1SJames Wright  *
95493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
96493642f1SJames Wright  * @param[in] path Path to the file
97493642f1SJames Wright  * @param[out] nrows Number of rows
98493642f1SJames Wright  */
99*2b916ea7SJeremy L Thompson static PetscErrorCode GetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
100493642f1SJames Wright   const PetscInt char_array_len = 512;
101493642f1SJames Wright   PetscInt       dims[2];
102493642f1SJames Wright   FILE          *fp;
103493642f1SJames Wright 
104493642f1SJames Wright   PetscFunctionBeginUser;
10534b5deb1SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
106493642f1SJames Wright   *nrows = dims[0];
10734b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
108*2b916ea7SJeremy L Thompson 
109493642f1SJames Wright   PetscFunctionReturn(0);
110493642f1SJames Wright }
111493642f1SJames Wright 
112493642f1SJames Wright /*
113493642f1SJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
114493642f1SJames Wright  *
115493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
116493642f1SJames Wright  * as the only two entries, separated by a single space.
117493642f1SJames Wright  * Assumes there are 14 columns in the file
118493642f1SJames Wright  *
119493642f1SJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
120493642f1SJames Wright  * profile found in the file
121493642f1SJames Wright  *
122493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
123493642f1SJames Wright  * @param[in] path Path to the STGInflow.dat file
124493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
125493642f1SJames Wright  */
126*2b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
127493642f1SJames Wright   PetscInt       ndims, dims[2];
128493642f1SJames Wright   FILE          *fp;
129493642f1SJames Wright   const PetscInt char_array_len = 512;
130493642f1SJames Wright   char           line[char_array_len];
131493642f1SJames Wright   char         **array;
132493642f1SJames Wright 
133493642f1SJames Wright   PetscFunctionBeginUser;
134493642f1SJames Wright 
13534b5deb1SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
136493642f1SJames Wright 
137493642f1SJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
138c77f3192SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
139493642f1SJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
140493642f1SJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
141*2b916ea7SJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
142493642f1SJames Wright 
143493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
14434b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
14534b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
146*2b916ea7SJeremy L Thompson     if (ndims < dims[1]) {
147*2b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
148*2b916ea7SJeremy L Thompson               ndims, dims[1]);
149*2b916ea7SJeremy L Thompson     }
150493642f1SJames Wright 
151c77f3192SJames Wright     wall_dist[i] = (CeedScalar)atof(array[0]);
152493642f1SJames Wright     ubar[0][i]   = (CeedScalar)atof(array[1]);
153493642f1SJames Wright     ubar[1][i]   = (CeedScalar)atof(array[2]);
154493642f1SJames Wright     ubar[2][i]   = (CeedScalar)atof(array[3]);
155493642f1SJames Wright     rij[0][i]    = (CeedScalar)atof(array[4]);
156493642f1SJames Wright     rij[1][i]    = (CeedScalar)atof(array[5]);
157493642f1SJames Wright     rij[2][i]    = (CeedScalar)atof(array[6]);
158493642f1SJames Wright     rij[3][i]    = (CeedScalar)atof(array[7]);
159493642f1SJames Wright     rij[4][i]    = (CeedScalar)atof(array[8]);
160493642f1SJames Wright     rij[5][i]    = (CeedScalar)atof(array[9]);
161493642f1SJames Wright     lt[i]        = (CeedScalar)atof(array[12]);
162493642f1SJames Wright     eps[i]       = (CeedScalar)atof(array[13]);
163493642f1SJames Wright 
164*2b916ea7SJeremy L Thompson     if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path);
165*2b916ea7SJeremy L Thompson     if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path);
166*2b916ea7SJeremy L Thompson     if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path);
167493642f1SJames Wright   }
168*2b916ea7SJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
16934b5deb1SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
17034b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
171*2b916ea7SJeremy L Thompson 
172493642f1SJames Wright   PetscFunctionReturn(0);
173493642f1SJames Wright }
174493642f1SJames Wright 
175493642f1SJames Wright /*
176493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
177493642f1SJames Wright  *
178493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
179493642f1SJames Wright  * as the only two entries, separated by a single space.
180493642f1SJames Wright  * Assumes there are 7 columns in the file
181493642f1SJames Wright  *
182493642f1SJames Wright  * @param[in]    comm    MPI_Comm for the program
183493642f1SJames Wright  * @param[in]    path    Path to the STGRand.dat file
184493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
185493642f1SJames Wright  */
186*2b916ea7SJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
187493642f1SJames Wright   PetscInt       ndims, dims[2];
188493642f1SJames Wright   FILE          *fp;
189493642f1SJames Wright   const PetscInt char_array_len = 512;
190493642f1SJames Wright   char           line[char_array_len];
191493642f1SJames Wright   char         **array;
192493642f1SJames Wright 
193493642f1SJames Wright   PetscFunctionBeginUser;
19434b5deb1SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
195493642f1SJames Wright 
196493642f1SJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
197*2b916ea7SJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
198*2b916ea7SJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
199493642f1SJames Wright 
200493642f1SJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
20134b5deb1SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
20234b5deb1SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
203*2b916ea7SJeremy L Thompson     if (ndims < dims[1]) {
204*2b916ea7SJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
205*2b916ea7SJeremy L Thompson               ndims, dims[1]);
206*2b916ea7SJeremy L Thompson     }
207493642f1SJames Wright 
208493642f1SJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
209493642f1SJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
210493642f1SJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
211493642f1SJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
212493642f1SJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
213493642f1SJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
214493642f1SJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
215493642f1SJames Wright   }
21634b5deb1SJames Wright   PetscCall(PetscFClose(comm, fp));
217*2b916ea7SJeremy L Thompson 
218493642f1SJames Wright   PetscFunctionReturn(0);
219493642f1SJames Wright }
220493642f1SJames Wright 
221493642f1SJames Wright /*
222493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
223493642f1SJames Wright  *
224493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
225493642f1SJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
226493642f1SJames Wright  * STGShur14Context instance.
227493642f1SJames Wright  *
228493642f1SJames Wright  * @param[in]    comm            MPI_Comm for the program
229493642f1SJames Wright  * @param[in]    dm              DM for the program
230493642f1SJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
231493642f1SJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
232493642f1SJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
233493642f1SJames Wright  */
234*2b916ea7SJeremy 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],
235*2b916ea7SJeremy L Thompson                                  STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) {
236493642f1SJames Wright   PetscInt         nmodes, nprofs;
237493642f1SJames Wright   STGShur14Context stg_ctx;
238493642f1SJames Wright   PetscFunctionBeginUser;
239493642f1SJames Wright 
240493642f1SJames Wright   // Get options
24134b5deb1SJames Wright   PetscCall(GetNRows(comm, stg_rand_path, &nmodes));
24234b5deb1SJames Wright   PetscCall(GetNRows(comm, stg_inflow_path, &nprofs));
243493642f1SJames Wright   if (nmodes > STG_NMODES_MAX)
244*2b916ea7SJeremy L Thompson     SETERRQ(comm, 1,
245*2b916ea7SJeremy L Thompson             "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT
246*2b916ea7SJeremy L Thompson             "). "
247*2b916ea7SJeremy L Thompson             "Change size of STG_NMODES_MAX and recompile",
248*2b916ea7SJeremy L Thompson             stg_rand_path, nmodes, STG_NMODES_MAX);
249493642f1SJames Wright 
250493642f1SJames Wright   {
251493642f1SJames Wright     STGShur14Context s;
25234b5deb1SJames Wright     PetscCall(PetscCalloc1(1, &s));
253493642f1SJames Wright     *s                         = **pstg_ctx;
254493642f1SJames Wright     s->nmodes                  = nmodes;
255493642f1SJames Wright     s->nprofs                  = nprofs;
256493642f1SJames Wright     s->offsets.sigma           = 0;
257493642f1SJames Wright     s->offsets.d               = nmodes * 3;
258493642f1SJames Wright     s->offsets.phi             = s->offsets.d + nmodes * 3;
259493642f1SJames Wright     s->offsets.kappa           = s->offsets.phi + nmodes;
260c77f3192SJames Wright     s->offsets.wall_dist       = s->offsets.kappa + nmodes;
261c77f3192SJames Wright     s->offsets.ubar            = s->offsets.wall_dist + nprofs;
262493642f1SJames Wright     s->offsets.cij             = s->offsets.ubar + nprofs * 3;
263493642f1SJames Wright     s->offsets.eps             = s->offsets.cij + nprofs * 6;
264493642f1SJames Wright     s->offsets.lt              = s->offsets.eps + nprofs;
265e6098bcdSJames Wright     s->offsets.ynodes          = s->offsets.lt + nprofs;
266e6098bcdSJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
267493642f1SJames Wright     s->total_bytes             = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]);
26834b5deb1SJames Wright     PetscCall(PetscMalloc(s->total_bytes, &stg_ctx));
269493642f1SJames Wright     *stg_ctx = *s;
27034b5deb1SJames Wright     PetscCall(PetscFree(s));
271493642f1SJames Wright   }
272493642f1SJames Wright 
27334b5deb1SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx));
27434b5deb1SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx));
275493642f1SJames Wright 
276e6098bcdSJames Wright   if (stg_ctx->nynodes > 0) {
277e6098bcdSJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
278e6098bcdSJames Wright     for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
279e6098bcdSJames Wright   }
280e6098bcdSJames Wright 
281493642f1SJames Wright   // -- Calculate kappa
282493642f1SJames Wright   {
283493642f1SJames Wright     CeedScalar *kappa     = &stg_ctx->data[stg_ctx->offsets.kappa];
284c77f3192SJames Wright     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
285493642f1SJames Wright     CeedScalar *lt        = &stg_ctx->data[stg_ctx->offsets.lt];
286493642f1SJames Wright     CeedScalar  le, le_max = 0;
287493642f1SJames Wright 
288*2b916ea7SJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
289c77f3192SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
290493642f1SJames Wright       if (le_max < le) le_max = le;
291493642f1SJames Wright     }
292493642f1SJames Wright     CeedScalar kmin = M_PI / le_max;
293493642f1SJames Wright 
294*2b916ea7SJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); }
295493642f1SJames Wright   }  // end calculate kappa
296493642f1SJames Wright 
29734b5deb1SJames Wright   PetscCall(PetscFree(*pstg_ctx));
298493642f1SJames Wright   *pstg_ctx = stg_ctx;
299493642f1SJames Wright   PetscFunctionReturn(0);
300493642f1SJames Wright }
301493642f1SJames Wright 
302*2b916ea7SJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
303*2b916ea7SJeremy L Thompson                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
304493642f1SJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
305493642f1SJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
306*2b916ea7SJeremy L Thompson   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
307*2b916ea7SJeremy L Thompson   CeedScalar               u0 = 0.0, alpha = 1.01;
308493642f1SJames Wright   CeedQFunctionContext     stg_context;
309493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
310493642f1SJames Wright   PetscFunctionBeginUser;
311493642f1SJames Wright 
312493642f1SJames Wright   // Get options
313493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
314*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
315*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
316*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
317*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
318*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
319*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
320*2b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
321*2b916ea7SJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
322493642f1SJames Wright   PetscOptionsEnd();
323493642f1SJames Wright 
32434b5deb1SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
3258085925cSJames Wright   global_stg_ctx->alpha              = alpha;
3268085925cSJames Wright   global_stg_ctx->u0                 = u0;
3278085925cSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
3288085925cSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
3298085925cSJames Wright   global_stg_ctx->mean_only          = mean_only;
330d4e0f297SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
3318085925cSJames Wright   global_stg_ctx->theta0             = theta0;
3328085925cSJames Wright   global_stg_ctx->P0                 = P0;
3338085925cSJames Wright   global_stg_ctx->nynodes            = nynodes;
334493642f1SJames Wright 
335493642f1SJames Wright   {
336493642f1SJames Wright     // Calculate dx assuming constant spacing
337493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
33834b5deb1SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
339493642f1SJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
340493642f1SJames Wright 
341493642f1SJames Wright     PetscInt nmax = 3, faces[3];
342*2b916ea7SJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
3438085925cSJames Wright     global_stg_ctx->dx = domain_size[0] / faces[0];
3448085925cSJames Wright     global_stg_ctx->dz = domain_size[2] / faces[2];
345493642f1SJames Wright   }
346493642f1SJames Wright 
347*2b916ea7SJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
3488085925cSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
349*2b916ea7SJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
350493642f1SJames Wright 
351*2b916ea7SJeremy L Thompson   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
352493642f1SJames Wright 
353493642f1SJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
354*2b916ea7SJeremy L Thompson   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
355*2b916ea7SJeremy L Thompson   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
356*2b916ea7SJeremy L Thompson   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
357493642f1SJames Wright 
35843bff553SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
35943bff553SJames Wright   problem->ics.qfunction         = ICsSTG;
36043bff553SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
36143bff553SJames Wright   problem->ics.qfunction_context = stg_context;
36243bff553SJames Wright 
363e6098bcdSJames Wright   if (use_stgstrong) {
3648085925cSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
3656d0190e2SJames Wright     problem->use_dirichlet_ceed = PETSC_TRUE;
366d7244455SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
367e6098bcdSJames Wright   } else {
368493642f1SJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
369493642f1SJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
370a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
371a6e8f989SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
372*2b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
373*2b916ea7SJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
374d7244455SJames Wright     problem->bc_from_ics = PETSC_TRUE;
375e6098bcdSJames Wright   }
376493642f1SJames Wright 
377493642f1SJames Wright   PetscFunctionReturn(0);
378493642f1SJames Wright }
379e6098bcdSJames Wright 
380*2b916ea7SJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
381e6098bcdSJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
382e6098bcdSJames Wright   // ^^assuming min(dy) is first element off the wall
383e6098bcdSJames Wright   PetscInt idx = -1;  // Index
384e6098bcdSJames Wright 
385e6098bcdSJames Wright   for (PetscInt i = 0; i < nynodes; i++) {
386e6098bcdSJames Wright     if (y < ynodes[i] + half_mindy) {
387*2b916ea7SJeremy L Thompson       idx = i;
388*2b916ea7SJeremy L Thompson       break;
389e6098bcdSJames Wright     }
390e6098bcdSJames Wright   }
391e6098bcdSJames Wright   if (idx == 0) return ynodes[1] - ynodes[0];
392e6098bcdSJames Wright   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
393e6098bcdSJames Wright   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
394e6098bcdSJames Wright }
395e6098bcdSJames Wright 
396e6098bcdSJames Wright // Function passed to DMAddBoundary
3976d0190e2SJames Wright // NOTE: Not used in favor of QFunction-based method
398*2b916ea7SJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
399e6098bcdSJames Wright   PetscFunctionBeginUser;
400*2b916ea7SJeremy L Thompson 
401e6098bcdSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
402e6098bcdSJames Wright   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
403e6098bcdSJames Wright   const bool             mean_only = stg_ctx->mean_only;
404e6098bcdSJames Wright   const PetscScalar      dx        = stg_ctx->dx;
405e6098bcdSJames Wright   const PetscScalar      dz        = stg_ctx->dz;
406e6098bcdSJames Wright   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
407e6098bcdSJames Wright   const PetscScalar      theta0    = stg_ctx->theta0;
408e6098bcdSJames Wright   const PetscScalar      P0        = stg_ctx->P0;
409e6098bcdSJames Wright   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
410e6098bcdSJames Wright   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
411e6098bcdSJames Wright   const PetscScalar      Rd        = cp - cv;
412e6098bcdSJames Wright 
413e6098bcdSJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
414e6098bcdSJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
415e6098bcdSJames Wright   if (!mean_only) {
416e6098bcdSJames Wright     const PetscInt     nynodes = stg_ctx->nynodes;
417e6098bcdSJames Wright     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
418e6098bcdSJames Wright     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
419e6098bcdSJames Wright     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
420e6098bcdSJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
421e6098bcdSJames Wright   } else {
422e6098bcdSJames Wright     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
423e6098bcdSJames Wright   }
424e6098bcdSJames Wright 
425e6098bcdSJames Wright   bcval[0] = rho;
426e6098bcdSJames Wright   bcval[1] = rho * u[0];
427e6098bcdSJames Wright   bcval[2] = rho * u[1];
428e6098bcdSJames Wright   bcval[3] = rho * u[2];
429e6098bcdSJames Wright   PetscFunctionReturn(0);
430e6098bcdSJames Wright }
431e6098bcdSJames Wright 
432*2b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
433e6098bcdSJames Wright   DMLabel label;
434e6098bcdSJames Wright   PetscFunctionBeginUser;
435e6098bcdSJames Wright 
436d374fb47SJames Wright   PetscInt comps[5], num_comps = 4;
4373636f6a4SJames Wright   switch (phys->state_var) {
4383636f6a4SJames Wright     case STATEVAR_CONSERVATIVE:
439d374fb47SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
440d374fb47SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
4413636f6a4SJames Wright       break;
4423636f6a4SJames Wright 
4433636f6a4SJames Wright     case STATEVAR_PRIMITIVE:
4443636f6a4SJames Wright       // {1,2,3,4} for u, v, w, T
4453636f6a4SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
4463636f6a4SJames Wright       break;
447d374fb47SJames Wright   }
448d374fb47SJames Wright 
44934b5deb1SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
450e6098bcdSJames Wright   // Set wall BCs
451e6098bcdSJames Wright   if (bc->num_inflow > 0) {
452*2b916ea7SJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
45334b5deb1SJames Wright                             NULL, global_stg_ctx, NULL));
454e6098bcdSJames Wright   }
455e6098bcdSJames Wright 
456e6098bcdSJames Wright   PetscFunctionReturn(0);
457e6098bcdSJames Wright }
4586d0190e2SJames Wright 
459*2b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
4609eeef72bSJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
4616d0190e2SJames Wright   CeedQFunction qf_strongbc;
4626d0190e2SJames Wright   PetscFunctionBeginUser;
463*2b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
464*2b916ea7SJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
4656d0190e2SJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
4666d0190e2SJames Wright   CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE);
4679eeef72bSJames Wright   CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4686d0190e2SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
4696d0190e2SJames Wright 
4706d0190e2SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
4716d0190e2SJames Wright   *pqf_strongbc = qf_strongbc;
4726d0190e2SJames Wright   PetscFunctionReturn(0);
4736d0190e2SJames Wright }
4749eeef72bSJames Wright 
475*2b916ea7SJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
4769eeef72bSJames Wright                                             CeedQFunction *pqf_strongbc) {
4779eeef72bSJames Wright   CeedQFunction qf_strongbc;
4789eeef72bSJames Wright   PetscFunctionBeginUser;
479*2b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc);
480*2b916ea7SJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
4819eeef72bSJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
4829eeef72bSJames Wright   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4839eeef72bSJames Wright 
4849eeef72bSJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
4859eeef72bSJames Wright   *pqf_strongbc = qf_strongbc;
4869eeef72bSJames Wright   PetscFunctionReturn(0);
4879eeef72bSJames Wright }
488