1a515125bSLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2a515125bSLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3a515125bSLeila Ghaffari // reserved. See files LICENSE and NOTICE for details. 4a515125bSLeila Ghaffari // 5a515125bSLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software 6a515125bSLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral 7a515125bSLeila Ghaffari // element discretizations for exascale applications. For more information and 8a515125bSLeila Ghaffari // source code availability see http://github.com/ceed. 9a515125bSLeila Ghaffari // 10a515125bSLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11a515125bSLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office 12a515125bSLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for 13a515125bSLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including 14a515125bSLeila Ghaffari // software, applications, hardware, advanced system engineering and early 15a515125bSLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative. 16a515125bSLeila Ghaffari 17a515125bSLeila Ghaffari /// @file 18a515125bSLeila Ghaffari /// Density current initial condition and operator for Navier-Stokes example using PETSc 19a515125bSLeila Ghaffari 20a515125bSLeila Ghaffari // Model from: 21a515125bSLeila Ghaffari // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 22a515125bSLeila Ghaffari // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 23a515125bSLeila Ghaffari 24a515125bSLeila Ghaffari #ifndef densitycurrent_h 25a515125bSLeila Ghaffari #define densitycurrent_h 26a515125bSLeila Ghaffari 27a515125bSLeila Ghaffari #include <math.h> 28*3a8779fbSJames Wright #include <ceed.h> 29a515125bSLeila Ghaffari 30a515125bSLeila Ghaffari #ifndef M_PI 31a515125bSLeila Ghaffari #define M_PI 3.14159265358979323846 32a515125bSLeila Ghaffari #endif 33a515125bSLeila Ghaffari 34a515125bSLeila Ghaffari #ifndef setup_context_struct 35a515125bSLeila Ghaffari #define setup_context_struct 36a515125bSLeila Ghaffari typedef struct SetupContext_ *SetupContext; 37a515125bSLeila Ghaffari struct SetupContext_ { 38a515125bSLeila Ghaffari CeedScalar theta0; 39a515125bSLeila Ghaffari CeedScalar thetaC; 40a515125bSLeila Ghaffari CeedScalar P0; 41a515125bSLeila Ghaffari CeedScalar N; 42a515125bSLeila Ghaffari CeedScalar cv; 43a515125bSLeila Ghaffari CeedScalar cp; 44a515125bSLeila Ghaffari CeedScalar g; 45a515125bSLeila Ghaffari CeedScalar rc; 46a515125bSLeila Ghaffari CeedScalar lx; 47a515125bSLeila Ghaffari CeedScalar ly; 48a515125bSLeila Ghaffari CeedScalar lz; 49a515125bSLeila Ghaffari CeedScalar center[3]; 50a515125bSLeila Ghaffari CeedScalar dc_axis[3]; 51a515125bSLeila Ghaffari CeedScalar wind[3]; 52a515125bSLeila Ghaffari CeedScalar time; 53a515125bSLeila Ghaffari int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 54a515125bSLeila Ghaffari int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 55a515125bSLeila Ghaffari int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 56a515125bSLeila Ghaffari }; 57a515125bSLeila Ghaffari #endif 58a515125bSLeila Ghaffari 59a515125bSLeila Ghaffari // ***************************************************************************** 60a515125bSLeila Ghaffari // This function sets the initial conditions and the boundary conditions 61a515125bSLeila Ghaffari // 62a515125bSLeila Ghaffari // These initial conditions are given in terms of potential temperature and 63a515125bSLeila Ghaffari // Exner pressure and then converted to density and total energy. 64a515125bSLeila Ghaffari // Initial momentum density is zero. 65a515125bSLeila Ghaffari // 66a515125bSLeila Ghaffari // Initial Conditions: 67a515125bSLeila Ghaffari // Potential Temperature: 68a515125bSLeila Ghaffari // theta = thetabar + delta_theta 69a515125bSLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 70a515125bSLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 71a515125bSLeila Ghaffari // r > rc : 0 72a515125bSLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 73a515125bSLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 74a515125bSLeila Ghaffari // Exner Pressure: 75a515125bSLeila Ghaffari // Pi = Pibar + deltaPi 76a515125bSLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 77a515125bSLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 78a515125bSLeila Ghaffari // Velocity/Momentum Density: 79a515125bSLeila Ghaffari // Ui = ui = 0 80a515125bSLeila Ghaffari // 81a515125bSLeila Ghaffari // Conversion to Conserved Variables: 82a515125bSLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 83a515125bSLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 84a515125bSLeila Ghaffari // 85a515125bSLeila Ghaffari // Boundary Conditions: 86a515125bSLeila Ghaffari // Mass Density: 87a515125bSLeila Ghaffari // 0.0 flux 88a515125bSLeila Ghaffari // Momentum Density: 89a515125bSLeila Ghaffari // 0.0 90a515125bSLeila Ghaffari // Energy Density: 91a515125bSLeila Ghaffari // 0.0 flux 92a515125bSLeila Ghaffari // 93a515125bSLeila Ghaffari // Constants: 94a515125bSLeila Ghaffari // theta0 , Potential temperature constant 95a515125bSLeila Ghaffari // thetaC , Potential temperature perturbation 96a515125bSLeila Ghaffari // P0 , Pressure at the surface 97a515125bSLeila Ghaffari // N , Brunt-Vaisala frequency 98a515125bSLeila Ghaffari // cv , Specific heat, constant volume 99a515125bSLeila Ghaffari // cp , Specific heat, constant pressure 100a515125bSLeila Ghaffari // Rd = cp - cv, Specific heat difference 101a515125bSLeila Ghaffari // g , Gravity 102a515125bSLeila Ghaffari // rc , Characteristic radius of thermal bubble 103a515125bSLeila Ghaffari // center , Location of bubble center 104a515125bSLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 105a515125bSLeila Ghaffari // ***************************************************************************** 106a515125bSLeila Ghaffari 107a515125bSLeila Ghaffari // ***************************************************************************** 108a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 109a515125bSLeila Ghaffari // (currently not implemented) and IC formulation for density current 110a515125bSLeila Ghaffari // ***************************************************************************** 111a515125bSLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_DC(CeedInt dim, CeedScalar time, 112a515125bSLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 113a515125bSLeila Ghaffari void *ctx) { 114a515125bSLeila Ghaffari // Context 115a515125bSLeila Ghaffari const SetupContext context = (SetupContext)ctx; 116a515125bSLeila Ghaffari const CeedScalar theta0 = context->theta0; 117a515125bSLeila Ghaffari const CeedScalar thetaC = context->thetaC; 118a515125bSLeila Ghaffari const CeedScalar P0 = context->P0; 119a515125bSLeila Ghaffari const CeedScalar N = context->N; 120a515125bSLeila Ghaffari const CeedScalar cv = context->cv; 121a515125bSLeila Ghaffari const CeedScalar cp = context->cp; 122a515125bSLeila Ghaffari const CeedScalar g = context->g; 123a515125bSLeila Ghaffari const CeedScalar rc = context->rc; 124a515125bSLeila Ghaffari const CeedScalar *center = context->center; 125a515125bSLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 126139613f2SLeila Ghaffari const CeedScalar Rd = cp - cv; 127a515125bSLeila Ghaffari 128a515125bSLeila Ghaffari // Setup 129a515125bSLeila Ghaffari // -- Coordinates 130a515125bSLeila Ghaffari const CeedScalar x = X[0]; 131a515125bSLeila Ghaffari const CeedScalar y = X[1]; 132a515125bSLeila Ghaffari const CeedScalar z = X[2]; 133a515125bSLeila Ghaffari 134a515125bSLeila Ghaffari // -- Potential temperature, density current 135a515125bSLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 136a515125bSLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 137a515125bSLeila Ghaffari for (CeedInt i=0; i<3; i++) 138a515125bSLeila Ghaffari rr[i] -= dc_axis[i] * 139a515125bSLeila Ghaffari (dc_axis[0]*rr[0] + dc_axis[1]*rr[1] + dc_axis[2]*rr[2]); 140a515125bSLeila Ghaffari const CeedScalar r = sqrt(rr[0]*rr[0] + rr[1]*rr[1] + rr[2]*rr[2]); 141a515125bSLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC*(1. + cos(M_PI*r/rc))/2. : 0.; 142a515125bSLeila Ghaffari const CeedScalar theta = theta0*exp(N*N*z/g) + delta_theta; 143a515125bSLeila Ghaffari 144a515125bSLeila Ghaffari // -- Exner pressure, hydrostatic balance 145a515125bSLeila Ghaffari const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N); 146a515125bSLeila Ghaffari // -- Density 147a515125bSLeila Ghaffari 148a515125bSLeila Ghaffari const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta); 149a515125bSLeila Ghaffari 150a515125bSLeila Ghaffari // Initial Conditions 151a515125bSLeila Ghaffari q[0] = rho; 152a515125bSLeila Ghaffari q[1] = 0.0; 153a515125bSLeila Ghaffari q[2] = 0.0; 154a515125bSLeila Ghaffari q[3] = 0.0; 155a515125bSLeila Ghaffari q[4] = rho * (cv*theta*Pi + g*z); 156a515125bSLeila Ghaffari 157a515125bSLeila Ghaffari return 0; 158a515125bSLeila Ghaffari } 159a515125bSLeila Ghaffari 160a515125bSLeila Ghaffari // ***************************************************************************** 161a515125bSLeila Ghaffari // This QFunction sets the initial conditions for density current 162a515125bSLeila Ghaffari // ***************************************************************************** 163a515125bSLeila Ghaffari CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, 164a515125bSLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 165a515125bSLeila Ghaffari // Inputs 166a515125bSLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 167a515125bSLeila Ghaffari 168a515125bSLeila Ghaffari // Outputs 169a515125bSLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 170a515125bSLeila Ghaffari 171a515125bSLeila Ghaffari CeedPragmaSIMD 172a515125bSLeila Ghaffari // Quadrature Point Loop 173a515125bSLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 174a515125bSLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 175139613f2SLeila Ghaffari CeedScalar q[5] = {0.}; 176a515125bSLeila Ghaffari 177a515125bSLeila Ghaffari Exact_DC(3, 0., x, 5, q, ctx); 178a515125bSLeila Ghaffari 179a515125bSLeila Ghaffari for (CeedInt j=0; j<5; j++) 180a515125bSLeila Ghaffari q0[j][i] = q[j]; 181a515125bSLeila Ghaffari } // End of Quadrature Point Loop 182a515125bSLeila Ghaffari 183a515125bSLeila Ghaffari return 0; 184a515125bSLeila Ghaffari } 185a515125bSLeila Ghaffari 186a515125bSLeila Ghaffari #endif // densitycurrent_h 187