xref: /honee/problems/stg_shur14.c (revision 493642f1e7bb5ccdccd1086ef1091462e675d35c)
1*493642f1SJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*493642f1SJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*493642f1SJames Wright //
4*493642f1SJames Wright // SPDX-License-Identifier: BSD-2-Clause
5*493642f1SJames Wright //
6*493642f1SJames Wright // This file is part of CEED:  http://github.com/ceed
7*493642f1SJames Wright 
8*493642f1SJames Wright /// @file
9*493642f1SJames Wright /// Implementation of the Synthetic Turbulence Generation (STG) algorithm
10*493642f1SJames Wright /// presented in Shur et al. 2014
11*493642f1SJames Wright 
12*493642f1SJames Wright #include <stdlib.h>
13*493642f1SJames Wright #include <math.h>
14*493642f1SJames Wright #include <petsc.h>
15*493642f1SJames Wright #include "../navierstokes.h"
16*493642f1SJames Wright #include "stg_shur14.h"
17*493642f1SJames Wright #include "../qfunctions/stg_shur14.h"
18*493642f1SJames Wright 
19*493642f1SJames Wright #ifndef M_PI
20*493642f1SJames Wright #define M_PI    3.14159265358979323846
21*493642f1SJames Wright #endif
22*493642f1SJames Wright 
23*493642f1SJames Wright /*
24*493642f1SJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
25*493642f1SJames Wright  *
26*493642f1SJames Wright  * This assumes the input matrices are in order [11,22,33,12,13,23]. This
27*493642f1SJames Wright  * format is also used for the output.
28*493642f1SJames Wright  *
29*493642f1SJames Wright  * @param[in]  comm   MPI_Comm
30*493642f1SJames Wright  * @param[in]  nprofs Number of matrices in Rij
31*493642f1SJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
32*493642f1SJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
33*493642f1SJames Wright  */
34*493642f1SJames Wright PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs,
35*493642f1SJames Wright                                   const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
36*493642f1SJames Wright 
37*493642f1SJames Wright   PetscFunctionBeginUser;
38*493642f1SJames Wright   for (PetscInt i=0; i<nprofs; i++) {
39*493642f1SJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
40*493642f1SJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
41*493642f1SJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2) );
42*493642f1SJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
43*493642f1SJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i]*Cij[4][i]) / Cij[1][i];
44*493642f1SJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
45*493642f1SJames Wright 
46*493642f1SJames Wright     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i]))
47*493642f1SJames Wright       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %d. "
48*493642f1SJames Wright               "Either STGInflow has non-SPD matrix or contains nan.", i+1);
49*493642f1SJames Wright   }
50*493642f1SJames Wright   PetscFunctionReturn(0);
51*493642f1SJames Wright }
52*493642f1SJames Wright 
53*493642f1SJames Wright 
54*493642f1SJames Wright /*
55*493642f1SJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
56*493642f1SJames Wright  *
57*493642f1SJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and
58*493642f1SJames Wright  * passes the file pointer in `fp`. It is not closed in this function, thus
59*493642f1SJames Wright  * `fp` must be closed sometime after this function has been called (using
60*493642f1SJames Wright  * `PetscFClose` for example).
61*493642f1SJames Wright  *
62*493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
63*493642f1SJames Wright  * as the only two entries, separated by a single space
64*493642f1SJames Wright  *
65*493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
66*493642f1SJames Wright  * @param[in] path Path to the file
67*493642f1SJames Wright  * @param[in] char_array_len Length of the character array that should contain each line
68*493642f1SJames Wright  * @param[out] dims Dimensions of the file, taken from the first line of the file
69*493642f1SJames Wright  * @param[out] fp File pointer to the opened file
70*493642f1SJames Wright  */
71*493642f1SJames Wright static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm,
72*493642f1SJames Wright                                         const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len,
73*493642f1SJames Wright                                         PetscInt dims[2], FILE **fp) {
74*493642f1SJames Wright   PetscErrorCode ierr;
75*493642f1SJames Wright   PetscInt ndims;
76*493642f1SJames Wright   char line[char_array_len];
77*493642f1SJames Wright   char **array;
78*493642f1SJames Wright 
79*493642f1SJames Wright   PetscFunctionBeginUser;
80*493642f1SJames Wright   ierr = PetscFOpen(comm, path, "r", fp); CHKERRQ(ierr);
81*493642f1SJames Wright   ierr = PetscSynchronizedFGets(comm, *fp, char_array_len, line); CHKERRQ(ierr);
82*493642f1SJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
83*493642f1SJames Wright   if (ndims != 2) SETERRQ(comm, -1,
84*493642f1SJames Wright                             "Found %d dimensions instead of 2 on the first line of %s",
85*493642f1SJames Wright                             ndims, path);
86*493642f1SJames Wright 
87*493642f1SJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
88*493642f1SJames Wright   ierr = PetscStrToArrayDestroy(ndims, array); CHKERRQ(ierr);
89*493642f1SJames Wright   PetscFunctionReturn(0);
90*493642f1SJames Wright }
91*493642f1SJames Wright 
92*493642f1SJames Wright /*
93*493642f1SJames Wright  * @brief Get the number of rows for the PHASTA file at path
94*493642f1SJames Wright  *
95*493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
96*493642f1SJames Wright  * as the only two entries, separated by a single space
97*493642f1SJames Wright  *
98*493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
99*493642f1SJames Wright  * @param[in] path Path to the file
100*493642f1SJames Wright  * @param[out] nrows Number of rows
101*493642f1SJames Wright  */
102*493642f1SJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm,
103*493642f1SJames Wright                                const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
104*493642f1SJames Wright   PetscErrorCode ierr;
105*493642f1SJames Wright   const PetscInt char_array_len = 512;
106*493642f1SJames Wright   PetscInt dims[2];
107*493642f1SJames Wright   FILE *fp;
108*493642f1SJames Wright 
109*493642f1SJames Wright   PetscFunctionBeginUser;
110*493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
111*493642f1SJames Wright   *nrows = dims[0];
112*493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
113*493642f1SJames Wright   PetscFunctionReturn(0);
114*493642f1SJames Wright }
115*493642f1SJames Wright 
116*493642f1SJames Wright /*
117*493642f1SJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
118*493642f1SJames Wright  *
119*493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
120*493642f1SJames Wright  * as the only two entries, separated by a single space.
121*493642f1SJames Wright  * Assumes there are 14 columns in the file
122*493642f1SJames Wright  *
123*493642f1SJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
124*493642f1SJames Wright  * profile found in the file
125*493642f1SJames Wright  *
126*493642f1SJames Wright  * @param[in] comm MPI_Comm for the program
127*493642f1SJames Wright  * @param[in] path Path to the STGInflow.dat file
128*493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
129*493642f1SJames Wright  */
130*493642f1SJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm,
131*493642f1SJames Wright                                     const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
132*493642f1SJames Wright   PetscErrorCode ierr;
133*493642f1SJames Wright   PetscInt ndims, dims[2];
134*493642f1SJames Wright   FILE *fp;
135*493642f1SJames Wright   const PetscInt char_array_len=512;
136*493642f1SJames Wright   char line[char_array_len];
137*493642f1SJames Wright   char **array;
138*493642f1SJames Wright 
139*493642f1SJames Wright   PetscFunctionBeginUser;
140*493642f1SJames Wright 
141*493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
142*493642f1SJames Wright 
143*493642f1SJames Wright   CeedScalar rij[6][stg_ctx->nprofs];
144*493642f1SJames Wright   CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw];
145*493642f1SJames Wright   CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps];
146*493642f1SJames Wright   CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
147*493642f1SJames Wright   CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs])
148*493642f1SJames Wright                                         &stg_ctx->data[stg_ctx->offsets.ubar];
149*493642f1SJames Wright 
150*493642f1SJames Wright   for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
151*493642f1SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
152*493642f1SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
153*493642f1SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
154*493642f1SJames Wright                                    "Line %d of %s does not contain enough columns (%d instead of %d)", i,
155*493642f1SJames Wright                                    path, ndims, dims[1]);
156*493642f1SJames Wright 
157*493642f1SJames Wright     prof_dw[i] = (CeedScalar) atof(array[0]);
158*493642f1SJames Wright     ubar[0][i] = (CeedScalar) atof(array[1]);
159*493642f1SJames Wright     ubar[1][i] = (CeedScalar) atof(array[2]);
160*493642f1SJames Wright     ubar[2][i] = (CeedScalar) atof(array[3]);
161*493642f1SJames Wright     rij[0][i]  = (CeedScalar) atof(array[4]);
162*493642f1SJames Wright     rij[1][i]  = (CeedScalar) atof(array[5]);
163*493642f1SJames Wright     rij[2][i]  = (CeedScalar) atof(array[6]);
164*493642f1SJames Wright     rij[3][i]  = (CeedScalar) atof(array[7]);
165*493642f1SJames Wright     rij[4][i]  = (CeedScalar) atof(array[8]);
166*493642f1SJames Wright     rij[5][i]  = (CeedScalar) atof(array[9]);
167*493642f1SJames Wright     lt[i]      = (CeedScalar) atof(array[12]);
168*493642f1SJames Wright     eps[i]     = (CeedScalar) atof(array[13]);
169*493642f1SJames Wright 
170*493642f1SJames Wright     if (prof_dw[i] < 0) SETERRQ(comm, -1,
171*493642f1SJames Wright                                   "Distance to wall in %s cannot be negative", path);
172*493642f1SJames Wright     if (lt[i] < 0) SETERRQ(comm, -1,
173*493642f1SJames Wright                              "Turbulent length scale in %s cannot be negative", path);
174*493642f1SJames Wright     if (eps[i] < 0) SETERRQ(comm, -1,
175*493642f1SJames Wright                               "Turbulent dissipation in %s cannot be negative", path);
176*493642f1SJames Wright 
177*493642f1SJames Wright   }
178*493642f1SJames Wright   CeedScalar (*cij)[stg_ctx->nprofs]  = (CeedScalar (*)[stg_ctx->nprofs])
179*493642f1SJames Wright                                         &stg_ctx->data[stg_ctx->offsets.cij];
180*493642f1SJames Wright   ierr = CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij); CHKERRQ(ierr);
181*493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
182*493642f1SJames Wright   PetscFunctionReturn(0);
183*493642f1SJames Wright }
184*493642f1SJames Wright 
185*493642f1SJames Wright /*
186*493642f1SJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
187*493642f1SJames Wright  *
188*493642f1SJames Wright  * Assumes that the first line of the file has the number of rows and columns
189*493642f1SJames Wright  * as the only two entries, separated by a single space.
190*493642f1SJames Wright  * Assumes there are 7 columns in the file
191*493642f1SJames Wright  *
192*493642f1SJames Wright  * @param[in]    comm    MPI_Comm for the program
193*493642f1SJames Wright  * @param[in]    path    Path to the STGRand.dat file
194*493642f1SJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
195*493642f1SJames Wright  */
196*493642f1SJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm,
197*493642f1SJames Wright                                   const char path[PETSC_MAX_PATH_LEN],
198*493642f1SJames Wright                                   STGShur14Context stg_ctx) {
199*493642f1SJames Wright 
200*493642f1SJames Wright   PetscErrorCode ierr;
201*493642f1SJames Wright   PetscInt ndims, dims[2];
202*493642f1SJames Wright   FILE *fp;
203*493642f1SJames Wright   const PetscInt char_array_len = 512;
204*493642f1SJames Wright   char line[char_array_len];
205*493642f1SJames Wright   char **array;
206*493642f1SJames Wright 
207*493642f1SJames Wright   PetscFunctionBeginUser;
208*493642f1SJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
209*493642f1SJames Wright 
210*493642f1SJames Wright   CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi];
211*493642f1SJames Wright   CeedScalar (*d)[stg_ctx->nmodes]     = (CeedScalar (*)[stg_ctx->nmodes])
212*493642f1SJames Wright                                          &stg_ctx->data[stg_ctx->offsets.d];
213*493642f1SJames Wright   CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes])
214*493642f1SJames Wright                                          &stg_ctx->data[stg_ctx->offsets.sigma];
215*493642f1SJames Wright 
216*493642f1SJames Wright   for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
217*493642f1SJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
218*493642f1SJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
219*493642f1SJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
220*493642f1SJames Wright                                    "Line %d of %s does not contain enough columns (%d instead of %d)", i,
221*493642f1SJames Wright                                    path, ndims, dims[1]);
222*493642f1SJames Wright 
223*493642f1SJames Wright     d[0][i]     = (CeedScalar) atof(array[0]);
224*493642f1SJames Wright     d[1][i]     = (CeedScalar) atof(array[1]);
225*493642f1SJames Wright     d[2][i]     = (CeedScalar) atof(array[2]);
226*493642f1SJames Wright     phi[i]      = (CeedScalar) atof(array[3]);
227*493642f1SJames Wright     sigma[0][i] = (CeedScalar) atof(array[4]);
228*493642f1SJames Wright     sigma[1][i] = (CeedScalar) atof(array[5]);
229*493642f1SJames Wright     sigma[2][i] = (CeedScalar) atof(array[6]);
230*493642f1SJames Wright   }
231*493642f1SJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
232*493642f1SJames Wright   PetscFunctionReturn(0);
233*493642f1SJames Wright }
234*493642f1SJames Wright 
235*493642f1SJames Wright /*
236*493642f1SJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
237*493642f1SJames Wright  *
238*493642f1SJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
239*493642f1SJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
240*493642f1SJames Wright  * STGShur14Context instance.
241*493642f1SJames Wright  *
242*493642f1SJames Wright  * @param[in]    comm            MPI_Comm for the program
243*493642f1SJames Wright  * @param[in]    dm              DM for the program
244*493642f1SJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
245*493642f1SJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
246*493642f1SJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
247*493642f1SJames Wright  */
248*493642f1SJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm,
249*493642f1SJames Wright                                  char stg_inflow_path[PETSC_MAX_PATH_LEN],
250*493642f1SJames Wright                                  char stg_rand_path[PETSC_MAX_PATH_LEN],
251*493642f1SJames Wright                                  STGShur14Context *pstg_ctx) {
252*493642f1SJames Wright   PetscErrorCode ierr;
253*493642f1SJames Wright   PetscInt nmodes, nprofs;
254*493642f1SJames Wright   STGShur14Context stg_ctx;
255*493642f1SJames Wright   PetscFunctionBeginUser;
256*493642f1SJames Wright 
257*493642f1SJames Wright   // Get options
258*493642f1SJames Wright   ierr = GetNRows(comm, stg_rand_path, &nmodes); CHKERRQ(ierr);
259*493642f1SJames Wright   ierr = GetNRows(comm, stg_inflow_path, &nprofs); CHKERRQ(ierr);
260*493642f1SJames Wright   if (nmodes > STG_NMODES_MAX)
261*493642f1SJames Wright     SETERRQ(comm, 1, "Number of wavemodes in %s (%d) exceeds STG_NMODES_MAX (%d). "
262*493642f1SJames Wright             "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes,
263*493642f1SJames Wright             STG_NMODES_MAX);
264*493642f1SJames Wright 
265*493642f1SJames Wright   {
266*493642f1SJames Wright     STGShur14Context s;
267*493642f1SJames Wright     ierr = PetscCalloc1(1, &s); CHKERRQ(ierr);
268*493642f1SJames Wright     *s = **pstg_ctx;
269*493642f1SJames Wright     s->nmodes = nmodes;
270*493642f1SJames Wright     s->nprofs = nprofs;
271*493642f1SJames Wright     s->offsets.sigma   = 0;
272*493642f1SJames Wright     s->offsets.d       = nmodes*3;
273*493642f1SJames Wright     s->offsets.phi     = s->offsets.d       + nmodes*3;
274*493642f1SJames Wright     s->offsets.kappa   = s->offsets.phi     + nmodes;
275*493642f1SJames Wright     s->offsets.prof_dw = s->offsets.kappa   + nmodes;
276*493642f1SJames Wright     s->offsets.ubar    = s->offsets.prof_dw + nprofs;
277*493642f1SJames Wright     s->offsets.cij     = s->offsets.ubar    + nprofs*3;
278*493642f1SJames Wright     s->offsets.eps     = s->offsets.cij     + nprofs*6;
279*493642f1SJames Wright     s->offsets.lt      = s->offsets.eps     + nprofs;
280*493642f1SJames Wright     PetscInt total_num_scalars = s->offsets.lt + nprofs;
281*493642f1SJames Wright     s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]);
282*493642f1SJames Wright     ierr = PetscMalloc(s->total_bytes, &stg_ctx); CHKERRQ(ierr);
283*493642f1SJames Wright     *stg_ctx = *s;
284*493642f1SJames Wright     ierr = PetscFree(s); CHKERRQ(ierr);
285*493642f1SJames Wright   }
286*493642f1SJames Wright 
287*493642f1SJames Wright   ierr = ReadSTGInflow(comm, stg_inflow_path, stg_ctx); CHKERRQ(ierr);
288*493642f1SJames Wright   ierr = ReadSTGRand(comm, stg_rand_path, stg_ctx); CHKERRQ(ierr);
289*493642f1SJames Wright 
290*493642f1SJames Wright   // -- Calculate kappa
291*493642f1SJames Wright   {
292*493642f1SJames Wright     CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
293*493642f1SJames Wright     CeedScalar *prof_dw = &stg_ctx->data[stg_ctx->offsets.prof_dw];
294*493642f1SJames Wright     CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
295*493642f1SJames Wright     CeedScalar le, le_max=0;
296*493642f1SJames Wright 
297*493642f1SJames Wright     CeedPragmaSIMD
298*493642f1SJames Wright     for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
299*493642f1SJames Wright       le = PetscMin(2*prof_dw[i], 3*lt[i]);
300*493642f1SJames Wright       if (le_max < le) le_max = le;
301*493642f1SJames Wright     }
302*493642f1SJames Wright     CeedScalar kmin = M_PI/le_max;
303*493642f1SJames Wright 
304*493642f1SJames Wright     CeedPragmaSIMD
305*493642f1SJames Wright     for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
306*493642f1SJames Wright       kappa[i] = kmin*pow(stg_ctx->alpha, i);
307*493642f1SJames Wright     }
308*493642f1SJames Wright   } //end calculate kappa
309*493642f1SJames Wright 
310*493642f1SJames Wright   *pstg_ctx = stg_ctx;
311*493642f1SJames Wright   PetscFunctionReturn(0);
312*493642f1SJames Wright }
313*493642f1SJames Wright 
314*493642f1SJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem,
315*493642f1SJames Wright                         User user, const bool prescribe_T,
316*493642f1SJames Wright                         const CeedScalar theta0, const CeedScalar P0) {
317*493642f1SJames Wright   PetscErrorCode ierr;
318*493642f1SJames Wright   char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
319*493642f1SJames Wright   char stg_rand_path[PETSC_MAX_PATH_LEN] = "./STGRand.dat";
320*493642f1SJames Wright   PetscBool mean_only = PETSC_FALSE;
321*493642f1SJames Wright   CeedScalar u0=0.0, alpha=1.01;
322*493642f1SJames Wright   STGShur14Context stg_ctx;
323*493642f1SJames Wright   CeedQFunctionContext stg_context;
324*493642f1SJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
325*493642f1SJames Wright   PetscFunctionBeginUser;
326*493642f1SJames Wright 
327*493642f1SJames Wright   // Get options
328*493642f1SJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
329*493642f1SJames Wright   ierr = PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL,
330*493642f1SJames Wright                             stg_inflow_path, stg_inflow_path,
331*493642f1SJames Wright                             sizeof(stg_inflow_path), NULL); CHKERRQ(ierr);
332*493642f1SJames Wright   ierr = PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL,
333*493642f1SJames Wright                             stg_rand_path,stg_rand_path,
334*493642f1SJames Wright                             sizeof(stg_rand_path), NULL); CHKERRQ(ierr);
335*493642f1SJames Wright   ierr = PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL,
336*493642f1SJames Wright                           alpha, &alpha, NULL); CHKERRQ(ierr);
337*493642f1SJames Wright   ierr = PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations",
338*493642f1SJames Wright                           NULL, u0, &u0, NULL); CHKERRQ(ierr);
339*493642f1SJames Wright   ierr = PetscOptionsBool("-stg_mean_only", "Only apply mean profile",
340*493642f1SJames Wright                           NULL, mean_only, &mean_only, NULL); CHKERRQ(ierr);
341*493642f1SJames Wright   PetscOptionsEnd();
342*493642f1SJames Wright 
343*493642f1SJames Wright   ierr = PetscCalloc1(1, &stg_ctx); CHKERRQ(ierr);
344*493642f1SJames Wright   stg_ctx->alpha         = alpha;
345*493642f1SJames Wright   stg_ctx->u0            = u0;
346*493642f1SJames Wright   stg_ctx->is_implicit   = user->phys->implicit;
347*493642f1SJames Wright   stg_ctx->prescribe_T   = prescribe_T;
348*493642f1SJames Wright   stg_ctx->mean_only     = mean_only;
349*493642f1SJames Wright   stg_ctx->theta0        = theta0;
350*493642f1SJames Wright   stg_ctx->P0            = P0;
351*493642f1SJames Wright 
352*493642f1SJames Wright   {
353*493642f1SJames Wright     // Calculate dx assuming constant spacing
354*493642f1SJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
355*493642f1SJames Wright     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
356*493642f1SJames Wright     for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
357*493642f1SJames Wright 
358*493642f1SJames Wright     PetscInt nmax = 3, faces[3];
359*493642f1SJames Wright     ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
360*493642f1SJames Wright                                    NULL); CHKERRQ(ierr);
361*493642f1SJames Wright     stg_ctx->dx = domain_size[0]/faces[0];
362*493642f1SJames Wright   }
363*493642f1SJames Wright 
364*493642f1SJames Wright   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
365*493642f1SJames Wright                               CEED_MEM_HOST, &newtonian_ig_ctx);
366*493642f1SJames Wright   stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
367*493642f1SJames Wright   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
368*493642f1SJames Wright                                   &newtonian_ig_ctx);
369*493642f1SJames Wright 
370*493642f1SJames Wright   ierr = GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path, &stg_ctx);
371*493642f1SJames Wright   CHKERRQ(ierr);
372*493642f1SJames Wright 
373*493642f1SJames Wright   CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context);
374*493642f1SJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
375*493642f1SJames Wright   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST,
376*493642f1SJames Wright                               CEED_USE_POINTER, stg_ctx->total_bytes, stg_ctx);
377*493642f1SJames Wright   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST,
378*493642f1SJames Wright                                      FreeContextPetsc);
379*493642f1SJames Wright   CeedQFunctionContextRegisterDouble(stg_context, "solution time",
380*493642f1SJames Wright                                      offsetof(struct STGShur14Context_, time), 1,
381*493642f1SJames Wright                                      "Phyiscal time of the solution");
382*493642f1SJames Wright 
383*493642f1SJames Wright   problem->apply_inflow.qfunction         = STGShur14_Inflow;
384*493642f1SJames Wright   problem->apply_inflow.qfunction_loc     = STGShur14_Inflow_loc;
385*493642f1SJames Wright   problem->apply_inflow.qfunction_context = stg_context;
386*493642f1SJames Wright 
387*493642f1SJames Wright   PetscFunctionReturn(0);
388*493642f1SJames Wright }
389