ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_DiffusionSrcForState_N.cpp File Reference
#include <ERF_Diffusion.H>
#include <ERF_EddyViscosity.H>
#include <ERF_PBLModels.H>
Include dependency graph for ERF_DiffusionSrcForState_N.cpp:

Functions

void DiffusionSrcForState_N (const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &exp_most, const Array4< const Real > &u, const Array4< const Real > &v, const Array4< const Real > &cell_data, const Array4< const Real > &cell_prim, const Array4< Real > &cell_rhs, const Array4< Real > &xflux, const Array4< Real > &yflux, const Array4< Real > &zflux, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_m, const Array4< const Real > &mf_u, const Array4< const Real > &mf_v, Array4< Real > &hfx_z, Array4< Real > &qfx1_z, Array4< Real > &qfx2_z, Array4< Real > &diss, const Array4< const Real > &mu_turb, const SolverChoice &solverChoice, const int level, const Array4< const Real > &tm_arr, const GpuArray< Real, AMREX_SPACEDIM > grav_gpu, const BCRec *bc_ptr, const bool use_most)
 

Function Documentation

◆ DiffusionSrcForState_N()

void DiffusionSrcForState_N ( const Box &  bx,
const Box &  domain,
int  start_comp,
int  num_comp,
const bool &  exp_most,
const Array4< const Real > &  u,
const Array4< const Real > &  v,
const Array4< const Real > &  cell_data,
const Array4< const Real > &  cell_prim,
const Array4< Real > &  cell_rhs,
const Array4< Real > &  xflux,
const Array4< Real > &  yflux,
const Array4< Real > &  zflux,
const GpuArray< Real, AMREX_SPACEDIM > &  cellSizeInv,
const Array4< const Real > &  SmnSmn_a,
const Array4< const Real > &  mf_m,
const Array4< const Real > &  mf_u,
const Array4< const Real > &  mf_v,
Array4< Real > &  hfx_z,
Array4< Real > &  qfx1_z,
Array4< Real > &  qfx2_z,
Array4< Real > &  diss,
const Array4< const Real > &  mu_turb,
const SolverChoice solverChoice,
const int  level,
const Array4< const Real > &  tm_arr,
const GpuArray< Real, AMREX_SPACEDIM >  grav_gpu,
const BCRec *  bc_ptr,
const bool  use_most 
)

Function for computing the scalar RHS for diffusion operator without terrain.

