ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InitCustomPert_Bomex.H
Go to the documentation of this file.
1 
2  ParmParse pp("prob");
3 
4  amrex::Real rho_0 = 0.0; pp.query("rho_0", rho_0);
5  amrex::Real T_0 = 0.0; pp.query("T_0" , T_0);
6  amrex::Real A_0 = 1.0; pp.query("A_0" , A_0);
7  amrex::Real KE_0 = 0.1; pp.query("KE_0" , KE_0);
8 
9  Real T_0_Pert_Mag = 0.0; pp.query("T_0_Pert_Mag", T_0_Pert_Mag);
10  Real qv_0_Pert_Mag = 0.0; pp.query("qv_0_Pert_Mag", qv_0_Pert_Mag);
11  Real pert_ref_height = 100.0; pp.query("pert_ref_height", pert_ref_height);
12  bool custom_TKE = false; pp.query("custom_TKE", custom_TKE);
13 
14  const bool use_moisture = (sc.moisture_type != MoistureType::None);
15 
16  const Real rdOcp = sc.rdOcp;
17 
18  ParallelForRNG(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k, const RandomEngine& engine) noexcept
19  {
20  // Geometry
21  const Real* prob_lo = geomdata.ProbLo();
22  const Real* prob_hi = geomdata.ProbHi();
23  const Real* dx = geomdata.CellSize();
24  const Real x = prob_lo[0] + (i + 0.5) * dx[0];
25  const Real y = prob_lo[1] + (j + 0.5) * dx[1];
26  const Real z = prob_lo[2] + (k + 0.5) * dx[2];
27 
28  // Define a point (xc,yc,zc) at the center of the domain
29  const Real xc = 0.5 * (prob_lo[0] + prob_hi[0]);
30  const Real yc = 0.5 * (prob_lo[1] + prob_hi[1]);
31  const Real zc = 0.5 * (prob_lo[2] + prob_hi[2]);
32 
33  const Real r = std::sqrt((x-xc)*(x-xc) + (y-yc)*(y-yc) + (z-zc)*(z-zc));
34  //
35  // Add temperature perturbations -- we want to keep pressure constant
36  // so these effectively end up as density perturbations
37  //
38  if ((z <= pert_ref_height) && (T_0_Pert_Mag != 0.0)) {
39 
40  Real rhotheta = state(i,j,k,RhoTheta_comp);
41  Real rho = state(i,j,k,Rho_comp);
42  Real qv = state(i,j,k,RhoQ1_comp) / rho;
43  Real Told = getTgivenRandRTh(rho,rhotheta,qv);
44  Real P = getPgivenRTh(rhotheta,qv);
45 
46  Real rand_double = amrex::Random(engine); // Between 0.0 and 1.0
47  Real Tpert = (rand_double*2.0 - 1.0)*T_0_Pert_Mag;
48  Real Tnew = Told + Tpert;
49 
50  Real theta_new = getThgivenTandP(Tnew,P,rdOcp);
51  Real rhonew = getRhogivenThetaPress(theta_new,P,rdOcp,qv);
52  state_pert(i, j, k, Rho_comp) = rhonew - rho;
53 
54  // Note we do not perturb this
55  state_pert(i, j, k, RhoTheta_comp) = 0.0;
56 
57  // Instead of perturbing (rho theta) we perturb T and hold (rho theta) fixed,
58  // which ends up being stored as a perturbation in rho
59  // state_pert(i, j, k, RhoTheta_comp) = (rand_double*2.0 - 1.0)*T_0_Pert_Mag;
60  }
61 
62  // Set scalar = A_0*exp(-10r^2), where r is distance from center of domain
63  state_pert(i, j, k, RhoScalar_comp) = A_0 * exp(-10.*r*r);
64 
65  // Set an initial value for KE
66  if (custom_TKE) {
67  state_pert(i, j, k, RhoKE_comp) = (1.0 - z/prob_hi[2]) * r_hse(i,j,k);
68  } else {
69  state_pert(i, j, k, RhoKE_comp) = KE_0;
70  }
71 
72  if (use_moisture) {
73  state_pert(i, j, k, RhoQ1_comp) = 0.0;
74  state_pert(i, j, k, RhoQ2_comp) = 0.0;
75  if ((z <= pert_ref_height) && (qv_0_Pert_Mag != 0.0))
76  {
77  Real rhoold = state(i,j,k,Rho_comp);
78  Real rhonew = rhoold + state_pert(i,j,k,Rho_comp);
79 
80  Real qvold = state(i,j,k,RhoQ1_comp) / rhoold;
81 
82  Real rand_double = amrex::Random(engine); // Between 0.0 and 1.0
83  Real qvnew = qvold + (rand_double*2.0 - 1.0)*qv_0_Pert_Mag;
84 
85  state_pert(i, j, k, RhoQ1_comp) = rhonew * qvnew - rhoold * qvold;
86  }
87  }
88  });
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 RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#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
amrex::Real A_0
Definition: ERF_InitCustomPert_Bomex.H:6
Real T_0_Pert_Mag
Definition: ERF_InitCustomPert_Bomex.H:9
Real pert_ref_height
Definition: ERF_InitCustomPert_Bomex.H:11
const bool use_moisture
Definition: ERF_InitCustomPert_Bomex.H:14
const Real rdOcp
Definition: ERF_InitCustomPert_Bomex.H:16
ParmParse pp("prob")
amrex::Real rho_0
Definition: ERF_InitCustomPert_Bomex.H:4
amrex::Real KE_0
Definition: ERF_InitCustomPert_Bomex.H:7
Real qv_0_Pert_Mag
Definition: ERF_InitCustomPert_Bomex.H:10
amrex::Real T_0
Definition: ERF_InitCustomPert_Bomex.H:5
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) { state_pert(i, j, k, RhoQ1_comp)=0.0;state_pert(i, j, k, RhoQ2_comp)=0.0;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;} } })
bool custom_TKE
Definition: ERF_InitCustomPert_Bomex.H:12
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
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ P
Definition: ERF_IndexDefines.H:148
@ qv
Definition: ERF_Kessler.H:28