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