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 State S_infty; 24 CeedScalar T_wall; // !< Temperature at the wall 25 CeedScalar x_inflow; // !< Location of inflow in x 26 CeedScalar n_cheb; // !< Number of Chebyshev terms 27 CeedScalar *X; // !< Chebyshev polynomial coordinate vector (CPU only) 28 CeedScalar eta_max; // !< Maximum eta in the domain 29 CeedScalar Tf_cheb[BLASIUS_MAX_N_CHEBYSHEV]; // !< Chebyshev coefficient for f 30 CeedScalar Th_cheb[BLASIUS_MAX_N_CHEBYSHEV - 1]; // !< Chebyshev coefficient for h 31 struct NewtonianIdealGasContext_ newtonian_ctx; 32 }; 33 34 // ***************************************************************************** 35 // This helper function evaluates Chebyshev polynomials with a set of coefficients with all their derivatives represented as a recurrence table. 36 // ***************************************************************************** 37 CEED_QFUNCTION_HELPER void ChebyshevEval(int N, const double *Tf, double x, double eta_max, double *f) { 38 double dX_deta = 2 / eta_max; 39 double table[4][3] = { 40 // Chebyshev polynomials T_0, T_1, T_2 of the first kind in (-1,1) 41 {1, x, 2 * x * x - 1}, 42 {0, 1, 4 * x }, 43 {0, 0, 4 }, 44 {0, 0, 0 } 45 }; 46 for (int i = 0; i < 4; i++) { 47 // i-th derivative of f 48 f[i] = table[i][0] * Tf[0] + table[i][1] * Tf[1] + table[i][2] * Tf[2]; 49 } 50 for (int i = 3; i < N; i++) { 51 // T_n(x) = 2xT_{n-1}(x) - T_{n-2}(x) 52 table[0][i % 3] = 2 * x * table[0][(i - 1) % 3] - table[0][(i - 2) % 3]; 53 // Differentiate Chebyshev polynomials with the recurrence relation 54 for (int j = 1; j < 4; j++) { 55 // T'_{n}(x)/n = 2T_{n-1}(x) + T'_{n-2}(x)/n-2 56 table[j][i % 3] = i * (2 * table[j - 1][(i - 1) % 3] + table[j][(i - 2) % 3] / (i - 2)); 57 } 58 for (int j = 0; j < 4; j++) { 59 f[j] += table[j][i % 3] * Tf[i]; 60 } 61 } 62 for (int i = 1; i < 4; i++) { 63 // Transform derivatives from Chebyshev [-1, 1] to [0, eta_max]. 64 for (int j = 0; j < i; j++) f[i] *= dX_deta; 65 } 66 } 67 68 // ***************************************************************************** 69 // This helper function computes the Blasius boundary layer solution. 70 // ***************************************************************************** 71 State CEED_QFUNCTION_HELPER(BlasiusSolution)(const BlasiusContext blasius, const CeedScalar x[3], const CeedScalar x0, const CeedScalar x_inflow, 72 const CeedScalar rho_infty, CeedScalar *t12) { 73 CeedInt N = blasius->n_cheb; 74 CeedScalar mu = blasius->newtonian_ctx.mu; 75 State S_infty = blasius->S_infty; 76 CeedScalar nu = mu / rho_infty; 77 CeedScalar U_infty = sqrt(Dot3(S_infty.Y.velocity, S_infty.Y.velocity)); 78 CeedScalar eta = x[1] * sqrt(U_infty / (nu * (x0 + x[0] - x_inflow))); 79 CeedScalar X = 2 * (eta / blasius->eta_max) - 1.; 80 CeedScalar Rd = GasConstant(&blasius->newtonian_ctx); 81 82 CeedScalar f[4], h[4]; 83 ChebyshevEval(N, blasius->Tf_cheb, X, blasius->eta_max, f); 84 ChebyshevEval(N - 1, blasius->Th_cheb, X, blasius->eta_max, h); 85 86 *t12 = mu * U_infty * f[2] * sqrt(U_infty / (nu * (x0 + x[0] - x_inflow))); 87 88 CeedScalar Y[5]; 89 Y[1] = U_infty * f[1]; 90 Y[2] = 0.5 * sqrt(nu * U_infty / (x0 + x[0] - x_inflow)) * (eta * f[1] - f[0]); 91 Y[3] = 0.; 92 Y[4] = S_infty.Y.temperature * h[0]; 93 Y[0] = rho_infty / h[0] * Rd * Y[4]; 94 return StateFromY(&blasius->newtonian_ctx, Y); 95 } 96 97 // ***************************************************************************** 98 // This QFunction sets a Blasius boundary layer for the initial condition 99 // ***************************************************************************** 100 CEED_QFUNCTION(ICsBlasius)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 101 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 102 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 103 104 const BlasiusContext context = (BlasiusContext)ctx; 105 const NewtonianIdealGasContext gas = &context->newtonian_ctx; 106 const CeedScalar mu = context->newtonian_ctx.mu; 107 const CeedScalar delta0 = context->delta0; 108 const CeedScalar x_inflow = context->x_inflow; 109 CeedScalar t12; 110 111 const State S_infty = context->S_infty; 112 const CeedScalar U_infty = sqrt(Dot3(S_infty.Y.velocity, S_infty.Y.velocity)); 113 114 const CeedScalar x0 = U_infty * S_infty.U.density / (mu * 25 / Square(delta0)); 115 116 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 117 const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 118 State s = BlasiusSolution(context, x, x0, x_inflow, S_infty.U.density, &t12); 119 CeedScalar q[5] = {0}; 120 121 StateToQ(gas, s, q, gas->state_var); 122 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 123 } 124 return 0; 125 } 126 127 // ***************************************************************************** 128 CEED_QFUNCTION(Blasius_Inflow)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 129 const BlasiusContext context = (BlasiusContext)ctx; 130 const CeedScalar(*q)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 131 const CeedScalar(*q_data_sur) = in[2]; 132 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 133 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 134 CeedScalar(*jac_data_sur) = context->newtonian_ctx.is_implicit ? out[1] : NULL; 135 136 const bool is_implicit = context->implicit; 137 const NewtonianIdealGasContext gas = &context->newtonian_ctx; 138 State S_infty = context->S_infty; 139 const CeedScalar rho_0 = S_infty.U.density; 140 const CeedScalar U_infty = sqrt(Dot3(S_infty.Y.velocity, S_infty.Y.velocity)); 141 const CeedScalar x0 = U_infty * rho_0 / (gas->mu * 25 / Square(context->delta0)); 142 const CeedScalar zeros[11] = {0.}; 143 144 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 145 CeedScalar wdetJb, norm[3]; 146 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 147 wdetJb *= is_implicit ? -1. : 1.; 148 149 // Calculate inflow values 150 const CeedScalar x[3] = {X[0][i], X[1][i], 0.}; 151 CeedScalar t12; 152 State s = BlasiusSolution(context, x, x0, context->x_inflow, rho_0, &t12); 153 CeedScalar qi[5]; 154 for (CeedInt j = 0; j < 5; j++) qi[j] = q[j][i]; 155 State s_int = StateFromU(gas, qi); 156 157 // enabling user to choose between weak T and weak rho inflow 158 if (context->weakT) { // density from the current solution 159 s.U.density = s_int.U.density; 160 s.Y = StatePrimitiveFromConservative(gas, s.U); 161 s.V = StateEntropyFromConservative(gas, s.U); 162 } else { // Total energy from current solution 163 s.U.E_total = s_int.U.E_total; 164 s.Y = StatePrimitiveFromConservative(gas, s.U); 165 s.V = StateEntropyFromConservative(gas, s.U); 166 } 167 168 StateConservative Flux_inviscid[3]; 169 FluxInviscid(&context->newtonian_ctx, s, Flux_inviscid); 170 171 const CeedScalar stress[3][3] = { 172 {0, t12, 0}, 173 {t12, 0, 0}, 174 {0, 0, 0} 175 }; 176 const CeedScalar Fe[3] = {0}; // TODO: viscous energy flux needs grad temperature 177 CeedScalar Flux[5]; 178 FluxTotal_Boundary(Flux_inviscid, stress, Fe, norm, Flux); 179 for (CeedInt j = 0; j < 5; j++) v[j][i] = -wdetJb * Flux[j]; 180 if (is_implicit) StoredValuesPack(Q, i, 0, 11, zeros, jac_data_sur); 181 } 182 return 0; 183 } 184 185 // ***************************************************************************** 186 CEED_QFUNCTION(Blasius_Inflow_Jacobian)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 187 const CeedScalar(*dq)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 188 const CeedScalar(*q_data_sur) = in[2]; 189 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[3]; 190 CeedScalar(*v)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 191 192 const BlasiusContext context = (BlasiusContext)ctx; 193 const NewtonianIdealGasContext gas = &context->newtonian_ctx; 194 const bool is_implicit = context->implicit; 195 const CeedScalar Rd = GasConstant(gas); 196 const CeedScalar gamma = HeatCapacityRatio(gas); 197 const State S_infty = context->S_infty; 198 const CeedScalar rho_0 = S_infty.U.density; 199 const CeedScalar U_infty = sqrt(Dot3(S_infty.Y.velocity, S_infty.Y.velocity)); 200 const CeedScalar x0 = U_infty * rho_0 / (gas->mu * 25 / Square(context->delta0)); 201 202 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 203 CeedScalar wdetJb, norm[3]; 204 QdataBoundaryUnpack_3D(Q, i, q_data_sur, &wdetJb, NULL, norm); 205 wdetJb *= is_implicit ? -1. : 1.; 206 207 // Calculate inflow values 208 const CeedScalar x[3] = {X[0][i], X[1][i], X[2][i]}; 209 CeedScalar t12; 210 State s = BlasiusSolution(context, x, x0, 0, rho_0, &t12); 211 212 // enabling user to choose between weak T and weak rho inflow 213 CeedScalar drho, dE, dP; 214 if (context->weakT) { 215 // rho should be from the current solution 216 drho = dq[0][i]; 217 CeedScalar dE_internal = drho * gas->cv * S_infty.Y.temperature; 218 CeedScalar dE_kinetic = .5 * drho * Dot3(s.Y.velocity, s.Y.velocity); 219 dE = dE_internal + dE_kinetic; 220 dP = drho * Rd * S_infty.Y.temperature; // interior rho with exterior T 221 } else { 222 // rho specified, E_internal from solution 223 drho = 0; 224 dE = dq[4][i]; 225 dP = dE * (gamma - 1.); 226 } 227 228 const CeedScalar u_normal = Dot3(norm, s.Y.velocity); 229 230 v[0][i] = -wdetJb * drho * u_normal; 231 for (int j = 0; j < 3; j++) { 232 v[j + 1][i] = -wdetJb * (drho * u_normal * s.Y.velocity[j] + norm[j] * dP); 233 } 234 v[4][i] = -wdetJb * u_normal * (dE + dP); 235 } 236 return 0; 237 } 238