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