xref: /libCEED/examples/fluids/qfunctions/densitycurrent.h (revision 7173c8dc81de0c321ecbfa65814b28460552c1f4)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Density current initial condition and operator for Navier-Stokes example using PETSc
10 
11 // Model from:
12 //   Semi-Implicit Formulations of the Navier-Stokes Equations: Application to
13 //   Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010).
14 
15 #ifndef densitycurrent_h
16 #define densitycurrent_h
17 
18 #include <math.h>
19 #include <ceed.h>
20 #include "newtonian_types.h"
21 
22 #ifndef M_PI
23 #define M_PI    3.14159265358979323846
24 #endif
25 
26 // *****************************************************************************
27 // This function sets the initial conditions and the boundary conditions
28 //
29 // These initial conditions are given in terms of potential temperature and
30 //   Exner pressure and then converted to density and total energy.
31 //   Initial momentum density is zero.
32 //
33 // Initial Conditions:
34 //   Potential Temperature:
35 //     theta = thetabar + delta_theta
36 //       thetabar   = theta0 exp( N**2 z / g )
37 //       delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2
38 //                     r > rc : 0
39 //         r        = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 )
40 //         with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble
41 //   Exner Pressure:
42 //     Pi = Pibar + deltaPi
43 //       Pibar      = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2)
44 //       deltaPi    = 0 (hydrostatic balance)
45 //   Velocity/Momentum Density:
46 //     Ui = ui = 0
47 //
48 // Conversion to Conserved Variables:
49 //   rho = P0 Pi**(cv/Rd) / (Rd theta)
50 //   E   = rho (cv T + (u u)/2 + g z)
51 //
52 //  Boundary Conditions:
53 //    Mass Density:
54 //      0.0 flux
55 //    Momentum Density:
56 //      0.0
57 //    Energy Density:
58 //      0.0 flux
59 //
60 // Constants:
61 //   theta0          ,  Potential temperature constant
62 //   thetaC          ,  Potential temperature perturbation
63 //   P0              ,  Pressure at the surface
64 //   N               ,  Brunt-Vaisala frequency
65 //   cv              ,  Specific heat, constant volume
66 //   cp              ,  Specific heat, constant pressure
67 //   Rd     = cp - cv,  Specific heat difference
68 //   g               ,  Gravity
69 //   rc              ,  Characteristic radius of thermal bubble
70 //   center          ,  Location of bubble center
71 //   dc_axis         ,  Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric
72 // *****************************************************************************
73 
74 // *****************************************************************************
75 // This helper function provides support for the exact, time-dependent solution
76 //   (currently not implemented) and IC formulation for density current
77 // *****************************************************************************
78 CEED_QFUNCTION_HELPER int Exact_DC(CeedInt dim, CeedScalar time,
79                                    const CeedScalar X[], CeedInt Nf, CeedScalar q[],
80                                    void *ctx) {
81   // Context
82   const SetupContext context = (SetupContext)ctx;
83   const CeedScalar theta0   = context->theta0;
84   const CeedScalar thetaC   = context->thetaC;
85   const CeedScalar P0       = context->P0;
86   const CeedScalar N        = context->N;
87   const CeedScalar cv       = context->cv;
88   const CeedScalar cp       = context->cp;
89   const CeedScalar *g_vec   = context->g;
90   const CeedScalar rc       = context->rc;
91   const CeedScalar *center  = context->center;
92   const CeedScalar *dc_axis = context->dc_axis;
93   const CeedScalar Rd       = cp - cv;
94   const CeedScalar g = -g_vec[2];
95 
96   // Setup
97   // -- Coordinates
98   const CeedScalar x = X[0];
99   const CeedScalar y = X[1];
100   const CeedScalar z = X[2];
101 
102   // -- Potential temperature, density current
103   CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]};
104   // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector)
105   for (CeedInt i=0; i<3; i++)
106     rr[i] -= dc_axis[i] *
107              (dc_axis[0]*rr[0] + dc_axis[1]*rr[1] + dc_axis[2]*rr[2]);
108   const CeedScalar r = sqrt(rr[0]*rr[0] + rr[1]*rr[1] + rr[2]*rr[2]);
109   const CeedScalar delta_theta = r <= rc ? thetaC*(1. + cos(M_PI*r/rc))/2. : 0.;
110   const CeedScalar theta = theta0*exp(N*N*z/g) + delta_theta;
111 
112   // -- Exner pressure, hydrostatic balance
113   const CeedScalar Pi = 1. + g*g*(exp(-N*N*z/g) - 1.) / (cp*theta0*N*N);
114   // -- Density
115 
116   const CeedScalar rho = P0 * pow(Pi, cv/Rd) / (Rd*theta);
117 
118   // Initial Conditions
119   q[0] = rho;
120   q[1] = 0.0;
121   q[2] = 0.0;
122   q[3] = 0.0;
123   q[4] = rho * (cv*theta*Pi + g*z);
124 
125   return 0;
126 }
127 
128 // *****************************************************************************
129 // This QFunction sets the initial conditions for density current
130 // *****************************************************************************
131 CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q,
132                       const CeedScalar *const *in, CeedScalar *const *out) {
133   // Inputs
134   const CeedScalar (*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
135 
136   // Outputs
137   CeedScalar (*q0)[CEED_Q_VLA] = (CeedScalar(*)[CEED_Q_VLA])out[0];
138 
139   CeedPragmaSIMD
140   // Quadrature Point Loop
141   for (CeedInt i=0; i<Q; i++) {
142     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
143     CeedScalar q[5] = {0.};
144 
145     Exact_DC(3, 0., x, 5, q, ctx);
146 
147     for (CeedInt j=0; j<5; j++)
148       q0[j][i] = q[j];
149   } // End of Quadrature Point Loop
150 
151   return 0;
152 }
153 
154 #endif // densitycurrent_h
155