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