13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Density current initial condition and operator for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 1177841947SLeila Ghaffari // Model from: 1277841947SLeila Ghaffari // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 1377841947SLeila Ghaffari // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 1477841947SLeila Ghaffari 1577841947SLeila Ghaffari #ifndef densitycurrent_h 1677841947SLeila Ghaffari #define densitycurrent_h 1777841947SLeila Ghaffari 1888b783a1SJames Wright #include <ceed.h> 19c9c2c079SJeremy L Thompson #include <math.h> 202b730f8bSJeremy L Thompson 21dc805cc4SLeila Ghaffari #include "newtonian_state.h" 22c9c2c079SJeremy L Thompson #include "newtonian_types.h" 2313fa47b2SJames Wright #include "utils.h" 2477841947SLeila Ghaffari 25dc805cc4SLeila Ghaffari typedef struct DensityCurrentContext_ *DensityCurrentContext; 26dc805cc4SLeila Ghaffari struct DensityCurrentContext_ { 27dc805cc4SLeila Ghaffari CeedScalar theta0; 28dc805cc4SLeila Ghaffari CeedScalar thetaC; 29dc805cc4SLeila Ghaffari CeedScalar P0; 30dc805cc4SLeila Ghaffari CeedScalar N; 31dc805cc4SLeila Ghaffari CeedScalar rc; 32dc805cc4SLeila Ghaffari CeedScalar center[3]; 33dc805cc4SLeila Ghaffari CeedScalar dc_axis[3]; 34dc805cc4SLeila Ghaffari struct NewtonianIdealGasContext_ newtonian_ctx; 35dc805cc4SLeila Ghaffari }; 36dc805cc4SLeila Ghaffari 3777841947SLeila Ghaffari // ***************************************************************************** 3877841947SLeila Ghaffari // This function sets the initial conditions and the boundary conditions 3977841947SLeila Ghaffari // 4077841947SLeila Ghaffari // These initial conditions are given in terms of potential temperature and 4177841947SLeila Ghaffari // Exner pressure and then converted to density and total energy. 4277841947SLeila Ghaffari // Initial momentum density is zero. 4377841947SLeila Ghaffari // 4477841947SLeila Ghaffari // Initial Conditions: 4577841947SLeila Ghaffari // Potential Temperature: 4677841947SLeila Ghaffari // theta = thetabar + delta_theta 4777841947SLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 4877841947SLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 4977841947SLeila Ghaffari // r > rc : 0 5077841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 5177841947SLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 5277841947SLeila Ghaffari // Exner Pressure: 5377841947SLeila Ghaffari // Pi = Pibar + deltaPi 5477841947SLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 5577841947SLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 5677841947SLeila Ghaffari // Velocity/Momentum Density: 5777841947SLeila Ghaffari // Ui = ui = 0 5877841947SLeila Ghaffari // 5977841947SLeila Ghaffari // Conversion to Conserved Variables: 6077841947SLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 6177841947SLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 6277841947SLeila Ghaffari // 6377841947SLeila Ghaffari // Boundary Conditions: 6477841947SLeila Ghaffari // Mass Density: 6577841947SLeila Ghaffari // 0.0 flux 6677841947SLeila Ghaffari // Momentum Density: 6777841947SLeila Ghaffari // 0.0 6877841947SLeila Ghaffari // Energy Density: 6977841947SLeila Ghaffari // 0.0 flux 7077841947SLeila Ghaffari // 7177841947SLeila Ghaffari // Constants: 7277841947SLeila Ghaffari // theta0 , Potential temperature constant 7377841947SLeila Ghaffari // thetaC , Potential temperature perturbation 7477841947SLeila Ghaffari // P0 , Pressure at the surface 7577841947SLeila Ghaffari // N , Brunt-Vaisala frequency 7677841947SLeila Ghaffari // cv , Specific heat, constant volume 7777841947SLeila Ghaffari // cp , Specific heat, constant pressure 7877841947SLeila Ghaffari // Rd = cp - cv, Specific heat difference 7977841947SLeila Ghaffari // g , Gravity 8077841947SLeila Ghaffari // rc , Characteristic radius of thermal bubble 8177841947SLeila Ghaffari // center , Location of bubble center 8277841947SLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 8377841947SLeila Ghaffari // ***************************************************************************** 8477841947SLeila Ghaffari 8577841947SLeila Ghaffari // ***************************************************************************** 8677841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 8777841947SLeila Ghaffari // (currently not implemented) and IC formulation for density current 8877841947SLeila Ghaffari // ***************************************************************************** 892b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) { 9077841947SLeila Ghaffari // Context 91dc805cc4SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 9277841947SLeila Ghaffari const CeedScalar theta0 = context->theta0; 9377841947SLeila Ghaffari const CeedScalar thetaC = context->thetaC; 9477841947SLeila Ghaffari const CeedScalar P0 = context->P0; 9577841947SLeila Ghaffari const CeedScalar N = context->N; 9677841947SLeila Ghaffari const CeedScalar rc = context->rc; 9777841947SLeila Ghaffari const CeedScalar *center = context->center; 9877841947SLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 99dc805cc4SLeila Ghaffari NewtonianIdealGasContext gas = &context->newtonian_ctx; 100dc805cc4SLeila Ghaffari const CeedScalar cp = gas->cp; 101dc805cc4SLeila Ghaffari const CeedScalar cv = gas->cv; 102e6225c47SLeila Ghaffari const CeedScalar Rd = cp - cv; 103dc805cc4SLeila Ghaffari const CeedScalar *g_vec = gas->g; 10488626eedSJames Wright const CeedScalar g = -g_vec[2]; 10577841947SLeila Ghaffari 10677841947SLeila Ghaffari // Setup 10777841947SLeila Ghaffari // -- Coordinates 10877841947SLeila Ghaffari const CeedScalar x = X[0]; 10977841947SLeila Ghaffari const CeedScalar y = X[1]; 11077841947SLeila Ghaffari const CeedScalar z = X[2]; 11177841947SLeila Ghaffari 11277841947SLeila Ghaffari // -- Potential temperature, density current 11377841947SLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 11477841947SLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 1152b730f8bSJeremy L Thompson for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr); 1162b89d87eSLeila Ghaffari const CeedScalar r = sqrt(Dot3(rr, rr)); 11777841947SLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.; 1182b89d87eSLeila Ghaffari const CeedScalar theta = theta0 * exp(Square(N) * z / g) + delta_theta; 11977841947SLeila Ghaffari 12077841947SLeila Ghaffari // -- Exner pressure, hydrostatic balance 1212b730f8bSJeremy L Thompson const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N)); 12277841947SLeila Ghaffari 12377841947SLeila Ghaffari // Initial Conditions 124dc805cc4SLeila Ghaffari CeedScalar Y[5] = {0.}; 125dc805cc4SLeila Ghaffari Y[0] = P0 * pow(Pi, cp / Rd); 126dc805cc4SLeila Ghaffari Y[1] = 0.0; 127dc805cc4SLeila Ghaffari Y[2] = 0.0; 128dc805cc4SLeila Ghaffari Y[3] = 0.0; 129dc805cc4SLeila Ghaffari Y[4] = Pi * theta; 13077841947SLeila Ghaffari 131dc805cc4SLeila Ghaffari return StateFromY(gas, Y, X); 13277841947SLeila Ghaffari } 13377841947SLeila Ghaffari 13477841947SLeila Ghaffari // ***************************************************************************** 13577841947SLeila Ghaffari // This QFunction sets the initial conditions for density current 13677841947SLeila Ghaffari // ***************************************************************************** 1372b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 13877841947SLeila Ghaffari // Inputs 13977841947SLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 14077841947SLeila Ghaffari 14177841947SLeila Ghaffari // Outputs 14277841947SLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 14377841947SLeila Ghaffari 144dc805cc4SLeila Ghaffari // Context 145dc805cc4SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 146dc805cc4SLeila Ghaffari 14777841947SLeila Ghaffari // Quadrature Point Loop 148*46603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 14977841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 150dc805cc4SLeila Ghaffari State s = Exact_DC(3, 0., x, 5, ctx); 1512b89d87eSLeila Ghaffari CeedScalar q[5] = {0}; 15297baf651SJames Wright switch (context->newtonian_ctx.state_var) { 15397baf651SJames Wright case STATEVAR_CONSERVATIVE: 1542b89d87eSLeila Ghaffari UnpackState_U(s.U, q); 15597baf651SJames Wright break; 15697baf651SJames Wright case STATEVAR_PRIMITIVE: 15797baf651SJames Wright UnpackState_Y(s.Y, q); 15897baf651SJames Wright break; 15997baf651SJames Wright } 1602b89d87eSLeila Ghaffari 1612b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 1622b89d87eSLeila Ghaffari 16377841947SLeila Ghaffari } // End of Quadrature Point Loop 16477841947SLeila Ghaffari 16577841947SLeila Ghaffari return 0; 16677841947SLeila Ghaffari } 16777841947SLeila Ghaffari 168dc805cc4SLeila Ghaffari // ***************************************************************************** 169dc805cc4SLeila Ghaffari 17077841947SLeila Ghaffari #endif // densitycurrent_h 171