xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision 175f00a62957e4f3f8aaa73c4289f1b68f0117f5)
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 
12ba6664aeSJames Wright #include <stdlib.h>
13ba6664aeSJames Wright #include <math.h>
14ba6664aeSJames Wright #include <petsc.h>
15ba6664aeSJames Wright #include "../navierstokes.h"
16ba6664aeSJames Wright #include "stg_shur14.h"
17ba6664aeSJames Wright #include "../qfunctions/stg_shur14.h"
18ba6664aeSJames Wright 
1965dd5cafSJames Wright STGShur14Context global_stg_ctx;
2065dd5cafSJames Wright 
21ba6664aeSJames Wright /*
22ba6664aeSJames Wright  * @brief Perform Cholesky decomposition on array of symmetric 3x3 matrices
23ba6664aeSJames Wright  *
24ba6664aeSJames Wright  * This assumes the input matrices are in order [11,22,33,12,13,23]. This
25ba6664aeSJames Wright  * format is also used for the output.
26ba6664aeSJames Wright  *
27ba6664aeSJames Wright  * @param[in]  comm   MPI_Comm
28ba6664aeSJames Wright  * @param[in]  nprofs Number of matrices in Rij
29ba6664aeSJames Wright  * @param[in]  Rij    Array of the symmetric matrices [6,nprofs]
30ba6664aeSJames Wright  * @param[out] Cij    Array of the Cholesky Decomposition matrices, [6,nprofs]
31ba6664aeSJames Wright  */
32ba6664aeSJames Wright PetscErrorCode CalcCholeskyDecomp(MPI_Comm comm, PetscInt nprofs,
33ba6664aeSJames Wright                                   const CeedScalar Rij[6][nprofs], CeedScalar Cij[6][nprofs]) {
34ba6664aeSJames Wright   PetscFunctionBeginUser;
35ba6664aeSJames Wright   for (PetscInt i=0; i<nprofs; i++) {
36ba6664aeSJames Wright     Cij[0][i] = sqrt(Rij[0][i]);
37ba6664aeSJames Wright     Cij[3][i] = Rij[3][i] / Cij[0][i];
38ba6664aeSJames Wright     Cij[1][i] = sqrt(Rij[1][i] - pow(Cij[3][i], 2) );
39ba6664aeSJames Wright     Cij[4][i] = Rij[4][i] / Cij[0][i];
40ba6664aeSJames Wright     Cij[5][i] = (Rij[5][i] - Cij[3][i]*Cij[4][i]) / Cij[1][i];
41ba6664aeSJames Wright     Cij[2][i] = sqrt(Rij[2][i] - pow(Cij[4][i], 2) - pow(Cij[5][i], 2));
42ba6664aeSJames Wright 
43ba6664aeSJames Wright     if (isnan(Cij[0][i]) || isnan(Cij[1][i]) || isnan(Cij[2][i]))
44990fdeb6SJeremy L Thompson       SETERRQ(comm, -1, "Cholesky decomposition failed at profile point %"
45990fdeb6SJeremy L Thompson               PetscInt_FMT ". Either STGInflow has non-SPD matrix or contains nan.", i+1);
46ba6664aeSJames Wright   }
47ba6664aeSJames Wright   PetscFunctionReturn(0);
48ba6664aeSJames Wright }
49ba6664aeSJames Wright 
50ba6664aeSJames Wright 
51ba6664aeSJames Wright /*
52ba6664aeSJames Wright  * @brief Open a PHASTA *.dat file, grabbing dimensions and file pointer
53ba6664aeSJames Wright  *
54ba6664aeSJames Wright  * This function opens the file specified by `path` using `PetscFOpen` and
55ba6664aeSJames Wright  * passes the file pointer in `fp`. It is not closed in this function, thus
56ba6664aeSJames Wright  * `fp` must be closed sometime after this function has been called (using
57ba6664aeSJames Wright  * `PetscFClose` for example).
58ba6664aeSJames Wright  *
59ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
60ba6664aeSJames Wright  * as the only two entries, separated by a single space
61ba6664aeSJames Wright  *
62ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
63ba6664aeSJames Wright  * @param[in] path Path to the file
64ba6664aeSJames Wright  * @param[in] char_array_len Length of the character array that should contain each line
65ba6664aeSJames Wright  * @param[out] dims Dimensions of the file, taken from the first line of the file
66ba6664aeSJames Wright  * @param[out] fp File pointer to the opened file
67ba6664aeSJames Wright  */
68ba6664aeSJames Wright static PetscErrorCode OpenPHASTADatFile(const MPI_Comm comm,
69ba6664aeSJames Wright                                         const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len,
70ba6664aeSJames Wright                                         PetscInt dims[2], FILE **fp) {
71ba6664aeSJames Wright   PetscErrorCode ierr;
72ba6664aeSJames Wright   PetscInt ndims;
73ba6664aeSJames Wright   char line[char_array_len];
74ba6664aeSJames Wright   char **array;
75ba6664aeSJames Wright 
76ba6664aeSJames Wright   PetscFunctionBeginUser;
77ba6664aeSJames Wright   ierr = PetscFOpen(comm, path, "r", fp); CHKERRQ(ierr);
78ba6664aeSJames Wright   ierr = PetscSynchronizedFGets(comm, *fp, char_array_len, line); CHKERRQ(ierr);
79ba6664aeSJames Wright   ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
80ba6664aeSJames Wright   if (ndims != 2) SETERRQ(comm, -1,
81990fdeb6SJeremy L Thompson                             "Found %" PetscInt_FMT" dimensions instead of 2 on the first line of %s",
82ba6664aeSJames Wright                             ndims, path);
83ba6664aeSJames Wright 
84ba6664aeSJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
85ba6664aeSJames Wright   ierr = PetscStrToArrayDestroy(ndims, array); CHKERRQ(ierr);
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  */
99ba6664aeSJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm,
100ba6664aeSJames Wright                                const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
101ba6664aeSJames Wright   PetscErrorCode ierr;
102ba6664aeSJames Wright   const PetscInt char_array_len = 512;
103ba6664aeSJames Wright   PetscInt dims[2];
104ba6664aeSJames Wright   FILE *fp;
105ba6664aeSJames Wright 
106ba6664aeSJames Wright   PetscFunctionBeginUser;
107ba6664aeSJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
108ba6664aeSJames Wright   *nrows = dims[0];
109ba6664aeSJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
110ba6664aeSJames Wright   PetscFunctionReturn(0);
111ba6664aeSJames Wright }
112ba6664aeSJames Wright 
113ba6664aeSJames Wright /*
114ba6664aeSJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
115ba6664aeSJames Wright  *
116ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
117ba6664aeSJames Wright  * as the only two entries, separated by a single space.
118ba6664aeSJames Wright  * Assumes there are 14 columns in the file
119ba6664aeSJames Wright  *
120ba6664aeSJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
121ba6664aeSJames Wright  * profile found in the file
122ba6664aeSJames Wright  *
123ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
124ba6664aeSJames Wright  * @param[in] path Path to the STGInflow.dat file
125ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
126ba6664aeSJames Wright  */
127ba6664aeSJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm,
128ba6664aeSJames Wright                                     const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
129ba6664aeSJames Wright   PetscErrorCode ierr;
130ba6664aeSJames Wright   PetscInt ndims, dims[2];
131ba6664aeSJames Wright   FILE *fp;
132ba6664aeSJames Wright   const PetscInt char_array_len=512;
133ba6664aeSJames Wright   char line[char_array_len];
134ba6664aeSJames Wright   char **array;
135ba6664aeSJames Wright 
136ba6664aeSJames Wright   PetscFunctionBeginUser;
137ba6664aeSJames Wright 
138ba6664aeSJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
139ba6664aeSJames Wright 
140ba6664aeSJames Wright   CeedScalar rij[6][stg_ctx->nprofs];
141*175f00a6SJames Wright   CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
142ba6664aeSJames Wright   CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps];
143ba6664aeSJames Wright   CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
144ba6664aeSJames Wright   CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs])
145ba6664aeSJames Wright                                         &stg_ctx->data[stg_ctx->offsets.ubar];
146ba6664aeSJames Wright 
147ba6664aeSJames Wright   for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
148ba6664aeSJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
149ba6664aeSJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
150ba6664aeSJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
151990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
152990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT")", i,
153ba6664aeSJames Wright                                    path, ndims, dims[1]);
154ba6664aeSJames Wright 
155*175f00a6SJames Wright     wall_dist[i] = (CeedScalar) atof(array[0]);
156ba6664aeSJames Wright     ubar[0][i]   = (CeedScalar) atof(array[1]);
157ba6664aeSJames Wright     ubar[1][i]   = (CeedScalar) atof(array[2]);
158ba6664aeSJames Wright     ubar[2][i]   = (CeedScalar) atof(array[3]);
159ba6664aeSJames Wright     rij[0][i]    = (CeedScalar) atof(array[4]);
160ba6664aeSJames Wright     rij[1][i]    = (CeedScalar) atof(array[5]);
161ba6664aeSJames Wright     rij[2][i]    = (CeedScalar) atof(array[6]);
162ba6664aeSJames Wright     rij[3][i]    = (CeedScalar) atof(array[7]);
163ba6664aeSJames Wright     rij[4][i]    = (CeedScalar) atof(array[8]);
164ba6664aeSJames Wright     rij[5][i]    = (CeedScalar) atof(array[9]);
165ba6664aeSJames Wright     lt[i]        = (CeedScalar) atof(array[12]);
166ba6664aeSJames Wright     eps[i]       = (CeedScalar) atof(array[13]);
167ba6664aeSJames Wright 
168*175f00a6SJames Wright     if (wall_dist[i] < 0) SETERRQ(comm, -1,
169ba6664aeSJames Wright                                     "Distance to wall in %s cannot be negative", path);
170ba6664aeSJames Wright     if (lt[i] < 0) SETERRQ(comm, -1,
171ba6664aeSJames Wright                              "Turbulent length scale in %s cannot be negative", path);
172ba6664aeSJames Wright     if (eps[i] < 0) SETERRQ(comm, -1,
173ba6664aeSJames Wright                               "Turbulent dissipation in %s cannot be negative", path);
174ba6664aeSJames Wright 
175ba6664aeSJames Wright   }
176ba6664aeSJames Wright   CeedScalar (*cij)[stg_ctx->nprofs]  = (CeedScalar (*)[stg_ctx->nprofs])
177ba6664aeSJames Wright                                         &stg_ctx->data[stg_ctx->offsets.cij];
178ba6664aeSJames Wright   ierr = CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij); CHKERRQ(ierr);
179ba6664aeSJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
180ba6664aeSJames Wright   PetscFunctionReturn(0);
181ba6664aeSJames Wright }
182ba6664aeSJames Wright 
183ba6664aeSJames Wright /*
184ba6664aeSJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
185ba6664aeSJames Wright  *
186ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
187ba6664aeSJames Wright  * as the only two entries, separated by a single space.
188ba6664aeSJames Wright  * Assumes there are 7 columns in the file
189ba6664aeSJames Wright  *
190ba6664aeSJames Wright  * @param[in]    comm    MPI_Comm for the program
191ba6664aeSJames Wright  * @param[in]    path    Path to the STGRand.dat file
192ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
193ba6664aeSJames Wright  */
194ba6664aeSJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm,
195ba6664aeSJames Wright                                   const char path[PETSC_MAX_PATH_LEN],
196ba6664aeSJames Wright                                   STGShur14Context stg_ctx) {
197ba6664aeSJames Wright   PetscErrorCode ierr;
198ba6664aeSJames Wright   PetscInt ndims, dims[2];
199ba6664aeSJames Wright   FILE *fp;
200ba6664aeSJames Wright   const PetscInt char_array_len = 512;
201ba6664aeSJames Wright   char line[char_array_len];
202ba6664aeSJames Wright   char **array;
203ba6664aeSJames Wright 
204ba6664aeSJames Wright   PetscFunctionBeginUser;
205ba6664aeSJames Wright   ierr = OpenPHASTADatFile(comm, path, char_array_len, dims, &fp); CHKERRQ(ierr);
206ba6664aeSJames Wright 
207ba6664aeSJames Wright   CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi];
208ba6664aeSJames Wright   CeedScalar (*d)[stg_ctx->nmodes]     = (CeedScalar (*)[stg_ctx->nmodes])
209ba6664aeSJames Wright                                          &stg_ctx->data[stg_ctx->offsets.d];
210ba6664aeSJames Wright   CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes])
211ba6664aeSJames Wright                                          &stg_ctx->data[stg_ctx->offsets.sigma];
212ba6664aeSJames Wright 
213ba6664aeSJames Wright   for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
214ba6664aeSJames Wright     ierr = PetscSynchronizedFGets(comm, fp, char_array_len, line); CHKERRQ(ierr);
215ba6664aeSJames Wright     ierr = PetscStrToArray(line, ' ', &ndims, &array); CHKERRQ(ierr);
216ba6664aeSJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
217990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
218990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT")", i,
219ba6664aeSJames Wright                                    path, ndims, dims[1]);
220ba6664aeSJames Wright 
221ba6664aeSJames Wright     d[0][i]     = (CeedScalar) atof(array[0]);
222ba6664aeSJames Wright     d[1][i]     = (CeedScalar) atof(array[1]);
223ba6664aeSJames Wright     d[2][i]     = (CeedScalar) atof(array[2]);
224ba6664aeSJames Wright     phi[i]      = (CeedScalar) atof(array[3]);
225ba6664aeSJames Wright     sigma[0][i] = (CeedScalar) atof(array[4]);
226ba6664aeSJames Wright     sigma[1][i] = (CeedScalar) atof(array[5]);
227ba6664aeSJames Wright     sigma[2][i] = (CeedScalar) atof(array[6]);
228ba6664aeSJames Wright   }
229ba6664aeSJames Wright   ierr = PetscFClose(comm, fp); CHKERRQ(ierr);
230ba6664aeSJames Wright   PetscFunctionReturn(0);
231ba6664aeSJames Wright }
232ba6664aeSJames Wright 
233ba6664aeSJames Wright /*
234ba6664aeSJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
235ba6664aeSJames Wright  *
236ba6664aeSJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
237ba6664aeSJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
238ba6664aeSJames Wright  * STGShur14Context instance.
239ba6664aeSJames Wright  *
240ba6664aeSJames Wright  * @param[in]    comm            MPI_Comm for the program
241ba6664aeSJames Wright  * @param[in]    dm              DM for the program
242ba6664aeSJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
243ba6664aeSJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
244ba6664aeSJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
245ba6664aeSJames Wright  */
246ba6664aeSJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm,
247ba6664aeSJames Wright                                  char stg_inflow_path[PETSC_MAX_PATH_LEN],
248ba6664aeSJames Wright                                  char stg_rand_path[PETSC_MAX_PATH_LEN],
2494e139266SJames Wright                                  STGShur14Context *pstg_ctx,
2504e139266SJames Wright                                  const CeedScalar ynodes[]) {
251ba6664aeSJames Wright   PetscErrorCode ierr;
252ba6664aeSJames Wright   PetscInt nmodes, nprofs;
253ba6664aeSJames Wright   STGShur14Context stg_ctx;
254ba6664aeSJames Wright   PetscFunctionBeginUser;
255ba6664aeSJames Wright 
256ba6664aeSJames Wright   // Get options
257ba6664aeSJames Wright   ierr = GetNRows(comm, stg_rand_path, &nmodes); CHKERRQ(ierr);
258ba6664aeSJames Wright   ierr = GetNRows(comm, stg_inflow_path, &nprofs); CHKERRQ(ierr);
259ba6664aeSJames Wright   if (nmodes > STG_NMODES_MAX)
260990fdeb6SJeremy L Thompson     SETERRQ(comm, 1, "Number of wavemodes in %s (%"
261990fdeb6SJeremy L Thompson             PetscInt_FMT") exceeds STG_NMODES_MAX (%" PetscInt_FMT"). "
262ba6664aeSJames Wright             "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes,
263ba6664aeSJames Wright             STG_NMODES_MAX);
264ba6664aeSJames Wright 
265ba6664aeSJames Wright   {
266ba6664aeSJames Wright     STGShur14Context s;
267ba6664aeSJames Wright     ierr = PetscCalloc1(1, &s); CHKERRQ(ierr);
268ba6664aeSJames Wright     *s = **pstg_ctx;
269ba6664aeSJames Wright     s->nmodes = nmodes;
270ba6664aeSJames Wright     s->nprofs = nprofs;
271ba6664aeSJames Wright     s->offsets.sigma   = 0;
272ba6664aeSJames Wright     s->offsets.d       = nmodes*3;
273ba6664aeSJames Wright     s->offsets.phi       = s->offsets.d         + nmodes*3;
274ba6664aeSJames Wright     s->offsets.kappa     = s->offsets.phi       + nmodes;
275*175f00a6SJames Wright     s->offsets.wall_dist = s->offsets.kappa     + nmodes;
276*175f00a6SJames Wright     s->offsets.ubar      = s->offsets.wall_dist + nprofs;
277ba6664aeSJames Wright     s->offsets.cij       = s->offsets.ubar      + nprofs*3;
278ba6664aeSJames Wright     s->offsets.eps       = s->offsets.cij       + nprofs*6;
279ba6664aeSJames Wright     s->offsets.lt        = s->offsets.eps       + nprofs;
2804e139266SJames Wright     s->offsets.ynodes    = s->offsets.lt        + nprofs;
2814e139266SJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
282ba6664aeSJames Wright     s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]);
283ba6664aeSJames Wright     ierr = PetscMalloc(s->total_bytes, &stg_ctx); CHKERRQ(ierr);
284ba6664aeSJames Wright     *stg_ctx = *s;
285ba6664aeSJames Wright     ierr = PetscFree(s); CHKERRQ(ierr);
286ba6664aeSJames Wright   }
287ba6664aeSJames Wright 
288ba6664aeSJames Wright   ierr = ReadSTGInflow(comm, stg_inflow_path, stg_ctx); CHKERRQ(ierr);
289ba6664aeSJames Wright   ierr = ReadSTGRand(comm, stg_rand_path, stg_ctx); CHKERRQ(ierr);
290ba6664aeSJames Wright 
2914e139266SJames Wright   if (stg_ctx->nynodes > 0) {
2924e139266SJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
2934e139266SJames Wright     for (PetscInt i=0; i<stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
2944e139266SJames Wright   }
2954e139266SJames Wright 
296ba6664aeSJames Wright   // -- Calculate kappa
297ba6664aeSJames Wright   {
298ba6664aeSJames Wright     CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
299*175f00a6SJames Wright     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
300ba6664aeSJames Wright     CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
301ba6664aeSJames Wright     CeedScalar le, le_max=0;
302ba6664aeSJames Wright 
303ba6664aeSJames Wright     CeedPragmaSIMD
304ba6664aeSJames Wright     for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
305*175f00a6SJames Wright       le = PetscMin(2*wall_dist[i], 3*lt[i]);
306ba6664aeSJames Wright       if (le_max < le) le_max = le;
307ba6664aeSJames Wright     }
308ba6664aeSJames Wright     CeedScalar kmin = M_PI/le_max;
309ba6664aeSJames Wright 
310ba6664aeSJames Wright     CeedPragmaSIMD
311ba6664aeSJames Wright     for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
312ba6664aeSJames Wright       kappa[i] = kmin*pow(stg_ctx->alpha, i);
313ba6664aeSJames Wright     }
314ba6664aeSJames Wright   } //end calculate kappa
315ba6664aeSJames Wright 
316d95c92dbSJames Wright   ierr = PetscFree(*pstg_ctx); CHKERRQ(ierr);
317ba6664aeSJames Wright   *pstg_ctx = stg_ctx;
318ba6664aeSJames Wright   PetscFunctionReturn(0);
319ba6664aeSJames Wright }
320ba6664aeSJames Wright 
321ba6664aeSJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem,
322ba6664aeSJames Wright                         User user, const bool prescribe_T,
3234e139266SJames Wright                         const CeedScalar theta0, const CeedScalar P0,
3244e139266SJames Wright                         const CeedScalar ynodes[], const CeedInt nynodes) {
325ba6664aeSJames Wright   PetscErrorCode ierr;
326ba6664aeSJames Wright   char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
327ba6664aeSJames Wright   char stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
3284e139266SJames Wright   PetscBool  mean_only     = PETSC_FALSE,
3294e139266SJames Wright              use_stgstrong = PETSC_FALSE;
3304e139266SJames Wright   CeedScalar u0            = 0.0,
3314e139266SJames Wright              alpha         = 1.01;
332ba6664aeSJames Wright   CeedQFunctionContext stg_context;
333ba6664aeSJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
334ba6664aeSJames Wright   PetscFunctionBeginUser;
335ba6664aeSJames Wright 
336ba6664aeSJames Wright   // Get options
337ba6664aeSJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
338ba6664aeSJames Wright   ierr = PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL,
339ba6664aeSJames Wright                             stg_inflow_path, stg_inflow_path,
340ba6664aeSJames Wright                             sizeof(stg_inflow_path), NULL); CHKERRQ(ierr);
341ba6664aeSJames Wright   ierr = PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL,
342ba6664aeSJames Wright                             stg_rand_path,stg_rand_path,
343ba6664aeSJames Wright                             sizeof(stg_rand_path), NULL); CHKERRQ(ierr);
344ba6664aeSJames Wright   ierr = PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL,
345ba6664aeSJames Wright                           alpha, &alpha, NULL); CHKERRQ(ierr);
346ba6664aeSJames Wright   ierr = PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations",
347ba6664aeSJames Wright                           NULL, u0, &u0, NULL); CHKERRQ(ierr);
348ba6664aeSJames Wright   ierr = PetscOptionsBool("-stg_mean_only", "Only apply mean profile",
349ba6664aeSJames Wright                           NULL, mean_only, &mean_only, NULL); CHKERRQ(ierr);
3504e139266SJames Wright   ierr = PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly",
3514e139266SJames Wright                           NULL, use_stgstrong, &use_stgstrong, NULL); CHKERRQ(ierr);
352ba6664aeSJames Wright   PetscOptionsEnd();
353ba6664aeSJames Wright 
35465dd5cafSJames Wright   ierr = PetscCalloc1(1, &global_stg_ctx); CHKERRQ(ierr);
35565dd5cafSJames Wright   global_stg_ctx->alpha         = alpha;
35665dd5cafSJames Wright   global_stg_ctx->u0            = u0;
35765dd5cafSJames Wright   global_stg_ctx->is_implicit   = user->phys->implicit;
35865dd5cafSJames Wright   global_stg_ctx->prescribe_T   = prescribe_T;
35965dd5cafSJames Wright   global_stg_ctx->mean_only     = mean_only;
36065dd5cafSJames Wright   global_stg_ctx->theta0        = theta0;
36165dd5cafSJames Wright   global_stg_ctx->P0            = P0;
36265dd5cafSJames Wright   global_stg_ctx->nynodes       = nynodes;
363ba6664aeSJames Wright 
364ba6664aeSJames Wright   {
365ba6664aeSJames Wright     // Calculate dx assuming constant spacing
366ba6664aeSJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
367ba6664aeSJames Wright     ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
368ba6664aeSJames Wright     for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
369ba6664aeSJames Wright 
370ba6664aeSJames Wright     PetscInt nmax = 3, faces[3];
371ba6664aeSJames Wright     ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &nmax,
372ba6664aeSJames Wright                                    NULL); CHKERRQ(ierr);
37365dd5cafSJames Wright     global_stg_ctx->dx = domain_size[0]/faces[0];
37465dd5cafSJames Wright     global_stg_ctx->dz = domain_size[2]/faces[2];
375ba6664aeSJames Wright   }
376ba6664aeSJames Wright 
377ba6664aeSJames Wright   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
378ba6664aeSJames Wright                               CEED_MEM_HOST, &newtonian_ig_ctx);
37965dd5cafSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
380ba6664aeSJames Wright   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
381ba6664aeSJames Wright                                   &newtonian_ig_ctx);
382ba6664aeSJames Wright 
38365dd5cafSJames Wright   ierr = GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path,
38465dd5cafSJames Wright                            &global_stg_ctx, ynodes); CHKERRQ(ierr);
385ba6664aeSJames Wright 
386ba6664aeSJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
387ba6664aeSJames Wright   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST,
38865dd5cafSJames Wright                               CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
389ba6664aeSJames Wright   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST,
390ba6664aeSJames Wright                                      FreeContextPetsc);
391ba6664aeSJames Wright   CeedQFunctionContextRegisterDouble(stg_context, "solution time",
392ba6664aeSJames Wright                                      offsetof(struct STGShur14Context_, time), 1,
39386c7fc99SJed Brown                                      "Physical time of the solution");
394ba6664aeSJames Wright 
395b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
396b77c53c9SJames Wright   problem->ics.qfunction         = ICsSTG;
397b77c53c9SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
398b77c53c9SJames Wright   problem->ics.qfunction_context = stg_context;
399b77c53c9SJames Wright 
4004e139266SJames Wright   if (use_stgstrong) {
40165dd5cafSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
402dada6cc0SJames Wright     problem->use_dirichlet_ceed = PETSC_TRUE;
40336b31e27SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
4044e139266SJames Wright   } else {
405ba6664aeSJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
406ba6664aeSJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
4074dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
4084dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
4094dbab5e5SJames Wright     CeedQFunctionContextReferenceCopy(stg_context,
4104dbab5e5SJames Wright                                       &problem->apply_inflow.qfunction_context);
4114dbab5e5SJames Wright     CeedQFunctionContextReferenceCopy(stg_context,
4124dbab5e5SJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
41336b31e27SJames Wright     problem->bc_from_ics = PETSC_TRUE;
4144e139266SJames Wright   }
415ba6664aeSJames Wright 
416ba6664aeSJames Wright   PetscFunctionReturn(0);
417ba6664aeSJames Wright }
4184e139266SJames Wright 
4194e139266SJames Wright static inline PetscScalar FindDy(const PetscScalar ynodes[],
4204e139266SJames Wright                                  const PetscInt nynodes, const PetscScalar y) {
4216ef2784eSLeila Ghaffari 
4224e139266SJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
4234e139266SJames Wright   // ^^assuming min(dy) is first element off the wall
4244e139266SJames Wright   PetscInt idx = -1; // Index
4254e139266SJames Wright 
4264e139266SJames Wright   for (PetscInt i=0; i<nynodes; i++) {
4274e139266SJames Wright     if (y < ynodes[i] + half_mindy) {
4284e139266SJames Wright       idx = i; break;
4294e139266SJames Wright     }
4304e139266SJames Wright   }
4314e139266SJames Wright   if      (idx == 0)          return ynodes[1] - ynodes[0];
4324e139266SJames Wright   else if (idx == nynodes-1)  return ynodes[nynodes-2] - ynodes[nynodes-1];
4334e139266SJames Wright   else                        return 0.5 * (ynodes[idx+1] - ynodes[idx-1]);
4344e139266SJames Wright }
4354e139266SJames Wright 
4364e139266SJames Wright // Function passed to DMAddBoundary
437dada6cc0SJames Wright // NOTE: Not used in favor of QFunction-based method
4384e139266SJames Wright PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time,
4394e139266SJames Wright                                const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
4404e139266SJames Wright   PetscFunctionBeginUser;
4414e139266SJames Wright 
4424e139266SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
4434e139266SJames Wright   PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
4444e139266SJames Wright   const bool mean_only      = stg_ctx->mean_only;
4454e139266SJames Wright   const PetscScalar dx      = stg_ctx->dx;
4464e139266SJames Wright   const PetscScalar dz      = stg_ctx->dz;
4474e139266SJames Wright   const PetscScalar mu      = stg_ctx->newtonian_ctx.mu;
4484e139266SJames Wright   const PetscScalar theta0  = stg_ctx->theta0;
4494e139266SJames Wright   const PetscScalar P0      = stg_ctx->P0;
4504e139266SJames Wright   const PetscScalar cv      = stg_ctx->newtonian_ctx.cv;
4514e139266SJames Wright   const PetscScalar cp      = stg_ctx->newtonian_ctx.cp;
4524e139266SJames Wright   const PetscScalar Rd      = cp - cv;
4534e139266SJames Wright 
4544e139266SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
4554e139266SJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
4564e139266SJames Wright   if (!mean_only) {
4574e139266SJames Wright     const PetscInt    nynodes = stg_ctx->nynodes;
4584e139266SJames Wright     const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes];
4594e139266SJames Wright     const PetscScalar h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
4604e139266SJames Wright     CalcSpectrum(x[1], eps, lt, h, mu/rho, qn, stg_ctx);
4614e139266SJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
4624e139266SJames Wright   } else {
4634e139266SJames Wright     for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
4644e139266SJames Wright   }
4654e139266SJames Wright 
4664e139266SJames Wright   bcval[0] = rho;
4674e139266SJames Wright   bcval[1] = rho * u[0];
4684e139266SJames Wright   bcval[2] = rho * u[1];
4694e139266SJames Wright   bcval[3] = rho * u[2];
4704e139266SJames Wright   PetscFunctionReturn(0);
4714e139266SJames Wright }
4724e139266SJames Wright 
473192a7459SJames Wright PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem,
474192a7459SJames Wright                               Physics phys) {
4754e139266SJames Wright   PetscErrorCode ierr;
4764e139266SJames Wright   DMLabel label;
4774e139266SJames Wright   PetscFunctionBeginUser;
4784e139266SJames Wright 
479192a7459SJames Wright   PetscInt comps[5], num_comps=4;
48097baf651SJames Wright   switch (phys->state_var) {
48197baf651SJames Wright   case STATEVAR_CONSERVATIVE:
482192a7459SJames Wright     // {0,1,2,3} for rho, rho*u, rho*v, rho*w
483192a7459SJames Wright     for(int i=0; i<4; i++) comps[i] = i;
48497baf651SJames Wright     break;
48597baf651SJames Wright 
48697baf651SJames Wright   case STATEVAR_PRIMITIVE:
48797baf651SJames Wright     // {1,2,3,4} for u, v, w, T
48897baf651SJames Wright     for(int i=0; i<4; i++) comps[i] = i+1;
48997baf651SJames Wright     break;
490192a7459SJames Wright   }
491192a7459SJames Wright 
4924e139266SJames Wright   ierr = DMGetLabel(dm, "Face Sets", &label); CHKERRQ(ierr);
4934e139266SJames Wright   // Set wall BCs
4944e139266SJames Wright   if (bc->num_inflow > 0) {
4954e139266SJames Wright     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label,
4964e139266SJames Wright                          bc->num_inflow, bc->inflows, 0, num_comps,
4974e139266SJames Wright                          comps, (void(*)(void))StrongSTGbcFunc,
49865dd5cafSJames Wright                          NULL, global_stg_ctx, NULL);  CHKERRQ(ierr);
4994e139266SJames Wright   }
5004e139266SJames Wright 
5014e139266SJames Wright   PetscFunctionReturn(0);
5024e139266SJames Wright }
503dada6cc0SJames Wright 
504dada6cc0SJames Wright PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem,
5055dc40723SJames Wright                                  CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
5065dc40723SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
507dada6cc0SJames Wright 
508dada6cc0SJames Wright   CeedQFunction qf_strongbc;
509dada6cc0SJames Wright   PetscFunctionBeginUser;
510dada6cc0SJames Wright   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF,
511dada6cc0SJames Wright                               STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
512dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur,
513dada6cc0SJames Wright                         CEED_EVAL_NONE);
514dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "x",        num_comp_x,    CEED_EVAL_NONE);
515dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "scale",    1,             CEED_EVAL_NONE);
5165dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "stg data", stg_data_size, CEED_EVAL_NONE);
517dada6cc0SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "q",        num_comp_q,    CEED_EVAL_NONE);
518dada6cc0SJames Wright 
519dada6cc0SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
520dada6cc0SJames Wright   *pqf_strongbc = qf_strongbc;
521dada6cc0SJames Wright   PetscFunctionReturn(0);
522dada6cc0SJames Wright }
5235dc40723SJames Wright 
5245dc40723SJames Wright PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem,
5255dc40723SJames Wright     CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
5265dc40723SJames Wright     CeedQFunction *pqf_strongbc) {
5275dc40723SJames Wright 
5285dc40723SJames Wright   CeedQFunction qf_strongbc;
5295dc40723SJames Wright   PetscFunctionBeginUser;
5305dc40723SJames Wright   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14,
5315dc40723SJames Wright                               Preprocess_STGShur14_loc, &qf_strongbc);
5325dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur,
5335dc40723SJames Wright                         CEED_EVAL_NONE);
5345dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "x",        num_comp_x,    CEED_EVAL_NONE);
5355dc40723SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
5365dc40723SJames Wright 
5375dc40723SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
5385dc40723SJames Wright   *pqf_strongbc = qf_strongbc;
5395dc40723SJames Wright   PetscFunctionReturn(0);
5405dc40723SJames Wright }
541