xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
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 
12*2b730f8bSJeremy L Thompson #include "stg_shur14.h"
13*2b730f8bSJeremy L Thompson 
14ba6664aeSJames Wright #include <math.h>
15ba6664aeSJames Wright #include <petsc.h>
16*2b730f8bSJeremy L Thompson #include <stdlib.h>
17*2b730f8bSJeremy L Thompson 
18ba6664aeSJames Wright #include "../navierstokes.h"
19ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h"
20ba6664aeSJames Wright 
2165dd5cafSJames Wright STGShur14Context global_stg_ctx;
2265dd5cafSJames Wright 
23ba6664aeSJames Wright /*
24ba6664aeSJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
25ba6664aeSJames Wright  *
26ba6664aeSJames Wright  * This assumes the input matrices are in order [11,22,33,12,13,23]. This
27ba6664aeSJames Wright  * format is also used for the output.
28ba6664aeSJames Wright  *
29ba6664aeSJames Wright  * @param[in]  comm   MPI_Comm
30ba6664aeSJames Wright  * @param[in]  nprofs Number of matrices in Rij
31ba6664aeSJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
32ba6664aeSJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
33ba6664aeSJames Wright  */
34*2b730f8bSJeremy L Thompson PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs, const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
35ba6664aeSJames Wright   PetscFunctionBeginUser;
36ba6664aeSJames Wright   for (PetscInt i = 0; i < nprofs; i++) {
37ba6664aeSJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
38ba6664aeSJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
39ba6664aeSJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2));
40ba6664aeSJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
41ba6664aeSJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i] * Cij[4][i]) / Cij[1][i];
42ba6664aeSJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
43ba6664aeSJames Wright 
44*2b730f8bSJeremy L Thompson     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i])) {
45*2b730f8bSJeremy L Thompson       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %" PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.",
46*2b730f8bSJeremy L Thompson               i + 1);
47*2b730f8bSJeremy L Thompson     }
48ba6664aeSJames Wright   }
49ba6664aeSJames Wright   PetscFunctionReturn(0);
50ba6664aeSJames Wright }
51ba6664aeSJames Wright 
52ba6664aeSJames Wright /*
53ba6664aeSJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
54ba6664aeSJames Wright  *
55ba6664aeSJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and
56ba6664aeSJames Wright  * passes the file pointer in `fp`. It is not closed in this function, thus
57ba6664aeSJames Wright  * `fp` must be closed sometime after this function has been called (using
58ba6664aeSJames Wright  * `PetscFClose` for example).
59ba6664aeSJames Wright  *
60ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
61ba6664aeSJames Wright  * as the only two entries, separated by a single space
62ba6664aeSJames Wright  *
63ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
64ba6664aeSJames Wright  * @param[in] path Path to the file
65ba6664aeSJames Wright  * @param[in] char_array_len Length of the character array that should contain each line
66ba6664aeSJames Wright  * @param[out] dims Dimensions of the file, taken from the first line of the file
67ba6664aeSJames Wright  * @param[out] fp File pointer to the opened file
68ba6664aeSJames Wright  */
69*2b730f8bSJeremy 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*2b730f8bSJeremy L Thompson                                         FILE **fp) {
71ba6664aeSJames Wright   PetscInt ndims;
72ba6664aeSJames Wright   char     line[char_array_len];
73ba6664aeSJames Wright   char   **array;
74ba6664aeSJames Wright 
75ba6664aeSJames Wright   PetscFunctionBeginUser;
7624941a69SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
7724941a69SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
7824941a69SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
79*2b730f8bSJeremy L Thompson   if (ndims != 2) {
80*2b730f8bSJeremy L Thompson     SETERRQ(comm, -1, "Found %" PetscInt_FMT " dimensions instead of 2 on the first line of %s", ndims, path);
81*2b730f8bSJeremy L Thompson   }
82ba6664aeSJames Wright 
83ba6664aeSJames Wright   for (PetscInt i = 0; i < ndims; i++) dims[i] = atoi(array[i]);
8424941a69SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
85*2b730f8bSJeremy L Thompson 
86ba6664aeSJames Wright   PetscFunctionReturn(0);
87ba6664aeSJames Wright }
88ba6664aeSJames Wright 
89ba6664aeSJames Wright /*
90ba6664aeSJames Wright  * @brief Get the number of rows for the PHASTA file at path
91ba6664aeSJames Wright  *
92ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
93ba6664aeSJames Wright  * as the only two entries, separated by a single space
94ba6664aeSJames Wright  *
95ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
96ba6664aeSJames Wright  * @param[in] path Path to the file
97ba6664aeSJames Wright  * @param[out] nrows Number of rows
98ba6664aeSJames Wright  */
99*2b730f8bSJeremy L Thompson static PetscErrorCode GetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
100ba6664aeSJames Wright   const PetscInt char_array_len = 512;
101ba6664aeSJames Wright   PetscInt       dims[2];
102ba6664aeSJames Wright   FILE          *fp;
103ba6664aeSJames Wright 
104ba6664aeSJames Wright   PetscFunctionBeginUser;
10524941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
106ba6664aeSJames Wright   *nrows = dims[0];
10724941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
108*2b730f8bSJeremy L Thompson 
109ba6664aeSJames Wright   PetscFunctionReturn(0);
110ba6664aeSJames Wright }
111ba6664aeSJames Wright 
112ba6664aeSJames Wright /*
113ba6664aeSJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
114ba6664aeSJames Wright  *
115ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
116ba6664aeSJames Wright  * as the only two entries, separated by a single space.
117ba6664aeSJames Wright  * Assumes there are 14 columns in the file
118ba6664aeSJames Wright  *
119ba6664aeSJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
120ba6664aeSJames Wright  * profile found in the file
121ba6664aeSJames Wright  *
122ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
123ba6664aeSJames Wright  * @param[in] path Path to the STGInflow.dat file
124ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
125ba6664aeSJames Wright  */
126*2b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGInflow(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
127ba6664aeSJames Wright   PetscInt       ndims, dims[2];
128ba6664aeSJames Wright   FILE          *fp;
129ba6664aeSJames Wright   const PetscInt char_array_len = 512;
130ba6664aeSJames Wright   char           line[char_array_len];
131ba6664aeSJames Wright   char         **array;
132ba6664aeSJames Wright 
133ba6664aeSJames Wright   PetscFunctionBeginUser;
134ba6664aeSJames Wright 
13524941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
136ba6664aeSJames Wright 
137ba6664aeSJames Wright   CeedScalar  rij[6][stg_ctx->nprofs];
138175f00a6SJames Wright   CeedScalar *wall_dist              = &stg_ctx->data[stg_ctx->offsets.wall_dist];
139ba6664aeSJames Wright   CeedScalar *eps                    = &stg_ctx->data[stg_ctx->offsets.eps];
140ba6664aeSJames Wright   CeedScalar *lt                     = &stg_ctx->data[stg_ctx->offsets.lt];
141*2b730f8bSJeremy L Thompson   CeedScalar(*ubar)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.ubar];
142ba6664aeSJames Wright 
143ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
14424941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
14524941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
146*2b730f8bSJeremy L Thompson     if (ndims < dims[1]) {
147*2b730f8bSJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
148*2b730f8bSJeremy L Thompson               ndims, dims[1]);
149*2b730f8bSJeremy L Thompson     }
150ba6664aeSJames Wright 
151175f00a6SJames Wright     wall_dist[i] = (CeedScalar)atof(array[0]);
152ba6664aeSJames Wright     ubar[0][i]   = (CeedScalar)atof(array[1]);
153ba6664aeSJames Wright     ubar[1][i]   = (CeedScalar)atof(array[2]);
154ba6664aeSJames Wright     ubar[2][i]   = (CeedScalar)atof(array[3]);
155ba6664aeSJames Wright     rij[0][i]    = (CeedScalar)atof(array[4]);
156ba6664aeSJames Wright     rij[1][i]    = (CeedScalar)atof(array[5]);
157ba6664aeSJames Wright     rij[2][i]    = (CeedScalar)atof(array[6]);
158ba6664aeSJames Wright     rij[3][i]    = (CeedScalar)atof(array[7]);
159ba6664aeSJames Wright     rij[4][i]    = (CeedScalar)atof(array[8]);
160ba6664aeSJames Wright     rij[5][i]    = (CeedScalar)atof(array[9]);
161ba6664aeSJames Wright     lt[i]        = (CeedScalar)atof(array[12]);
162ba6664aeSJames Wright     eps[i]       = (CeedScalar)atof(array[13]);
163ba6664aeSJames Wright 
164*2b730f8bSJeremy L Thompson     if (wall_dist[i] < 0) SETERRQ(comm, -1, "Distance to wall in %s cannot be negative", path);
165*2b730f8bSJeremy L Thompson     if (lt[i] < 0) SETERRQ(comm, -1, "Turbulent length scale in %s cannot be negative", path);
166*2b730f8bSJeremy L Thompson     if (eps[i] < 0) SETERRQ(comm, -1, "Turbulent dissipation in %s cannot be negative", path);
167ba6664aeSJames Wright   }
168*2b730f8bSJeremy L Thompson   CeedScalar(*cij)[stg_ctx->nprofs] = (CeedScalar(*)[stg_ctx->nprofs]) & stg_ctx->data[stg_ctx->offsets.cij];
16924941a69SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
17024941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
171*2b730f8bSJeremy L Thompson 
172ba6664aeSJames Wright   PetscFunctionReturn(0);
173ba6664aeSJames Wright }
174ba6664aeSJames Wright 
175ba6664aeSJames Wright /*
176ba6664aeSJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
177ba6664aeSJames Wright  *
178ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
179ba6664aeSJames Wright  * as the only two entries, separated by a single space.
180ba6664aeSJames Wright  * Assumes there are 7 columns in the file
181ba6664aeSJames Wright  *
182ba6664aeSJames Wright  * @param[in]    comm    MPI_Comm for the program
183ba6664aeSJames Wright  * @param[in]    path    Path to the STGRand.dat file
184ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
185ba6664aeSJames Wright  */
186*2b730f8bSJeremy L Thompson static PetscErrorCode ReadSTGRand(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
187ba6664aeSJames Wright   PetscInt       ndims, dims[2];
188ba6664aeSJames Wright   FILE          *fp;
189ba6664aeSJames Wright   const PetscInt char_array_len = 512;
190ba6664aeSJames Wright   char           line[char_array_len];
191ba6664aeSJames Wright   char         **array;
192ba6664aeSJames Wright 
193ba6664aeSJames Wright   PetscFunctionBeginUser;
19424941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
195ba6664aeSJames Wright 
196ba6664aeSJames Wright   CeedScalar *phi                     = &stg_ctx->data[stg_ctx->offsets.phi];
197*2b730f8bSJeremy L Thompson   CeedScalar(*d)[stg_ctx->nmodes]     = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.d];
198*2b730f8bSJeremy L Thompson   CeedScalar(*sigma)[stg_ctx->nmodes] = (CeedScalar(*)[stg_ctx->nmodes]) & stg_ctx->data[stg_ctx->offsets.sigma];
199ba6664aeSJames Wright 
200ba6664aeSJames Wright   for (PetscInt i = 0; i < stg_ctx->nmodes; i++) {
20124941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
20224941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
203*2b730f8bSJeremy L Thompson     if (ndims < dims[1]) {
204*2b730f8bSJeremy L Thompson       SETERRQ(comm, -1, "Line %" PetscInt_FMT " of %s does not contain enough columns (%" PetscInt_FMT " instead of %" PetscInt_FMT ")", i, path,
205*2b730f8bSJeremy L Thompson               ndims, dims[1]);
206*2b730f8bSJeremy L Thompson     }
207ba6664aeSJames Wright 
208ba6664aeSJames Wright     d[0][i]     = (CeedScalar)atof(array[0]);
209ba6664aeSJames Wright     d[1][i]     = (CeedScalar)atof(array[1]);
210ba6664aeSJames Wright     d[2][i]     = (CeedScalar)atof(array[2]);
211ba6664aeSJames Wright     phi[i]      = (CeedScalar)atof(array[3]);
212ba6664aeSJames Wright     sigma[0][i] = (CeedScalar)atof(array[4]);
213ba6664aeSJames Wright     sigma[1][i] = (CeedScalar)atof(array[5]);
214ba6664aeSJames Wright     sigma[2][i] = (CeedScalar)atof(array[6]);
215ba6664aeSJames Wright   }
21624941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
217*2b730f8bSJeremy L Thompson 
218ba6664aeSJames Wright   PetscFunctionReturn(0);
219ba6664aeSJames Wright }
220ba6664aeSJames Wright 
221ba6664aeSJames Wright /*
222ba6664aeSJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
223ba6664aeSJames Wright  *
224ba6664aeSJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
225ba6664aeSJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
226ba6664aeSJames Wright  * STGShur14Context instance.
227ba6664aeSJames Wright  *
228ba6664aeSJames Wright  * @param[in]    comm            MPI_Comm for the program
229ba6664aeSJames Wright  * @param[in]    dm              DM for the program
230ba6664aeSJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
231ba6664aeSJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
232ba6664aeSJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
233ba6664aeSJames Wright  */
234*2b730f8bSJeremy 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*2b730f8bSJeremy L Thompson                                  STGShur14Context *pstg_ctx, const CeedScalar ynodes[]) {
236ba6664aeSJames Wright   PetscInt         nmodes, nprofs;
237ba6664aeSJames Wright   STGShur14Context stg_ctx;
238ba6664aeSJames Wright   PetscFunctionBeginUser;
239ba6664aeSJames Wright 
240ba6664aeSJames Wright   // Get options
24124941a69SJames Wright   PetscCall(GetNRows(comm, stg_rand_path, &nmodes));
24224941a69SJames Wright   PetscCall(GetNRows(comm, stg_inflow_path, &nprofs));
243ba6664aeSJames Wright   if (nmodes > STG_NMODES_MAX)
244*2b730f8bSJeremy L Thompson     SETERRQ(comm, 1,
245*2b730f8bSJeremy L Thompson             "Number of wavemodes in %s (%" PetscInt_FMT ") exceeds STG_NMODES_MAX (%" PetscInt_FMT
246*2b730f8bSJeremy L Thompson             "). "
247*2b730f8bSJeremy L Thompson             "Change size of STG_NMODES_MAX and recompile",
248*2b730f8bSJeremy L Thompson             stg_rand_path, nmodes, STG_NMODES_MAX);
249ba6664aeSJames Wright 
250ba6664aeSJames Wright   {
251ba6664aeSJames Wright     STGShur14Context s;
25224941a69SJames Wright     PetscCall(PetscCalloc1(1, &s));
253ba6664aeSJames Wright     *s                         = **pstg_ctx;
254ba6664aeSJames Wright     s->nmodes                  = nmodes;
255ba6664aeSJames Wright     s->nprofs                  = nprofs;
256ba6664aeSJames Wright     s->offsets.sigma           = 0;
257ba6664aeSJames Wright     s->offsets.d               = nmodes * 3;
258ba6664aeSJames Wright     s->offsets.phi             = s->offsets.d + nmodes * 3;
259ba6664aeSJames Wright     s->offsets.kappa           = s->offsets.phi + nmodes;
260175f00a6SJames Wright     s->offsets.wall_dist       = s->offsets.kappa + nmodes;
261175f00a6SJames Wright     s->offsets.ubar            = s->offsets.wall_dist + nprofs;
262ba6664aeSJames Wright     s->offsets.cij             = s->offsets.ubar + nprofs * 3;
263ba6664aeSJames Wright     s->offsets.eps             = s->offsets.cij + nprofs * 6;
264ba6664aeSJames Wright     s->offsets.lt              = s->offsets.eps + nprofs;
2654e139266SJames Wright     s->offsets.ynodes          = s->offsets.lt + nprofs;
2664e139266SJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
267ba6664aeSJames Wright     s->total_bytes             = sizeof(*stg_ctx) + total_num_scalars * sizeof(stg_ctx->data[0]);
26824941a69SJames Wright     PetscCall(PetscMalloc(s->total_bytes, &stg_ctx));
269ba6664aeSJames Wright     *stg_ctx = *s;
27024941a69SJames Wright     PetscCall(PetscFree(s));
271ba6664aeSJames Wright   }
272ba6664aeSJames Wright 
27324941a69SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx));
27424941a69SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx));
275ba6664aeSJames Wright 
2764e139266SJames Wright   if (stg_ctx->nynodes > 0) {
2774e139266SJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
2784e139266SJames Wright     for (PetscInt i = 0; i < stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
2794e139266SJames Wright   }
2804e139266SJames Wright 
281ba6664aeSJames Wright   // -- Calculate kappa
282ba6664aeSJames Wright   {
283ba6664aeSJames Wright     CeedScalar *kappa     = &stg_ctx->data[stg_ctx->offsets.kappa];
284175f00a6SJames Wright     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
285ba6664aeSJames Wright     CeedScalar *lt        = &stg_ctx->data[stg_ctx->offsets.lt];
286ba6664aeSJames Wright     CeedScalar  le, le_max = 0;
287ba6664aeSJames Wright 
288*2b730f8bSJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nprofs; i++) {
289175f00a6SJames Wright       le = PetscMin(2 * wall_dist[i], 3 * lt[i]);
290ba6664aeSJames Wright       if (le_max < le) le_max = le;
291ba6664aeSJames Wright     }
292ba6664aeSJames Wright     CeedScalar kmin = M_PI / le_max;
293ba6664aeSJames Wright 
294*2b730f8bSJeremy L Thompson     CeedPragmaSIMD for (PetscInt i = 0; i < stg_ctx->nmodes; i++) { kappa[i] = kmin * pow(stg_ctx->alpha, i); }
295ba6664aeSJames Wright   }  // end calculate kappa
296ba6664aeSJames Wright 
29724941a69SJames Wright   PetscCall(PetscFree(*pstg_ctx));
298ba6664aeSJames Wright   *pstg_ctx = stg_ctx;
299ba6664aeSJames Wright   PetscFunctionReturn(0);
300ba6664aeSJames Wright }
301ba6664aeSJames Wright 
302*2b730f8bSJeremy L Thompson PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem, User user, const bool prescribe_T, const CeedScalar theta0,
303*2b730f8bSJeremy L Thompson                         const CeedScalar P0, const CeedScalar ynodes[], const CeedInt nynodes) {
304ba6664aeSJames Wright   char                     stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
305ba6664aeSJames Wright   char                     stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
306*2b730f8bSJeremy L Thompson   PetscBool                mean_only = PETSC_FALSE, use_stgstrong = PETSC_FALSE, use_fluctuating_IC = PETSC_FALSE;
307*2b730f8bSJeremy L Thompson   CeedScalar               u0 = 0.0, alpha = 1.01;
308ba6664aeSJames Wright   CeedQFunctionContext     stg_context;
309ba6664aeSJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
310ba6664aeSJames Wright   PetscFunctionBeginUser;
311ba6664aeSJames Wright 
312ba6664aeSJames Wright   // Get options
313ba6664aeSJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
314*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL, stg_inflow_path, stg_inflow_path, sizeof(stg_inflow_path), NULL));
315*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL, stg_rand_path, stg_rand_path, sizeof(stg_rand_path), NULL));
316*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL, alpha, &alpha, NULL));
317*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations", NULL, u0, &u0, NULL));
318*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile", NULL, mean_only, &mean_only, NULL));
319*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly", NULL, use_stgstrong, &use_stgstrong, NULL));
320*2b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-stg_fluctuating_IC", "\"Extrude\" the fluctuations through the domain as an initial condition", NULL,
321*2b730f8bSJeremy L Thompson                              use_fluctuating_IC, &use_fluctuating_IC, NULL));
322ba6664aeSJames Wright   PetscOptionsEnd();
323ba6664aeSJames Wright 
32424941a69SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
32565dd5cafSJames Wright   global_stg_ctx->alpha              = alpha;
32665dd5cafSJames Wright   global_stg_ctx->u0                 = u0;
32765dd5cafSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
32865dd5cafSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
32965dd5cafSJames Wright   global_stg_ctx->mean_only          = mean_only;
33089060322SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
33165dd5cafSJames Wright   global_stg_ctx->theta0             = theta0;
33265dd5cafSJames Wright   global_stg_ctx->P0                 = P0;
33365dd5cafSJames Wright   global_stg_ctx->nynodes            = nynodes;
334ba6664aeSJames Wright 
335ba6664aeSJames Wright   {
336ba6664aeSJames Wright     // Calculate dx assuming constant spacing
337ba6664aeSJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
33824941a69SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
339ba6664aeSJames Wright     for (PetscInt i = 0; i < 3; i++) domain_size[i] = domain_max[i] - domain_min[i];
340ba6664aeSJames Wright 
341ba6664aeSJames Wright     PetscInt nmax = 3, faces[3];
342*2b730f8bSJeremy L Thompson     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax, NULL));
34365dd5cafSJames Wright     global_stg_ctx->dx = domain_size[0] / faces[0];
34465dd5cafSJames Wright     global_stg_ctx->dz = domain_size[2] / faces[2];
345ba6664aeSJames Wright   }
346ba6664aeSJames Wright 
347*2b730f8bSJeremy L Thompson   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
34865dd5cafSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
349*2b730f8bSJeremy L Thompson   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
350ba6664aeSJames Wright 
351*2b730f8bSJeremy L Thompson   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &global_stg_ctx, ynodes));
352ba6664aeSJames Wright 
353ba6664aeSJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
354*2b730f8bSJeremy L Thompson   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST, CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
355*2b730f8bSJeremy L Thompson   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST, FreeContextPetsc);
356*2b730f8bSJeremy L Thompson   CeedQFunctionContextRegisterDouble(stg_context, "solution time", offsetof(struct STGShur14Context_, time), 1, "Physical time of the solution");
357ba6664aeSJames Wright 
358b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
359b77c53c9SJames Wright   problem->ics.qfunction         = ICsSTG;
360b77c53c9SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
361b77c53c9SJames Wright   problem->ics.qfunction_context = stg_context;
362b77c53c9SJames Wright 
3634e139266SJames Wright   if (use_stgstrong) {
36465dd5cafSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
365dada6cc0SJames Wright     problem->use_dirichlet_ceed = PETSC_TRUE;
36636b31e27SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
3674e139266SJames Wright   } else {
368ba6664aeSJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
369ba6664aeSJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
3704dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
3714dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
372*2b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow.qfunction_context);
373*2b730f8bSJeremy L Thompson     CeedQFunctionContextReferenceCopy(stg_context, &problem->apply_inflow_jacobian.qfunction_context);
37436b31e27SJames Wright     problem->bc_from_ics = PETSC_TRUE;
3754e139266SJames Wright   }
376ba6664aeSJames Wright 
377ba6664aeSJames Wright   PetscFunctionReturn(0);
378ba6664aeSJames Wright }
3794e139266SJames Wright 
380*2b730f8bSJeremy L Thompson static inline PetscScalar FindDy(const PetscScalar ynodes[], const PetscInt nynodes, const PetscScalar y) {
3814e139266SJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
3824e139266SJames Wright   // ^^assuming min(dy) is first element off the wall
3834e139266SJames Wright   PetscInt idx = -1;  // Index
3844e139266SJames Wright 
3854e139266SJames Wright   for (PetscInt i = 0; i < nynodes; i++) {
3864e139266SJames Wright     if (y < ynodes[i] + half_mindy) {
387*2b730f8bSJeremy L Thompson       idx = i;
388*2b730f8bSJeremy L Thompson       break;
3894e139266SJames Wright     }
3904e139266SJames Wright   }
3914e139266SJames Wright   if (idx == 0) return ynodes[1] - ynodes[0];
3924e139266SJames Wright   else if (idx == nynodes - 1) return ynodes[nynodes - 2] - ynodes[nynodes - 1];
3934e139266SJames Wright   else return 0.5 * (ynodes[idx + 1] - ynodes[idx - 1]);
3944e139266SJames Wright }
3954e139266SJames Wright 
3964e139266SJames Wright // Function passed to DMAddBoundary
397dada6cc0SJames Wright // NOTE: Not used in favor of QFunction-based method
398*2b730f8bSJeremy L Thompson PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
3994e139266SJames Wright   PetscFunctionBeginUser;
400*2b730f8bSJeremy L Thompson 
4014e139266SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context)ctx;
4024e139266SJames Wright   PetscScalar            qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
4034e139266SJames Wright   const bool             mean_only = stg_ctx->mean_only;
4044e139266SJames Wright   const PetscScalar      dx        = stg_ctx->dx;
4054e139266SJames Wright   const PetscScalar      dz        = stg_ctx->dz;
4064e139266SJames Wright   const PetscScalar      mu        = stg_ctx->newtonian_ctx.mu;
4074e139266SJames Wright   const PetscScalar      theta0    = stg_ctx->theta0;
4084e139266SJames Wright   const PetscScalar      P0        = stg_ctx->P0;
4094e139266SJames Wright   const PetscScalar      cv        = stg_ctx->newtonian_ctx.cv;
4104e139266SJames Wright   const PetscScalar      cp        = stg_ctx->newtonian_ctx.cp;
4114e139266SJames Wright   const PetscScalar      Rd        = cp - cv;
4124e139266SJames Wright 
4134e139266SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
4144e139266SJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
4154e139266SJames Wright   if (!mean_only) {
4164e139266SJames Wright     const PetscInt     nynodes = stg_ctx->nynodes;
4174e139266SJames Wright     const PetscScalar *ynodes  = &stg_ctx->data[stg_ctx->offsets.ynodes];
4184e139266SJames Wright     const PetscScalar  h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
4194e139266SJames Wright     CalcSpectrum(x[1], eps, lt, h, mu / rho, qn, stg_ctx);
4204e139266SJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
4214e139266SJames Wright   } else {
4224e139266SJames Wright     for (CeedInt j = 0; j < 3; j++) u[j] = ubar[j];
4234e139266SJames Wright   }
4244e139266SJames Wright 
4254e139266SJames Wright   bcval[0] = rho;
4264e139266SJames Wright   bcval[1] = rho * u[0];
4274e139266SJames Wright   bcval[2] = rho * u[1];
4284e139266SJames Wright   bcval[3] = rho * u[2];
4294e139266SJames Wright   PetscFunctionReturn(0);
4304e139266SJames Wright }
4314e139266SJames Wright 
432*2b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem, Physics phys) {
4334e139266SJames Wright   DMLabel label;
4344e139266SJames Wright   PetscFunctionBeginUser;
4354e139266SJames Wright 
436192a7459SJames Wright   PetscInt comps[5], num_comps = 4;
43797baf651SJames Wright   switch (phys->state_var) {
43897baf651SJames Wright     case STATEVAR_CONSERVATIVE:
439192a7459SJames Wright       // {0,1,2,3} for rho, rho*u, rho*v, rho*w
440192a7459SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i;
44197baf651SJames Wright       break;
44297baf651SJames Wright 
44397baf651SJames Wright     case STATEVAR_PRIMITIVE:
44497baf651SJames Wright       // {1,2,3,4} for u, v, w, T
44597baf651SJames Wright       for (int i = 0; i < 4; i++) comps[i] = i + 1;
44697baf651SJames Wright       break;
447192a7459SJames Wright   }
448192a7459SJames Wright 
44924941a69SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
4504e139266SJames Wright   // Set wall BCs
4514e139266SJames Wright   if (bc->num_inflow > 0) {
452*2b730f8bSJeremy L Thompson     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label, bc->num_inflow, bc->inflows, 0, num_comps, comps, (void (*)(void))StrongSTGbcFunc,
45324941a69SJames Wright                             NULL, global_stg_ctx, NULL));
4544e139266SJames Wright   }
4554e139266SJames Wright 
4564e139266SJames Wright   PetscFunctionReturn(0);
4574e139266SJames Wright }
458dada6cc0SJames Wright 
459*2b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
4605dc40723SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
461dada6cc0SJames Wright   CeedQFunction qf_strongbc;
462dada6cc0SJames Wright   PetscFunctionBeginUser;
463*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF, STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
464*2b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
465dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
466dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc, "scale", 1, CEED_EVAL_NONE);
4675dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
468dada6cc0SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "q", num_comp_q, CEED_EVAL_NONE);
469dada6cc0SJames Wright 
470dada6cc0SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
471dada6cc0SJames Wright   *pqf_strongbc = qf_strongbc;
472dada6cc0SJames Wright   PetscFunctionReturn(0);
473dada6cc0SJames Wright }
4745dc40723SJames Wright 
475*2b730f8bSJeremy L Thompson PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem, CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
4765dc40723SJames Wright                                             CeedQFunction *pqf_strongbc) {
4775dc40723SJames Wright   CeedQFunction qf_strongbc;
4785dc40723SJames Wright   PetscFunctionBeginUser;
479*2b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14, Preprocess_STGShur14_loc, &qf_strongbc);
480*2b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
4815dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc, "x", num_comp_x, CEED_EVAL_NONE);
4825dc40723SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
4835dc40723SJames Wright 
4845dc40723SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
4855dc40723SJames Wright   *pqf_strongbc = qf_strongbc;
4865dc40723SJames Wright   PetscFunctionReturn(0);
4875dc40723SJames Wright }
488