ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
AdvanceEWP.cpp File Reference
#include <EWP.H>
#include <IndexDefines.H>
Include dependency graph for AdvanceEWP.cpp:

Functions

void ewp_advance (int lev, const Geometry &geom, const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, MultiFab &W_old, MultiFab &mf_vars_ewp, const amrex::MultiFab &mf_Nturb)
 
void ewp_update (const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, const MultiFab &mf_vars_ewp)
 
void ewp_source_terms_cellcentered (const Geometry &geom, const MultiFab &cons_in, const MultiFab &U_old, const MultiFab &V_old, const MultiFab &W_old, MultiFab &mf_vars_ewp, const amrex::MultiFab &mf_Nturb)
 

Variables

Real R_ewp = 30.0
 
Real z_c_ewp = 100.0
 
Real C_T_ewp = 0.8
 
Real C_TKE_ewp = 0.0
 

Function Documentation

◆ ewp_advance()

void ewp_advance ( int  lev,
const Geometry &  geom,
const Real &  dt_advance,
MultiFab &  cons_in,
MultiFab &  U_old,
MultiFab &  V_old,
MultiFab &  W_old,
MultiFab &  mf_vars_ewp,
const amrex::MultiFab &  mf_Nturb 
)
15 {
16  ewp_source_terms_cellcentered(geom, cons_in, U_old, V_old, W_old, mf_vars_ewp, mf_Nturb);
17  ewp_update(dt_advance, cons_in, U_old, V_old, mf_vars_ewp);
18 }
void ewp_source_terms_cellcentered(const Geometry &geom, const MultiFab &cons_in, const MultiFab &U_old, const MultiFab &V_old, const MultiFab &W_old, MultiFab &mf_vars_ewp, const amrex::MultiFab &mf_Nturb)
Definition: AdvanceEWP.cpp:54
void ewp_update(const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, const MultiFab &mf_vars_ewp)
Definition: AdvanceEWP.cpp:21

Referenced by ERF::Advance().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ewp_source_terms_cellcentered()

void ewp_source_terms_cellcentered ( const Geometry &  geom,
const MultiFab &  cons_in,
const MultiFab &  U_old,
const MultiFab &  V_old,
const MultiFab &  W_old,
MultiFab &  mf_vars_ewp,
const amrex::MultiFab &  mf_Nturb 
)
58 {
59 
60  auto dx = geom.CellSizeArray();
61  auto ProbHiArr = geom.ProbHiArray();
62  auto ProbLoArr = geom.ProbLoArray();
63 
64  // Domain valid box
65  const amrex::Box& domain = geom.Domain();
66  int domlo_x = domain.smallEnd(0);
67  int domhi_x = domain.bigEnd(0) + 1;
68  int domlo_y = domain.smallEnd(1);
69  int domhi_y = domain.bigEnd(1) + 1;
70  int domlo_z = domain.smallEnd(2);
71  int domhi_z = domain.bigEnd(2) + 1;
72 
73  // The order of variables are - Vabs dVabsdt, dudt, dvdt, dTKEdt
74  mf_vars_ewp.setVal(0.0);
75 
76  for ( MFIter mfi(cons_in,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
77 
78  const Box& bx = mfi.tilebox();
79  const Box& gbx = mfi.growntilebox(1);
80  auto cons_array = cons_in.array(mfi);
81  auto ewp_array = mf_vars_ewp.array(mfi);
82  auto Nturb_array = mf_Nturb.array(mfi);
83  auto u_vel = U_old.array(mfi);
84  auto v_vel = V_old.array(mfi);
85  auto w_vel = W_old.array(mfi);
86 
87  amrex::IntVect lo = bx.smallEnd();
88 
89  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
90  int ii = amrex::min(amrex::max(i, domlo_x), domhi_x);
91  int jj = amrex::min(amrex::max(j, domlo_y), domhi_y);
92  int kk = amrex::min(amrex::max(k, domlo_z), domhi_z);
93 
94 
95  Real x = ProbLoArr[0] + (ii+0.5) * dx[0];
96  Real y = ProbLoArr[1] + (jj+0.5) * dx[1];
97  Real z = ProbLoArr[2] + (kk+0.5) * dx[2];
98 
99  // Compute Fitch source terms
100 
101  Real Vabs = std::pow(u_vel(i,j,k)*u_vel(i,j,k) +
102  v_vel(i,j,k)*v_vel(i,j,k) +
103  w_vel(i,j,kk)*w_vel(i,j,kk), 0.5);
104 
105  Real L_wake = std::pow(dx[0]*dx[1],0.5);
106  Real K_turb = 100.0;
107  Real sigma_0 = 60.0;
108  Real sigma_e = Vabs/(3.0*K_turb*L_wake)*(std::pow(2.0*K_turb*L_wake/Vabs + std::pow(sigma_0,2),3.0/2.0) - std::pow(sigma_0,3));
109 
110  Real phi = std::atan2(v_vel(i,j,k),u_vel(i,j,k)); // Wind direction w.r.t the x-dreiction
111  Real fac = -std::pow(M_PI/8.0,0.5)*C_T_ewp*std::pow(R_ewp,2)*std::pow(Vabs,2)/(dx[0]*dx[1]*sigma_e)*std::exp(-0.5*std::pow((z - z_c_ewp)/sigma_e,2));
112  ewp_array(i,j,k,0) = fac*std::cos(phi)*Nturb_array(i,j,k);
113  ewp_array(i,j,k,1) = fac*std::sin(phi)*Nturb_array(i,j,k);
114  ewp_array(i,j,k,2) = 0.0;
115  });
116  }
117 }
Real C_T_ewp
Definition: AdvanceEWP.cpp:7
Real R_ewp
Definition: AdvanceEWP.cpp:5
Real z_c_ewp
Definition: AdvanceEWP.cpp:6

Referenced by ewp_advance().

Here is the caller graph for this function:

◆ ewp_update()

void ewp_update ( const Real &  dt_advance,
MultiFab &  cons_in,
MultiFab &  U_old,
MultiFab &  V_old,
const MultiFab &  mf_vars_ewp 
)
25 {
26 
27  for ( MFIter mfi(cons_in,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
28 
29  Box bx = mfi.tilebox();
30  Box tbx = mfi.nodaltilebox(0);
31  Box tby = mfi.nodaltilebox(1);
32 
33  auto cons_array = cons_in.array(mfi);
34  auto ewp_array = mf_vars_ewp.array(mfi);
35  auto u_vel = U_old.array(mfi);
36  auto v_vel = V_old.array(mfi);
37 
38  ParallelFor(tbx, tby, bx,
39  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
40  {
41  u_vel(i,j,k) = u_vel(i,j,k) + (ewp_array(i-1,j,k,0) + ewp_array(i,j,k,0))/2.0*dt_advance;
42  },
43  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
44  {
45  v_vel(i,j,k) = v_vel(i,j,k) + (ewp_array(i,j-1,k,1) + ewp_array(i,j,k,1))/2.0*dt_advance;
46  },
47  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
48  {
49  cons_array(i,j,k,RhoQKE_comp) = cons_array(i,j,k,RhoQKE_comp) + ewp_array(i,j,k,2)*dt_advance;
50  });
51  }
52 }
#define RhoQKE_comp
Definition: IndexDefines.H:15

Referenced by ewp_advance().

Here is the caller graph for this function:

Variable Documentation

◆ C_T_ewp

Real C_T_ewp = 0.8

◆ C_TKE_ewp

Real C_TKE_ewp = 0.0

◆ R_ewp

Real R_ewp = 30.0

◆ z_c_ewp

Real z_c_ewp = 100.0