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