xref: /libCEED/examples/fluids/qfunctions/densitycurrent.h (revision a0add3c91c6ddcd77c0d376840911ac920dd4230)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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).
1477841947SLeila Ghaffari 
1577841947SLeila Ghaffari #ifndef densitycurrent_h
1677841947SLeila Ghaffari #define densitycurrent_h
1777841947SLeila Ghaffari 
1877841947SLeila Ghaffari #include <math.h>
1988b783a1SJames Wright #include <ceed.h>
20*a0add3c9SJed Brown #include "newtonian_types.h"
2177841947SLeila Ghaffari 
2277841947SLeila Ghaffari #ifndef M_PI
2377841947SLeila Ghaffari #define M_PI    3.14159265358979323846
2477841947SLeila Ghaffari #endif
2577841947SLeila Ghaffari 
2677841947SLeila Ghaffari // *****************************************************************************
2777841947SLeila Ghaffari // This function sets the initial conditions and the boundary conditions
2877841947SLeila Ghaffari //
2977841947SLeila Ghaffari // These initial conditions are given in terms of potential temperature and
3077841947SLeila Ghaffari //   Exner pressure and then converted to density and total energy.
3177841947SLeila Ghaffari //   Initial momentum density is zero.
3277841947SLeila Ghaffari //
3377841947SLeila Ghaffari // Initial Conditions:
3477841947SLeila Ghaffari //   Potential Temperature:
3577841947SLeila Ghaffari //     theta = thetabar + delta_theta
3677841947SLeila Ghaffari //       thetabar   = theta0 exp( N**2 z / g )
3777841947SLeila Ghaffari //       delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2
3877841947SLeila Ghaffari //                     r > rc : 0
3977841947SLeila Ghaffari //         r        = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 )
4077841947SLeila Ghaffari //         with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble
4177841947SLeila Ghaffari //   Exner Pressure:
4277841947SLeila Ghaffari //     Pi = Pibar + deltaPi
4377841947SLeila Ghaffari //       Pibar      = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2)
4477841947SLeila Ghaffari //       deltaPi    = 0 (hydrostatic balance)
4577841947SLeila Ghaffari //   Velocity/Momentum Density:
4677841947SLeila Ghaffari //     Ui = ui = 0
4777841947SLeila Ghaffari //
4877841947SLeila Ghaffari // Conversion to Conserved Variables:
4977841947SLeila Ghaffari //   rho = P0 Pi**(cv/Rd) / (Rd theta)
5077841947SLeila Ghaffari //   E   = rho (cv T + (u u)/2 + g z)
5177841947SLeila Ghaffari //
5277841947SLeila Ghaffari //  Boundary Conditions:
5377841947SLeila Ghaffari //    Mass Density:
5477841947SLeila Ghaffari //      0.0 flux
5577841947SLeila Ghaffari //    Momentum Density:
5677841947SLeila Ghaffari //      0.0
5777841947SLeila Ghaffari //    Energy Density:
5877841947SLeila Ghaffari //      0.0 flux
5977841947SLeila Ghaffari //
6077841947SLeila Ghaffari // Constants:
6177841947SLeila Ghaffari //   theta0          ,  Potential temperature constant
6277841947SLeila Ghaffari //   thetaC          ,  Potential temperature perturbation
6377841947SLeila Ghaffari //   P0              ,  Pressure at the surface
6477841947SLeila Ghaffari //   N               ,  Brunt-Vaisala frequency
6577841947SLeila Ghaffari //   cv              ,  Specific heat, constant volume
6677841947SLeila Ghaffari //   cp              ,  Specific heat, constant pressure
6777841947SLeila Ghaffari //   Rd     = cp - cv,  Specific heat difference
6877841947SLeila Ghaffari //   g               ,  Gravity
6977841947SLeila Ghaffari //   rc              ,  Characteristic radius of thermal bubble
7077841947SLeila Ghaffari //   center          ,  Location of bubble center
7177841947SLeila Ghaffari //   dc_axis         ,  Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric
7277841947SLeila Ghaffari // *****************************************************************************
7377841947SLeila Ghaffari 
7477841947SLeila Ghaffari // *****************************************************************************
7577841947SLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
7677841947SLeila Ghaffari //   (currently not implemented) and IC formulation for density current
7777841947SLeila Ghaffari // *****************************************************************************
7877841947SLeila Ghaffari CEED_QFUNCTION_HELPER int Exact_DC(CeedInt dim, CeedScalar time,
7977841947SLeila Ghaffari                                    const CeedScalar X[], CeedInt Nf, CeedScalar q[],
8077841947SLeila Ghaffari                                    void *ctx) {
8177841947SLeila Ghaffari   // Context
8277841947SLeila Ghaffari   const SetupContext context = (SetupContext)ctx;
8377841947SLeila Ghaffari   const CeedScalar theta0   = context->theta0;
8477841947SLeila Ghaffari   const CeedScalar thetaC   = context->thetaC;
8577841947SLeila Ghaffari   const CeedScalar P0       = context->P0;
8677841947SLeila Ghaffari   const CeedScalar N        = context->N;
8777841947SLeila Ghaffari   const CeedScalar cv       = context->cv;
8877841947SLeila Ghaffari   const CeedScalar cp       = context->cp;
8988626eedSJames Wright   const CeedScalar *g_vec   = context->g;
9077841947SLeila Ghaffari   const CeedScalar rc       = context->rc;
9177841947SLeila Ghaffari   const CeedScalar *center  = context->center;
9277841947SLeila Ghaffari   const CeedScalar *dc_axis = context->dc_axis;
93e6225c47SLeila Ghaffari   const CeedScalar Rd       = cp - cv;
9488626eedSJames Wright   const CeedScalar g = -g_vec[2];
9577841947SLeila Ghaffari 
9677841947SLeila Ghaffari   // Setup
9777841947SLeila Ghaffari   // -- Coordinates
9877841947SLeila Ghaffari   const CeedScalar x = X[0];
9977841947SLeila Ghaffari   const CeedScalar y = X[1];
10077841947SLeila Ghaffari   const CeedScalar z = X[2];
10177841947SLeila Ghaffari 
10277841947SLeila Ghaffari   // -- Potential temperature, density current
10377841947SLeila Ghaffari   CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]};
10477841947SLeila Ghaffari   // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector)
10577841947SLeila Ghaffari   for (CeedInt i=0; i<3; i++)
10677841947SLeila Ghaffari     rr[i] -= dc_axis[i] *
10777841947SLeila Ghaffari              (dc_axis[0]*rr[0] + dc_axis[1]*rr[1] + dc_axis[2]*rr[2]);
10877841947SLeila Ghaffari   const CeedScalar r = sqrt(rr[0]*rr[0] + rr[1]*rr[1] + rr[2]*rr[2]);
10977841947SLeila Ghaffari   const CeedScalar delta_theta = r <= rc ? thetaC*(1. + cos(M_PI*r/rc))/2. : 0.;
11077841947SLeila Ghaffari   const CeedScalar theta = theta0*exp(N*N*z/g) + delta_theta;
11177841947SLeila Ghaffari 
11277841947SLeila Ghaffari   // -- Exner pressure, hydrostatic balance
11377841947SLeila Ghaffari   const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N);
11477841947SLeila Ghaffari   // -- Density
11577841947SLeila Ghaffari 
11677841947SLeila Ghaffari   const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta);
11777841947SLeila Ghaffari 
11877841947SLeila Ghaffari   // Initial Conditions
11977841947SLeila Ghaffari   q[0] = rho;
12077841947SLeila Ghaffari   q[1] = 0.0;
12177841947SLeila Ghaffari   q[2] = 0.0;
12277841947SLeila Ghaffari   q[3] = 0.0;
12377841947SLeila Ghaffari   q[4] = rho * (cv*theta*Pi + g*z);
12477841947SLeila Ghaffari 
12577841947SLeila Ghaffari   return 0;
12677841947SLeila Ghaffari }
12777841947SLeila Ghaffari 
12877841947SLeila Ghaffari // *****************************************************************************
12977841947SLeila Ghaffari // This QFunction sets the initial conditions for density current
13077841947SLeila Ghaffari // *****************************************************************************
13177841947SLeila Ghaffari CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q,
13277841947SLeila Ghaffari                       const CeedScalar *const *in, CeedScalar *const *out) {
13377841947SLeila Ghaffari   // Inputs
13477841947SLeila Ghaffari   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
13577841947SLeila Ghaffari 
13677841947SLeila Ghaffari   // Outputs
13777841947SLeila Ghaffari   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
13877841947SLeila Ghaffari 
13977841947SLeila Ghaffari   CeedPragmaSIMD
14077841947SLeila Ghaffari   // Quadrature Point Loop
14177841947SLeila Ghaffari   for (CeedInt i=0; i<Q; i++) {
14277841947SLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
143e6225c47SLeila Ghaffari     CeedScalar q[5] = {0.};
14477841947SLeila Ghaffari 
14577841947SLeila Ghaffari     Exact_DC(3, 0., x, 5, q, ctx);
14677841947SLeila Ghaffari 
14777841947SLeila Ghaffari     for (CeedInt j=0; j<5; j++)
14877841947SLeila Ghaffari       q0[j][i] = q[j];
14977841947SLeila Ghaffari   } // End of Quadrature Point Loop
15077841947SLeila Ghaffari 
15177841947SLeila Ghaffari   return 0;
15277841947SLeila Ghaffari }
15377841947SLeila Ghaffari 
15477841947SLeila Ghaffari #endif // densitycurrent_h
155