1 // Copyright (c) 2017-2022, 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 /// Structs and helper functions regarding the state of a newtonian simulation 10 11 #ifndef newtonian_state_h 12 #define newtonian_state_h 13 14 #include <ceed.h> 15 #include <math.h> 16 17 #include "newtonian_types.h" 18 #include "utils.h" 19 20 typedef struct { 21 CeedScalar density; 22 CeedScalar momentum[3]; 23 CeedScalar E_total; 24 } StateConservative; 25 26 typedef struct { 27 StateConservative U; 28 StatePrimitive Y; 29 } State; 30 31 CEED_QFUNCTION_HELPER void UnpackState_U(StateConservative s, CeedScalar U[5]) { 32 U[0] = s.density; 33 for (int i = 0; i < 3; i++) U[i + 1] = s.momentum[i]; 34 U[4] = s.E_total; 35 } 36 37 CEED_QFUNCTION_HELPER void UnpackState_Y(StatePrimitive s, CeedScalar Y[5]) { 38 Y[0] = s.pressure; 39 for (int i = 0; i < 3; i++) Y[i + 1] = s.velocity[i]; 40 Y[4] = s.temperature; 41 } 42 43 CEED_QFUNCTION_HELPER CeedScalar HeatCapacityRatio(NewtonianIdealGasContext gas) { return gas->cp / gas->cv; } 44 45 CEED_QFUNCTION_HELPER CeedScalar GasConstant(NewtonianIdealGasContext gas) { return gas->cp - gas->cv; } 46 47 CEED_QFUNCTION_HELPER CeedScalar Prandtl(NewtonianIdealGasContext gas) { return gas->cp * gas->mu / gas->k; } 48 49 CEED_QFUNCTION_HELPER CeedScalar SoundSpeed(NewtonianIdealGasContext gas, CeedScalar T) { return sqrt(gas->cp * (HeatCapacityRatio(gas) - 1.) * T); } 50 51 CEED_QFUNCTION_HELPER CeedScalar Mach(NewtonianIdealGasContext gas, CeedScalar T, CeedScalar u) { return u / SoundSpeed(gas, T); } 52 53 CEED_QFUNCTION_HELPER CeedScalar TotalSpecificEnthalpy(NewtonianIdealGasContext gas, const State s) { 54 // Ignoring potential energy 55 CeedScalar e_internal = gas->cv * s.Y.temperature; 56 CeedScalar e_kinetic = 0.5 * Dot3(s.Y.velocity, s.Y.velocity); 57 return e_internal + e_kinetic + s.Y.pressure / s.U.density; 58 } 59 60 CEED_QFUNCTION_HELPER CeedScalar TotalSpecificEnthalpy_fwd(NewtonianIdealGasContext gas, const State s, const State ds) { 61 // Ignoring potential energy 62 CeedScalar de_kinetic = Dot3(ds.Y.velocity, s.Y.velocity); 63 CeedScalar de_internal = gas->cv * ds.Y.temperature; 64 return de_internal + de_kinetic + ds.Y.pressure / s.U.density - s.Y.pressure / Square(s.U.density) * ds.U.density; 65 } 66 67 CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative(NewtonianIdealGasContext gas, StateConservative U) { 68 StatePrimitive Y; 69 for (CeedInt i = 0; i < 3; i++) Y.velocity[i] = U.momentum[i] / U.density; 70 CeedScalar e_kinetic = .5 * Dot3(Y.velocity, Y.velocity); 71 CeedScalar e_total = U.E_total / U.density; 72 CeedScalar e_internal = e_total - e_kinetic; 73 Y.temperature = e_internal / gas->cv; 74 Y.pressure = (HeatCapacityRatio(gas) - 1) * U.density * e_internal; 75 return Y; 76 } 77 78 CEED_QFUNCTION_HELPER StatePrimitive StatePrimitiveFromConservative_fwd(NewtonianIdealGasContext gas, State s, StateConservative dU) { 79 StatePrimitive dY; 80 for (CeedInt i = 0; i < 3; i++) { 81 dY.velocity[i] = (dU.momentum[i] - s.Y.velocity[i] * dU.density) / s.U.density; 82 } 83 CeedScalar e_kinetic = .5 * Dot3(s.Y.velocity, s.Y.velocity); 84 CeedScalar de_kinetic = Dot3(dY.velocity, s.Y.velocity); 85 CeedScalar e_total = s.U.E_total / s.U.density; 86 CeedScalar de_total = (dU.E_total - e_total * dU.density) / s.U.density; 87 CeedScalar e_internal = e_total - e_kinetic; 88 CeedScalar de_internal = de_total - de_kinetic; 89 dY.temperature = de_internal / gas->cv; 90 dY.pressure = (HeatCapacityRatio(gas) - 1) * (dU.density * e_internal + s.U.density * de_internal); 91 return dY; 92 } 93 94 CEED_QFUNCTION_HELPER StateConservative StateConservativeFromPrimitive(NewtonianIdealGasContext gas, StatePrimitive Y) { 95 StateConservative U; 96 U.density = Y.pressure / (GasConstant(gas) * Y.temperature); 97 for (int i = 0; i < 3; i++) U.momentum[i] = U.density * Y.velocity[i]; 98 CeedScalar e_internal = gas->cv * Y.temperature; 99 CeedScalar e_kinetic = .5 * Dot3(Y.velocity, Y.velocity); 100 CeedScalar e_total = e_internal + e_kinetic; 101 U.E_total = U.density * e_total; 102 return U; 103 } 104 105 CEED_QFUNCTION_HELPER StateConservative StateConservativeFromPrimitive_fwd(NewtonianIdealGasContext gas, State s, StatePrimitive dY) { 106 StateConservative dU; 107 dU.density = (dY.pressure * s.Y.temperature - s.Y.pressure * dY.temperature) / (GasConstant(gas) * s.Y.temperature * s.Y.temperature); 108 for (int i = 0; i < 3; i++) { 109 dU.momentum[i] = dU.density * s.Y.velocity[i] + s.U.density * dY.velocity[i]; 110 } 111 CeedScalar e_kinetic = .5 * Dot3(s.Y.velocity, s.Y.velocity); 112 CeedScalar de_kinetic = Dot3(dY.velocity, s.Y.velocity); 113 CeedScalar e_internal = gas->cv * s.Y.temperature; 114 CeedScalar de_internal = gas->cv * dY.temperature; 115 CeedScalar e_total = e_internal + e_kinetic; 116 CeedScalar de_total = de_internal + de_kinetic; 117 dU.E_total = dU.density * e_total + s.U.density * de_total; 118 return dU; 119 } 120 121 CEED_QFUNCTION_HELPER State StateFromPrimitive(NewtonianIdealGasContext gas, StatePrimitive Y) { 122 StateConservative U = StateConservativeFromPrimitive(gas, Y); 123 State s; 124 s.U = U; 125 s.Y = Y; 126 return s; 127 } 128 129 CEED_QFUNCTION_HELPER State StateFromPrimitive_fwd(NewtonianIdealGasContext gas, State s, StatePrimitive dY) { 130 StateConservative dU = StateConservativeFromPrimitive_fwd(gas, s, dY); 131 State ds; 132 ds.U = dU; 133 ds.Y = dY; 134 return ds; 135 } 136 137 // linear combination of n states 138 CEED_QFUNCTION_HELPER StateConservative StateConservativeMult(CeedInt n, const CeedScalar a[], const StateConservative X[]) { 139 StateConservative R = {0}; 140 for (CeedInt i = 0; i < n; i++) { 141 R.density += a[i] * X[i].density; 142 for (int j = 0; j < 3; j++) R.momentum[j] += a[i] * X[i].momentum[j]; 143 R.E_total += a[i] * X[i].E_total; 144 } 145 return R; 146 } 147 148 CEED_QFUNCTION_HELPER StateConservative StateConservativeAXPBYPCZ(CeedScalar a, StateConservative X, CeedScalar b, StateConservative Y, CeedScalar c, 149 StateConservative Z) { 150 StateConservative R; 151 R.density = a * X.density + b * Y.density + c * Z.density; 152 for (int i = 0; i < 3; i++) R.momentum[i] = a * X.momentum[i] + b * Y.momentum[i] + c * Z.momentum[i]; 153 R.E_total = a * X.E_total + b * Y.E_total + c * Z.E_total; 154 return R; 155 } 156 157 CEED_QFUNCTION_HELPER void StateToU(NewtonianIdealGasContext gas, const State input, CeedScalar U[5]) { UnpackState_U(input.U, U); } 158 159 CEED_QFUNCTION_HELPER void StateToY(NewtonianIdealGasContext gas, const State input, CeedScalar Y[5]) { UnpackState_Y(input.Y, Y); } 160 161 CEED_QFUNCTION_HELPER void StateToQ(NewtonianIdealGasContext gas, const State input, CeedScalar Q[5], StateVariable state_var) { 162 switch (state_var) { 163 case STATEVAR_CONSERVATIVE: 164 StateToU(gas, input, Q); 165 break; 166 case STATEVAR_PRIMITIVE: 167 StateToY(gas, input, Q); 168 break; 169 } 170 } 171 172 CEED_QFUNCTION_HELPER State StateFromU(NewtonianIdealGasContext gas, const CeedScalar U[5]) { 173 State s; 174 s.U.density = U[0]; 175 s.U.momentum[0] = U[1]; 176 s.U.momentum[1] = U[2]; 177 s.U.momentum[2] = U[3]; 178 s.U.E_total = U[4]; 179 s.Y = StatePrimitiveFromConservative(gas, s.U); 180 return s; 181 } 182 183 CEED_QFUNCTION_HELPER State StateFromU_fwd(NewtonianIdealGasContext gas, State s, const CeedScalar dU[5]) { 184 State ds; 185 ds.U.density = dU[0]; 186 ds.U.momentum[0] = dU[1]; 187 ds.U.momentum[1] = dU[2]; 188 ds.U.momentum[2] = dU[3]; 189 ds.U.E_total = dU[4]; 190 ds.Y = StatePrimitiveFromConservative_fwd(gas, s, ds.U); 191 return ds; 192 } 193 194 CEED_QFUNCTION_HELPER State StateFromY(NewtonianIdealGasContext gas, const CeedScalar Y[5]) { 195 State s; 196 s.Y.pressure = Y[0]; 197 s.Y.velocity[0] = Y[1]; 198 s.Y.velocity[1] = Y[2]; 199 s.Y.velocity[2] = Y[3]; 200 s.Y.temperature = Y[4]; 201 s.U = StateConservativeFromPrimitive(gas, s.Y); 202 return s; 203 } 204 205 CEED_QFUNCTION_HELPER State StateFromY_fwd(NewtonianIdealGasContext gas, State s, const CeedScalar dY[5]) { 206 State ds; 207 ds.Y.pressure = dY[0]; 208 ds.Y.velocity[0] = dY[1]; 209 ds.Y.velocity[1] = dY[2]; 210 ds.Y.velocity[2] = dY[3]; 211 ds.Y.temperature = dY[4]; 212 ds.U = StateConservativeFromPrimitive_fwd(gas, s, ds.Y); 213 return ds; 214 } 215 216 CEED_QFUNCTION_HELPER State StateFromQ(NewtonianIdealGasContext gas, const CeedScalar Q[5], StateVariable state_var) { 217 State s; 218 switch (state_var) { 219 case STATEVAR_CONSERVATIVE: 220 s = StateFromU(gas, Q); 221 break; 222 case STATEVAR_PRIMITIVE: 223 s = StateFromY(gas, Q); 224 break; 225 } 226 return s; 227 } 228 229 CEED_QFUNCTION_HELPER State StateFromQ_fwd(NewtonianIdealGasContext gas, State s, const CeedScalar dQ[5], StateVariable state_var) { 230 State ds; 231 switch (state_var) { 232 case STATEVAR_CONSERVATIVE: 233 ds = StateFromU_fwd(gas, s, dQ); 234 break; 235 case STATEVAR_PRIMITIVE: 236 ds = StateFromY_fwd(gas, s, dQ); 237 break; 238 } 239 return ds; 240 } 241 242 CEED_QFUNCTION_HELPER void FluxInviscid(NewtonianIdealGasContext gas, State s, StateConservative Flux[3]) { 243 for (CeedInt i = 0; i < 3; i++) { 244 Flux[i].density = s.U.momentum[i]; 245 for (CeedInt j = 0; j < 3; j++) Flux[i].momentum[j] = s.U.momentum[i] * s.Y.velocity[j] + s.Y.pressure * (i == j); 246 Flux[i].E_total = (s.U.E_total + s.Y.pressure) * s.Y.velocity[i]; 247 } 248 } 249 250 CEED_QFUNCTION_HELPER void FluxInviscid_fwd(NewtonianIdealGasContext gas, State s, State ds, StateConservative dFlux[3]) { 251 for (CeedInt i = 0; i < 3; i++) { 252 dFlux[i].density = ds.U.momentum[i]; 253 for (CeedInt j = 0; j < 3; j++) { 254 dFlux[i].momentum[j] = ds.U.momentum[i] * s.Y.velocity[j] + s.U.momentum[i] * ds.Y.velocity[j] + ds.Y.pressure * (i == j); 255 } 256 dFlux[i].E_total = (ds.U.E_total + ds.Y.pressure) * s.Y.velocity[i] + (s.U.E_total + s.Y.pressure) * ds.Y.velocity[i]; 257 } 258 } 259 260 CEED_QFUNCTION_HELPER StateConservative FluxInviscidDotNormal(NewtonianIdealGasContext gas, State s, const CeedScalar normal[3]) { 261 StateConservative Flux[3], Flux_dot_n = {0}; 262 FluxInviscid(gas, s, Flux); 263 for (CeedInt i = 0; i < 3; i++) { 264 Flux_dot_n.density += Flux[i].density * normal[i]; 265 for (CeedInt j = 0; j < 3; j++) Flux_dot_n.momentum[j] += Flux[i].momentum[j] * normal[i]; 266 Flux_dot_n.E_total += Flux[i].E_total * normal[i]; 267 } 268 return Flux_dot_n; 269 } 270 271 CEED_QFUNCTION_HELPER StateConservative FluxInviscidDotNormal_fwd(NewtonianIdealGasContext gas, State s, State ds, const CeedScalar normal[3]) { 272 StateConservative dFlux[3], Flux_dot_n = {0}; 273 FluxInviscid_fwd(gas, s, ds, dFlux); 274 for (CeedInt i = 0; i < 3; i++) { 275 Flux_dot_n.density += dFlux[i].density * normal[i]; 276 for (CeedInt j = 0; j < 3; j++) Flux_dot_n.momentum[j] += dFlux[i].momentum[j] * normal[i]; 277 Flux_dot_n.E_total += dFlux[i].E_total * normal[i]; 278 } 279 return Flux_dot_n; 280 } 281 282 CEED_QFUNCTION_HELPER void FluxInviscidStrong(NewtonianIdealGasContext gas, State s, State ds[3], CeedScalar strong_conv[5]) { 283 for (CeedInt i = 0; i < 5; i++) strong_conv[i] = 0; 284 for (CeedInt i = 0; i < 3; i++) { 285 StateConservative dF[3]; 286 FluxInviscid_fwd(gas, s, ds[i], dF); 287 CeedScalar dF_i[5]; 288 UnpackState_U(dF[i], dF_i); 289 for (CeedInt j = 0; j < 5; j++) strong_conv[j] += dF_i[j]; 290 } 291 } 292 293 CEED_QFUNCTION_HELPER void FluxTotal(const StateConservative F_inviscid[3], CeedScalar stress[3][3], CeedScalar Fe[3], CeedScalar Flux[5][3]) { 294 for (CeedInt j = 0; j < 3; j++) { 295 Flux[0][j] = F_inviscid[j].density; 296 for (CeedInt k = 0; k < 3; k++) Flux[k + 1][j] = F_inviscid[j].momentum[k] - stress[k][j]; 297 Flux[4][j] = F_inviscid[j].E_total + Fe[j]; 298 } 299 } 300 301 CEED_QFUNCTION_HELPER void FluxTotal_Boundary(const StateConservative F_inviscid[3], const CeedScalar stress[3][3], const CeedScalar Fe[3], 302 const CeedScalar normal[3], CeedScalar Flux[5]) { 303 for (CeedInt j = 0; j < 5; j++) Flux[j] = 0.; 304 for (CeedInt j = 0; j < 3; j++) { 305 Flux[0] += F_inviscid[j].density * normal[j]; 306 for (CeedInt k = 0; k < 3; k++) { 307 Flux[k + 1] += (F_inviscid[j].momentum[k] - stress[k][j]) * normal[j]; 308 } 309 Flux[4] += (F_inviscid[j].E_total + Fe[j]) * normal[j]; 310 } 311 } 312 313 CEED_QFUNCTION_HELPER void FluxTotal_RiemannBoundary(const StateConservative F_inviscid_normal, const CeedScalar stress[3][3], const CeedScalar Fe[3], 314 const CeedScalar normal[3], CeedScalar Flux[5]) { 315 Flux[0] = F_inviscid_normal.density; 316 for (CeedInt k = 0; k < 3; k++) Flux[k + 1] = F_inviscid_normal.momentum[k]; 317 Flux[4] = F_inviscid_normal.E_total; 318 for (CeedInt j = 0; j < 3; j++) { 319 for (CeedInt k = 0; k < 3; k++) { 320 Flux[k + 1] -= stress[k][j] * normal[j]; 321 } 322 Flux[4] += Fe[j] * normal[j]; 323 } 324 } 325 326 CEED_QFUNCTION_HELPER void VelocityGradient(const State grad_s[3], CeedScalar grad_velocity[3][3]) { 327 grad_velocity[0][0] = grad_s[0].Y.velocity[0]; 328 grad_velocity[0][1] = grad_s[1].Y.velocity[0]; 329 grad_velocity[0][2] = grad_s[2].Y.velocity[0]; 330 grad_velocity[1][0] = grad_s[0].Y.velocity[1]; 331 grad_velocity[1][1] = grad_s[1].Y.velocity[1]; 332 grad_velocity[1][2] = grad_s[2].Y.velocity[1]; 333 grad_velocity[2][0] = grad_s[0].Y.velocity[2]; 334 grad_velocity[2][1] = grad_s[1].Y.velocity[2]; 335 grad_velocity[2][2] = grad_s[2].Y.velocity[2]; 336 } 337 338 CEED_QFUNCTION_HELPER void KMStrainRate(const CeedScalar grad_velocity[3][3], CeedScalar strain_rate[6]) { 339 const CeedScalar weight = 1 / sqrt(2.); // Really sqrt(2.) / 2 340 strain_rate[0] = grad_velocity[0][0]; 341 strain_rate[1] = grad_velocity[1][1]; 342 strain_rate[2] = grad_velocity[2][2]; 343 strain_rate[3] = weight * (grad_velocity[1][2] + grad_velocity[2][1]); 344 strain_rate[4] = weight * (grad_velocity[0][2] + grad_velocity[2][0]); 345 strain_rate[5] = weight * (grad_velocity[0][1] + grad_velocity[1][0]); 346 } 347 348 // Kelvin-Mandel notation 349 CEED_QFUNCTION_HELPER void KMStrainRate_State(const State grad_s[3], CeedScalar strain_rate[6]) { 350 CeedScalar grad_velocity[3][3]; 351 VelocityGradient(grad_s, grad_velocity); 352 KMStrainRate(grad_velocity, strain_rate); 353 } 354 355 CEED_QFUNCTION_HELPER void RotationRate(const CeedScalar grad_velocity[3][3], CeedScalar rotation_rate[3][3]) { 356 rotation_rate[0][0] = 0; 357 rotation_rate[1][1] = 0; 358 rotation_rate[2][2] = 0; 359 rotation_rate[1][2] = 0.5 * (grad_velocity[1][2] - grad_velocity[2][1]); 360 rotation_rate[0][2] = 0.5 * (grad_velocity[0][2] - grad_velocity[2][0]); 361 rotation_rate[0][1] = 0.5 * (grad_velocity[0][1] - grad_velocity[1][0]); 362 rotation_rate[2][1] = -rotation_rate[1][2]; 363 rotation_rate[2][0] = -rotation_rate[0][2]; 364 rotation_rate[1][0] = -rotation_rate[0][1]; 365 } 366 367 CEED_QFUNCTION_HELPER void NewtonianStress(NewtonianIdealGasContext gas, const CeedScalar strain_rate[6], CeedScalar stress[6]) { 368 CeedScalar div_u = strain_rate[0] + strain_rate[1] + strain_rate[2]; 369 for (CeedInt i = 0; i < 6; i++) { 370 stress[i] = gas->mu * (2 * strain_rate[i] + gas->lambda * div_u * (i < 3)); 371 } 372 } 373 374 CEED_QFUNCTION_HELPER void ViscousEnergyFlux(NewtonianIdealGasContext gas, StatePrimitive Y, const State grad_s[3], const CeedScalar stress[3][3], 375 CeedScalar Fe[3]) { 376 for (CeedInt i = 0; i < 3; i++) { 377 Fe[i] = -Y.velocity[0] * stress[0][i] - Y.velocity[1] * stress[1][i] - Y.velocity[2] * stress[2][i] - gas->k * grad_s[i].Y.temperature; 378 } 379 } 380 381 CEED_QFUNCTION_HELPER void ViscousEnergyFlux_fwd(NewtonianIdealGasContext gas, StatePrimitive Y, StatePrimitive dY, const State grad_ds[3], 382 const CeedScalar stress[3][3], const CeedScalar dstress[3][3], CeedScalar dFe[3]) { 383 for (CeedInt i = 0; i < 3; i++) { 384 dFe[i] = -Y.velocity[0] * dstress[0][i] - dY.velocity[0] * stress[0][i] - Y.velocity[1] * dstress[1][i] - dY.velocity[1] * stress[1][i] - 385 Y.velocity[2] * dstress[2][i] - dY.velocity[2] * stress[2][i] - gas->k * grad_ds[i].Y.temperature; 386 } 387 } 388 389 CEED_QFUNCTION_HELPER void Vorticity(const State grad_s[3], CeedScalar vorticity[3]) { 390 CeedScalar grad_velocity[3][3]; 391 VelocityGradient(grad_s, grad_velocity); 392 Curl3(grad_velocity, vorticity); 393 } 394 395 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference(CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, StateVariable state_var, 396 const CeedScalar *grad_q, const CeedScalar dXdx[3][3], State grad_s[3]) { 397 for (CeedInt k = 0; k < 3; k++) { 398 CeedScalar dqi[5]; 399 for (CeedInt j = 0; j < 5; j++) { 400 dqi[j] = 401 grad_q[(Q * 5) * 0 + Q * j + i] * dXdx[0][k] + grad_q[(Q * 5) * 1 + Q * j + i] * dXdx[1][k] + grad_q[(Q * 5) * 2 + Q * j + i] * dXdx[2][k]; 402 } 403 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 404 } 405 } 406 407 CEED_QFUNCTION_HELPER void StatePhysicalGradientFromReference_Boundary(CeedInt Q, CeedInt i, NewtonianIdealGasContext gas, State s, 408 StateVariable state_var, const CeedScalar *grad_q, const CeedScalar dXdx[2][3], 409 State grad_s[3]) { 410 for (CeedInt k = 0; k < 3; k++) { 411 CeedScalar dqi[5]; 412 for (CeedInt j = 0; j < 5; j++) { 413 dqi[j] = grad_q[(Q * 5) * 0 + Q * j + i] * dXdx[0][k] + grad_q[(Q * 5) * 1 + Q * j + i] * dXdx[1][k]; 414 } 415 grad_s[k] = StateFromQ_fwd(gas, s, dqi, state_var); 416 } 417 } 418 419 #endif // newtonian_state_h 420