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