xref: /libCEED/examples/fluids/problems/stg_shur14.c (revision 24941a69eaf54c83b32846ca05428878af166e84)
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   PetscInt ndims;
72ba6664aeSJames Wright   char line[char_array_len];
73ba6664aeSJames Wright   char **array;
74ba6664aeSJames Wright 
75ba6664aeSJames Wright   PetscFunctionBeginUser;
76*24941a69SJames Wright   PetscCall(PetscFOpen(comm, path, "r", fp));
77*24941a69SJames Wright   PetscCall(PetscSynchronizedFGets(comm, *fp, char_array_len, line));
78*24941a69SJames Wright   PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
79ba6664aeSJames Wright   if (ndims != 2) SETERRQ(comm, -1,
80990fdeb6SJeremy L Thompson                             "Found %" PetscInt_FMT" dimensions instead of 2 on the first line of %s",
81ba6664aeSJames Wright                             ndims, path);
82ba6664aeSJames Wright 
83ba6664aeSJames Wright   for (PetscInt i=0; i<ndims; i++)  dims[i] = atoi(array[i]);
84*24941a69SJames Wright   PetscCall(PetscStrToArrayDestroy(ndims, array));
85ba6664aeSJames Wright   PetscFunctionReturn(0);
86ba6664aeSJames Wright }
87ba6664aeSJames Wright 
88ba6664aeSJames Wright /*
89ba6664aeSJames Wright  * @brief Get the number of rows for the PHASTA file at path
90ba6664aeSJames Wright  *
91ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
92ba6664aeSJames Wright  * as the only two entries, separated by a single space
93ba6664aeSJames Wright  *
94ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
95ba6664aeSJames Wright  * @param[in] path Path to the file
96ba6664aeSJames Wright  * @param[out] nrows Number of rows
97ba6664aeSJames Wright  */
98ba6664aeSJames Wright static PetscErrorCode GetNRows(const MPI_Comm comm,
99ba6664aeSJames Wright                                const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows) {
100ba6664aeSJames Wright   const PetscInt char_array_len = 512;
101ba6664aeSJames Wright   PetscInt dims[2];
102ba6664aeSJames Wright   FILE *fp;
103ba6664aeSJames Wright 
104ba6664aeSJames Wright   PetscFunctionBeginUser;
105*24941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
106ba6664aeSJames Wright   *nrows = dims[0];
107*24941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
108ba6664aeSJames Wright   PetscFunctionReturn(0);
109ba6664aeSJames Wright }
110ba6664aeSJames Wright 
111ba6664aeSJames Wright /*
112ba6664aeSJames Wright  * @brief Read the STGInflow file and load the contents into stg_ctx
113ba6664aeSJames Wright  *
114ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
115ba6664aeSJames Wright  * as the only two entries, separated by a single space.
116ba6664aeSJames Wright  * Assumes there are 14 columns in the file
117ba6664aeSJames Wright  *
118ba6664aeSJames Wright  * Function calculates the Cholesky decomposition from the Reynolds stress
119ba6664aeSJames Wright  * profile found in the file
120ba6664aeSJames Wright  *
121ba6664aeSJames Wright  * @param[in] comm MPI_Comm for the program
122ba6664aeSJames Wright  * @param[in] path Path to the STGInflow.dat file
123ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
124ba6664aeSJames Wright  */
125ba6664aeSJames Wright static PetscErrorCode ReadSTGInflow(const MPI_Comm comm,
126ba6664aeSJames Wright                                     const char path[PETSC_MAX_PATH_LEN], STGShur14Context stg_ctx) {
127ba6664aeSJames Wright   PetscInt ndims, dims[2];
128ba6664aeSJames Wright   FILE *fp;
129ba6664aeSJames Wright   const PetscInt char_array_len=512;
130ba6664aeSJames Wright   char line[char_array_len];
131ba6664aeSJames Wright   char **array;
132ba6664aeSJames Wright 
133ba6664aeSJames Wright   PetscFunctionBeginUser;
134ba6664aeSJames Wright 
135*24941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
136ba6664aeSJames Wright 
137ba6664aeSJames Wright   CeedScalar rij[6][stg_ctx->nprofs];
138175f00a6SJames Wright   CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
139ba6664aeSJames Wright   CeedScalar *eps = &stg_ctx->data[stg_ctx->offsets.eps];
140ba6664aeSJames Wright   CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
141ba6664aeSJames Wright   CeedScalar (*ubar)[stg_ctx->nprofs] = (CeedScalar (*)[stg_ctx->nprofs])
142ba6664aeSJames Wright                                         &stg_ctx->data[stg_ctx->offsets.ubar];
143ba6664aeSJames Wright 
144ba6664aeSJames Wright   for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
145*24941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
146*24941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
147ba6664aeSJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
148990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
149990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT")", i,
150ba6664aeSJames Wright                                    path, ndims, dims[1]);
151ba6664aeSJames Wright 
152175f00a6SJames Wright     wall_dist[i] = (CeedScalar) atof(array[0]);
153ba6664aeSJames Wright     ubar[0][i]   = (CeedScalar) atof(array[1]);
154ba6664aeSJames Wright     ubar[1][i]   = (CeedScalar) atof(array[2]);
155ba6664aeSJames Wright     ubar[2][i]   = (CeedScalar) atof(array[3]);
156ba6664aeSJames Wright     rij[0][i]    = (CeedScalar) atof(array[4]);
157ba6664aeSJames Wright     rij[1][i]    = (CeedScalar) atof(array[5]);
158ba6664aeSJames Wright     rij[2][i]    = (CeedScalar) atof(array[6]);
159ba6664aeSJames Wright     rij[3][i]    = (CeedScalar) atof(array[7]);
160ba6664aeSJames Wright     rij[4][i]    = (CeedScalar) atof(array[8]);
161ba6664aeSJames Wright     rij[5][i]    = (CeedScalar) atof(array[9]);
162ba6664aeSJames Wright     lt[i]        = (CeedScalar) atof(array[12]);
163ba6664aeSJames Wright     eps[i]       = (CeedScalar) atof(array[13]);
164ba6664aeSJames Wright 
165175f00a6SJames Wright     if (wall_dist[i] < 0) SETERRQ(comm, -1,
166ba6664aeSJames Wright                                     "Distance to wall in %s cannot be negative", path);
167ba6664aeSJames Wright     if (lt[i] < 0) SETERRQ(comm, -1,
168ba6664aeSJames Wright                              "Turbulent length scale in %s cannot be negative", path);
169ba6664aeSJames Wright     if (eps[i] < 0) SETERRQ(comm, -1,
170ba6664aeSJames Wright                               "Turbulent dissipation in %s cannot be negative", path);
171ba6664aeSJames Wright 
172ba6664aeSJames Wright   }
173ba6664aeSJames Wright   CeedScalar (*cij)[stg_ctx->nprofs]  = (CeedScalar (*)[stg_ctx->nprofs])
174ba6664aeSJames Wright                                         &stg_ctx->data[stg_ctx->offsets.cij];
175*24941a69SJames Wright   PetscCall(CalcCholeskyDecomp(comm, stg_ctx->nprofs, rij, cij));
176*24941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
177ba6664aeSJames Wright   PetscFunctionReturn(0);
178ba6664aeSJames Wright }
179ba6664aeSJames Wright 
180ba6664aeSJames Wright /*
181ba6664aeSJames Wright  * @brief Read the STGRand file and load the contents into stg_ctx
182ba6664aeSJames Wright  *
183ba6664aeSJames Wright  * Assumes that the first line of the file has the number of rows and columns
184ba6664aeSJames Wright  * as the only two entries, separated by a single space.
185ba6664aeSJames Wright  * Assumes there are 7 columns in the file
186ba6664aeSJames Wright  *
187ba6664aeSJames Wright  * @param[in]    comm    MPI_Comm for the program
188ba6664aeSJames Wright  * @param[in]    path    Path to the STGRand.dat file
189ba6664aeSJames Wright  * @param[inout] stg_ctx STGShur14Context where the data will be loaded into
190ba6664aeSJames Wright  */
191ba6664aeSJames Wright static PetscErrorCode ReadSTGRand(const MPI_Comm comm,
192ba6664aeSJames Wright                                   const char path[PETSC_MAX_PATH_LEN],
193ba6664aeSJames Wright                                   STGShur14Context stg_ctx) {
194ba6664aeSJames Wright   PetscInt ndims, dims[2];
195ba6664aeSJames Wright   FILE *fp;
196ba6664aeSJames Wright   const PetscInt char_array_len = 512;
197ba6664aeSJames Wright   char line[char_array_len];
198ba6664aeSJames Wright   char **array;
199ba6664aeSJames Wright 
200ba6664aeSJames Wright   PetscFunctionBeginUser;
201*24941a69SJames Wright   PetscCall(OpenPHASTADatFile(comm, path, char_array_len, dims, &fp));
202ba6664aeSJames Wright 
203ba6664aeSJames Wright   CeedScalar *phi = &stg_ctx->data[stg_ctx->offsets.phi];
204ba6664aeSJames Wright   CeedScalar (*d)[stg_ctx->nmodes]     = (CeedScalar (*)[stg_ctx->nmodes])
205ba6664aeSJames Wright                                          &stg_ctx->data[stg_ctx->offsets.d];
206ba6664aeSJames Wright   CeedScalar (*sigma)[stg_ctx->nmodes] = (CeedScalar (*)[stg_ctx->nmodes])
207ba6664aeSJames Wright                                          &stg_ctx->data[stg_ctx->offsets.sigma];
208ba6664aeSJames Wright 
209ba6664aeSJames Wright   for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
210*24941a69SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, char_array_len, line));
211*24941a69SJames Wright     PetscCall(PetscStrToArray(line, ' ', &ndims, &array));
212ba6664aeSJames Wright     if (ndims < dims[1]) SETERRQ(comm, -1,
213990fdeb6SJeremy L Thompson                                    "Line %" PetscInt_FMT" of %s does not contain enough columns (%"
214990fdeb6SJeremy L Thompson                                    PetscInt_FMT" instead of %" PetscInt_FMT")", i,
215ba6664aeSJames Wright                                    path, ndims, dims[1]);
216ba6664aeSJames Wright 
217ba6664aeSJames Wright     d[0][i]     = (CeedScalar) atof(array[0]);
218ba6664aeSJames Wright     d[1][i]     = (CeedScalar) atof(array[1]);
219ba6664aeSJames Wright     d[2][i]     = (CeedScalar) atof(array[2]);
220ba6664aeSJames Wright     phi[i]      = (CeedScalar) atof(array[3]);
221ba6664aeSJames Wright     sigma[0][i] = (CeedScalar) atof(array[4]);
222ba6664aeSJames Wright     sigma[1][i] = (CeedScalar) atof(array[5]);
223ba6664aeSJames Wright     sigma[2][i] = (CeedScalar) atof(array[6]);
224ba6664aeSJames Wright   }
225*24941a69SJames Wright   PetscCall(PetscFClose(comm, fp));
226ba6664aeSJames Wright   PetscFunctionReturn(0);
227ba6664aeSJames Wright }
228ba6664aeSJames Wright 
229ba6664aeSJames Wright /*
230ba6664aeSJames Wright  * @brief Read STG data from input paths and put in STGShur14Context
231ba6664aeSJames Wright  *
232ba6664aeSJames Wright  * Reads data from input paths and puts them into a STGShur14Context object.
233ba6664aeSJames Wright  * Data stored initially in `*pstg_ctx` will be copied over to the new
234ba6664aeSJames Wright  * STGShur14Context instance.
235ba6664aeSJames Wright  *
236ba6664aeSJames Wright  * @param[in]    comm            MPI_Comm for the program
237ba6664aeSJames Wright  * @param[in]    dm              DM for the program
238ba6664aeSJames Wright  * @param[in]    stg_inflow_path Path to STGInflow.dat file
239ba6664aeSJames Wright  * @param[in]    stg_rand_path   Path to STGRand.dat file
240ba6664aeSJames Wright  * @param[inout] pstg_ctx        Pointer to STGShur14Context where the data will be loaded into
241ba6664aeSJames Wright  */
242ba6664aeSJames Wright PetscErrorCode GetSTGContextData(const MPI_Comm comm, const DM dm,
243ba6664aeSJames Wright                                  char stg_inflow_path[PETSC_MAX_PATH_LEN],
244ba6664aeSJames Wright                                  char stg_rand_path[PETSC_MAX_PATH_LEN],
2454e139266SJames Wright                                  STGShur14Context *pstg_ctx,
2464e139266SJames Wright                                  const CeedScalar ynodes[]) {
247ba6664aeSJames Wright   PetscInt nmodes, nprofs;
248ba6664aeSJames Wright   STGShur14Context stg_ctx;
249ba6664aeSJames Wright   PetscFunctionBeginUser;
250ba6664aeSJames Wright 
251ba6664aeSJames Wright   // Get options
252*24941a69SJames Wright   PetscCall(GetNRows(comm, stg_rand_path, &nmodes));
253*24941a69SJames Wright   PetscCall(GetNRows(comm, stg_inflow_path, &nprofs));
254ba6664aeSJames Wright   if (nmodes > STG_NMODES_MAX)
255990fdeb6SJeremy L Thompson     SETERRQ(comm, 1, "Number of wavemodes in %s (%"
256990fdeb6SJeremy L Thompson             PetscInt_FMT") exceeds STG_NMODES_MAX (%" PetscInt_FMT"). "
257ba6664aeSJames Wright             "Change size of STG_NMODES_MAX and recompile", stg_rand_path, nmodes,
258ba6664aeSJames Wright             STG_NMODES_MAX);
259ba6664aeSJames Wright 
260ba6664aeSJames Wright   {
261ba6664aeSJames Wright     STGShur14Context s;
262*24941a69SJames Wright     PetscCall(PetscCalloc1(1, &s));
263ba6664aeSJames Wright     *s = **pstg_ctx;
264ba6664aeSJames Wright     s->nmodes = nmodes;
265ba6664aeSJames Wright     s->nprofs = nprofs;
266ba6664aeSJames Wright     s->offsets.sigma   = 0;
267ba6664aeSJames Wright     s->offsets.d       = nmodes*3;
268ba6664aeSJames Wright     s->offsets.phi       = s->offsets.d         + nmodes*3;
269ba6664aeSJames Wright     s->offsets.kappa     = s->offsets.phi       + nmodes;
270175f00a6SJames Wright     s->offsets.wall_dist = s->offsets.kappa     + nmodes;
271175f00a6SJames Wright     s->offsets.ubar      = s->offsets.wall_dist + nprofs;
272ba6664aeSJames Wright     s->offsets.cij       = s->offsets.ubar      + nprofs*3;
273ba6664aeSJames Wright     s->offsets.eps       = s->offsets.cij       + nprofs*6;
274ba6664aeSJames Wright     s->offsets.lt        = s->offsets.eps       + nprofs;
2754e139266SJames Wright     s->offsets.ynodes    = s->offsets.lt        + nprofs;
2764e139266SJames Wright     PetscInt total_num_scalars = s->offsets.ynodes + s->nynodes;
277ba6664aeSJames Wright     s->total_bytes = sizeof(*stg_ctx) + total_num_scalars*sizeof(stg_ctx->data[0]);
278*24941a69SJames Wright     PetscCall(PetscMalloc(s->total_bytes, &stg_ctx));
279ba6664aeSJames Wright     *stg_ctx = *s;
280*24941a69SJames Wright     PetscCall(PetscFree(s));
281ba6664aeSJames Wright   }
282ba6664aeSJames Wright 
283*24941a69SJames Wright   PetscCall(ReadSTGInflow(comm, stg_inflow_path, stg_ctx));
284*24941a69SJames Wright   PetscCall(ReadSTGRand(comm, stg_rand_path, stg_ctx));
285ba6664aeSJames Wright 
2864e139266SJames Wright   if (stg_ctx->nynodes > 0) {
2874e139266SJames Wright     CeedScalar *ynodes_ctx = &stg_ctx->data[stg_ctx->offsets.ynodes];
2884e139266SJames Wright     for (PetscInt i=0; i<stg_ctx->nynodes; i++) ynodes_ctx[i] = ynodes[i];
2894e139266SJames Wright   }
2904e139266SJames Wright 
291ba6664aeSJames Wright   // -- Calculate kappa
292ba6664aeSJames Wright   {
293ba6664aeSJames Wright     CeedScalar *kappa = &stg_ctx->data[stg_ctx->offsets.kappa];
294175f00a6SJames Wright     CeedScalar *wall_dist = &stg_ctx->data[stg_ctx->offsets.wall_dist];
295ba6664aeSJames Wright     CeedScalar *lt = &stg_ctx->data[stg_ctx->offsets.lt];
296ba6664aeSJames Wright     CeedScalar le, le_max=0;
297ba6664aeSJames Wright 
298ba6664aeSJames Wright     CeedPragmaSIMD
299ba6664aeSJames Wright     for (PetscInt i=0; i<stg_ctx->nprofs; i++) {
300175f00a6SJames Wright       le = PetscMin(2*wall_dist[i], 3*lt[i]);
301ba6664aeSJames Wright       if (le_max < le) le_max = le;
302ba6664aeSJames Wright     }
303ba6664aeSJames Wright     CeedScalar kmin = M_PI/le_max;
304ba6664aeSJames Wright 
305ba6664aeSJames Wright     CeedPragmaSIMD
306ba6664aeSJames Wright     for (PetscInt i=0; i<stg_ctx->nmodes; i++) {
307ba6664aeSJames Wright       kappa[i] = kmin*pow(stg_ctx->alpha, i);
308ba6664aeSJames Wright     }
309ba6664aeSJames Wright   } //end calculate kappa
310ba6664aeSJames Wright 
311*24941a69SJames Wright   PetscCall(PetscFree(*pstg_ctx));
312ba6664aeSJames Wright   *pstg_ctx = stg_ctx;
313ba6664aeSJames Wright   PetscFunctionReturn(0);
314ba6664aeSJames Wright }
315ba6664aeSJames Wright 
316ba6664aeSJames Wright PetscErrorCode SetupSTG(const MPI_Comm comm, const DM dm, ProblemData *problem,
317ba6664aeSJames Wright                         User user, const bool prescribe_T,
3184e139266SJames Wright                         const CeedScalar theta0, const CeedScalar P0,
3194e139266SJames Wright                         const CeedScalar ynodes[], const CeedInt nynodes) {
320ba6664aeSJames Wright   char stg_inflow_path[PETSC_MAX_PATH_LEN] = "./STGInflow.dat";
321ba6664aeSJames Wright   char stg_rand_path[PETSC_MAX_PATH_LEN]   = "./STGRand.dat";
3224e139266SJames Wright   PetscBool  mean_only          = PETSC_FALSE,
32389060322SJames Wright              use_stgstrong      = PETSC_FALSE,
32489060322SJames Wright              use_fluctuating_IC = PETSC_FALSE;
3254e139266SJames Wright   CeedScalar u0            = 0.0,
3264e139266SJames Wright              alpha         = 1.01;
327ba6664aeSJames Wright   CeedQFunctionContext stg_context;
328ba6664aeSJames Wright   NewtonianIdealGasContext newtonian_ig_ctx;
329ba6664aeSJames Wright   PetscFunctionBeginUser;
330ba6664aeSJames Wright 
331ba6664aeSJames Wright   // Get options
332ba6664aeSJames Wright   PetscOptionsBegin(comm, NULL, "STG Boundary Condition Options", NULL);
333*24941a69SJames Wright   PetscCall(PetscOptionsString("-stg_inflow_path", "Path to STGInflow.dat", NULL,
334ba6664aeSJames Wright                                stg_inflow_path, stg_inflow_path,
335*24941a69SJames Wright                                sizeof(stg_inflow_path), NULL));
336*24941a69SJames Wright   PetscCall(PetscOptionsString("-stg_rand_path", "Path to STGInflow.dat", NULL,
337ba6664aeSJames Wright                                stg_rand_path,stg_rand_path,
338*24941a69SJames Wright                                sizeof(stg_rand_path), NULL));
339*24941a69SJames Wright   PetscCall(PetscOptionsReal("-stg_alpha", "Growth rate of the wavemodes", NULL,
340*24941a69SJames Wright                              alpha, &alpha, NULL));
341*24941a69SJames Wright   PetscCall(PetscOptionsReal("-stg_u0", "Advective velocity for the fluctuations",
342*24941a69SJames Wright                              NULL, u0, &u0, NULL));
343*24941a69SJames Wright   PetscCall(PetscOptionsBool("-stg_mean_only", "Only apply mean profile",
344*24941a69SJames Wright                              NULL, mean_only, &mean_only, NULL));
345*24941a69SJames Wright   PetscCall(PetscOptionsBool("-stg_strong", "Enforce STG inflow strongly",
346*24941a69SJames Wright                              NULL, use_stgstrong, &use_stgstrong, NULL));
34789060322SJames Wright   PetscCall(PetscOptionsBool("-stg_fluctuating_IC",
34889060322SJames Wright                              "\"Extrude\" the fluctuations through the domain as an initial condition",
34989060322SJames Wright                              NULL, use_fluctuating_IC, &use_fluctuating_IC, NULL));
350ba6664aeSJames Wright   PetscOptionsEnd();
351ba6664aeSJames Wright 
352*24941a69SJames Wright   PetscCall(PetscCalloc1(1, &global_stg_ctx));
35365dd5cafSJames Wright   global_stg_ctx->alpha              = alpha;
35465dd5cafSJames Wright   global_stg_ctx->u0                 = u0;
35565dd5cafSJames Wright   global_stg_ctx->is_implicit        = user->phys->implicit;
35665dd5cafSJames Wright   global_stg_ctx->prescribe_T        = prescribe_T;
35765dd5cafSJames Wright   global_stg_ctx->mean_only          = mean_only;
35889060322SJames Wright   global_stg_ctx->use_fluctuating_IC = use_fluctuating_IC;
35965dd5cafSJames Wright   global_stg_ctx->theta0             = theta0;
36065dd5cafSJames Wright   global_stg_ctx->P0                 = P0;
36165dd5cafSJames Wright   global_stg_ctx->nynodes            = nynodes;
362ba6664aeSJames Wright 
363ba6664aeSJames Wright   {
364ba6664aeSJames Wright     // Calculate dx assuming constant spacing
365ba6664aeSJames Wright     PetscReal domain_min[3], domain_max[3], domain_size[3];
366*24941a69SJames Wright     PetscCall(DMGetBoundingBox(dm, domain_min, domain_max));
367ba6664aeSJames Wright     for (PetscInt i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
368ba6664aeSJames Wright 
369ba6664aeSJames Wright     PetscInt nmax = 3, faces[3];
370*24941a69SJames Wright     PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces,
371*24941a69SJames Wright                                       &nmax, NULL));
37265dd5cafSJames Wright     global_stg_ctx->dx = domain_size[0]/faces[0];
37365dd5cafSJames Wright     global_stg_ctx->dz = domain_size[2]/faces[2];
374ba6664aeSJames Wright   }
375ba6664aeSJames Wright 
376ba6664aeSJames Wright   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context,
377ba6664aeSJames Wright                               CEED_MEM_HOST, &newtonian_ig_ctx);
37865dd5cafSJames Wright   global_stg_ctx->newtonian_ctx = *newtonian_ig_ctx;
379ba6664aeSJames Wright   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context,
380ba6664aeSJames Wright                                   &newtonian_ig_ctx);
381ba6664aeSJames Wright 
382*24941a69SJames Wright   PetscCall(GetSTGContextData(comm, dm, stg_inflow_path, stg_rand_path,
383*24941a69SJames Wright                               &global_stg_ctx, ynodes));
384ba6664aeSJames Wright 
385ba6664aeSJames Wright   CeedQFunctionContextCreate(user->ceed, &stg_context);
386ba6664aeSJames Wright   CeedQFunctionContextSetData(stg_context, CEED_MEM_HOST,
38765dd5cafSJames Wright                               CEED_USE_POINTER, global_stg_ctx->total_bytes, global_stg_ctx);
388ba6664aeSJames Wright   CeedQFunctionContextSetDataDestroy(stg_context, CEED_MEM_HOST,
389ba6664aeSJames Wright                                      FreeContextPetsc);
390ba6664aeSJames Wright   CeedQFunctionContextRegisterDouble(stg_context, "solution time",
391ba6664aeSJames Wright                                      offsetof(struct STGShur14Context_, time), 1,
39286c7fc99SJed Brown                                      "Physical time of the solution");
393ba6664aeSJames Wright 
394b77c53c9SJames Wright   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
395b77c53c9SJames Wright   problem->ics.qfunction         = ICsSTG;
396b77c53c9SJames Wright   problem->ics.qfunction_loc     = ICsSTG_loc;
397b77c53c9SJames Wright   problem->ics.qfunction_context = stg_context;
398b77c53c9SJames Wright 
3994e139266SJames Wright   if (use_stgstrong) {
40065dd5cafSJames Wright     // Use default boundary integral QF (BoundaryIntegral) in newtonian.h
401dada6cc0SJames Wright     problem->use_dirichlet_ceed = PETSC_TRUE;
40236b31e27SJames Wright     problem->bc_from_ics        = PETSC_FALSE;
4034e139266SJames Wright   } else {
404ba6664aeSJames Wright     problem->apply_inflow.qfunction              = STGShur14_Inflow;
405ba6664aeSJames Wright     problem->apply_inflow.qfunction_loc          = STGShur14_Inflow_loc;
4064dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction     = STGShur14_Inflow_Jacobian;
4074dbab5e5SJames Wright     problem->apply_inflow_jacobian.qfunction_loc = STGShur14_Inflow_Jacobian_loc;
4084dbab5e5SJames Wright     CeedQFunctionContextReferenceCopy(stg_context,
4094dbab5e5SJames Wright                                       &problem->apply_inflow.qfunction_context);
4104dbab5e5SJames Wright     CeedQFunctionContextReferenceCopy(stg_context,
4114dbab5e5SJames Wright                                       &problem->apply_inflow_jacobian.qfunction_context);
41236b31e27SJames Wright     problem->bc_from_ics = PETSC_TRUE;
4134e139266SJames Wright   }
414ba6664aeSJames Wright 
415ba6664aeSJames Wright   PetscFunctionReturn(0);
416ba6664aeSJames Wright }
4174e139266SJames Wright 
4184e139266SJames Wright static inline PetscScalar FindDy(const PetscScalar ynodes[],
4194e139266SJames Wright                                  const PetscInt nynodes, const PetscScalar y) {
4204e139266SJames Wright   const PetscScalar half_mindy = 0.5 * (ynodes[1] - ynodes[0]);
4214e139266SJames Wright   // ^^assuming min(dy) is first element off the wall
4224e139266SJames Wright   PetscInt idx = -1; // Index
4234e139266SJames Wright 
4244e139266SJames Wright   for (PetscInt i=0; i<nynodes; i++) {
4254e139266SJames Wright     if (y < ynodes[i] + half_mindy) {
4264e139266SJames Wright       idx = i; break;
4274e139266SJames Wright     }
4284e139266SJames Wright   }
4294e139266SJames Wright   if      (idx == 0)          return ynodes[1] - ynodes[0];
4304e139266SJames Wright   else if (idx == nynodes-1)  return ynodes[nynodes-2] - ynodes[nynodes-1];
4314e139266SJames Wright   else                        return 0.5 * (ynodes[idx+1] - ynodes[idx-1]);
4324e139266SJames Wright }
4334e139266SJames Wright 
4344e139266SJames Wright // Function passed to DMAddBoundary
435dada6cc0SJames Wright // NOTE: Not used in favor of QFunction-based method
4364e139266SJames Wright PetscErrorCode StrongSTGbcFunc(PetscInt dim, PetscReal time,
4374e139266SJames Wright                                const PetscReal x[], PetscInt Nc, PetscScalar bcval[], void *ctx) {
4384e139266SJames Wright   PetscFunctionBeginUser;
4394e139266SJames Wright   const STGShur14Context stg_ctx = (STGShur14Context) ctx;
4404e139266SJames Wright   PetscScalar qn[stg_ctx->nmodes], u[3], ubar[3], cij[6], eps, lt;
4414e139266SJames Wright   const bool mean_only      = stg_ctx->mean_only;
4424e139266SJames Wright   const PetscScalar dx      = stg_ctx->dx;
4434e139266SJames Wright   const PetscScalar dz      = stg_ctx->dz;
4444e139266SJames Wright   const PetscScalar mu      = stg_ctx->newtonian_ctx.mu;
4454e139266SJames Wright   const PetscScalar theta0  = stg_ctx->theta0;
4464e139266SJames Wright   const PetscScalar P0      = stg_ctx->P0;
4474e139266SJames Wright   const PetscScalar cv      = stg_ctx->newtonian_ctx.cv;
4484e139266SJames Wright   const PetscScalar cp      = stg_ctx->newtonian_ctx.cp;
4494e139266SJames Wright   const PetscScalar Rd      = cp - cv;
4504e139266SJames Wright 
4514e139266SJames Wright   const CeedScalar rho = P0 / (Rd * theta0);
4524e139266SJames Wright   InterpolateProfile(x[1], ubar, cij, &eps, &lt, stg_ctx);
4534e139266SJames Wright   if (!mean_only) {
4544e139266SJames Wright     const PetscInt    nynodes = stg_ctx->nynodes;
4554e139266SJames Wright     const PetscScalar *ynodes = &stg_ctx->data[stg_ctx->offsets.ynodes];
4564e139266SJames Wright     const PetscScalar h[3]    = {dx, FindDy(ynodes, nynodes, x[1]), dz};
4574e139266SJames Wright     CalcSpectrum(x[1], eps, lt, h, mu/rho, qn, stg_ctx);
4584e139266SJames Wright     STGShur14_Calc(x, time, ubar, cij, qn, u, stg_ctx);
4594e139266SJames Wright   } else {
4604e139266SJames Wright     for (CeedInt j=0; j<3; j++) u[j] = ubar[j];
4614e139266SJames Wright   }
4624e139266SJames Wright 
4634e139266SJames Wright   bcval[0] = rho;
4644e139266SJames Wright   bcval[1] = rho * u[0];
4654e139266SJames Wright   bcval[2] = rho * u[1];
4664e139266SJames Wright   bcval[3] = rho * u[2];
4674e139266SJames Wright   PetscFunctionReturn(0);
4684e139266SJames Wright }
4694e139266SJames Wright 
470192a7459SJames Wright PetscErrorCode SetupStrongSTG(DM dm, SimpleBC bc, ProblemData *problem,
471192a7459SJames Wright                               Physics phys) {
4724e139266SJames Wright   DMLabel label;
4734e139266SJames Wright   PetscFunctionBeginUser;
4744e139266SJames Wright 
475192a7459SJames Wright   PetscInt comps[5], num_comps=4;
47697baf651SJames Wright   switch (phys->state_var) {
47797baf651SJames Wright   case STATEVAR_CONSERVATIVE:
478192a7459SJames Wright     // {0,1,2,3} for rho, rho*u, rho*v, rho*w
479192a7459SJames Wright     for(int i=0; i<4; i++) comps[i] = i;
48097baf651SJames Wright     break;
48197baf651SJames Wright 
48297baf651SJames Wright   case STATEVAR_PRIMITIVE:
48397baf651SJames Wright     // {1,2,3,4} for u, v, w, T
48497baf651SJames Wright     for(int i=0; i<4; i++) comps[i] = i+1;
48597baf651SJames Wright     break;
486192a7459SJames Wright   }
487192a7459SJames Wright 
488*24941a69SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &label));
4894e139266SJames Wright   // Set wall BCs
4904e139266SJames Wright   if (bc->num_inflow > 0) {
491*24941a69SJames Wright     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "STG", label,
4924e139266SJames Wright                             bc->num_inflow, bc->inflows, 0, num_comps,
4934e139266SJames Wright                             comps, (void(*)(void))StrongSTGbcFunc,
494*24941a69SJames Wright                             NULL, global_stg_ctx, NULL));
4954e139266SJames Wright   }
4964e139266SJames Wright 
4974e139266SJames Wright   PetscFunctionReturn(0);
4984e139266SJames Wright }
499dada6cc0SJames Wright 
500dada6cc0SJames Wright PetscErrorCode SetupStrongSTG_QF(Ceed ceed, ProblemData *problem,
5015dc40723SJames Wright                                  CeedInt num_comp_x, CeedInt num_comp_q, CeedInt stg_data_size,
5025dc40723SJames Wright                                  CeedInt q_data_size_sur, CeedQFunction *pqf_strongbc) {
503dada6cc0SJames Wright 
504dada6cc0SJames Wright   CeedQFunction qf_strongbc;
505dada6cc0SJames Wright   PetscFunctionBeginUser;
506dada6cc0SJames Wright   CeedQFunctionCreateInterior(ceed, 1, STGShur14_Inflow_StrongQF,
507dada6cc0SJames Wright                               STGShur14_Inflow_StrongQF_loc, &qf_strongbc);
508dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur,
509dada6cc0SJames Wright                         CEED_EVAL_NONE);
510dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "x",        num_comp_x,    CEED_EVAL_NONE);
511dada6cc0SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "scale",    1,             CEED_EVAL_NONE);
5125dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "stg data", stg_data_size, CEED_EVAL_NONE);
513dada6cc0SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "q",        num_comp_q,    CEED_EVAL_NONE);
514dada6cc0SJames Wright 
515dada6cc0SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
516dada6cc0SJames Wright   *pqf_strongbc = qf_strongbc;
517dada6cc0SJames Wright   PetscFunctionReturn(0);
518dada6cc0SJames Wright }
5195dc40723SJames Wright 
5205dc40723SJames Wright PetscErrorCode SetupStrongSTG_PreProcessing(Ceed ceed, ProblemData *problem,
5215dc40723SJames Wright     CeedInt num_comp_x, CeedInt stg_data_size, CeedInt q_data_size_sur,
5225dc40723SJames Wright     CeedQFunction *pqf_strongbc) {
5235dc40723SJames Wright 
5245dc40723SJames Wright   CeedQFunction qf_strongbc;
5255dc40723SJames Wright   PetscFunctionBeginUser;
5265dc40723SJames Wright   CeedQFunctionCreateInterior(ceed, 1, Preprocess_STGShur14,
5275dc40723SJames Wright                               Preprocess_STGShur14_loc, &qf_strongbc);
5285dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc, "surface qdata", q_data_size_sur,
5295dc40723SJames Wright                         CEED_EVAL_NONE);
5305dc40723SJames Wright   CeedQFunctionAddInput(qf_strongbc,  "x",        num_comp_x,    CEED_EVAL_NONE);
5315dc40723SJames Wright   CeedQFunctionAddOutput(qf_strongbc, "stg data", stg_data_size, CEED_EVAL_NONE);
5325dc40723SJames Wright 
5335dc40723SJames Wright   CeedQFunctionSetContext(qf_strongbc, problem->ics.qfunction_context);
5345dc40723SJames Wright   *pqf_strongbc = qf_strongbc;
5355dc40723SJames Wright   PetscFunctionReturn(0);
5365dc40723SJames Wright }
537