ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InitCustomPert_RICO.H
Go to the documentation of this file.
1  ParmParse pp_prob("prob");
2 
3  Real A_0 = 1.0; pp_prob.query("A_0", A_0);
4  Real KE_0 = 0.1; pp_prob.query("KE_0", KE_0);
5 
6  // sinusoidal perturbations for regression tests
7  Real pert_deltaT = 0.0; pp_prob.query("pert_deltaT", pert_deltaT);
8  Real pert_deltaQV = 0.0; pp_prob.query("pert_deltaQV", pert_deltaQV);
9  Real pert_periods_T = 0.0; pp_prob.query("pert_periods_T", pert_periods_T);
10  Real pert_periods_QV = 0.0; pp_prob.query("pert_periods_QV", pert_periods_QV);
11 
12  Real T_0_Pert_Mag = 0.0; pp_prob.query("T_0_Pert_Mag", T_0_Pert_Mag);
13  Real qv_0_Pert_Mag = 0.0; pp_prob.query("qv_0_Pert_Mag", qv_0_Pert_Mag);
14 
15  Real wbar_sub_max = -0.65; pp_prob.query("wbar_sub_max", wbar_sub_max);
16 
17  Real pert_ref_height = 100.0; pp_prob.query("pert_ref_height", pert_ref_height);
18 
19  bool custom_TKE = false; pp_prob.query("custom_TKE", custom_TKE);
20 
21  const Real prob_lo_x = geomdata.ProbLo()[0];
22  const Real prob_hi_x = geomdata.ProbHi()[0];
23 
26 
27  Real tfac = pert_deltaT * std::exp(0.5) / pert_ref_height;
28  Real qvfac = pert_deltaQV * std::exp(0.5) / pert_ref_height;
29 
30  const bool use_moisture = (sc.moisture_type != MoistureType::None);
31 
32  const Real rdOcp = sc.rdOcp;
33 
34  ParallelForRNG(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const RandomEngine& engine) noexcept
35  {
36  // Geometry
37  const Real* prob_lo = geomdata.ProbLo();
38  const Real* prob_hi = geomdata.ProbHi();
39  const Real* dx = geomdata.CellSize();
40  const Real x = prob_lo[0] + (i + 0.5) * dx[0];
41  const Real y = prob_lo[1] + (j + 0.5) * dx[1];
42  const Real z = prob_lo[2] + (k + 0.5) * dx[2];
43 
44  // Define a point (xc,yc,zc) at the center of the domain
45  const Real xc = 0.5 * (prob_lo[0] + prob_hi[0]);
46  const Real yc = 0.5 * (prob_lo[1] + prob_hi[1]);
47  const Real zc = 0.5 * (prob_lo[2] + prob_hi[2]);
48 
49  const Real r = std::sqrt((x-xc)*(x-xc) + (y-yc)*(y-yc) + (z-zc)*(z-zc));
50  //
51  // Add temperature perturbations -- we want to keep pressure constant
52  // so these effectively end up as density perturbations
53  //
54  if ((z <= pert_ref_height) && (T_0_Pert_Mag != 0.0)) {
55 
56  Real rhotheta = state(i,j,k,RhoTheta_comp);
57  Real rho = state(i,j,k,Rho_comp);
58  Real qv = state(i,j,k,RhoQ1_comp) / rho;
59  Real Told = getTgivenRandRTh(rho,rhotheta,qv);
60  Real P = getPgivenRTh(rhotheta,qv);
61 
62  Real rand_double = amrex::Random(engine); // Between 0.0 and 1.0
63  Real Tpert = (rand_double*2.0 - 1.0)*T_0_Pert_Mag;
64  Real Tnew = Told + Tpert;
65 
66  Real theta_new = getThgivenTandP(Tnew,P,rdOcp);
67  Real rhonew = getRhogivenThetaPress(theta_new,P,rdOcp,qv);
68  state_pert(i, j, k, Rho_comp) = rhonew - rho;
69 
70  // Note we do not perturb this
71  state_pert(i, j, k, RhoTheta_comp) = 0.0;
72 
73  // Instead of perturbing (rho theta) we perturb T and hold (rho theta) fixed,
74  // which ends up being stored as a perturbation in rho
75  // state_pert(i, j, k, RhoTheta_comp) = (rand_double*2.0 - 1.0)*T_0_Pert_Mag;
76  }
77 
78  // Set scalar = A_0*exp(-10r^2), where r is distance from center of domain
79  state_pert(i, j, k, RhoScalar_comp) = A_0 * exp(-10.*r*r);
80 
81  // Set an initial value for KE
82  if (custom_TKE) {
83  state_pert(i, j, k, RhoKE_comp) = (1.0 - z/prob_hi[2]) * r_hse(i,j,k);
84  } else {
85  state_pert(i, j, k, RhoKE_comp) = KE_0;
86  }
87 
88  if (use_moisture) {
89  if ((z <= pert_ref_height) && (qv_0_Pert_Mag != 0.0))
90  {
91  Real rhoold = state(i,j,k,Rho_comp);
92  Real rhonew = rhoold + state_pert(i,j,k,Rho_comp);
93 
94  Real qvold = state(i,j,k,RhoQ1_comp) / rhoold;
95 
96  Real rand_double = amrex::Random(engine); // Between 0.0 and 1.0
97  Real qvnew = qvold + (rand_double*2.0 - 1.0)*qv_0_Pert_Mag;
98 
99  state_pert(i, j, k, RhoQ1_comp) = rhonew * qvnew - rhoold * qvold;
100  }
101  }
102 
103  // sinusoidal variation of theta and qv for regression tests
104  if (pert_deltaT != 0.0)
105  {
106  const Real zl = z / pert_ref_height;
107  const Real damp = std::exp(-0.5 * zl * zl);
108 
109  Real rhotheta = state(i,j,k,RhoTheta_comp);
110  Real rho = state(i,j,k,Rho_comp);
111  Real qv = state(i,j,k,RhoQ1_comp) / rho;
112  Real Told = getTgivenRandRTh(rho,rhotheta,qv);
113  Real P = getPgivenRTh(rhotheta,qv);
114 
115  Real Tpert = tfac * damp * z * std::cos(cval * (x - xc));
116  Real Tnew = Told + Tpert;
117 
118  Real theta_new = getThgivenTandP(Tnew,P,rdOcp);
119  Real rhonew = getRhogivenThetaPress(theta_new,P,rdOcp,qv);
120 
121  state_pert(i, j, k, Rho_comp) = rhonew - rho;
122  state_pert(i, j, k, RhoTheta_comp) = 0.0;
123  }
124 
125  if (use_moisture && pert_deltaQV != 0.0)
126  {
127  const Real zl = z / pert_ref_height;
128  const Real damp = std::exp(-0.5 * zl * zl);
129 
130  Real rhoold = state(i,j,k,Rho_comp);
131  Real rhonew = rhoold + state_pert(i,j,k,Rho_comp);
132 
133  Real qvold = state(i,j,k,RhoQ1_comp) / rhoold;
134  Real qvnew = qvold + qvfac * damp * z * std::cos(dval * (x - xc));
135 
136  state_pert(i, j, k, RhoQ1_comp) = rhonew * qvnew - rhoold * qvold;
137  }
138  });
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getRhogivenThetaPress(const amrex::Real th, const amrex::Real p, const amrex::Real rdOcp, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:96
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getThgivenTandP(const amrex::Real T, const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:18
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
const Real zc
Definition: ERF_InitCustomPert_ABL.H:35
const Real yc
Definition: ERF_InitCustomPert_ABL.H:34
const Real xc
Definition: ERF_InitCustomPert_ABL.H:33
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
state_pert(i, j, k, RhoTheta_comp)
rho
Definition: ERF_InitCustomPert_Bubble.H:106
const amrex::Real * prob_lo
Definition: ERF_InitCustomPert_IsentropicVortex.H:16
const amrex::Real * prob_hi
Definition: ERF_InitCustomPert_IsentropicVortex.H:17
Real pert_periods_QV
Definition: ERF_InitCustomPert_RICO.H:10
Real T_0_Pert_Mag
Definition: ERF_InitCustomPert_RICO.H:12
Real pert_ref_height
Definition: ERF_InitCustomPert_RICO.H:17
const bool use_moisture
Definition: ERF_InitCustomPert_RICO.H:30
Real KE_0
Definition: ERF_InitCustomPert_RICO.H:4
Real dval
Definition: ERF_InitCustomPert_RICO.H:25
const Real rdOcp
Definition: ERF_InitCustomPert_RICO.H:32
ParmParse pp_prob("prob")
ParallelForRNG(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const RandomEngine &engine) noexcept { const Real *prob_lo=geomdata.ProbLo();const Real *prob_hi=geomdata.ProbHi();const Real *dx=geomdata.CellSize();const Real x=prob_lo[0]+(i+0.5) *dx[0];const Real y=prob_lo[1]+(j+0.5) *dx[1];const Real z=prob_lo[2]+(k+0.5) *dx[2];const Real xc=0.5 *(prob_lo[0]+prob_hi[0]);const Real yc=0.5 *(prob_lo[1]+prob_hi[1]);const Real zc=0.5 *(prob_lo[2]+prob_hi[2]);const Real r=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc)+(z-zc) *(z-zc));if((z<=pert_ref_height) &&(T_0_Pert_Mag !=0.0)) { Real rhotheta=state(i, j, k, RhoTheta_comp);Real rho=state(i, j, k, Rho_comp);Real qv=state(i, j, k, RhoQ1_comp)/rho;Real Told=getTgivenRandRTh(rho, rhotheta, qv);Real P=getPgivenRTh(rhotheta, qv);Real rand_double=amrex::Random(engine);Real Tpert=(rand_double *2.0 - 1.0) *T_0_Pert_Mag;Real Tnew=Told+Tpert;Real theta_new=getThgivenTandP(Tnew, P, rdOcp);Real rhonew=getRhogivenThetaPress(theta_new, P, rdOcp, qv);state_pert(i, j, k, Rho_comp)=rhonew - rho;state_pert(i, j, k, RhoTheta_comp)=0.0;} state_pert(i, j, k, RhoScalar_comp)=A_0 *exp(-10.*r *r);if(custom_TKE) { state_pert(i, j, k, RhoKE_comp)=(1.0 - z/prob_hi[2]) *r_hse(i, j, k);} else { state_pert(i, j, k, RhoKE_comp)=KE_0;} if(use_moisture) { if((z<=pert_ref_height) &&(qv_0_Pert_Mag !=0.0)) { Real rhoold=state(i, j, k, Rho_comp);Real rhonew=rhoold+state_pert(i, j, k, Rho_comp);Real qvold=state(i, j, k, RhoQ1_comp)/rhoold;Real rand_double=amrex::Random(engine);Real qvnew=qvold+(rand_double *2.0 - 1.0) *qv_0_Pert_Mag;state_pert(i, j, k, RhoQ1_comp)=rhonew *qvnew - rhoold *qvold;} } if(pert_deltaT !=0.0) { const Real zl=z/pert_ref_height;const Real damp=std::exp(-0.5 *zl *zl);Real rhotheta=state(i, j, k, RhoTheta_comp);Real rho=state(i, j, k, Rho_comp);Real qv=state(i, j, k, RhoQ1_comp)/rho;Real Told=getTgivenRandRTh(rho, rhotheta, qv);Real P=getPgivenRTh(rhotheta, qv);Real Tpert=tfac *damp *z *std::cos(cval *(x - xc));Real Tnew=Told+Tpert;Real theta_new=getThgivenTandP(Tnew, P, rdOcp);Real rhonew=getRhogivenThetaPress(theta_new, P, rdOcp, qv);state_pert(i, j, k, Rho_comp)=rhonew - rho;state_pert(i, j, k, RhoTheta_comp)=0.0;} if(use_moisture &&pert_deltaQV !=0.0) { const Real zl=z/pert_ref_height;const Real damp=std::exp(-0.5 *zl *zl);Real rhoold=state(i, j, k, Rho_comp);Real rhonew=rhoold+state_pert(i, j, k, Rho_comp);Real qvold=state(i, j, k, RhoQ1_comp)/rhoold;Real qvnew=qvold+qvfac *damp *z *std::cos(dval *(x - xc));state_pert(i, j, k, RhoQ1_comp)=rhonew *qvnew - rhoold *qvold;} })
const Real prob_lo_x
Definition: ERF_InitCustomPert_RICO.H:21
Real tfac
Definition: ERF_InitCustomPert_RICO.H:27
const Real prob_hi_x
Definition: ERF_InitCustomPert_RICO.H:22
Real wbar_sub_max
Definition: ERF_InitCustomPert_RICO.H:15
Real pert_deltaQV
Definition: ERF_InitCustomPert_RICO.H:8
Real pert_periods_T
Definition: ERF_InitCustomPert_RICO.H:9
Real qvfac
Definition: ERF_InitCustomPert_RICO.H:28
Real A_0
Definition: ERF_InitCustomPert_RICO.H:3
Real qv_0_Pert_Mag
Definition: ERF_InitCustomPert_RICO.H:13
Real pert_deltaT
Definition: ERF_InitCustomPert_RICO.H:7
Real cval
Definition: ERF_InitCustomPert_RICO.H:24
bool custom_TKE
Definition: ERF_InitCustomPert_RICO.H:19
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ P
Definition: ERF_IndexDefines.H:148
@ qv
Definition: ERF_Kessler.H:28