xref: /libCEED/examples/fluids/qfunctions/blasius.h (revision 07d5b98a8feba68a643190b8ea9bcdac5c3e6570)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Operator for Navier-Stokes example using PETSc
10 #include <ceed.h>
11 
12 #include "newtonian_state.h"
13 #include "newtonian_types.h"
14 #include "utils.h"
15 
16 #define BLASIUS_MAX_N_CHEBYSHEV 50
17 
18 typedef struct BlasiusContext_ *BlasiusContext;
19 struct BlasiusContext_ {
20   bool                             implicit;                              // !< Using implicit timesteping or not
21   bool                             weakT;                                 // !< flag to set Temperature weakly at inflow
22   CeedScalar                       delta0;                                // !< Boundary layer height at inflow
23   CeedScalar                       U_inf;                                 // !< Velocity at boundary layer edge
24   CeedScalar                       T_inf;                                 // !< Temperature at boundary layer edge
25   CeedScalar                       T_wall;                                // !< Temperature at the wall
26   CeedScalar                       P0;                                    // !< Pressure at outflow
27   CeedScalar                       x_inflow;                              // !< Location of inflow in x
28   CeedScalar                       n_cheb;                                // !< Number of Chebyshev terms
29   CeedScalar                      *X;                                     // !< Chebyshev polynomial coordinate vector (CPU only)
30   CeedScalar                       eta_max;                               // !< Maximum eta in the domain
31   CeedScalar                       Tf_cheb[BLASIUS_MAX_N_CHEBYSHEV];      // !< Chebyshev coefficient for f
32   CeedScalar                       Th_cheb[BLASIUS_MAX_N_CHEBYSHEV - 1];  // !< Chebyshev coefficient for h
33   struct NewtonianIdealGasContext_ newtonian_ctx;
34 };
35 
36 // *****************************************************************************
37 // This helper function evaluates Chebyshev polynomials with a set of coefficients with all their derivatives represented as a recurrence table.
38 // *****************************************************************************
39 CEED_QFUNCTION_HELPER void ChebyshevEval(int N, const double *Tf, double x, double eta_max, double *f) {
40   double dX_deta     = 2 / eta_max;
41   double table[4][3] = {
42   // Chebyshev polynomials T_0, T_1, T_2 of the first kind in (-1,1)
43       {1, x, 2 * x * x - 1},
44       {0, 1, 4 * x        },
45       {0, 0, 4            },
46       {0, 0, 0            }
47   };
48   for (int i = 0; i < 4; i++) {
49     // i-th derivative of f
50     f[i] = table[i][0] * Tf[0] + table[i][1] * Tf[1] + table[i][2] * Tf[2];
51   }
52   for (int i = 3; i < N; i++) {
53     // T_n(x) = 2xT_{n-1}(x) - T_{n-2}(x)
54     table[0][i % 3] = 2 * x * table[0][(i - 1) % 3] - table[0][(i - 2) % 3];
55     // Differentiate Chebyshev polynomials with the recurrence relation
56     for (int j = 1; j < 4; j++) {
57       // T'_{n}(x)/n = 2T_{n-1}(x) + T'_{n-2}(x)/n-2
58       table[j][i % 3] = i * (2 * table[j - 1][(i - 1) % 3] + table[j][(i - 2) % 3] / (i - 2));
59     }
60     for (int j = 0; j < 4; j++) {
61       f[j] += table[j][i % 3] * Tf[i];
62     }
63   }
64   for (int i = 1; i < 4; i++) {
65     // Transform derivatives from Chebyshev [-1, 1] to [0, eta_max].
66     for (int j = 0; j < i; j++) f[i] *= dX_deta;
67   }
68 }
69 
70 // *****************************************************************************
71 // This helper function computes the Blasius boundary layer solution.
72 // *****************************************************************************
73 State CEED_QFUNCTION_HELPER(BlasiusSolution)(const BlasiusContext blasius, const CeedScalar x[3], const CeedScalar x0, const CeedScalar x_inflow,
74                                              const CeedScalar rho_infty, CeedScalar *t12) {
75   CeedInt    N     = blasius->n_cheb;
76   CeedScalar mu    = blasius->newtonian_ctx.mu;
77   CeedScalar nu    = mu / rho_infty;
78   CeedScalar eta   = x[1] * sqrt(blasius->U_inf / (nu * (x0 + x[0] - x_inflow)));
79   CeedScalar X     = 2 * (eta / blasius->eta_max) - 1.;
80   CeedScalar U_inf = blasius->U_inf;
81   CeedScalar Rd    = GasConstant(&blasius->newtonian_ctx);
82 
83   CeedScalar f[4], h[4];
84   ChebyshevEval(N, blasius->Tf_cheb, X, blasius->eta_max, f);
85   ChebyshevEval(N - 1, blasius->Th_cheb, X, blasius->eta_max, h);
86 
87   *t12 = mu * U_inf * f[2] * sqrt(U_inf / (nu * (x0 + x[0] - x_inflow)));
88 
89   CeedScalar Y[5];
90   Y[1] = U_inf * f[1];
91   Y[2] = 0.5 * sqrt(nu * U_inf / (x0 + x[0] - x_inflow)) * (eta * f[1] - f[0]);
92   Y[3] = 0.;
93   Y[4] = blasius->T_inf * h[0];
94   Y[0] = rho_infty / h[0] * Rd * Y[4];
95   return StateFromY(&blasius->newtonian_ctx, Y);
96 }
97 
98 // *****************************************************************************
99 // This QFunction sets a Blasius boundary layer for the initial condition
100 // *****************************************************************************
101 CEED_QFUNCTION(ICsBlasius)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
102   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
103   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
104 
105   const BlasiusContext           context  = (BlasiusContext)ctx;
106   const NewtonianIdealGasContext gas      = &context->newtonian_ctx;
107   const CeedScalar               mu       = context->newtonian_ctx.mu;
108   const CeedScalar               delta0   = context->delta0;
109   const CeedScalar               x_inflow = context->x_inflow;
110   CeedScalar                     t12;
111 
112   const CeedScalar Y_inf[5] = {context->P0, context->U_inf, 0, 0, context->T_inf};
113   const State      s_inf    = StateFromY(gas, Y_inf);
114 
115   const CeedScalar x0 = context->U_inf * s_inf.U.density / (mu * 25 / Square(delta0));
116 
117   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
118     const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]};
119     State            s    = BlasiusSolution(context, x, x0, x_inflow, s_inf.U.density, &t12);
120     CeedScalar       q[5] = {0};
121 
122     switch (gas->state_var) {
123       case STATEVAR_CONSERVATIVE:
124         UnpackState_U(s.U, q);
125         break;
126       case STATEVAR_PRIMITIVE:
127         UnpackState_Y(s.Y, q);
128         break;
129     }
130     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
131   }
132   return 0;
133 }
134 
135 // *****************************************************************************
136 CEED_QFUNCTION(Blasius_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
137   // Inputs
138   const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
139   const CeedScalar(*q_data_sur)    = in[2];
140   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3];
141 
142   // Outputs
143   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
144   CeedScalar(*jac_data_sur)  = out[1];
145 
146   const BlasiusContext     context     = (BlasiusContext)ctx;
147   const bool               is_implicit = context->implicit;
148   NewtonianIdealGasContext gas         = &context->newtonian_ctx;
149   const CeedScalar         mu          = context->newtonian_ctx.mu;
150   const CeedScalar         Rd          = GasConstant(&context->newtonian_ctx);
151   const CeedScalar         T_inf       = context->T_inf;
152   const CeedScalar         P0          = context->P0;
153   const CeedScalar         delta0      = context->delta0;
154   const CeedScalar         U_inf       = context->U_inf;
155   const CeedScalar         x_inflow    = context->x_inflow;
156   const bool               weakT       = context->weakT;
157   const CeedScalar         rho_0       = P0 / (Rd * T_inf);
158   const CeedScalar         x0          = U_inf * rho_0 / (mu * 25 / Square(delta0));
159   const CeedScalar         zeros[11]   = {0.};
160 
161   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
162     CeedScalar wdetJb, norm[3];
163     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
164     wdetJb *= is_implicit ? -1. : 1.;
165 
166     // Calculate inflow values
167     const CeedScalar x[3] = {X[0][i], X[1][i], 0.};
168     CeedScalar       t12;
169     State            s = BlasiusSolution(context, x, x0, x_inflow, rho_0, &t12);
170     CeedScalar       qi[5];
171     for (CeedInt j = 0; j < 5; j++) qi[j] = q[j][i];
172     State s_int = StateFromU(gas, qi);
173 
174     // enabling user to choose between weak T and weak rho inflow
175     if (weakT) {  // density from the current solution
176       s.U.density = s_int.U.density;
177       s.Y         = StatePrimitiveFromConservative(gas, s.U);
178     } else {  // Total energy from current solution
179       s.U.E_total = s_int.U.E_total;
180       s.Y         = StatePrimitiveFromConservative(gas, s.U);
181     }
182 
183     StateConservative Flux_inviscid[3];
184     FluxInviscid(&context->newtonian_ctx, s, Flux_inviscid);
185 
186     const CeedScalar stress[3][3] = {
187         {0,   t12, 0},
188         {t12, 0,   0},
189         {0,   0,   0}
190     };
191     const CeedScalar Fe[3] = {0};  // TODO: viscous energy flux needs grad temperature
192     CeedScalar       Flux[5];
193     FluxTotal_Boundary(Flux_inviscid, stress, Fe, norm, Flux);
194     for (CeedInt j = 0; j < 5; j++) v[j][i] = -wdetJb * Flux[j];
195     StoredValuesPack(Q, i, 0, 11, zeros, jac_data_sur);
196   }
197   return 0;
198 }
199 
200 // *****************************************************************************
201 CEED_QFUNCTION(Blasius_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
202   // Inputs
203   const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
204   const CeedScalar(*q_data_sur)     = in[2];
205   const CeedScalar(*X)[CEED_Q_VLA]  = (const CeedScalar(*)[CEED_Q_VLA])in[3];
206 
207   // Outputs
208   CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
209 
210   const BlasiusContext context     = (BlasiusContext)ctx;
211   const bool           is_implicit = context->implicit;
212   const CeedScalar     mu          = context->newtonian_ctx.mu;
213   const CeedScalar     cv          = context->newtonian_ctx.cv;
214   const CeedScalar     Rd          = GasConstant(&context->newtonian_ctx);
215   const CeedScalar     gamma       = HeatCapacityRatio(&context->newtonian_ctx);
216   const CeedScalar     T_inf       = context->T_inf;
217   const CeedScalar     P0          = context->P0;
218   const CeedScalar     delta0      = context->delta0;
219   const CeedScalar     U_inf       = context->U_inf;
220   const bool           weakT       = context->weakT;
221   const CeedScalar     rho_0       = P0 / (Rd * T_inf);
222   const CeedScalar     x0          = U_inf * rho_0 / (mu * 25 / (delta0 * delta0));
223 
224   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
225     CeedScalar wdetJb, norm[3];
226     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
227     wdetJb *= is_implicit ? -1. : 1.;
228 
229     // Calculate inflow values
230     const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]};
231     CeedScalar       t12;
232     State            s = BlasiusSolution(context, x, x0, 0, rho_0, &t12);
233 
234     // enabling user to choose between weak T and weak rho inflow
235     CeedScalar drho, dE, dP;
236     if (weakT) {
237       // rho should be from the current solution
238       drho                   = dq[0][i];
239       CeedScalar dE_internal = drho * cv * T_inf;
240       CeedScalar dE_kinetic  = .5 * drho * Dot3(s.Y.velocity, s.Y.velocity);
241       dE                     = dE_internal + dE_kinetic;
242       dP                     = drho * Rd * T_inf;  // interior rho with exterior T
243     } else {                                       // rho specified, E_internal from solution
244       drho = 0;
245       dE   = dq[4][i];
246       dP   = dE * (gamma - 1.);
247     }
248 
249     const CeedScalar u_normal = Dot3(norm, s.Y.velocity);
250 
251     v[0][i] = -wdetJb * drho * u_normal;
252     for (int j = 0; j < 3; j++) {
253       v[j + 1][i] = -wdetJb * (drho * u_normal * s.Y.velocity[j] + norm[j] * dP);
254     }
255     v[4][i] = -wdetJb * u_normal * (dE + dP);
256   }
257   return 0;
258 }
259