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 /// Density current initial condition and operator for Navier-Stokes example using PETSc 10 11 // Model from: 12 // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 13 // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 14 15 #ifndef densitycurrent_h 16 #define densitycurrent_h 17 18 #include <ceed.h> 19 #include <math.h> 20 21 #include "newtonian_state.h" 22 #include "newtonian_types.h" 23 #include "utils.h" 24 25 typedef struct DensityCurrentContext_ *DensityCurrentContext; 26 struct DensityCurrentContext_ { 27 CeedScalar theta0; 28 CeedScalar thetaC; 29 CeedScalar P0; 30 CeedScalar N; 31 CeedScalar rc; 32 CeedScalar center[3]; 33 CeedScalar dc_axis[3]; 34 struct NewtonianIdealGasContext_ newtonian_ctx; 35 }; 36 37 // ***************************************************************************** 38 // This function sets the initial conditions and the boundary conditions 39 // 40 // These initial conditions are given in terms of potential temperature and Exner pressure and then converted to density and total energy. 41 // Initial momentum density is zero. 42 // 43 // Initial Conditions: 44 // Potential Temperature: 45 // theta = thetabar + delta_theta 46 // thetabar = theta0 exp( N**2 z / g ) 47 // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 48 // r > rc : 0 49 // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 50 // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 51 // Exner Pressure: 52 // Pi = Pibar + deltaPi 53 // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 54 // deltaPi = 0 (hydrostatic balance) 55 // Velocity/Momentum Density: 56 // Ui = ui = 0 57 // 58 // Conversion to Conserved Variables: 59 // rho = P0 Pi**(cv/Rd) / (Rd theta) 60 // E = rho (cv T + (u u)/2 + g z) 61 // 62 // Boundary Conditions: 63 // Mass Density: 64 // 0.0 flux 65 // Momentum Density: 66 // 0.0 67 // Energy Density: 68 // 0.0 flux 69 // 70 // Constants: 71 // theta0 , Potential temperature constant 72 // thetaC , Potential temperature perturbation 73 // P0 , Pressure at the surface 74 // N , Brunt-Vaisala frequency 75 // cv , Specific heat, constant volume 76 // cp , Specific heat, constant pressure 77 // Rd = cp - cv, Specific heat difference 78 // g , Gravity 79 // rc , Characteristic radius of thermal bubble 80 // center , Location of bubble center 81 // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 82 // ***************************************************************************** 83 84 // ***************************************************************************** 85 // This helper function provides support for the exact, time-dependent solution 86 // (currently not implemented) and IC formulation for density current 87 // ***************************************************************************** 88 CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) { 89 // Context 90 const DensityCurrentContext context = (DensityCurrentContext)ctx; 91 const CeedScalar theta0 = context->theta0; 92 const CeedScalar thetaC = context->thetaC; 93 const CeedScalar P0 = context->P0; 94 const CeedScalar N = context->N; 95 const CeedScalar rc = context->rc; 96 const CeedScalar *center = context->center; 97 const CeedScalar *dc_axis = context->dc_axis; 98 NewtonianIdealGasContext gas = &context->newtonian_ctx; 99 const CeedScalar cp = gas->cp; 100 const CeedScalar cv = gas->cv; 101 const CeedScalar Rd = cp - cv; 102 const CeedScalar *g_vec = gas->g; 103 const CeedScalar g = -g_vec[2]; 104 105 // Setup 106 // -- Coordinates 107 const CeedScalar x = X[0]; 108 const CeedScalar y = X[1]; 109 const CeedScalar z = X[2]; 110 111 // -- Potential temperature, density current 112 CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 113 // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 114 for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr); 115 const CeedScalar r = sqrt(Dot3(rr, rr)); 116 const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.; 117 const CeedScalar theta = theta0 * exp(Square(N) * z / g) + delta_theta; 118 119 // -- Exner pressure, hydrostatic balance 120 const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N)); 121 122 // Initial Conditions 123 CeedScalar Y[5] = {0.}; 124 Y[0] = P0 * pow(Pi, cp / Rd); 125 Y[1] = 0.0; 126 Y[2] = 0.0; 127 Y[3] = 0.0; 128 Y[4] = Pi * theta; 129 130 return StateFromY(gas, Y, X); 131 } 132 133 // ***************************************************************************** 134 // This QFunction sets the initial conditions for density current 135 // ***************************************************************************** 136 CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 137 // Inputs 138 const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 139 140 // Outputs 141 CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 142 143 // Context 144 const DensityCurrentContext context = (DensityCurrentContext)ctx; 145 146 // Quadrature Point Loop 147 CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 148 const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 149 State s = Exact_DC(3, 0., x, 5, ctx); 150 CeedScalar q[5] = {0}; 151 switch (context->newtonian_ctx.state_var) { 152 case STATEVAR_CONSERVATIVE: 153 UnpackState_U(s.U, q); 154 break; 155 case STATEVAR_PRIMITIVE: 156 UnpackState_Y(s.Y, q); 157 break; 158 } 159 160 for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 161 162 } // End of Quadrature Point Loop 163 164 return 0; 165 } 166 167 // ***************************************************************************** 168 169 #endif // densitycurrent_h 170