xref: /libCEED/examples/fluids/qfunctions/blasius.h (revision fb133d4b0da93cd2b1f52e39cd1305dc4c4450c3)
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   const BlasiusContext context     = (BlasiusContext)ctx;
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   CeedScalar(*v)[CEED_Q_VLA]       = (CeedScalar(*)[CEED_Q_VLA])out[0];
142   CeedScalar(*jac_data_sur)        = context->newtonian_ctx.is_implicit ? out[1] : NULL;
143 
144   const bool                     is_implicit = context->implicit;
145   const NewtonianIdealGasContext gas         = &context->newtonian_ctx;
146   const CeedScalar               rho_0       = context->P0 / (GasConstant(gas) * context->T_inf);
147   const CeedScalar               x0          = context->U_inf * rho_0 / (gas->mu * 25 / Square(context->delta0));
148   const CeedScalar               zeros[11]   = {0.};
149 
150   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
151     CeedScalar wdetJb, norm[3];
152     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
153     wdetJb *= is_implicit ? -1. : 1.;
154 
155     // Calculate inflow values
156     const CeedScalar x[3] = {X[0][i], X[1][i], 0.};
157     CeedScalar       t12;
158     State            s = BlasiusSolution(context, x, x0, context->x_inflow, rho_0, &t12);
159     CeedScalar       qi[5];
160     for (CeedInt j = 0; j < 5; j++) qi[j] = q[j][i];
161     State s_int = StateFromU(gas, qi);
162 
163     // enabling user to choose between weak T and weak rho inflow
164     if (context->weakT) {  // density from the current solution
165       s.U.density = s_int.U.density;
166       s.Y         = StatePrimitiveFromConservative(gas, s.U);
167     } else {  // Total energy from current solution
168       s.U.E_total = s_int.U.E_total;
169       s.Y         = StatePrimitiveFromConservative(gas, s.U);
170     }
171 
172     StateConservative Flux_inviscid[3];
173     FluxInviscid(&context->newtonian_ctx, s, Flux_inviscid);
174 
175     const CeedScalar stress[3][3] = {
176         {0,   t12, 0},
177         {t12, 0,   0},
178         {0,   0,   0}
179     };
180     const CeedScalar Fe[3] = {0};  // TODO: viscous energy flux needs grad temperature
181     CeedScalar       Flux[5];
182     FluxTotal_Boundary(Flux_inviscid, stress, Fe, norm, Flux);
183     for (CeedInt j = 0; j < 5; j++) v[j][i] = -wdetJb * Flux[j];
184     if (is_implicit) StoredValuesPack(Q, i, 0, 11, zeros, jac_data_sur);
185   }
186   return 0;
187 }
188 
189 // *****************************************************************************
190 CEED_QFUNCTION(Blasius_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
191   const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
192   const CeedScalar(*q_data_sur)     = in[2];
193   const CeedScalar(*X)[CEED_Q_VLA]  = (const CeedScalar(*)[CEED_Q_VLA])in[3];
194   CeedScalar(*v)[CEED_Q_VLA]        = (CeedScalar(*)[CEED_Q_VLA])out[0];
195 
196   const BlasiusContext           context     = (BlasiusContext)ctx;
197   const NewtonianIdealGasContext gas         = &context->newtonian_ctx;
198   const bool                     is_implicit = context->implicit;
199   const CeedScalar               Rd          = GasConstant(gas);
200   const CeedScalar               gamma       = HeatCapacityRatio(gas);
201   const CeedScalar               rho_0       = context->P0 / (Rd * context->T_inf);
202   const CeedScalar               x0          = context->U_inf * rho_0 / (gas->mu * 25 / (Square(context->delta0)));
203 
204   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
205     CeedScalar wdetJb, norm[3];
206     QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm);
207     wdetJb *= is_implicit ? -1. : 1.;
208 
209     // Calculate inflow values
210     const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]};
211     CeedScalar       t12;
212     State            s = BlasiusSolution(context, x, x0, 0, rho_0, &t12);
213 
214     // enabling user to choose between weak T and weak rho inflow
215     CeedScalar drho, dE, dP;
216     if (context->weakT) {
217       // rho should be from the current solution
218       drho                   = dq[0][i];
219       CeedScalar dE_internal = drho * gas->cv * context->T_inf;
220       CeedScalar dE_kinetic  = .5 * drho * Dot3(s.Y.velocity, s.Y.velocity);
221       dE                     = dE_internal + dE_kinetic;
222       dP                     = drho * Rd * context->T_inf;  // interior rho with exterior T
223     } else {                                                // rho specified, E_internal from solution
224       drho = 0;
225       dE   = dq[4][i];
226       dP   = dE * (gamma - 1.);
227     }
228 
229     const CeedScalar u_normal = Dot3(norm, s.Y.velocity);
230 
231     v[0][i] = -wdetJb * drho * u_normal;
232     for (int j = 0; j < 3; j++) {
233       v[j + 1][i] = -wdetJb * (drho * u_normal * s.Y.velocity[j] + norm[j] * dP);
234     }
235     v[4][i] = -wdetJb * u_normal * (dE + dP);
236   }
237   return 0;
238 }
239