15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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). 1488b783a1SJames Wright #include <ceed.h> 15c9c2c079SJeremy L Thompson #include <math.h> 162b730f8bSJeremy L Thompson 17dc805cc4SLeila Ghaffari #include "newtonian_state.h" 18c9c2c079SJeremy L Thompson #include "newtonian_types.h" 1913fa47b2SJames Wright #include "utils.h" 2077841947SLeila Ghaffari 21dc805cc4SLeila Ghaffari typedef struct DensityCurrentContext_ *DensityCurrentContext; 22dc805cc4SLeila Ghaffari struct DensityCurrentContext_ { 23dc805cc4SLeila Ghaffari CeedScalar theta0; 24dc805cc4SLeila Ghaffari CeedScalar thetaC; 25dc805cc4SLeila Ghaffari CeedScalar P0; 26dc805cc4SLeila Ghaffari CeedScalar N; 27dc805cc4SLeila Ghaffari CeedScalar rc; 28dc805cc4SLeila Ghaffari CeedScalar center[3]; 29dc805cc4SLeila Ghaffari CeedScalar dc_axis[3]; 30dc805cc4SLeila Ghaffari struct NewtonianIdealGasContext_ newtonian_ctx; 31dc805cc4SLeila Ghaffari }; 32dc805cc4SLeila Ghaffari 3377841947SLeila Ghaffari // ***************************************************************************** 3477841947SLeila Ghaffari // This function sets the initial conditions and the boundary conditions 3577841947SLeila Ghaffari // 36ea61e9acSJeremy L Thompson // These initial conditions are given in terms of potential temperature and Exner pressure and then converted to density and total energy. 3777841947SLeila Ghaffari // Initial momentum density is zero. 3877841947SLeila Ghaffari // 3977841947SLeila Ghaffari // Initial Conditions: 4077841947SLeila Ghaffari // Potential Temperature: 4177841947SLeila Ghaffari // theta = thetabar + delta_theta 4277841947SLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 4377841947SLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 4477841947SLeila Ghaffari // r > rc : 0 4577841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 4677841947SLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 4777841947SLeila Ghaffari // Exner Pressure: 4877841947SLeila Ghaffari // Pi = Pibar + deltaPi 4977841947SLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 5077841947SLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 5177841947SLeila Ghaffari // Velocity/Momentum Density: 5277841947SLeila Ghaffari // Ui = ui = 0 5377841947SLeila Ghaffari // 5477841947SLeila Ghaffari // Conversion to Conserved Variables: 5577841947SLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 5677841947SLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 5777841947SLeila Ghaffari // 5877841947SLeila Ghaffari // Boundary Conditions: 5977841947SLeila Ghaffari // Mass Density: 6077841947SLeila Ghaffari // 0.0 flux 6177841947SLeila Ghaffari // Momentum Density: 6277841947SLeila Ghaffari // 0.0 6377841947SLeila Ghaffari // Energy Density: 6477841947SLeila Ghaffari // 0.0 flux 6577841947SLeila Ghaffari // 6677841947SLeila Ghaffari // Constants: 6777841947SLeila Ghaffari // theta0 , Potential temperature constant 6877841947SLeila Ghaffari // thetaC , Potential temperature perturbation 6977841947SLeila Ghaffari // P0 , Pressure at the surface 7077841947SLeila Ghaffari // N , Brunt-Vaisala frequency 7177841947SLeila Ghaffari // cv , Specific heat, constant volume 7277841947SLeila Ghaffari // cp , Specific heat, constant pressure 7377841947SLeila Ghaffari // Rd = cp - cv, Specific heat difference 7477841947SLeila Ghaffari // g , Gravity 7577841947SLeila Ghaffari // rc , Characteristic radius of thermal bubble 7677841947SLeila Ghaffari // center , Location of bubble center 7777841947SLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 7877841947SLeila Ghaffari // ***************************************************************************** 7977841947SLeila Ghaffari 8077841947SLeila Ghaffari // ***************************************************************************** 8177841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 8277841947SLeila Ghaffari // (currently not implemented) and IC formulation for density current 8377841947SLeila Ghaffari // ***************************************************************************** 842b730f8bSJeremy L Thompson CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) { 8577841947SLeila Ghaffari // Context 86dc805cc4SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 8777841947SLeila Ghaffari const CeedScalar theta0 = context->theta0; 8877841947SLeila Ghaffari const CeedScalar thetaC = context->thetaC; 8977841947SLeila Ghaffari const CeedScalar P0 = context->P0; 9077841947SLeila Ghaffari const CeedScalar N = context->N; 9177841947SLeila Ghaffari const CeedScalar rc = context->rc; 9277841947SLeila Ghaffari const CeedScalar *center = context->center; 9377841947SLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 94dc805cc4SLeila Ghaffari NewtonianIdealGasContext gas = &context->newtonian_ctx; 95dc805cc4SLeila Ghaffari const CeedScalar cp = gas->cp; 96dc805cc4SLeila Ghaffari const CeedScalar cv = gas->cv; 97e6225c47SLeila Ghaffari const CeedScalar Rd = cp - cv; 98dc805cc4SLeila Ghaffari const CeedScalar *g_vec = gas->g; 9988626eedSJames Wright const CeedScalar g = -g_vec[2]; 10077841947SLeila Ghaffari 10177841947SLeila Ghaffari // Setup 10277841947SLeila Ghaffari // -- Coordinates 10377841947SLeila Ghaffari const CeedScalar x = X[0]; 10477841947SLeila Ghaffari const CeedScalar y = X[1]; 10577841947SLeila Ghaffari const CeedScalar z = X[2]; 10677841947SLeila Ghaffari 10777841947SLeila Ghaffari // -- Potential temperature, density current 10877841947SLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 10977841947SLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 1102b730f8bSJeremy L Thompson for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr); 1112b89d87eSLeila Ghaffari const CeedScalar r = sqrt(Dot3(rr, rr)); 11277841947SLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.; 1132b89d87eSLeila Ghaffari const CeedScalar theta = theta0 * exp(Square(N) * z / g) + delta_theta; 11477841947SLeila Ghaffari 11577841947SLeila Ghaffari // -- Exner pressure, hydrostatic balance 1162b730f8bSJeremy L Thompson const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N)); 11777841947SLeila Ghaffari 11877841947SLeila Ghaffari // Initial Conditions 119dc805cc4SLeila Ghaffari CeedScalar Y[5] = {0.}; 120dc805cc4SLeila Ghaffari Y[0] = P0 * pow(Pi, cp / Rd); 121dc805cc4SLeila Ghaffari Y[1] = 0.0; 122dc805cc4SLeila Ghaffari Y[2] = 0.0; 123dc805cc4SLeila Ghaffari Y[3] = 0.0; 1240d674129SKenneth E. Jansen Y[4] = Pi * theta; 12577841947SLeila Ghaffari 1263bd61617SKenneth E. Jansen return StateFromY(gas, Y); 12777841947SLeila Ghaffari } 12877841947SLeila Ghaffari 12977841947SLeila Ghaffari // ***************************************************************************** 13077841947SLeila Ghaffari // This QFunction sets the initial conditions for density current 13177841947SLeila Ghaffari // ***************************************************************************** 1322b730f8bSJeremy L Thompson CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 13377841947SLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 13477841947SLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 13577841947SLeila Ghaffari 136dc805cc4SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 137*a2d72b6fSJames Wright const NewtonianIdealGasContext gas = &context->newtonian_ctx; 138dc805cc4SLeila Ghaffari 13946603fc5SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 14077841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 141dc805cc4SLeila Ghaffari State s = Exact_DC(3, 0., x, 5, ctx); 1422b89d87eSLeila Ghaffari CeedScalar q[5] = {0}; 143*a2d72b6fSJames Wright StateToQ(gas, s, q, gas->state_var); 1442b730f8bSJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 145f0b01153SJames Wright } 14677841947SLeila Ghaffari return 0; 14777841947SLeila Ghaffari } 148