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