1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, 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). 143a8779fbSJames Wright #include <ceed.h> 15d0cce58aSJeremy L Thompson #include <math.h> 162b916ea7SJeremy L Thompson 17cbe60e31SLeila Ghaffari #include "newtonian_state.h" 18d0cce58aSJeremy L Thompson #include "newtonian_types.h" 19704b8bbeSJames Wright #include "utils.h" 20a515125bSLeila Ghaffari 21cbe60e31SLeila Ghaffari typedef struct DensityCurrentContext_ *DensityCurrentContext; 22cbe60e31SLeila Ghaffari struct DensityCurrentContext_ { 23cbe60e31SLeila Ghaffari CeedScalar theta0; 24cbe60e31SLeila Ghaffari CeedScalar thetaC; 25cbe60e31SLeila Ghaffari CeedScalar P0; 26cbe60e31SLeila Ghaffari CeedScalar N; 27cbe60e31SLeila Ghaffari CeedScalar rc; 28cbe60e31SLeila Ghaffari CeedScalar center[3]; 29cbe60e31SLeila Ghaffari CeedScalar dc_axis[3]; 30cbe60e31SLeila Ghaffari struct NewtonianIdealGasContext_ newtonian_ctx; 31cbe60e31SLeila Ghaffari }; 32cbe60e31SLeila Ghaffari 33a515125bSLeila Ghaffari // ***************************************************************************** 34a515125bSLeila Ghaffari // This function sets the initial conditions and the boundary conditions 35a515125bSLeila Ghaffari // 3604e40bb6SJeremy L Thompson // These initial conditions are given in terms of potential temperature and Exner pressure and then converted to density and total energy. 37a515125bSLeila Ghaffari // Initial momentum density is zero. 38a515125bSLeila Ghaffari // 39a515125bSLeila Ghaffari // Initial Conditions: 40a515125bSLeila Ghaffari // Potential Temperature: 41a515125bSLeila Ghaffari // theta = thetabar + delta_theta 42a515125bSLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 43a515125bSLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 44a515125bSLeila Ghaffari // r > rc : 0 45a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 46a515125bSLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 47a515125bSLeila Ghaffari // Exner Pressure: 48a515125bSLeila Ghaffari // Pi = Pibar + deltaPi 49a515125bSLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 50a515125bSLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 51a515125bSLeila Ghaffari // Velocity/Momentum Density: 52a515125bSLeila Ghaffari // Ui = ui = 0 53a515125bSLeila Ghaffari // 54a515125bSLeila Ghaffari // Conversion to Conserved Variables: 55a515125bSLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 56a515125bSLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 57a515125bSLeila Ghaffari // 58a515125bSLeila Ghaffari // Boundary Conditions: 59a515125bSLeila Ghaffari // Mass Density: 60a515125bSLeila Ghaffari // 0.0 flux 61a515125bSLeila Ghaffari // Momentum Density: 62a515125bSLeila Ghaffari // 0.0 63a515125bSLeila Ghaffari // Energy Density: 64a515125bSLeila Ghaffari // 0.0 flux 65a515125bSLeila Ghaffari // 66a515125bSLeila Ghaffari // Constants: 67a515125bSLeila Ghaffari // theta0 , Potential temperature constant 68a515125bSLeila Ghaffari // thetaC , Potential temperature perturbation 69a515125bSLeila Ghaffari // P0 , Pressure at the surface 70a515125bSLeila Ghaffari // N , Brunt-Vaisala frequency 71a515125bSLeila Ghaffari // cv , Specific heat, constant volume 72a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 73a515125bSLeila Ghaffari // Rd = cp - cv, Specific heat difference 74a515125bSLeila Ghaffari // g , Gravity 75a515125bSLeila Ghaffari // rc , Characteristic radius of thermal bubble 76a515125bSLeila Ghaffari // center , Location of bubble center 77a515125bSLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 78a515125bSLeila Ghaffari // ***************************************************************************** 79a515125bSLeila Ghaffari 80a515125bSLeila Ghaffari // ***************************************************************************** 81a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 82a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for density current 83a515125bSLeila Ghaffari // ***************************************************************************** 842b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) { 85a515125bSLeila Ghaffari // Context 86cbe60e31SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 87a515125bSLeila Ghaffari const CeedScalar theta0 = context->theta0; 88a515125bSLeila Ghaffari const CeedScalar thetaC = context->thetaC; 89a515125bSLeila Ghaffari const CeedScalar P0 = context->P0; 90a515125bSLeila Ghaffari const CeedScalar N = context->N; 91a515125bSLeila Ghaffari const CeedScalar rc = context->rc; 92a515125bSLeila Ghaffari const CeedScalar *center = context->center; 93a515125bSLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 94cbe60e31SLeila Ghaffari NewtonianIdealGasContext gas = &context->newtonian_ctx; 95cbe60e31SLeila Ghaffari const CeedScalar cp = gas->cp; 96cbe60e31SLeila Ghaffari const CeedScalar cv = gas->cv; 97139613f2SLeila Ghaffari const CeedScalar Rd = cp - cv; 98cbe60e31SLeila Ghaffari const CeedScalar *g_vec = gas->g; 99bb8a0c61SJames Wright const CeedScalar g = -g_vec[2]; 100a515125bSLeila Ghaffari 101a515125bSLeila Ghaffari // Setup 102a515125bSLeila Ghaffari // -- Coordinates 103a515125bSLeila Ghaffari const CeedScalar x = X[0]; 104a515125bSLeila Ghaffari const CeedScalar y = X[1]; 105a515125bSLeila Ghaffari const CeedScalar z = X[2]; 106a515125bSLeila Ghaffari 107a515125bSLeila Ghaffari // -- Potential temperature, density current 108a515125bSLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 109a515125bSLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 1102b916ea7SJeremy L Thompson for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr); 111d1b9ef12SLeila Ghaffari const CeedScalar r = sqrt(Dot3(rr, rr)); 112a515125bSLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.; 113d1b9ef12SLeila Ghaffari const CeedScalar theta = theta0 * exp(Square(N) * z / g) + delta_theta; 114a515125bSLeila Ghaffari 115a515125bSLeila Ghaffari // -- Exner pressure, hydrostatic balance 1162b916ea7SJeremy L Thompson const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N)); 117a515125bSLeila Ghaffari 118a515125bSLeila Ghaffari // Initial Conditions 119cbe60e31SLeila Ghaffari CeedScalar Y[5] = {0.}; 120cbe60e31SLeila Ghaffari Y[0] = P0 * pow(Pi, cp / Rd); 121cbe60e31SLeila Ghaffari Y[1] = 0.0; 122cbe60e31SLeila Ghaffari Y[2] = 0.0; 123cbe60e31SLeila Ghaffari Y[3] = 0.0; 1242ca21690SKenneth E. Jansen Y[4] = Pi * theta; 125a515125bSLeila Ghaffari 126edcfef1bSKenneth E. Jansen return StateFromY(gas, Y); 127a515125bSLeila Ghaffari } 128a515125bSLeila Ghaffari 129a515125bSLeila Ghaffari // ***************************************************************************** 130a515125bSLeila Ghaffari // This QFunction sets the initial conditions for density current 131a515125bSLeila Ghaffari // ***************************************************************************** 1322b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) { 133a515125bSLeila Ghaffari const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 134a515125bSLeila Ghaffari CeedScalar(*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 135a515125bSLeila Ghaffari 136cbe60e31SLeila Ghaffari const DensityCurrentContext context = (DensityCurrentContext)ctx; 137cbe60e31SLeila Ghaffari 1383d65b166SJames Wright CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) { 139a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 140cbe60e31SLeila Ghaffari State s = Exact_DC(3, 0., x, 5, ctx); 141d1b9ef12SLeila Ghaffari CeedScalar q[5] = {0}; 1423636f6a4SJames Wright switch (context->newtonian_ctx.state_var) { 1433636f6a4SJames Wright case STATEVAR_CONSERVATIVE: 144d1b9ef12SLeila Ghaffari UnpackState_U(s.U, q); 1453636f6a4SJames Wright break; 1463636f6a4SJames Wright case STATEVAR_PRIMITIVE: 1473636f6a4SJames Wright UnpackState_Y(s.Y, q); 1483636f6a4SJames Wright break; 1493636f6a4SJames Wright } 150d1b9ef12SLeila Ghaffari 1512b916ea7SJeremy L Thompson for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j]; 152*b193fadcSJames Wright } 153a515125bSLeila Ghaffari return 0; 154a515125bSLeila Ghaffari } 155