Parameters
[in]bxcell center box to loop over
[in]domainbox of the whole domain
[in]start_compstarting component index
[in]num_compnumber of components
[in]uvelocity in x-dir
[in]vvelocity in y-dir
[in]cell_dataconserved cell center vars
[in]cell_primprimitive cell center vars
[out]cell_rhsRHS for cell center vars
[in]xfluxflux in x-dir
[in]yfluxflux in y-dir
[in]zfluxflux in z-dir
[in]cellSizeInvinverse cell size array
[in]SmnSmn_astrain rate magnitude
[in]mf_mmap factor at cell center
[in]mf_umap factor at x-face
[in]mf_vmap factor at y-face
[in,out]hfx_zheat flux in z-dir
[in,out]qfx1_zheat flux in z-dir
[out]qfx2_zheat flux in z-dir
[in]dissdissipation of TKE
[in]mu_turbturbulent viscosity
[in]diffChoicecontainer of diffusion parameters
[in]turbChoicecontainer of turbulence parameters
[in]tm_arrtheta mean array
[in]grav_gpugravity vector
[in]bc_ptrcontainer with boundary conditions
[in]use_mostwhether we have turned on MOST BCs
67 {
68  BL_PROFILE_VAR("DiffusionSrcForState_N()",DiffusionSrcForState_N);
69 
70  DiffChoice diffChoice = solverChoice.diffChoice;
71  TurbChoice turbChoice = solverChoice.turbChoice[level];
72 
73  amrex::ignore_unused(use_most);
74 
75  const Real dx_inv = cellSizeInv[0];
76  const Real dy_inv = cellSizeInv[1];
77  const Real dz_inv = cellSizeInv[2];
78 
79  const auto& dom_lo = lbound(domain);
80  const auto& dom_hi = ubound(domain);
81 
82  Real l_inv_theta0 = 1.0 / turbChoice.theta_ref;
83  Real l_abs_g = std::abs(grav_gpu[2]);
84 
85  bool l_use_ddorf = (turbChoice.les_type == LESType::Deardorff);
86  bool l_use_mynn = (turbChoice.pbl_type == PBLType::MYNN25);
87 
88  bool l_consA = (diffChoice.molec_diff_type == MolecDiffType::ConstantAlpha);
89  bool l_turb = ( (turbChoice.les_type == LESType::Smagorinsky) ||
90  (turbChoice.les_type == LESType::Deardorff ) ||
91  (turbChoice.pbl_type == PBLType::MYNN25 ) ||
92  (turbChoice.pbl_type == PBLType::YSU ) );
93 
94  const Box xbx = surroundingNodes(bx,0);
95  const Box ybx = surroundingNodes(bx,1);
96  const Box zbx = surroundingNodes(bx,2);
97 
98  const int end_comp = start_comp + num_comp - 1;
99 
100  // Theta, KE, Scalar
101  Vector<Real> alpha_eff(NPRIMVAR_max, 0.0);
102  if (l_consA) {
103  for (int i = 0; i < NPRIMVAR_max; ++i) {
104  switch (i) {
105  case PrimTheta_comp:
106  alpha_eff[PrimTheta_comp] = diffChoice.alpha_T;
107  break;
108  case PrimScalar_comp:
109  alpha_eff[PrimScalar_comp] = diffChoice.alpha_C;
110  break;
111  case PrimQ1_comp:
112  alpha_eff[PrimQ1_comp] = diffChoice.alpha_C;
113  break;
114  case PrimQ2_comp:
115  alpha_eff[PrimQ2_comp] = diffChoice.alpha_C;
116  break;
117  case PrimQ3_comp:
118  alpha_eff[PrimQ3_comp] = diffChoice.alpha_C;
119  break;
120  case PrimQ4_comp:
121  alpha_eff[PrimQ4_comp] = diffChoice.alpha_C;
122  break;
123  case PrimQ5_comp:
124  alpha_eff[PrimQ5_comp] = diffChoice.alpha_C;
125  break;
126  case PrimQ6_comp:
127  alpha_eff[PrimQ6_comp] = diffChoice.alpha_C;
128  break;
129  default:
130  alpha_eff[i] = 0.0;
131  break;
132  }
133  }
134  } else {
135  for (int i = 0; i < NPRIMVAR_max; ++i) {
136  switch (i) {
137  case PrimTheta_comp:
138  alpha_eff[PrimTheta_comp] = diffChoice.rhoAlpha_T;
139  break;
140  case PrimScalar_comp:
141  alpha_eff[PrimScalar_comp] = diffChoice.rhoAlpha_C;
142  break;
143  case PrimQ1_comp:
144  alpha_eff[PrimQ1_comp] = diffChoice.rhoAlpha_C;
145  break;
146  case PrimQ2_comp:
147  alpha_eff[PrimQ2_comp] = diffChoice.rhoAlpha_C;
148  break;
149  case PrimQ3_comp:
150  alpha_eff[PrimQ3_comp] = diffChoice.rhoAlpha_C;
151  break;
152  case PrimQ4_comp:
153  alpha_eff[PrimQ4_comp] = diffChoice.rhoAlpha_C;
154  break;
155  case PrimQ5_comp:
156  alpha_eff[PrimQ5_comp] = diffChoice.rhoAlpha_C;
157  break;
158  case PrimQ6_comp:
159  alpha_eff[PrimQ6_comp] = diffChoice.rhoAlpha_C;
160  break;
161  default:
162  alpha_eff[i] = 0.0;
163  break;
164  }
165  }
166  }
167 
168  Vector<int> eddy_diff_idx{EddyDiff::Theta_h, EddyDiff::KE_h, EddyDiff::Scalar_h,
171  Vector<int> eddy_diff_idy{EddyDiff::Theta_h, EddyDiff::KE_h, EddyDiff::Scalar_h,
174  Vector<int> eddy_diff_idz{EddyDiff::Theta_v, EddyDiff::KE_v, EddyDiff::Scalar_v,
177 
178  // Device vectors
179  Gpu::AsyncVector<Real> alpha_eff_d;
180  Gpu::AsyncVector<int> eddy_diff_idx_d,eddy_diff_idy_d,eddy_diff_idz_d;
181  alpha_eff_d.resize(alpha_eff.size());
182  eddy_diff_idx_d.resize(eddy_diff_idx.size());
183  eddy_diff_idy_d.resize(eddy_diff_idy.size());
184  eddy_diff_idz_d.resize(eddy_diff_idz.size());
185 
186  Gpu::copy(Gpu::hostToDevice, alpha_eff.begin() , alpha_eff.end() , alpha_eff_d.begin());
187  Gpu::copy(Gpu::hostToDevice, eddy_diff_idx.begin(), eddy_diff_idx.end(), eddy_diff_idx_d.begin());
188  Gpu::copy(Gpu::hostToDevice, eddy_diff_idy.begin(), eddy_diff_idy.end(), eddy_diff_idy_d.begin());
189  Gpu::copy(Gpu::hostToDevice, eddy_diff_idz.begin(), eddy_diff_idz.end(), eddy_diff_idz_d.begin());
190 
191  // Capture pointers for device code
192  Real* d_alpha_eff = alpha_eff_d.data();
193  int* d_eddy_diff_idx = eddy_diff_idx_d.data();
194  int* d_eddy_diff_idy = eddy_diff_idy_d.data();
195  int* d_eddy_diff_idz = eddy_diff_idz_d.data();
196 
197  // Compute fluxes at each face
198  if (l_consA && l_turb) {
199  ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
200  {
201  const int qty_index = start_comp + n;
202  const int prim_index = qty_index - 1;
203  // const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
204  const int prim_scal_index = prim_index;
205 
206  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
207  BCVars::RhoScalar_bc_comp : qty_index;
208  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
209  bool ext_dir_on_xlo = ( ((bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir) ||
210  (bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir_prim))
211  && i == dom_lo.x);
212  bool ext_dir_on_xhi = ( ((bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir) ||
213  (bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir_prim))
214  && i == dom_hi.x+1);
215 
216  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i-1, j, k, Rho_comp) );
217  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
218  rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_scal_index])
219  + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_scal_index]) );
220 
221  if (ext_dir_on_xlo) {
222  xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
223  + 3. * cell_prim(i , j, k, prim_index)
224  - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
225  } else if (ext_dir_on_xhi) {
226  xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
227  - 3. * cell_prim(i-1, j, k, prim_index)
228  + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
229  } else {
230  xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) -
231  cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
232  }
233  });
234  ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
235  {
236  const int qty_index = start_comp + n;
237  const int prim_index = qty_index - 1;
238  // const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
239  const int prim_scal_index = prim_index;
240 
241  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
242  BCVars::RhoScalar_bc_comp : qty_index;
243  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
244  bool ext_dir_on_ylo = ( ((bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir) ||
245  (bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir_prim))
246  && j == dom_lo.y);
247  bool ext_dir_on_yhi = ( ((bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir) ||
248  (bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir_prim))
249  && j == dom_hi.y+1);
250 
251  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
252  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
253  rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_scal_index])
254  + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_scal_index]) );
255 
256  if (ext_dir_on_ylo) {
257  yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
258  + 3. * cell_prim(i, j , k, prim_index)
259  - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
260  } else if (ext_dir_on_yhi) {
261  yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
262  - 3. * cell_prim(i, j-1, k, prim_index)
263  + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
264  } else {
265  yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
266  }
267  });
268  ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
269  {
270  const int qty_index = start_comp + n;
271  const int prim_index = qty_index - 1;
272  // const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
273  const int prim_scal_index = prim_index;
274 
275  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
276  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
277  rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_scal_index])
278  + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_scal_index]) );
279 
280  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
281  BCVars::RhoScalar_bc_comp : qty_index;
282  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
283  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
284  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
285  && k == dom_lo.z);
286  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
287  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
288  && k == dom_hi.z+1);
289  bool most_on_zlo = ( use_most && exp_most &&
290  (bc_ptr[BCVars::cons_bc+qty_index].lo(2) == ERFBCType::foextrap) &&
291  k == dom_lo.z);
292 
293  if (ext_dir_on_zlo) {
294  zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
295  + 3. * cell_prim(i, j, k , prim_index)
296  - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
297  } else if (ext_dir_on_zhi) {
298  zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
299  - 3. * cell_prim(i, j, k-1, prim_index)
300  + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
301  } else if (most_on_zlo && (qty_index == RhoTheta_comp)) {
302  zflux(i,j,k,qty_index) = hfx_z(i,j,0);
303  } else if (most_on_zlo && (qty_index == RhoQ1_comp)) {
304  zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
305  } else {
306  zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
307  }
308 
309  if (qty_index == RhoTheta_comp) {
310  if (!most_on_zlo) {
311  hfx_z(i,j,k) = zflux(i,j,k,qty_index);
312  }
313  } else if (qty_index == RhoQ1_comp) {
314  if (!most_on_zlo) {
315  qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
316  }
317  } else if (qty_index == RhoQ2_comp) {
318  qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
319  }
320  });
321  } else if (l_turb) {
322  // with MolecDiffType::Constant or None
323  ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
324  {
325  const int qty_index = start_comp + n;
326  const int prim_index = qty_index - 1;
327 
328  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
329  BCVars::RhoScalar_bc_comp : qty_index;
330  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
331  bool ext_dir_on_xlo = ( ((bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir) ||
332  (bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir_prim))
333  && i == dom_lo.x);
334  bool ext_dir_on_xhi = ( ((bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir) ||
335  (bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir_prim))
336  && i == dom_hi.x+1);
337 
338  Real rhoAlpha = d_alpha_eff[prim_index];
339  rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_index])
340  + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_index]) );
341 
342  if (ext_dir_on_xlo) {
343  xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
344  + 3. * cell_prim(i , j, k, prim_index)
345  - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
346  } else if (ext_dir_on_xhi) {
347  xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
348  - 3. * cell_prim(i-1, j, k, prim_index)
349  + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
350  } else {
351  xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
352  }
353  });
354  ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
355  {
356  const int qty_index = start_comp + n;
357  const int prim_index = qty_index - 1;
358 
359  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
360  BCVars::RhoScalar_bc_comp : qty_index;
361  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
362  bool ext_dir_on_ylo = ( ((bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir) ||
363  (bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir_prim))
364  && j == dom_lo.y);
365  bool ext_dir_on_yhi = ( ((bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir) ||
366  (bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir_prim))
367  && j == dom_hi.y+1);
368 
369  Real rhoAlpha = d_alpha_eff[prim_index];
370  rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_index])
371  + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_index]) );
372 
373  if (ext_dir_on_ylo) {
374  yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
375  + 3. * cell_prim(i, j , k, prim_index)
376  - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
377  } else if (ext_dir_on_yhi) {
378  yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
379  - 3. * cell_prim(i, j-1, k, prim_index)
380  + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
381  } else {
382  yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
383  }
384  });
385  ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
386  {
387  const int qty_index = start_comp + n;
388  const int prim_index = qty_index - 1;
389 
390  Real rhoAlpha = d_alpha_eff[prim_index];
391  rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_index])
392  + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_index]) );
393 
394  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
395  BCVars::RhoScalar_bc_comp : qty_index;
396  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
397  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
398  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
399  && k == dom_lo.z);
400  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
401  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
402  && k == dom_hi.z+1);
403  bool most_on_zlo = ( use_most && exp_most &&
404  (bc_ptr[bc_comp].lo(2) == ERFBCType::foextrap) &&
405  k == dom_lo.z);
406 
407  if (ext_dir_on_zlo) {
408  zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
409  + 3. * cell_prim(i, j, k , prim_index)
410  - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
411  } else if (ext_dir_on_zhi) {
412  zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
413  - 3. * cell_prim(i, j, k-1, prim_index)
414  + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
415  } else if (most_on_zlo && (qty_index == RhoTheta_comp)) {
416  zflux(i,j,k,qty_index) = hfx_z(i,j,0);
417  } else if (most_on_zlo && (qty_index == RhoQ1_comp)) {
418  zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
419  } else {
420  zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
421  }
422 
423  if (qty_index == RhoTheta_comp) {
424  if (!most_on_zlo) {
425  hfx_z(i,j,k) = zflux(i,j,k,qty_index);
426  }
427  } else if (qty_index == RhoQ1_comp) {
428  if (!most_on_zlo) {
429  qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
430  }
431  } else if (qty_index == RhoQ2_comp) {
432  qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
433  }
434  });
435  } else if(l_consA) {
436  // without an LES/PBL model
437  ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
438  {
439  const int qty_index = start_comp + n;
440  const int prim_index = qty_index - 1;
441 
442  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
443  BCVars::RhoScalar_bc_comp : qty_index;
444  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
445  bool ext_dir_on_xlo = ( ((bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir) ||
446  (bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir_prim))
447  && i == dom_lo.x);
448  bool ext_dir_on_xhi = ( ((bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir) ||
449  (bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir_prim))
450  && i == dom_hi.x+1);
451 
452  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i-1, j, k, Rho_comp) );
453  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
454 
455  if (ext_dir_on_xlo) {
456  xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
457  + 3. * cell_prim(i , j, k, prim_index)
458  - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
459  } else if (ext_dir_on_xhi) {
460  xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
461  - 3. * cell_prim(i-1, j, k, prim_index)
462  + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
463  } else {
464  xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
465  }
466  });
467  ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
468  {
469  const int qty_index = start_comp + n;
470  const int prim_index = qty_index - 1;
471 
472  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
473  BCVars::RhoScalar_bc_comp : qty_index;
474  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
475  bool ext_dir_on_ylo = ( ((bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir) ||
476  (bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir_prim))
477  && j == dom_lo.y);
478  bool ext_dir_on_yhi = ( ((bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir) ||
479  (bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir_prim))
480  && j == dom_hi.y+1);
481 
482  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
483  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
484 
485  if (ext_dir_on_ylo) {
486  yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
487  + 3. * cell_prim(i, j , k, prim_index)
488  - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
489  } else if (ext_dir_on_yhi) {
490  yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
491  - 3. * cell_prim(i, j-1, k, prim_index)
492  + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
493  } else {
494  yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
495  }
496  });
497  ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
498  {
499  const int qty_index = start_comp + n;
500  const int prim_index = qty_index - 1;
501 
502  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
503  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
504 
505  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
506  BCVars::RhoScalar_bc_comp : qty_index;
507  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
508  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
509  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
510  && k == dom_lo.z);
511  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
512  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
513  && k == dom_hi.z+1);
514  bool most_on_zlo = ( use_most && exp_most &&
515  (bc_ptr[bc_comp].lo(2) == ERFBCType::foextrap) &&
516  k == dom_lo.z);
517 
518  if (ext_dir_on_zlo) {
519  zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
520  + 3. * cell_prim(i, j, k , prim_index)
521  - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
522  } else if (ext_dir_on_zhi) {
523  zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
524  - 3. * cell_prim(i, j, k-1, prim_index)
525  + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
526  } else if (most_on_zlo && (qty_index == RhoTheta_comp)) {
527  zflux(i,j,k,qty_index) = hfx_z(i,j,0);
528  } else if (most_on_zlo && (qty_index == RhoQ1_comp)) {
529  zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
530  } else {
531  zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
532  }
533 
534  if (qty_index == RhoTheta_comp) {
535  if (!most_on_zlo) {
536  hfx_z(i,j,k) = zflux(i,j,k,qty_index);
537  }
538  } else if (qty_index == RhoQ1_comp) {
539  if (!most_on_zlo) {
540  qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
541  }
542  } else if (qty_index == RhoQ2_comp) {
543  qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
544  }
545  });
546  } else {
547  // with MolecDiffType::Constant or None
548  // without an LES/PBL model
549  ParallelFor(xbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
550  {
551  const int qty_index = start_comp + n;
552  const int prim_index = qty_index - 1;
553 
554  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
555  BCVars::RhoScalar_bc_comp : qty_index;
556  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
557  bool ext_dir_on_xlo = ( ((bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir) ||
558  (bc_ptr[bc_comp].lo(0) == ERFBCType::ext_dir_prim))
559  && i == dom_lo.x);
560  bool ext_dir_on_xhi = ( ((bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir) ||
561  (bc_ptr[bc_comp].hi(0) == ERFBCType::ext_dir_prim))
562  && i == dom_hi.x+1);
563 
564  Real rhoAlpha = d_alpha_eff[prim_index];
565 
566  if (ext_dir_on_xlo) {
567  xflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i-1, j, k, prim_index)
568  + 3. * cell_prim(i , j, k, prim_index)
569  - (1./3.) * cell_prim(i+1, j, k, prim_index) ) * dx_inv;
570  } else if (ext_dir_on_xhi) {
571  xflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i , j, k, prim_index)
572  - 3. * cell_prim(i-1, j, k, prim_index)
573  + (1./3.) * cell_prim(i-2, j, k, prim_index) ) * dx_inv;
574  } else {
575  xflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i-1, j, k, prim_index)) * dx_inv * mf_u(i,j,0);
576  }
577  });
578  ParallelFor(ybx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
579  {
580  const int qty_index = start_comp + n;
581  const int prim_index = qty_index - 1;
582 
583  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
584  BCVars::RhoScalar_bc_comp : qty_index;
585  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
586  bool ext_dir_on_ylo = ( ((bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir) ||
587  (bc_ptr[bc_comp].lo(1) == ERFBCType::ext_dir_prim))
588  && j == dom_lo.y);
589  bool ext_dir_on_yhi = ( ((bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir) ||
590  (bc_ptr[bc_comp].hi(1) == ERFBCType::ext_dir_prim))
591  && j == dom_hi.y+1);
592 
593  Real rhoAlpha = d_alpha_eff[prim_index];
594 
595  if (ext_dir_on_ylo) {
596  yflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j-1, k, prim_index)
597  + 3. * cell_prim(i, j , k, prim_index)
598  - (1./3.) * cell_prim(i, j+1, k, prim_index) ) * dy_inv;
599  } else if (ext_dir_on_yhi) {
600  yflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j , k, prim_index)
601  - 3. * cell_prim(i, j-1, k, prim_index)
602  + (1./3.) * cell_prim(i, j-2, k, prim_index) ) * dy_inv;
603  } else {
604  yflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j-1, k, prim_index)) * dy_inv * mf_v(i,j,0);
605  }
606  });
607  ParallelFor(zbx, num_comp,[=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
608  {
609  const int qty_index = start_comp + n;
610  const int prim_index = qty_index - 1;
611 
612  Real rhoAlpha = d_alpha_eff[prim_index];
613 
614  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
615  BCVars::RhoScalar_bc_comp : qty_index;
616  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
617  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
618  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
619  && k == dom_lo.z);
620  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
621  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
622  && k == dom_hi.z+1);
623  bool most_on_zlo = ( use_most && exp_most &&
624  (bc_ptr[BCVars::cons_bc+qty_index].lo(2) == ERFBCType::foextrap) &&
625  k == dom_lo.z);
626 
627  if (ext_dir_on_zlo) {
628  zflux(i,j,k,qty_index) = -rhoAlpha * ( -(8./3.) * cell_prim(i, j, k-1, prim_index)
629  + 3. * cell_prim(i, j, k , prim_index)
630  - (1./3.) * cell_prim(i, j, k+1, prim_index) ) * dz_inv;
631  } else if (ext_dir_on_zhi) {
632  zflux(i,j,k,qty_index) = -rhoAlpha * ( (8./3.) * cell_prim(i, j, k , prim_index)
633  - 3. * cell_prim(i, j, k-1, prim_index)
634  + (1./3.) * cell_prim(i, j, k-2, prim_index) ) * dz_inv;
635  } else if (most_on_zlo && (qty_index == RhoTheta_comp)) {
636  zflux(i,j,k,qty_index) = hfx_z(i,j,0);
637  } else if (most_on_zlo && (qty_index == RhoQ1_comp)) {
638  zflux(i,j,k,qty_index) = qfx1_z(i,j,0);
639  } else {
640  zflux(i,j,k,qty_index) = -rhoAlpha * (cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index)) * dz_inv;
641  }
642 
643  if (qty_index == RhoTheta_comp) {
644  if (!most_on_zlo) {
645  hfx_z(i,j,k) = zflux(i,j,k,qty_index);
646  }
647  } else if (qty_index == RhoQ1_comp) {
648  if (!most_on_zlo) {
649  qfx1_z(i,j,k) = zflux(i,j,k,qty_index);
650  }
651  } else if (qty_index == RhoQ2_comp) {
652  qfx2_z(i,j,k) = zflux(i,j,k,qty_index);
653  }
654  });
655  }
656 
657  // Use fluxes to compute RHS
658  for (int n(0); n < num_comp; ++n)
659  {
660  int qty_index = start_comp + n;
661  ParallelFor(bx,[=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
662  {
663 
664  cell_rhs(i,j,k,qty_index) -= (xflux(i+1,j ,k ,qty_index) - xflux(i, j, k, qty_index)) * dx_inv * mf_m(i,j,0) // Diffusive flux in x-dir
665  +(yflux(i ,j+1,k ,qty_index) - yflux(i, j, k, qty_index)) * dy_inv * mf_m(i,j,0) // Diffusive flux in y-dir
666  +(zflux(i ,j ,k+1,qty_index) - zflux(i, j, k, qty_index)) * dz_inv; // Diffusive flux in z-dir
667  });
668  }
669 
670  // Using Deardorff (see Sullivan et al 1994)
671  //
672  // Note: At this point, the thermal diffusivity ("Khv" field in ERF), the
673  // subgrid heat flux ("hfx_z" here), and the subgrid dissipation
674  // ("diss" here) have been updated by ComputeTurbulentViscosityLES --
675  // at the beginning of each timestep.
676  // The strain rate magnitude is updated at the beginning of the first
677  // RK stage only, therefore the shear production term also does not
678  // change between RK stages.
679  // The surface heat flux hfx_z(i,j,-1) is updated in MOSTStress at
680  // each RK stage if using the ERF_EXPLICIT_MOST_STRESS path, but that
681  // does not change the buoyancy production term here.
682  if (l_use_ddorf && (start_comp <= RhoKE_comp) && (end_comp >=RhoKE_comp)) {
683  int qty_index = RhoKE_comp;
684  ParallelFor(bx,[=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
685  {
686  // Add Buoyancy Source
687  // where the SGS buoyancy flux tau_{theta,i} = -KH * dtheta/dx_i,
688  // such that for dtheta/dz < 0, there is a positive (upward) heat
689  // flux; the TKE buoyancy production is then
690  // B = g/theta_0 * tau_{theta,w}
691  // for a dry atmosphere.
692  // TODO: To account for moisture, the Brunt-Vaisala frequency,
693  // N^2 = g[1/theta * dtheta/dz + ...]
694  // **should** be a function of the water vapor and total water
695  // mixing ratios, depending on whether conditions are saturated or
696  // not (see the WRF model description, Skamarock et al 2019).
697  cell_rhs(i,j,k,qty_index) += l_abs_g * l_inv_theta0 * hfx_z(i,j,k);
698 
699  // TKE shear production
700  // P = -tau_ij * S_ij = 2 * mu_turb * S_ij * S_ij
701  // Note: This assumes that the horizontal and vertical diffusivities
702  // of momentum are equal
703  cell_rhs(i,j,k,qty_index) += 2.0*mu_turb(i,j,k,EddyDiff::Mom_v) * SmnSmn_a(i,j,k);
704 
705  // TKE dissipation
706  cell_rhs(i,j,k,qty_index) -= diss(i,j,k);
707  });
708  }
709 
710  // Using PBL
711  if (l_use_mynn && start_comp <= RhoKE_comp && end_comp >=RhoKE_comp) {
712  int qty_index = RhoKE_comp;
713  auto pbl_mynn_B1_l = turbChoice.pbl_mynn.B1;
714 
715  const int rhoqv_comp = solverChoice.RhoQv_comp;
716  const int rhoqr_comp = solverChoice.RhoQr_comp;
717 
718  ParallelFor(bx,[=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
719  {
720  bool c_ext_dir_on_zlo = ( (bc_ptr[BCVars::cons_bc].lo(2) == ERFBCType::ext_dir) );
721  bool c_ext_dir_on_zhi = ( (bc_ptr[BCVars::cons_bc].lo(5) == ERFBCType::ext_dir) );
722  bool u_ext_dir_on_zlo = ( (bc_ptr[BCVars::xvel_bc].lo(2) == ERFBCType::ext_dir) );
723  bool u_ext_dir_on_zhi = ( (bc_ptr[BCVars::xvel_bc].lo(5) == ERFBCType::ext_dir) );
724  bool v_ext_dir_on_zlo = ( (bc_ptr[BCVars::yvel_bc].lo(2) == ERFBCType::ext_dir) );
725  bool v_ext_dir_on_zhi = ( (bc_ptr[BCVars::yvel_bc].lo(5) == ERFBCType::ext_dir) );
726 
727  // This computes shear production, buoyancy production, and dissipation terms only.
728  cell_rhs(i, j, k, qty_index) += ComputeQKESourceTerms(i,j,k,u,v,cell_data,cell_prim,
729  mu_turb,cellSizeInv,domain,
730  pbl_mynn_B1_l,tm_arr(i,j,0),
731  rhoqv_comp, rhoqr_comp,
732  c_ext_dir_on_zlo, c_ext_dir_on_zhi,
733  u_ext_dir_on_zlo, u_ext_dir_on_zhi,
734  v_ext_dir_on_zlo, v_ext_dir_on_zhi,
735  use_most);
736  });
737  }
738 }
void DiffusionSrcForState_N(const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &exp_most, const Array4< const Real > &u, const Array4< const Real > &v, const Array4< const Real > &cell_data, const Array4< const Real > &cell_prim, const Array4< Real > &cell_rhs, const Array4< Real > &xflux, const Array4< Real > &yflux, const Array4< Real > &zflux, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_m, const Array4< const Real > &mf_u, const Array4< const Real > &mf_v, Array4< Real > &hfx_z, Array4< Real > &qfx1_z, Array4< Real > &qfx2_z, Array4< Real > &diss, const Array4< const Real > &mu_turb, const SolverChoice &solverChoice, const int level, const Array4< const Real > &tm_arr, const GpuArray< Real, AMREX_SPACEDIM > grav_gpu, const BCRec *bc_ptr, const bool use_most)
Definition: ERF_DiffusionSrcForState_N.cpp:40
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define PrimQ1_comp
Definition: ERF_IndexDefines.H:53
#define NPRIMVAR_max
Definition: ERF_IndexDefines.H:33
#define PrimQ2_comp
Definition: ERF_IndexDefines.H:54
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define PrimQ4_comp
Definition: ERF_IndexDefines.H:56
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define PrimQ6_comp
Definition: ERF_IndexDefines.H:58
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define PrimQ5_comp
Definition: ERF_IndexDefines.H:57
#define PrimTheta_comp
Definition: ERF_IndexDefines.H:50
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define PrimScalar_comp
Definition: ERF_IndexDefines.H:52
#define PrimQ3_comp
Definition: ERF_IndexDefines.H:55
#define RhoKE_comp
Definition: ERF_IndexDefines.H:38
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real ComputeQKESourceTerms(int i, int j, int k, const amrex::Array4< const amrex::Real > &uvel, const amrex::Array4< const amrex::Real > &vvel, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &K_turb, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Box &domain, amrex::Real pbl_mynn_B1_l, const amrex::Real theta_mean, const int RhoQv_comp, const int RhoQr_comp, bool c_ext_dir_on_zlo, bool c_ext_dir_on_zhi, bool u_ext_dir_on_zlo, bool u_ext_dir_on_zhi, bool v_ext_dir_on_zlo, bool v_ext_dir_on_zhi, const bool use_most=false, const amrex::Real met_h_zeta=1.0)
Definition: ERF_PBLModels.H:157
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:80
@ yvel_bc
Definition: ERF_IndexDefines.H:89
@ cons_bc
Definition: ERF_IndexDefines.H:76
@ xvel_bc
Definition: ERF_IndexDefines.H:88
@ foextrap
Definition: ERF_IndexDefines.H:179
@ ext_dir
Definition: ERF_IndexDefines.H:180
@ ext_dir_prim
Definition: ERF_IndexDefines.H:182
@ Theta_v
Definition: ERF_IndexDefines.H:157
@ Scalar_v
Definition: ERF_IndexDefines.H:159
@ Q_v
Definition: ERF_IndexDefines.H:160
@ Q_h
Definition: ERF_IndexDefines.H:155
@ Mom_v
Definition: ERF_IndexDefines.H:156
@ Scalar_h
Definition: ERF_IndexDefines.H:154
@ KE_v
Definition: ERF_IndexDefines.H:158
@ Theta_h
Definition: ERF_IndexDefines.H:152
@ KE_h
Definition: ERF_IndexDefines.H:153
Definition: ERF_DiffStruct.H:19
amrex::Real rhoAlpha_C
Definition: ERF_DiffStruct.H:87
amrex::Real rhoAlpha_T
Definition: ERF_DiffStruct.H:86
amrex::Real alpha_T
Definition: ERF_DiffStruct.H:79
amrex::Real alpha_C
Definition: ERF_DiffStruct.H:80
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:76
amrex::Real B1
Definition: ERF_MYNNStruct.H:43
int RhoQr_comp
Definition: ERF_DataStruct.H:649
DiffChoice diffChoice
Definition: ERF_DataStruct.H:538
int RhoQv_comp
Definition: ERF_DataStruct.H:643
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:540
Definition: ERF_TurbStruct.H:29
MYNNLevel25 pbl_mynn
Definition: ERF_TurbStruct.H:194
PBLType pbl_type
Definition: ERF_TurbStruct.H:192
LESType les_type
Definition: ERF_TurbStruct.H:169
amrex::Real theta_ref
Definition: ERF_TurbStruct.H:189
Here is the call graph for this function: