1*ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2*ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3a515125bSLeila Ghaffari 4a515125bSLeila Ghaffari /// @file 5a515125bSLeila Ghaffari /// Density current initial condition and operator for Navier-Stokes example using PETSc 6a515125bSLeila Ghaffari 7a515125bSLeila Ghaffari // Model from: 8a515125bSLeila Ghaffari // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 9a515125bSLeila Ghaffari // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 103a8779fbSJames Wright #include <ceed.h> 11d0cce58aSJeremy L Thompson #include <math.h> 122b916ea7SJeremy L Thompson 13cbe60e31SLeila Ghaffari #include "newtonian_state.h" 14d0cce58aSJeremy L Thompson #include "newtonian_types.h" 15704b8bbeSJames Wright #include "utils.h" 16a515125bSLeila Ghaffari 17cbe60e31SLeila Ghaffari typedef struct DensityCurrentContext_ *DensityCurrentContext; 18cbe60e31SLeila Ghaffari struct DensityCurrentContext_ { 19cbe60e31SLeila Ghaffari CeedScalar theta0; 20cbe60e31SLeila Ghaffari CeedScalar thetaC; 21cbe60e31SLeila Ghaffari CeedScalar P0; 22cbe60e31SLeila Ghaffari CeedScalar N; 23cbe60e31SLeila Ghaffari CeedScalar rc; 24cbe60e31SLeila Ghaffari CeedScalar center[3]; 25cbe60e31SLeila Ghaffari CeedScalar dc_axis[3]; 26cbe60e31SLeila Ghaffari struct NewtonianIdealGasContext_ newtonian_ctx; 27cbe60e31SLeila Ghaffari }; 28cbe60e31SLeila Ghaffari 29a515125bSLeila Ghaffari // ***************************************************************************** 30a515125bSLeila Ghaffari // This function sets the initial conditions and the boundary conditions 31a515125bSLeila Ghaffari // 3204e40bb6SJeremy L Thompson // These initial conditions are given in terms of potential temperature and Exner pressure and then converted to density and total energy. 33a515125bSLeila Ghaffari // Initial momentum density is zero. 34a515125bSLeila Ghaffari // 35a515125bSLeila Ghaffari // Initial Conditions: 36a515125bSLeila Ghaffari // Potential Temperature: 37a515125bSLeila Ghaffari // theta = thetabar + delta_theta 38a515125bSLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 39a515125bSLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 40a515125bSLeila Ghaffari // r > rc : 0 41a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 42a515125bSLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 43a515125bSLeila Ghaffari // Exner Pressure: 44a515125bSLeila Ghaffari // Pi = Pibar + deltaPi 45a515125bSLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 46a515125bSLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 47a515125bSLeila Ghaffari // Velocity/Momentum Density: 48a515125bSLeila Ghaffari // Ui = ui = 0 49a515125bSLeila Ghaffari // 50a515125bSLeila Ghaffari // Conversion to Conserved Variables: 51a515125bSLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 52a515125bSLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 53a515125bSLeila Ghaffari // 54a515125bSLeila Ghaffari // Boundary Conditions: 55a515125bSLeila Ghaffari // Mass Density: 56a515125bSLeila Ghaffari // 0.0 flux 57a515125bSLeila Ghaffari // Momentum Density: 58a515125bSLeila Ghaffari // 0.0 59a515125bSLeila Ghaffari // Energy Density: 60a515125bSLeila Ghaffari // 0.0 flux 61a515125bSLeila Ghaffari // 62a515125bSLeila Ghaffari // Constants: 63a515125bSLeila Ghaffari // theta0 , Potential temperature constant 64a515125bSLeila Ghaffari // thetaC , Potential temperature perturbation 65a515125bSLeila Ghaffari // P0 , Pressure at the surface 66a515125bSLeila Ghaffari // N , Brunt-Vaisala frequency 67a515125bSLeila Ghaffari // cv , Specific heat, constant volume 68a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 69a515125bSLeila Ghaffari // Rd = cp - cv, Specific heat difference 70a515125bSLeila Ghaffari // g , Gravity 71a515125bSLeila Ghaffari // rc , Characteristic radius of thermal bubble 72a515125bSLeila Ghaffari // center , Location of bubble center 73a515125bSLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 74a515125bSLeila Ghaffari // ***************************************************************************** 75a515125bSLeila Ghaffari 76a515125bSLeila Ghaffari // ***************************************************************************** 77a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 78a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for density current 79a515125bSLeila Ghaffari // ***************************************************************************** 802b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) { 81a515125bSLeila Ghaffari // Context 82cbe60e31SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 83a515125bSLeila Ghaffari const CeedScalar theta0 = context->theta0; 84a515125bSLeila Ghaffari const CeedScalar thetaC = context->thetaC; 85a515125bSLeila Ghaffari const CeedScalar P0 = context->P0; 86a515125bSLeila Ghaffari const CeedScalar N = context->N; 87a515125bSLeila Ghaffari const CeedScalar rc = context->rc; 88a515125bSLeila Ghaffari const CeedScalar *center = context->center; 89a515125bSLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 90cbe60e31SLeila Ghaffari NewtonianIdealGasContext gas = &context->newtonian_ctx; 91cbe60e31SLeila Ghaffari const CeedScalar cp = gas->cp; 92cbe60e31SLeila Ghaffari const CeedScalar cv = gas->cv; 93139613f2SLeila Ghaffari const CeedScalar Rd = cp - cv; 94cbe60e31SLeila Ghaffari const CeedScalar *g_vec = gas->g; 95bb8a0c61SJames Wright const CeedScalar g = -g_vec[2]; 96a515125bSLeila Ghaffari 97a515125bSLeila Ghaffari // Setup 98a515125bSLeila Ghaffari // -- Coordinates 99a515125bSLeila Ghaffari const CeedScalar x = X[0]; 100a515125bSLeila Ghaffari const CeedScalar y = X[1]; 101a515125bSLeila Ghaffari const CeedScalar z = X[2]; 102a515125bSLeila Ghaffari 103a515125bSLeila Ghaffari // -- Potential temperature, density current 104a515125bSLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 105a515125bSLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 1062b916ea7SJeremy L Thompson for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr); 10764667825SJames Wright const CeedScalar r = Norm3(rr); 108a515125bSLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.; 109d1b9ef12SLeila Ghaffari const CeedScalar theta = theta0 * exp(Square(N) * z / g) + delta_theta; 110a515125bSLeila Ghaffari 111a515125bSLeila Ghaffari // -- Exner pressure, hydrostatic balance 1122b916ea7SJeremy L Thompson const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N)); 113a515125bSLeila Ghaffari 114a515125bSLeila Ghaffari // Initial Conditions 115cbe60e31SLeila Ghaffari CeedScalar Y[5] = {0.}; 116cbe60e31SLeila Ghaffari Y[0] = P0 * pow(Pi, cp / Rd); 117cbe60e31SLeila Ghaffari Y[1] = 0.0; 118cbe60e31SLeila Ghaffari Y[2] = 0.0; 119cbe60e31SLeila Ghaffari Y[3] = 0.0; 1202ca21690SKenneth E. Jansen Y[4] = Pi * theta; 121a515125bSLeila Ghaffari 122edcfef1bSKenneth E. Jansen return StateFromY(gas, Y); 123a515125bSLeila Ghaffari } 124a515125bSLeila Ghaffari 125a515125bSLeila Ghaffari // ***************************************************************************** 126a515125bSLeila Ghaffari // This QFunction sets the initial conditions for density current 127a515125bSLeila Ghaffari // ***************************************************************************** 1282b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 129a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 130a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 131a515125bSLeila Ghaffari 132cbe60e31SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 1339b103f75SJames Wright const NewtonianIdealGasContext gas = &context->newtonian_ctx; 134cbe60e31SLeila Ghaffari 1353d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 136a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 137cbe60e31SLeila Ghaffari State s = Exact_DC(3, 0., x, 5, ctx); 138a541e550SJames Wright CeedScalar q[5]; 1399b103f75SJames Wright StateToQ(gas, s, q, gas->state_var); 1402b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 141b193fadcSJames Wright } 142a515125bSLeila Ghaffari return 0; 143a515125bSLeila Ghaffari } 144