xref: /honee/problems/stg_shur14.c (revision 6dd99bedd556801b3601df82fac20dedc5308a81)
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 
12493642f1SJames Wright #include <stdlib.h>
13493642f1SJames Wright #include <math.h>
14493642f1SJames Wright #include <petsc.h>
15493642f1SJames Wright #include "../navierstokes.h"
16493642f1SJames Wright #include "stg_shur14.h"
17493642f1SJames Wright #include "../qfunctions/stg_shur14.h"
18493642f1SJames Wright 
19493642f1SJames Wright #ifndef M_PI
20493642f1SJames Wright #define M_PI    3.14159265358979323846
21493642f1SJames Wright #endif
22493642f1SJames 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  */
34493642f1SJames Wright PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs,
35493642f1SJames Wright                                   const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
36493642f1SJames Wright   PetscFunctionBeginUser;
37493642f1SJames Wright   for (PetscInt i=0; i<nprofs; i++) {
38493642f1SJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
39493642f1SJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
40493642f1SJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2) );
41493642f1SJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
42493642f1SJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i]*Cij[4][i]) / Cij[1][i];
43493642f1SJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
44493642f1SJames Wright 
45493642f1SJames Wright     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i]))
46493642f1SJames Wright       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %d. "
47493642f1SJames Wright               "Either STGInflow has non-SPD matrix or contains nan.", i+1);
48493642f1SJames Wright   }
49493642f1SJames Wright   PetscFunctionReturn(0);
50493642f1SJames Wright }
51493642f1SJames Wright 
52493642f1SJames Wright 
53493642f1SJames Wright /*
54493642f1SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
55493642f1SJames Wright  *
56493642f1SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and
57493642f1SJames Wright  * passes the file pointer in `fp`. It is not closed in this function, thus
58493642f1SJames Wright  * `fp` must be closed sometime after this function has been called (using
59493642f1SJames Wright  * `PetscFClose` for example).
60493642f1SJames Wright  *
61493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
62493642f1SJames Wright  * as the only two entries, separated by a single space
63493642f1SJames Wright  *
64493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
65493642f1SJames Wright  * @param[in] path Path to the file
66493642f1SJames Wright  * @param[in] char_array_len Length of the character array that should contain each line
67493642f1SJames Wright  * @param[out] dims Dimensions of the file, taken from the first line of the file
68493642f1SJames Wright  * @param[out] fp File pointer to the opened file
69493642f1SJames Wright  */
70493642f1SJames Wright static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm,
71493642f1SJames Wright                                         const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len,
72493642f1SJames Wright                                         PetscInt dims[2], FILE **fp) {
73493642f1SJames Wright   PetscErrorCode ierr;
74493642f1SJames Wright   PetscInt ndims;
75493642f1SJames Wright   char line[char_array_len];
76493642f1SJames Wright   char **array;
77493642f1SJames Wright 
78493642f1SJames Wright   PetscFunctionBeginUser;
79493642f1SJames Wright   ierr = PetscFOpen(comm, path, "r", fp); CHKERRQ(ierr);
80493642f1SJames Wright   ierr = PetscSynchronizedFGets(comm, *fp, char_array_len, line); CHKERRQ(ierr);
81493642f1SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
82493642f1SJames Wright   if (ndims != 2) SETERRQ(comm, -1,
83493642f1SJames Wright                             "Found %d dimensions instead of 2 on the first line of %s",
84493642f1SJames Wright                             ndims, path);
85493642f1SJames Wright 
86493642f1SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
87493642f1SJames Wright   ierr = PetscStrToArrayDestroy(ndims, array); CHKERRQ(ierr);
88493642f1SJames Wright   PetscFunctionReturn(0);
89493642f1SJames Wright }
90493642f1SJames Wright 
91493642f1SJames Wright /*
92493642f1SJames Wright  * @brief Get the number of rows for the PHASTA file at path
93493642f1SJames Wright  *
94493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
95493642f1SJames Wright  * as the only two entries, separated by a single space
96493642f1SJames Wright  *
97493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
98493642f1SJames Wright  * @param[in] path Path to the file
99493642f1SJames Wright  * @param[out] nrows Number of rows
100493642f1SJames Wright  */
101493642f1SJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm,
102493642f1SJames Wright                                const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
103493642f1SJames Wright   PetscErrorCode ierr;
104493642f1SJames Wright   const PetscInt char_array_len = 512;
105493642f1SJames Wright   PetscInt dims[2];
106493642f1SJames Wright   FILE *fp;
107493642f1SJames Wright 
108493642f1SJames Wright   PetscFunctionBeginUser;
109493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
110493642f1SJames Wright   *nrows = dims[0];
111493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
112493642f1SJames Wright   PetscFunctionReturn(0);
113493642f1SJames Wright }
114493642f1SJames Wright 
115493642f1SJames Wright /*
116493642f1SJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
117493642f1SJames Wright  *
118493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
119493642f1SJames Wright  * as the only two entries, separated by a single space.
120493642f1SJames Wright  * Assumes there are 14 columns in the file
121493642f1SJames Wright  *
122493642f1SJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
123493642f1SJames Wright  * profile found in the file
124493642f1SJames Wright  *
125493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
126493642f1SJames Wright  * @param[in] path Path to the STGInflow.dat file
127493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
128493642f1SJames Wright  */
129493642f1SJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm,
130493642f1SJames Wright                                     const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
131493642f1SJames Wright   PetscErrorCode ierr;
132493642f1SJames Wright   PetscInt ndims, dims[2];
133493642f1SJames Wright   FILE *fp;
134493642f1SJames Wright   const PetscInt char_array_len=512;
135493642f1SJames Wright   char line[char_array_len];
136493642f1SJames Wright   char **array;
137493642f1SJames Wright 
138493642f1SJames Wright   PetscFunctionBeginUser;
139493642f1SJames Wright 
140493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
141493642f1SJames Wright 
142493642f1SJames Wright   CeedScalar rij[6][stg_ctx->nprofs];
143493642f1SJames Wright   CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw];
144493642f1SJames Wright   CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps];
145493642f1SJames Wright   CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
146493642f1SJames Wright   CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs])
147493642f1SJames Wright                                         &stg_ctx->data[stg_ctx->offsets.ubar];
148493642f1SJames Wright 
149493642f1SJames Wright   for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
150493642f1SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
151493642f1SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
152493642f1SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
153493642f1SJames Wright                                    "Line %d of %s does not contain enough columns (%d instead of %d)", i,
154493642f1SJames Wright                                    path, ndims, dims[1]);
155493642f1SJames Wright 
156493642f1SJames Wright     prof_dw[i] = (CeedScalar) atof(array[0]);
157493642f1SJames Wright     ubar[0][i] = (CeedScalar) atof(array[1]);
158493642f1SJames Wright     ubar[1][i] = (CeedScalar) atof(array[2]);
159493642f1SJames Wright     ubar[2][i] = (CeedScalar) atof(array[3]);
160493642f1SJames Wright     rij[0][i]  = (CeedScalar) atof(array[4]);
161493642f1SJames Wright     rij[1][i]  = (CeedScalar) atof(array[5]);
162493642f1SJames Wright     rij[2][i]  = (CeedScalar) atof(array[6]);
163493642f1SJames Wright     rij[3][i]  = (CeedScalar) atof(array[7]);
164493642f1SJames Wright     rij[4][i]  = (CeedScalar) atof(array[8]);
165493642f1SJames Wright     rij[5][i]  = (CeedScalar) atof(array[9]);
166493642f1SJames Wright     lt[i]      = (CeedScalar) atof(array[12]);
167493642f1SJames Wright     eps[i]     = (CeedScalar) atof(array[13]);
168493642f1SJames Wright 
169493642f1SJames Wright     if (prof_dw[i] < 0) SETERRQ(comm, -1,
170493642f1SJames Wright                                   "Distance to wall in %s cannot be negative", path);
171493642f1SJames Wright     if (lt[i] < 0) SETERRQ(comm, -1,
172493642f1SJames Wright                              "Turbulent length scale in %s cannot be negative", path);
173493642f1SJames Wright     if (eps[i] < 0) SETERRQ(comm, -1,
174493642f1SJames Wright                               "Turbulent dissipation in %s cannot be negative", path);
175493642f1SJames Wright 
176493642f1SJames Wright   }
177493642f1SJames Wright   CeedScalar (*cij)[stg_ctx->nprofs]  = (CeedScalar (*)[stg_ctx->nprofs])
178493642f1SJames Wright                                         &stg_ctx->data[stg_ctx->offsets.cij];
179493642f1SJames Wright   ierr = CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij); CHKERRQ(ierr);
180493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
181493642f1SJames Wright   PetscFunctionReturn(0);
182493642f1SJames Wright }
183493642f1SJames Wright 
184493642f1SJames Wright /*
185493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
186493642f1SJames Wright  *
187493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
188493642f1SJames Wright  * as the only two entries, separated by a single space.
189493642f1SJames Wright  * Assumes there are 7 columns in the file
190493642f1SJames Wright  *
191493642f1SJames Wright  * @param[in]    comm    MPI_Comm for the program
192493642f1SJames Wright  * @param[in]    path    Path to the STGRand.dat file
193493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
194493642f1SJames Wright  */
195493642f1SJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm,
196493642f1SJames Wright                                   const char path[PETSC_MAX_PATH_LEN],
197493642f1SJames Wright                                   STGShur14Context stg_ctx) {
198493642f1SJames Wright   PetscErrorCode ierr;
199493642f1SJames Wright   PetscInt ndims, dims[2];
200493642f1SJames Wright   FILE *fp;
201493642f1SJames Wright   const PetscInt char_array_len = 512;
202493642f1SJames Wright   char line[char_array_len];
203493642f1SJames Wright   char **array;
204493642f1SJames Wright 
205493642f1SJames Wright   PetscFunctionBeginUser;
206493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
207493642f1SJames Wright 
208493642f1SJames Wright   CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi];
209493642f1SJames Wright   CeedScalar (*d)[stg_ctx->nmodes]     = (CeedScalar (*)[stg_ctx->nmodes])
210493642f1SJames Wright                                          &stg_ctx->data[stg_ctx->offsets.d];
211493642f1SJames Wright   CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes])
212493642f1SJames Wright                                          &stg_ctx->data[stg_ctx->offsets.sigma];
213493642f1SJames Wright 
214493642f1SJames Wright   for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
215493642f1SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
216493642f1SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
217493642f1SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
218493642f1SJames Wright                                    "Line %d of %s does not contain enough columns (%d instead of %d)", i,
219493642f1SJames Wright                                    path, ndims, dims[1]);
220493642f1SJames Wright 
221493642f1SJames Wright     d[0][i]     = (CeedScalar) atof(array[0]);
222493642f1SJames Wright     d[1][i]     = (CeedScalar) atof(array[1]);
223493642f1SJames Wright     d[2][i]     = (CeedScalar) atof(array[2]);
224493642f1SJames Wright     phi[i]      = (CeedScalar) atof(array[3]);
225493642f1SJames Wright     sigma[0][i] = (CeedScalar) atof(array[4]);
226493642f1SJames Wright     sigma[1][i] = (CeedScalar) atof(array[5]);
227493642f1SJames Wright     sigma[2][i] = (CeedScalar) atof(array[6]);
228493642f1SJames Wright   }
229493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
230493642f1SJames Wright   PetscFunctionReturn(0);
231493642f1SJames Wright }
232493642f1SJames Wright 
233493642f1SJames Wright /*
234493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
235493642f1SJames Wright  *
236493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
237493642f1SJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
238493642f1SJames Wright  * STGShur14Context instance.
239493642f1SJames Wright  *
240493642f1SJames Wright  * @param[in]    comm            MPI_Comm for the program
241493642f1SJames Wright  * @param[in]    dm              DM for the program
242493642f1SJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
243493642f1SJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
244493642f1SJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
245493642f1SJames Wright  */
246493642f1SJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm,
247493642f1SJames Wright                                  char stg_inflow_path[PETSC_MAX_PATH_LEN],
248493642f1SJames Wright                                  char stg_rand_path[PETSC_MAX_PATH_LEN],
249e6098bcdSJames Wright                                  STGShur14Context *pstg_ctx,
250e6098bcdSJames Wright                                  const CeedScalar ynodes[]) {
251493642f1SJames Wright   PetscErrorCode ierr;
252493642f1SJames Wright   PetscInt nmodes, nprofs;
253493642f1SJames Wright   STGShur14Context stg_ctx;
254493642f1SJames Wright   PetscFunctionBeginUser;
255493642f1SJames Wright 
256493642f1SJames Wright   // Get options
257493642f1SJames Wright   ierr = GetNRows(comm, stg_rand_path, &nmodes); CHKERRQ(ierr);
258493642f1SJames Wright   ierr = GetNRows(comm, stg_inflow_path, &nprofs); CHKERRQ(ierr);
259493642f1SJames Wright   if (nmodes > STG_NMODES_MAX)
260493642f1SJames Wright     SETERRQ(comm, 1, "Number of wavemodes in %s (%d) exceeds STG_NMODES_MAX (%d). "
261493642f1SJames Wright             "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes,
262493642f1SJames Wright             STG_NMODES_MAX);
263493642f1SJames Wright 
264493642f1SJames Wright   {
265493642f1SJames Wright     STGShur14Context s;
266493642f1SJames Wright     ierr = PetscCalloc1(1, &s); CHKERRQ(ierr);
267493642f1SJames Wright     *s = **pstg_ctx;
268493642f1SJames Wright     s->nmodes = nmodes;
269493642f1SJames Wright     s->nprofs = nprofs;
270493642f1SJames Wright     s->offsets.sigma   = 0;
271493642f1SJames Wright     s->offsets.d       = nmodes*3;
272493642f1SJames Wright     s->offsets.phi     = s->offsets.d       + nmodes*3;
273493642f1SJames Wright     s->offsets.kappa   = s->offsets.phi     + nmodes;
274493642f1SJames Wright     s->offsets.prof_dw = s->offsets.kappa   + nmodes;
275493642f1SJames Wright     s->offsets.ubar    = s->offsets.prof_dw + nprofs;
276493642f1SJames Wright     s->offsets.cij     = s->offsets.ubar    + nprofs*3;
277493642f1SJames Wright     s->offsets.eps     = s->offsets.cij     + nprofs*6;
278493642f1SJames Wright     s->offsets.lt      = s->offsets.eps     + nprofs;
279e6098bcdSJames Wright     s->offsets.ynodes  = s->offsets.lt      + nprofs;
280e6098bcdSJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
281493642f1SJames Wright     s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]);
282493642f1SJames Wright     ierr = PetscMalloc(s->total_bytes, &stg_ctx); CHKERRQ(ierr);
283493642f1SJames Wright     *stg_ctx = *s;
284493642f1SJames Wright     ierr = PetscFree(s); CHKERRQ(ierr);
285493642f1SJames Wright   }
286493642f1SJames Wright 
287493642f1SJames Wright   ierr = ReadSTGInflow(comm, stg_inflow_path, stg_ctx); CHKERRQ(ierr);
288493642f1SJames Wright   ierr = ReadSTGRand(comm, stg_rand_path, stg_ctx); CHKERRQ(ierr);
289493642f1SJames Wright 
290e6098bcdSJames Wright   if (stg_ctx->nynodes > 0) {
291e6098bcdSJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
292e6098bcdSJames Wright     for (PetscInt i=0; i<stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
293e6098bcdSJames Wright   }
294e6098bcdSJames Wright 
295493642f1SJames Wright   // -- Calculate kappa
296493642f1SJames Wright   {
297493642f1SJames Wright     CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
298493642f1SJames Wright     CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw];
299493642f1SJames Wright     CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
300493642f1SJames Wright     CeedScalar le, le_max=0;
301493642f1SJames Wright 
302493642f1SJames Wright     CeedPragmaSIMD
303493642f1SJames Wright     for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
304493642f1SJames Wright       le = PetscMin(2*prof_dw[i], 3*lt[i]);
305493642f1SJames Wright       if (le_max < le) le_max = le;
306493642f1SJames Wright     }
307493642f1SJames Wright     CeedScalar kmin = M_PI/le_max;
308493642f1SJames Wright 
309493642f1SJames Wright     CeedPragmaSIMD
310493642f1SJames Wright     for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
311493642f1SJames Wright       kappa[i] = kmin*pow(stg_ctx->alpha, i);
312493642f1SJames Wright     }
313493642f1SJames Wright   } //end calculate kappa
314493642f1SJames Wright 
31576a21e77SJames Wright   ierr = PetscFree(*pstg_ctx); CHKERRQ(ierr);
316493642f1SJames Wright   *pstg_ctx = stg_ctx;
317493642f1SJames Wright   PetscFunctionReturn(0);
318493642f1SJames Wright }
319493642f1SJames Wright 
320493642f1SJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem,
321493642f1SJames Wright                         User user, const bool prescribe_T,
322e6098bcdSJames Wright                         const CeedScalar theta0, const CeedScalar P0,
323e6098bcdSJames Wright                         const CeedScalar ynodes[], const CeedInt nynodes) {
324493642f1SJames Wright   PetscErrorCode ierr;
325493642f1SJames Wright   char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
326493642f1SJames Wright   char stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
327e6098bcdSJames Wright   PetscBool  mean_only     = PETSC_FALSE,
328e6098bcdSJames Wright              use_stgstrong = PETSC_FALSE;
329e6098bcdSJames Wright   CeedScalar u0            = 0.0,
330e6098bcdSJames Wright              alpha         = 1.01;
331493642f1SJames Wright   STGShur14Context stg_ctx;
332493642f1SJames Wright   CeedQFunctionContext stg_context;
333493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
334493642f1SJames Wright   PetscFunctionBeginUser;
335493642f1SJames Wright 
336493642f1SJames Wright   // Get options
337493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
338493642f1SJames Wright   ierr = PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL,
339493642f1SJames Wright                             stg_inflow_path, stg_inflow_path,
340493642f1SJames Wright                             sizeof(stg_inflow_path), NULL); CHKERRQ(ierr);
341493642f1SJames Wright   ierr = PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL,
342493642f1SJames Wright                             stg_rand_path,stg_rand_path,
343493642f1SJames Wright                             sizeof(stg_rand_path), NULL); CHKERRQ(ierr);
344493642f1SJames Wright   ierr = PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL,
345493642f1SJames Wright                           alpha, &alpha, NULL); CHKERRQ(ierr);
346493642f1SJames Wright   ierr = PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations",
347493642f1SJames Wright                           NULL, u0, &u0, NULL); CHKERRQ(ierr);
348493642f1SJames Wright   ierr = PetscOptionsBool("-stg_mean_only", "Only apply mean profile",
349493642f1SJames Wright                           NULL, mean_only, &mean_only, NULL); CHKERRQ(ierr);
350e6098bcdSJames Wright   ierr = PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly",
351e6098bcdSJames Wright                           NULL, use_stgstrong, &use_stgstrong, NULL); CHKERRQ(ierr);
352493642f1SJames Wright   PetscOptionsEnd();
353493642f1SJames Wright 
354493642f1SJames Wright   ierr = PetscCalloc1(1, &stg_ctx); CHKERRQ(ierr);
355493642f1SJames Wright   stg_ctx->alpha         = alpha;
356493642f1SJames Wright   stg_ctx->u0            = u0;
357493642f1SJames Wright   stg_ctx->is_implicit   = user->phys->implicit;
358493642f1SJames Wright   stg_ctx->prescribe_T   = prescribe_T;
359493642f1SJames Wright   stg_ctx->mean_only     = mean_only;
360493642f1SJames Wright   stg_ctx->theta0        = theta0;
361493642f1SJames Wright   stg_ctx->P0            = P0;
362e6098bcdSJames Wright   stg_ctx->nynodes       = nynodes;
363493642f1SJames Wright 
364493642f1SJames Wright   {
365493642f1SJames Wright     // Calculate dx assuming constant spacing
366493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
367493642f1SJames Wright     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
368493642f1SJames Wright     for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
369493642f1SJames Wright 
370493642f1SJames Wright     PetscInt nmax = 3, faces[3];
371493642f1SJames Wright     ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
372493642f1SJames Wright                                    NULL); CHKERRQ(ierr);
373493642f1SJames Wright     stg_ctx->dx = domain_size[0]/faces[0];
374e6098bcdSJames Wright     stg_ctx->dz = domain_size[2]/faces[2];
375493642f1SJames Wright   }
376493642f1SJames Wright 
377493642f1SJames Wright   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
378493642f1SJames Wright                               CEED_MEM_HOST, &newtonian_ig_ctx);
379493642f1SJames Wright   stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
380493642f1SJames Wright   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
381493642f1SJames Wright                                   &newtonian_ig_ctx);
382493642f1SJames Wright 
383e6098bcdSJames Wright   ierr = GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &stg_ctx,
384e6098bcdSJames Wright                            ynodes); CHKERRQ(ierr);
385493642f1SJames Wright 
386493642f1SJames Wright   CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context);
387493642f1SJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
388493642f1SJames Wright   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST,
389493642f1SJames Wright                               CEED_USE_POINTER, stg_ctx->total_bytes, stg_ctx);
390493642f1SJames Wright   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST,
391493642f1SJames Wright                                      FreeContextPetsc);
392493642f1SJames Wright   CeedQFunctionContextRegisterDouble(stg_context, "solution time",
393493642f1SJames Wright                                      offsetof(struct STGShur14Context_, time), 1,
394493642f1SJames Wright                                      "Phyiscal time of the solution");
395493642f1SJames Wright 
396e6098bcdSJames Wright   if (use_stgstrong) {
397e6098bcdSJames Wright     problem->apply_inflow.qfunction     = STGShur14_Inflow_Strong;
398e6098bcdSJames Wright     problem->apply_inflow.qfunction_loc = STGShur14_Inflow_Strong_loc;
399d7244455SJames Wright     problem->bc_from_ics                = PETSC_FALSE;
400e6098bcdSJames Wright   } else {
401493642f1SJames Wright     problem->apply_inflow.qfunction     = STGShur14_Inflow;
402493642f1SJames Wright     problem->apply_inflow.qfunction_loc = STGShur14_Inflow_loc;
403d7244455SJames Wright     problem->bc_from_ics                = PETSC_TRUE;
404e6098bcdSJames Wright   }
405493642f1SJames Wright   problem->apply_inflow.qfunction_context = stg_context;
406493642f1SJames Wright 
407493642f1SJames Wright   PetscFunctionReturn(0);
408493642f1SJames Wright }
409e6098bcdSJames Wright 
410e6098bcdSJames Wright static inline PetscScalar FindDy(const PetscScalar ynodes[],
411e6098bcdSJames Wright                                  const PetscInt nynodes, const PetscScalar y) {
412*6dd99bedSLeila Ghaffari 
413e6098bcdSJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
414e6098bcdSJames Wright   // ^^assuming min(dy) is first element off the wall
415e6098bcdSJames Wright   PetscInt idx = -1; // Index
416e6098bcdSJames Wright 
417e6098bcdSJames Wright   for (PetscInt i=0; i<nynodes; i++) {
418e6098bcdSJames Wright     if (y < ynodes[i] + half_mindy) {
419e6098bcdSJames Wright       idx = i; break;
420e6098bcdSJames Wright     }
421e6098bcdSJames Wright   }
422e6098bcdSJames Wright   if      (idx == 0)          return ynodes[1] - ynodes[0];
423e6098bcdSJames Wright   else if (idx == nynodes-1)  return ynodes[nynodes-2] - ynodes[nynodes-1];
424e6098bcdSJames Wright   else                        return 0.5 * (ynodes[idx+1] - ynodes[idx-1]);
425e6098bcdSJames Wright }
426e6098bcdSJames Wright 
427e6098bcdSJames Wright // Function passed to DMAddBoundary
428e6098bcdSJames Wright PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time,
429e6098bcdSJames Wright                                const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
430e6098bcdSJames Wright   PetscFunctionBeginUser;
431e6098bcdSJames Wright 
432e6098bcdSJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
433e6098bcdSJames Wright   PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
434e6098bcdSJames Wright   const bool mean_only      = stg_ctx->mean_only;
435e6098bcdSJames Wright   const PetscScalar dx      = stg_ctx->dx;
436e6098bcdSJames Wright   const PetscScalar dz      = stg_ctx->dz;
437e6098bcdSJames Wright   const PetscScalar mu      = stg_ctx->newtonian_ctx.mu;
438e6098bcdSJames Wright   const PetscScalar theta0  = stg_ctx->theta0;
439e6098bcdSJames Wright   const PetscScalar P0      = stg_ctx->P0;
440e6098bcdSJames Wright   const PetscScalar cv      = stg_ctx->newtonian_ctx.cv;
441e6098bcdSJames Wright   const PetscScalar cp      = stg_ctx->newtonian_ctx.cp;
442e6098bcdSJames Wright   const PetscScalar Rd      = cp - cv;
443e6098bcdSJames Wright 
444e6098bcdSJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
445e6098bcdSJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
446e6098bcdSJames Wright   if (!mean_only) {
447e6098bcdSJames Wright     const PetscInt    nynodes = stg_ctx->nynodes;
448e6098bcdSJames Wright     const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes];
449e6098bcdSJames Wright     const PetscScalar h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
450e6098bcdSJames Wright     CalcSpectrum(x[1], eps, lt, h, mu/rho, qn, stg_ctx);
451e6098bcdSJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
452e6098bcdSJames Wright   } else {
453e6098bcdSJames Wright     for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
454e6098bcdSJames Wright   }
455e6098bcdSJames Wright 
456e6098bcdSJames Wright   bcval[0] = rho;
457e6098bcdSJames Wright   bcval[1] = rho * u[0];
458e6098bcdSJames Wright   bcval[2] = rho * u[1];
459e6098bcdSJames Wright   bcval[3] = rho * u[2];
460e6098bcdSJames Wright   PetscFunctionReturn(0);
461e6098bcdSJames Wright }
462e6098bcdSJames Wright 
463e6098bcdSJames Wright PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem,
464e6098bcdSJames Wright                               STGShur14Context stg_ctx) {
465e6098bcdSJames Wright   PetscErrorCode ierr;
466e6098bcdSJames Wright   DMLabel label;
467e6098bcdSJames Wright   const PetscInt comps[] = {0, 1, 2, 3};
468e6098bcdSJames Wright   const PetscInt num_comps = 4;
469e6098bcdSJames Wright   PetscFunctionBeginUser;
470e6098bcdSJames Wright 
471e6098bcdSJames Wright   ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
472e6098bcdSJames Wright   // Set wall BCs
473e6098bcdSJames Wright   if (bc->num_inflow > 0) {
474e6098bcdSJames Wright     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label,
475e6098bcdSJames Wright                          bc->num_inflow, bc->inflows, 0, num_comps,
476e6098bcdSJames Wright                          comps, (void(*)(void))StrongSTGbcFunc,
477e6098bcdSJames Wright                          NULL, stg_ctx, NULL);  CHKERRQ(ierr);
478e6098bcdSJames Wright   }
479e6098bcdSJames Wright 
480e6098bcdSJames Wright   PetscFunctionReturn(0);
481e6098bcdSJames Wright }
482