177841947SLeila Ghaffari // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 277841947SLeila Ghaffari // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 377841947SLeila Ghaffari // reserved. See files LICENSE and NOTICE for details. 477841947SLeila Ghaffari // 577841947SLeila Ghaffari // This file is part of CEED, a collection of benchmarks, miniapps, software 677841947SLeila Ghaffari // libraries and APIs for efficient high-order finite element and spectral 777841947SLeila Ghaffari // element discretizations for exascale applications. For more information and 877841947SLeila Ghaffari // source code availability see http://github.com/ceed. 977841947SLeila Ghaffari // 1077841947SLeila Ghaffari // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 1177841947SLeila Ghaffari // a collaborative effort of two U.S. Department of Energy organizations (Office 1277841947SLeila Ghaffari // of Science and the National Nuclear Security Administration) responsible for 1377841947SLeila Ghaffari // the planning and preparation of a capable exascale ecosystem, including 1477841947SLeila Ghaffari // software, applications, hardware, advanced system engineering and early 1577841947SLeila Ghaffari // testbed platforms, in support of the nation's exascale computing imperative. 1677841947SLeila Ghaffari 1777841947SLeila Ghaffari /// @file 1877841947SLeila Ghaffari /// Density current initial condition and operator for Navier-Stokes example using PETSc 1977841947SLeila Ghaffari 2077841947SLeila Ghaffari // Model from: 2177841947SLeila Ghaffari // Semi-Implicit Formulations of the Navier-Stokes Equations: Application to 2277841947SLeila Ghaffari // Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010). 2377841947SLeila Ghaffari 2477841947SLeila Ghaffari #ifndef densitycurrent_h 2577841947SLeila Ghaffari #define densitycurrent_h 2677841947SLeila Ghaffari 2777841947SLeila Ghaffari #include <math.h> 28*88b783a1SJames Wright #include <ceed.h> 2977841947SLeila Ghaffari 3077841947SLeila Ghaffari #ifndef M_PI 3177841947SLeila Ghaffari #define M_PI 3.14159265358979323846 3277841947SLeila Ghaffari #endif 3377841947SLeila Ghaffari 3477841947SLeila Ghaffari #ifndef setup_context_struct 3577841947SLeila Ghaffari #define setup_context_struct 3677841947SLeila Ghaffari typedef struct SetupContext_ *SetupContext; 3777841947SLeila Ghaffari struct SetupContext_ { 3877841947SLeila Ghaffari CeedScalar theta0; 3977841947SLeila Ghaffari CeedScalar thetaC; 4077841947SLeila Ghaffari CeedScalar P0; 4177841947SLeila Ghaffari CeedScalar N; 4277841947SLeila Ghaffari CeedScalar cv; 4377841947SLeila Ghaffari CeedScalar cp; 4477841947SLeila Ghaffari CeedScalar g; 4577841947SLeila Ghaffari CeedScalar rc; 4677841947SLeila Ghaffari CeedScalar lx; 4777841947SLeila Ghaffari CeedScalar ly; 4877841947SLeila Ghaffari CeedScalar lz; 4977841947SLeila Ghaffari CeedScalar center[3]; 5077841947SLeila Ghaffari CeedScalar dc_axis[3]; 5177841947SLeila Ghaffari CeedScalar wind[3]; 5277841947SLeila Ghaffari CeedScalar time; 5377841947SLeila Ghaffari int wind_type; // See WindType: 0=ROTATION, 1=TRANSLATION 5477841947SLeila Ghaffari int bubble_type; // See BubbleType: 0=SPHERE, 1=CYLINDER 5577841947SLeila Ghaffari int bubble_continuity_type; // See BubbleContinuityType: 0=SMOOTH, 1=BACK_SHARP 2=THICK 5677841947SLeila Ghaffari }; 5777841947SLeila Ghaffari #endif 5877841947SLeila Ghaffari 5977841947SLeila Ghaffari // ***************************************************************************** 6077841947SLeila Ghaffari // This function sets the initial conditions and the boundary conditions 6177841947SLeila Ghaffari // 6277841947SLeila Ghaffari // These initial conditions are given in terms of potential temperature and 6377841947SLeila Ghaffari // Exner pressure and then converted to density and total energy. 6477841947SLeila Ghaffari // Initial momentum density is zero. 6577841947SLeila Ghaffari // 6677841947SLeila Ghaffari // Initial Conditions: 6777841947SLeila Ghaffari // Potential Temperature: 6877841947SLeila Ghaffari // theta = thetabar + delta_theta 6977841947SLeila Ghaffari // thetabar = theta0 exp( N**2 z / g ) 7077841947SLeila Ghaffari // delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2 7177841947SLeila Ghaffari // r > rc : 0 7277841947SLeila Ghaffari // r = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 ) 7377841947SLeila Ghaffari // with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble 7477841947SLeila Ghaffari // Exner Pressure: 7577841947SLeila Ghaffari // Pi = Pibar + deltaPi 7677841947SLeila Ghaffari // Pibar = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2) 7777841947SLeila Ghaffari // deltaPi = 0 (hydrostatic balance) 7877841947SLeila Ghaffari // Velocity/Momentum Density: 7977841947SLeila Ghaffari // Ui = ui = 0 8077841947SLeila Ghaffari // 8177841947SLeila Ghaffari // Conversion to Conserved Variables: 8277841947SLeila Ghaffari // rho = P0 Pi**(cv/Rd) / (Rd theta) 8377841947SLeila Ghaffari // E = rho (cv T + (u u)/2 + g z) 8477841947SLeila Ghaffari // 8577841947SLeila Ghaffari // Boundary Conditions: 8677841947SLeila Ghaffari // Mass Density: 8777841947SLeila Ghaffari // 0.0 flux 8877841947SLeila Ghaffari // Momentum Density: 8977841947SLeila Ghaffari // 0.0 9077841947SLeila Ghaffari // Energy Density: 9177841947SLeila Ghaffari // 0.0 flux 9277841947SLeila Ghaffari // 9377841947SLeila Ghaffari // Constants: 9477841947SLeila Ghaffari // theta0 , Potential temperature constant 9577841947SLeila Ghaffari // thetaC , Potential temperature perturbation 9677841947SLeila Ghaffari // P0 , Pressure at the surface 9777841947SLeila Ghaffari // N , Brunt-Vaisala frequency 9877841947SLeila Ghaffari // cv , Specific heat, constant volume 9977841947SLeila Ghaffari // cp , Specific heat, constant pressure 10077841947SLeila Ghaffari // Rd = cp - cv, Specific heat difference 10177841947SLeila Ghaffari // g , Gravity 10277841947SLeila Ghaffari // rc , Characteristic radius of thermal bubble 10377841947SLeila Ghaffari // center , Location of bubble center 10477841947SLeila Ghaffari // dc_axis , Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric 10577841947SLeila Ghaffari // ***************************************************************************** 10677841947SLeila Ghaffari 10777841947SLeila Ghaffari // ***************************************************************************** 10877841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution 10977841947SLeila Ghaffari // (currently not implemented) and IC formulation for density current 11077841947SLeila Ghaffari // ***************************************************************************** 11177841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_DC(CeedInt dim, CeedScalar time, 11277841947SLeila Ghaffari const CeedScalar X[], CeedInt Nf, CeedScalar q[], 11377841947SLeila Ghaffari void *ctx) { 11477841947SLeila Ghaffari // Context 11577841947SLeila Ghaffari const SetupContext context = (SetupContext)ctx; 11677841947SLeila Ghaffari const CeedScalar theta0 = context->theta0; 11777841947SLeila Ghaffari const CeedScalar thetaC = context->thetaC; 11877841947SLeila Ghaffari const CeedScalar P0 = context->P0; 11977841947SLeila Ghaffari const CeedScalar N = context->N; 12077841947SLeila Ghaffari const CeedScalar cv = context->cv; 12177841947SLeila Ghaffari const CeedScalar cp = context->cp; 12277841947SLeila Ghaffari const CeedScalar g = context->g; 12377841947SLeila Ghaffari const CeedScalar rc = context->rc; 12477841947SLeila Ghaffari const CeedScalar *center = context->center; 12577841947SLeila Ghaffari const CeedScalar *dc_axis = context->dc_axis; 126e6225c47SLeila Ghaffari const CeedScalar Rd = cp - cv; 12777841947SLeila Ghaffari 12877841947SLeila Ghaffari // Setup 12977841947SLeila Ghaffari // -- Coordinates 13077841947SLeila Ghaffari const CeedScalar x = X[0]; 13177841947SLeila Ghaffari const CeedScalar y = X[1]; 13277841947SLeila Ghaffari const CeedScalar z = X[2]; 13377841947SLeila Ghaffari 13477841947SLeila Ghaffari // -- Potential temperature, density current 13577841947SLeila Ghaffari CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]}; 13677841947SLeila Ghaffari // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector) 13777841947SLeila Ghaffari for (CeedInt i=0; i<3; i++) 13877841947SLeila Ghaffari rr[i] -= dc_axis[i] * 13977841947SLeila Ghaffari (dc_axis[0]*rr[0] + dc_axis[1]*rr[1] + dc_axis[2]*rr[2]); 14077841947SLeila Ghaffari const CeedScalar r = sqrt(rr[0]*rr[0] + rr[1]*rr[1] + rr[2]*rr[2]); 14177841947SLeila Ghaffari const CeedScalar delta_theta = r <= rc ? thetaC*(1. + cos(M_PI*r/rc))/2. : 0.; 14277841947SLeila Ghaffari const CeedScalar theta = theta0*exp(N*N*z/g) + delta_theta; 14377841947SLeila Ghaffari 14477841947SLeila Ghaffari // -- Exner pressure, hydrostatic balance 14577841947SLeila Ghaffari const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N); 14677841947SLeila Ghaffari // -- Density 14777841947SLeila Ghaffari 14877841947SLeila Ghaffari const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta); 14977841947SLeila Ghaffari 15077841947SLeila Ghaffari // Initial Conditions 15177841947SLeila Ghaffari q[0] = rho; 15277841947SLeila Ghaffari q[1] = 0.0; 15377841947SLeila Ghaffari q[2] = 0.0; 15477841947SLeila Ghaffari q[3] = 0.0; 15577841947SLeila Ghaffari q[4] = rho * (cv*theta*Pi + g*z); 15677841947SLeila Ghaffari 15777841947SLeila Ghaffari return 0; 15877841947SLeila Ghaffari } 15977841947SLeila Ghaffari 16077841947SLeila Ghaffari // ***************************************************************************** 16177841947SLeila Ghaffari // This QFunction sets the initial conditions for density current 16277841947SLeila Ghaffari // ***************************************************************************** 16377841947SLeila Ghaffari CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, 16477841947SLeila Ghaffari const CeedScalar *const *in, CeedScalar *const *out) { 16577841947SLeila Ghaffari // Inputs 16677841947SLeila Ghaffari const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0]; 16777841947SLeila Ghaffari 16877841947SLeila Ghaffari // Outputs 16977841947SLeila Ghaffari CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0]; 17077841947SLeila Ghaffari 17177841947SLeila Ghaffari CeedPragmaSIMD 17277841947SLeila Ghaffari // Quadrature Point Loop 17377841947SLeila Ghaffari for (CeedInt i=0; i<Q; i++) { 17477841947SLeila Ghaffari const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]}; 175e6225c47SLeila Ghaffari CeedScalar q[5] = {0.}; 17677841947SLeila Ghaffari 17777841947SLeila Ghaffari Exact_DC(3, 0., x, 5, q, ctx); 17877841947SLeila Ghaffari 17977841947SLeila Ghaffari for (CeedInt j=0; j<5; j++) 18077841947SLeila Ghaffari q0[j][i] = q[j]; 18177841947SLeila Ghaffari } // End of Quadrature Point Loop 18277841947SLeila Ghaffari 18377841947SLeila Ghaffari return 0; 18477841947SLeila Ghaffari } 18577841947SLeila Ghaffari 18677841947SLeila Ghaffari #endif // densitycurrent_h 187