xref: /honee/qfunctions/densitycurrent.h (revision 475f0cac5d40259768f4556cf888e8f2448554cb)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5ea615d4cSJames Wright /// Density current initial condition and operator for HONEE
6a515125bSLeila Ghaffari 
7a515125bSLeila Ghaffari // Model from:
8a515125bSLeila Ghaffari //   Semi-Implicit Formulations of the Navier-Stokes Equations: Application to
9a515125bSLeila Ghaffari //   Nonhydrostatic Atmospheric Modeling, Giraldo, Restelli, and Lauter (2010).
103e17a7a1SJames Wright #include <ceed/types.h>
112b916ea7SJeremy L Thompson 
12cbe60e31SLeila Ghaffari #include "newtonian_state.h"
13d0cce58aSJeremy L Thompson #include "newtonian_types.h"
14704b8bbeSJames Wright #include "utils.h"
15a515125bSLeila Ghaffari 
16cbe60e31SLeila Ghaffari typedef struct DensityCurrentContext_ *DensityCurrentContext;
17cbe60e31SLeila Ghaffari struct DensityCurrentContext_ {
18cbe60e31SLeila Ghaffari   CeedScalar                       theta0;
19cbe60e31SLeila Ghaffari   CeedScalar                       thetaC;
20cbe60e31SLeila Ghaffari   CeedScalar                       P0;
21cbe60e31SLeila Ghaffari   CeedScalar                       N;
22cbe60e31SLeila Ghaffari   CeedScalar                       rc;
23cbe60e31SLeila Ghaffari   CeedScalar                       center[3];
24cbe60e31SLeila Ghaffari   CeedScalar                       dc_axis[3];
25*cde3d787SJames Wright   struct NewtonianIdealGasContext_ newt_ctx;
26cbe60e31SLeila Ghaffari };
27cbe60e31SLeila Ghaffari 
28a515125bSLeila Ghaffari // *****************************************************************************
29a515125bSLeila Ghaffari // This function sets the initial conditions and the boundary conditions
30a515125bSLeila Ghaffari //
3104e40bb6SJeremy L Thompson // These initial conditions are given in terms of potential temperature and Exner pressure and then converted to density and total energy.
32a515125bSLeila Ghaffari //   Initial momentum density is zero.
33a515125bSLeila Ghaffari //
34a515125bSLeila Ghaffari // Initial Conditions:
35a515125bSLeila Ghaffari //   Potential Temperature:
36a515125bSLeila Ghaffari //     theta = thetabar + delta_theta
37a515125bSLeila Ghaffari //       thetabar   = theta0 exp( N**2 z / g )
38a515125bSLeila Ghaffari //       delta_theta = r <= rc : thetaC(1 + cos(pi r/rc)) / 2
39a515125bSLeila Ghaffari //                     r > rc : 0
40a515125bSLeila Ghaffari //         r        = sqrt( (x - xc)**2 + (y - yc)**2 + (z - zc)**2 )
41a515125bSLeila Ghaffari //         with (xc,yc,zc) center of domain, rc characteristic radius of thermal bubble
42a515125bSLeila Ghaffari //   Exner Pressure:
43a515125bSLeila Ghaffari //     Pi = Pibar + deltaPi
44a515125bSLeila Ghaffari //       Pibar      = 1. + g**2 (exp( - N**2 z / g ) - 1) / (cp theta0 N**2)
45a515125bSLeila Ghaffari //       deltaPi    = 0 (hydrostatic balance)
46a515125bSLeila Ghaffari //   Velocity/Momentum Density:
47a515125bSLeila Ghaffari //     Ui = ui = 0
48a515125bSLeila Ghaffari //
49a515125bSLeila Ghaffari // Conversion to Conserved Variables:
50a515125bSLeila Ghaffari //   rho = P0 Pi**(cv/Rd) / (Rd theta)
51a515125bSLeila Ghaffari //   E   = rho (cv T + (u u)/2 + g z)
52a515125bSLeila Ghaffari //
53a515125bSLeila Ghaffari //  Boundary Conditions:
54a515125bSLeila Ghaffari //    Mass Density:
55a515125bSLeila Ghaffari //      0.0 flux
56a515125bSLeila Ghaffari //    Momentum Density:
57a515125bSLeila Ghaffari //      0.0
58a515125bSLeila Ghaffari //    Energy Density:
59a515125bSLeila Ghaffari //      0.0 flux
60a515125bSLeila Ghaffari //
61a515125bSLeila Ghaffari // Constants:
62a515125bSLeila Ghaffari //   theta0          ,  Potential temperature constant
63a515125bSLeila Ghaffari //   thetaC          ,  Potential temperature perturbation
64a515125bSLeila Ghaffari //   P0              ,  Pressure at the surface
65a515125bSLeila Ghaffari //   N               ,  Brunt-Vaisala frequency
66a515125bSLeila Ghaffari //   cv              ,  Specific heat, constant volume
67a515125bSLeila Ghaffari //   cp              ,  Specific heat, constant pressure
68a515125bSLeila Ghaffari //   Rd     = cp - cv,  Specific heat difference
69a515125bSLeila Ghaffari //   g               ,  Gravity
70a515125bSLeila Ghaffari //   rc              ,  Characteristic radius of thermal bubble
71a515125bSLeila Ghaffari //   center          ,  Location of bubble center
72a515125bSLeila Ghaffari //   dc_axis         ,  Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric
73a515125bSLeila Ghaffari // *****************************************************************************
74a515125bSLeila Ghaffari 
75a515125bSLeila Ghaffari // *****************************************************************************
76a515125bSLeila Ghaffari // This helper function provides support for the exact, time-dependent solution
77a515125bSLeila Ghaffari //   (currently not implemented) and IC formulation for density current
78a515125bSLeila Ghaffari // *****************************************************************************
Exact_DC(CeedInt dim,CeedScalar time,const CeedScalar X[],CeedInt Nf,void * ctx)792b916ea7SJeremy L Thompson CEED_QFUNCTION_HELPER State Exact_DC(CeedInt dim, CeedScalar time, const CeedScalar X[], CeedInt Nf, void *ctx) {
80a515125bSLeila Ghaffari   // Context
81cbe60e31SLeila Ghaffari   const DensityCurrentContext context = (DensityCurrentContext)ctx;
82a515125bSLeila Ghaffari   const CeedScalar            theta0  = context->theta0;
83a515125bSLeila Ghaffari   const CeedScalar            thetaC  = context->thetaC;
84a515125bSLeila Ghaffari   const CeedScalar            P0      = context->P0;
85a515125bSLeila Ghaffari   const CeedScalar            N       = context->N;
86a515125bSLeila Ghaffari   const CeedScalar            rc      = context->rc;
87a515125bSLeila Ghaffari   const CeedScalar           *center  = context->center;
88a515125bSLeila Ghaffari   const CeedScalar           *dc_axis = context->dc_axis;
89*cde3d787SJames Wright   NewtonianIGProperties       gas     = context->newt_ctx.gas;
90*cde3d787SJames Wright   const CeedScalar            cp      = gas.cp;
91*cde3d787SJames Wright   const CeedScalar            Rd      = GasConstant(gas);
92*cde3d787SJames Wright   const CeedScalar           *g_vec   = context->newt_ctx.g;
93bb8a0c61SJames Wright   const CeedScalar            g       = -g_vec[2];
94a515125bSLeila Ghaffari 
95a515125bSLeila Ghaffari   // Setup
96a515125bSLeila Ghaffari   // -- Coordinates
97a515125bSLeila Ghaffari   const CeedScalar x = X[0];
98a515125bSLeila Ghaffari   const CeedScalar y = X[1];
99a515125bSLeila Ghaffari   const CeedScalar z = X[2];
100a515125bSLeila Ghaffari 
101a515125bSLeila Ghaffari   // -- Potential temperature, density current
102a515125bSLeila Ghaffari   CeedScalar rr[3] = {x - center[0], y - center[1], z - center[2]};
103a515125bSLeila Ghaffari   // (I - q q^T) r: distance from dc_axis (or from center if dc_axis is the zero vector)
1042b916ea7SJeremy L Thompson   for (CeedInt i = 0; i < 3; i++) rr[i] -= dc_axis[i] * Dot3(dc_axis, rr);
10564667825SJames Wright   const CeedScalar r           = Norm3(rr);
106a515125bSLeila Ghaffari   const CeedScalar delta_theta = r <= rc ? thetaC * (1. + cos(M_PI * r / rc)) / 2. : 0.;
107d1b9ef12SLeila Ghaffari   const CeedScalar theta       = theta0 * exp(Square(N) * z / g) + delta_theta;
108a515125bSLeila Ghaffari 
109a515125bSLeila Ghaffari   // -- Exner pressure, hydrostatic balance
1102b916ea7SJeremy L Thompson   const CeedScalar Pi = 1. + Square(g) * (exp(-Square(N) * z / g) - 1.) / (cp * theta0 * Square(N));
111a515125bSLeila Ghaffari 
112a515125bSLeila Ghaffari   // Initial Conditions
113cbe60e31SLeila Ghaffari   CeedScalar Y[5] = {0.};
114cbe60e31SLeila Ghaffari   Y[0]            = P0 * pow(Pi, cp / Rd);
115cbe60e31SLeila Ghaffari   Y[1]            = 0.0;
116cbe60e31SLeila Ghaffari   Y[2]            = 0.0;
117cbe60e31SLeila Ghaffari   Y[3]            = 0.0;
1182ca21690SKenneth E. Jansen   Y[4]            = Pi * theta;
119a515125bSLeila Ghaffari 
120edcfef1bSKenneth E. Jansen   return StateFromY(gas, Y);
121a515125bSLeila Ghaffari }
122a515125bSLeila Ghaffari 
123a515125bSLeila Ghaffari // *****************************************************************************
124a515125bSLeila Ghaffari // This QFunction sets the initial conditions for density current
125a515125bSLeila Ghaffari // *****************************************************************************
ICsDC(void * ctx,CeedInt Q,const CeedScalar * const * in,CeedScalar * const * out)1262b916ea7SJeremy L Thompson CEED_QFUNCTION(ICsDC)(void *ctx, CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) {
127a515125bSLeila Ghaffari   const CeedScalar(*X)[CEED_Q_VLA] = (const CeedScalar(*)[CEED_Q_VLA])in[0];
128a515125bSLeila Ghaffari   CeedScalar(*q0)[CEED_Q_VLA]      = (CeedScalar(*)[CEED_Q_VLA])out[0];
129a515125bSLeila Ghaffari 
130cbe60e31SLeila Ghaffari   const DensityCurrentContext context = (DensityCurrentContext)ctx;
131*cde3d787SJames Wright   NewtonianIGProperties       gas     = context->newt_ctx.gas;
132cbe60e31SLeila Ghaffari 
1333d65b166SJames Wright   CeedPragmaSIMD for (CeedInt i = 0; i < Q; i++) {
134a515125bSLeila Ghaffari     const CeedScalar x[] = {X[0][i], X[1][i], X[2][i]};
135cbe60e31SLeila Ghaffari     State            s   = Exact_DC(3, 0., x, 5, ctx);
136a541e550SJames Wright     CeedScalar       q[5];
137*cde3d787SJames Wright     StateToQ(gas, s, q, context->newt_ctx.state_var);
1382b916ea7SJeremy L Thompson     for (CeedInt j = 0; j < 5; j++) q0[j][i] = q[j];
139b193fadcSJames Wright   }
140a515125bSLeila Ghaffari   return 0;
141a515125bSLeila Ghaffari }
142