ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_DiffusionSrcForState_T.cpp File Reference
#include "ERF_Diffusion.H"
#include "ERF_EddyViscosity.H"
#include "ERF_TerrainMetrics.H"
#include "ERF_PBLModels.H"
#include "ERF_SetupDiff.H"
#include "ERF_AddTKESources.H"
#include "ERF_AddQKESources.H"
Include dependency graph for ERF_DiffusionSrcForState_T.cpp:

Functions

void DiffusionSrcForState_T (const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &rotate, 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 Array4< const Real > &z_nd, const Array4< const Real > &z_cc, const Array4< const Real > &ax, const Array4< const Real > &ay, const Array4< const Real > &, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_ux, const Array4< const Real > &mf_vx, const Array4< const Real > &mf_my, const Array4< const Real > &mf_uy, const Array4< const Real > &mf_vy, Array4< Real > &hfx_x, Array4< Real > &hfx_y, Array4< Real > &hfx_z, Array4< Real > &qfx1_x, Array4< Real > &qfx1_y, 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_SurfLayer, const Real implicit_fac)
 

Function Documentation

◆ DiffusionSrcForState_T()

void DiffusionSrcForState_T ( const Box &  bx,
const Box &  domain,
int  start_comp,
int  num_comp,
const bool &  rotate,
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 Array4< const Real > &  z_nd,
const Array4< const Real > &  z_cc,
const Array4< const Real > &  ax,
const Array4< const Real > &  ay,
const Array4< const Real > &  ,
const Array4< const Real > &  detJ,
const GpuArray< Real, AMREX_SPACEDIM > &  cellSizeInv,
const Array4< const Real > &  SmnSmn_a,
const Array4< const Real > &  mf_mx,
const Array4< const Real > &  mf_ux,
const Array4< const Real > &  mf_vx,
const Array4< const Real > &  mf_my,
const Array4< const Real > &  mf_uy,
const Array4< const Real > &  mf_vy,
Array4< Real > &  hfx_x,
Array4< Real > &  hfx_y,
Array4< Real > &  hfx_z,
Array4< Real > &  qfx1_x,
Array4< Real > &  qfx1_y,
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_SurfLayer,
const Real  implicit_fac 
)

Function for computing the scalar RHS for diffusion operator with terrain-fitted coordinates.

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]z_ndphysical z height
[in]detJJacobian determinant
[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_SurfLayerwhether we have turned on subgrid diffusion
[in]implicit_fac– factor of implicitness for vertical differences only
85 {
86  BL_PROFILE_VAR("DiffusionSrcForState_T()",DiffusionSrcForState_T);
87 
88  const Real explicit_fac = one - implicit_fac;
89 
90 #include "ERF_SetupDiff.H"
91  Real l_abs_g = std::abs(grav_gpu[2]);
92 
93  const Real dz_inv = cellSizeInv[2];
94 
95  // We need to grow these boxes in the vertical direction when tiling so that we can access xflux and yflux
96  // to modify zflux
97  Box xbx_g1(xbx); Box ybx_g1(ybx);
98  if (xbx_g1.smallEnd(2) != dom_lo.z) xbx_g1.growLo(2,1);
99  if (ybx_g1.smallEnd(2) != dom_lo.z) ybx_g1.growLo(2,1);
100  if (xbx_g1.bigEnd(2) != dom_hi.z) xbx_g1.growHi(2,1);
101  if (ybx_g1.bigEnd(2) != dom_hi.z) ybx_g1.growHi(2,1);
102 
103  for (int n(0); n<num_comp; ++n) {
104  const int qty_index = start_comp + n;
105 
106  // Constant alpha & Turb model
107  if (l_consA && l_turb) {
108  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
109  {
110  const int prim_index = qty_index - 1;
111  const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
112 
113  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i-1, j, k, Rho_comp) );
114  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
115  rhoAlpha += myhalf * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_scal_index])
116  + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_scal_index]) );
117 
118  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
119 
120  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
121  BCVars::RhoScalar_bc_comp : qty_index;
122  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
123 
124  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
125 
126  Real idz_hi = one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
127  Real idz_lo = one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
128  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
129  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
130  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
131 
132  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
133  xflux(i,j,k) = hfx_x(i,j,0);
134  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
135  xflux(i,j,k) = qfx1_x(i,j,0);
136  } else {
137  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
138  }
139  });
140  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
141  {
142  const int prim_index = qty_index - 1;
143  const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
144 
145  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
146  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
147  rhoAlpha += myhalf * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_scal_index])
148  + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_scal_index]) );
149 
150  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
151 
152  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
153  BCVars::RhoScalar_bc_comp : qty_index;
154  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
155  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
156 
157  Real idz_hi = one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
158  Real idz_lo = one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
159  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
160  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
161  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
162 
163  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
164  yflux(i,j,k) = hfx_y(i,j,0);
165  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
166  yflux(i,j,k) = qfx1_y(i,j,0);
167  } else {
168  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
169  }
170  });
171  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
172  {
173  const int prim_index = qty_index - 1;
174  const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
175 
176  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
177  Real rhoAlpha = rhoFace * d_alpha_eff[prim_scal_index];
178  rhoAlpha += myhalf * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_scal_index])
179  + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_scal_index]) );
180 
181  Real GradCz;
182  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
183  BCVars::RhoScalar_bc_comp : qty_index;
184  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
185  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
186  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
187  && k == dom_lo.z);
188  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
189  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim) )
190  && k == dom_hi.z+1);
191  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
192 
193  if (ext_dir_on_zlo) {
194  // Third order stencil with variable dz
195  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
196  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
197  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
198  Real idz0 = one / dz0;
199  Real f = (dz1 / dz0) + two;
200  Real f2 = f*f;
201  Real c3 = two / (f - f2);
202  Real c2 = -f2*c3;
203  Real c1 = -(one-f2)*c3;
204  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
205  + c2 * cell_prim(i, j, k , prim_index)
206  + c3 * cell_prim(i, j, k+1, prim_index) );
207  } else if (ext_dir_on_zhi) {
208  // Third order stencil with variable dz
209  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
210  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
211  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
212  Real idz0 = one / dz0;
213  Real f = (dz1 / dz0) + two;
214  Real f2 = f*f;
215  Real c3 = two / (f - f2);
216  Real c2 = -f2*c3;
217  Real c1 = -(one-f2)*c3;
218  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
219  + c2 * cell_prim(i, j, k-1, prim_index)
220  + c3 * cell_prim(i, j, k-2, prim_index) ) );
221  } else {
222  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
223  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
224  }
225 
226  if (SurfLayer_on_zlo) {
227  if (qty_index == RhoTheta_comp) {
228  zflux(i,j,k) = hfx_z(i,j,0);
229  } else if (qty_index == RhoQ1_comp) {
230  zflux(i,j,k) = qfx1_z(i,j,0);
231  } else {
232  zflux(i,j,k) = zero;
233  }
234  } else {
235  zflux(i,j,k) = -rhoAlpha * GradCz;
236  }
237 
238  if (qty_index == RhoTheta_comp) {
239  if (!SurfLayer_on_zlo) {
240  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
241  }
242  } else if (qty_index == RhoQ1_comp) {
243  if (!SurfLayer_on_zlo) {
244  qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
245  }
246  } else if (qty_index == RhoQ2_comp) {
247  qfx2_z(i,j,k) = zflux(i,j,k);
248  }
249  });
250  // Constant rho*alpha & Turb model
251  } else if (l_turb) {
252  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
253  {
254  const int prim_index = qty_index - 1;
255 
256  Real rhoAlpha = d_alpha_eff[prim_index];
257  rhoAlpha += myhalf * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_index])
258  + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_index]) );
259 
260  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
261 
262  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
263  BCVars::RhoScalar_bc_comp : qty_index;
264  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
265  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
266 
267  Real idz_hi = one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
268  Real idz_lo = one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
269  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
270  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
271  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
272 
273  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
274  xflux(i,j,k) = hfx_x(i,j,0);
275  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
276  xflux(i,j,k) = qfx1_x(i,j,0);
277  } else {
278  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
279  }
280  });
281  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
282  {
283  const int prim_index = qty_index - 1;
284 
285  Real rhoAlpha = d_alpha_eff[prim_index];
286  rhoAlpha += myhalf * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_index])
287  + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_index]) );
288 
289  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
290 
291  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
292  BCVars::RhoScalar_bc_comp : qty_index;
293  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
294  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
295 
296  Real idz_hi = one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
297  Real idz_lo = one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
298  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
299  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
300  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
301 
302  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
303  yflux(i,j,k) = hfx_y(i,j,0);
304  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
305  yflux(i,j,k) = qfx1_y(i,j,0);
306  } else {
307  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
308  }
309  });
310  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
311  {
312  const int prim_index = qty_index - 1;
313 
314  Real rhoAlpha = d_alpha_eff[prim_index];
315  rhoAlpha += myhalf * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_index])
316  + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_index]) );
317 
318  Real GradCz;
319  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
320  BCVars::RhoScalar_bc_comp : qty_index;
321  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
322  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
323  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
324  && k == dom_lo.z);
325  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
326  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
327  && k == dom_hi.z+1);
328  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
329 
330  if (ext_dir_on_zlo) {
331  // Third order stencil with variable dz
332  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
333  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
334  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
335  Real idz0 = one / dz0;
336  Real f = (dz1 / dz0) + two;
337  Real f2 = f*f;
338  Real c3 = two / (f - f2);
339  Real c2 = -f2*c3;
340  Real c1 = -(one-f2)*c3;
341  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
342  + c2 * cell_prim(i, j, k , prim_index)
343  + c3 * cell_prim(i, j, k+1, prim_index) );
344  } else if (ext_dir_on_zhi) {
345  // Third order stencil with variable dz
346  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
347  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
348  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
349  Real idz0 = one / dz0;
350  Real f = (dz1 / dz0) + two;
351  Real f2 = f*f;
352  Real c3 = two / (f - f2);
353  Real c2 = -f2*c3;
354  Real c1 = -(one-f2)*c3;
355  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
356  + c2 * cell_prim(i, j, k-1, prim_index)
357  + c3 * cell_prim(i, j, k-2, prim_index) ) );
358  } else {
359  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
360  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
361  }
362 
363  if (SurfLayer_on_zlo) {
364  if (qty_index == RhoTheta_comp) {
365  zflux(i,j,k) = hfx_z(i,j,0);
366  } else if (qty_index == RhoQ1_comp) {
367  zflux(i,j,k) = qfx1_z(i,j,0);
368  } else {
369  zflux(i,j,k) = zero;
370  }
371  } else {
372  zflux(i,j,k) = -rhoAlpha * GradCz;
373  }
374 
375  if (qty_index == RhoTheta_comp) {
376  if (!SurfLayer_on_zlo) {
377  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
378  }
379  } else if (qty_index == RhoQ1_comp) {
380  if (!SurfLayer_on_zlo) {
381  qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
382  }
383  } else if (qty_index == RhoQ2_comp) {
384  qfx2_z(i,j,k) = zflux(i,j,k);
385  }
386 
387  });
388  // Constant alpha & no LES/PBL model
389  } else if(l_consA) {
390  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
391  {
392  const int prim_index = qty_index - 1;
393 
394  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i-1, j, k, Rho_comp) );
395  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
396 
397  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
398 
399  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
400  BCVars::RhoScalar_bc_comp : qty_index;
401  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
402  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
403 
404  Real idz_hi = one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
405  Real idz_lo = one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
406  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
407  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
408  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
409 
410  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
411  xflux(i,j,k) = hfx_x(i,j,0);
412  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
413  xflux(i,j,k) = qfx1_x(i,j,0);
414  } else {
415  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
416  }
417  });
418  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
419  {
420  const int prim_index = qty_index - 1;
421 
422  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
423  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
424 
425  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
426 
427  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
428  BCVars::RhoScalar_bc_comp : qty_index;
429  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
430  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
431 
432  Real idz_hi = one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
433  Real idz_lo = one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
434  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
435  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
436  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
437 
438  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
439  yflux(i,j,k) = hfx_y(i,j,0);
440  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
441  yflux(i,j,k) = qfx1_y(i,j,0);
442  } else {
443  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
444  }
445  });
446  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
447  {
448  const int prim_index = qty_index - 1;
449 
450  Real rhoFace = myhalf * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
451  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
452 
453  Real GradCz;
454  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
455  BCVars::RhoScalar_bc_comp : qty_index;
456  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
457  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
458  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
459  && k == dom_lo.z);
460  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
461  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
462  && k == dom_hi.z+1);
463  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
464 
465  if (ext_dir_on_zlo) {
466  // Third order stencil with variable dz
467  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
468  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
469  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
470  Real idz0 = one / dz0;
471  Real f = (dz1 / dz0) + two;
472  Real f2 = f*f;
473  Real c3 = two / (f - f2);
474  Real c2 = -f2*c3;
475  Real c1 = -(one-f2)*c3;
476  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
477  + c2 * cell_prim(i, j, k , prim_index)
478  + c3 * cell_prim(i, j, k+1, prim_index) );
479  } else if (ext_dir_on_zhi) {
480  // Third order stencil with variable dz
481  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
482  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
483  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
484  Real idz0 = one / dz0;
485  Real f = (dz1 / dz0) + two;
486  Real f2 = f*f;
487  Real c3 = two / (f - f2);
488  Real c2 = -f2*c3;
489  Real c1 = -(one-f2)*c3;
490  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
491  + c2 * cell_prim(i, j, k-1, prim_index)
492  + c3 * cell_prim(i, j, k-2, prim_index) ) );
493  } else {
494  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
495  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
496  }
497 
498  if (SurfLayer_on_zlo) {
499  if (qty_index == RhoTheta_comp) {
500  zflux(i,j,k) = hfx_z(i,j,0);
501  } else if (qty_index == RhoQ1_comp) {
502  zflux(i,j,k) = qfx1_z(i,j,0);
503  } else {
504  zflux(i,j,k) = zero;
505  }
506  } else {
507  zflux(i,j,k) = -rhoAlpha * GradCz;
508  }
509 
510  if (qty_index == RhoTheta_comp) {
511  if (!SurfLayer_on_zlo) {
512  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
513  }
514  } else if (qty_index == RhoQ1_comp) {
515  if (!SurfLayer_on_zlo) {
516  qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
517  }
518  } else if (qty_index == RhoQ2_comp) {
519  qfx2_z(i,j,k) = zflux(i,j,k);
520  }
521  });
522  // Constant rho*alpha & no LES/PBL model
523  } else {
524  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
525  {
526  const int prim_index = qty_index - 1;
527 
528  Real rhoAlpha = d_alpha_eff[prim_index];
529 
530  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
531 
532  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
533  BCVars::RhoScalar_bc_comp : qty_index;
534  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
535  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
536 
537  Real idz_hi = one / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
538  Real idz_lo = one / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
539  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
540  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
541  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
542 
543  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
544  xflux(i,j,k) = hfx_x(i,j,0);
545  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
546  xflux(i,j,k) = qfx1_x(i,j,0);
547  } else {
548  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
549  }
550  });
551  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
552  {
553  const int prim_index = qty_index - 1;
554 
555  Real rhoAlpha = d_alpha_eff[prim_index];
556 
557  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
558 
559  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
560  BCVars::RhoScalar_bc_comp : qty_index;
561  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
562  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
563 
564  Real idz_hi = one / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
565  Real idz_lo = one / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
566  Real GradCz = myhalf * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
567  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
568  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
569 
570  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
571  yflux(i,j,k) = hfx_y(i,j,0);
572  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
573  yflux(i,j,k) = qfx1_y(i,j,0);
574  } else {
575  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
576  }
577  });
578  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
579  {
580  const int prim_index = qty_index - 1;
581 
582  Real rhoAlpha = d_alpha_eff[prim_index];
583 
584 
585  Real GradCz;
586  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
587  BCVars::RhoScalar_bc_comp : qty_index;
588  if (bc_comp > BCVars::RhoScalar_bc_comp) bc_comp -= (NSCALARS-1);
589  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
590  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
591  && k == dom_lo.z);
592  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir) ||
593  (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim))
594  && k == dom_hi.z+1);
595  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
596 
597  if (ext_dir_on_zlo) {
598  // Third order stencil with variable dz
599  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
600  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
601  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
602  Real idz0 = one / dz0;
603  Real f = (dz1 / dz0) + two;
604  Real f2 = f*f;
605  Real c3 = two / (f - f2);
606  Real c2 = -f2*c3;
607  Real c1 = -(one-f2)*c3;
608  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
609  + c2 * cell_prim(i, j, k , prim_index)
610  + c3 * cell_prim(i, j, k+1, prim_index) );
611  } else if (ext_dir_on_zhi) {
612  // Third order stencil with variable dz
613  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
614  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
615  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
616  Real idz0 = one / dz0;
617  Real f = (dz1 / dz0) + two;
618  Real f2 = f*f;
619  Real c3 = two / (f - f2);
620  Real c2 = -f2*c3;
621  Real c1 = -(one-f2)*c3;
622  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
623  + c2 * cell_prim(i, j, k-1, prim_index)
624  + c3 * cell_prim(i, j, k-2, prim_index) ) );
625  } else {
626  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
627  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
628  }
629 
630  if (SurfLayer_on_zlo) {
631  if (qty_index == RhoTheta_comp) {
632  zflux(i,j,k) = hfx_z(i,j,0);
633  } else if (qty_index == RhoQ1_comp) {
634  zflux(i,j,k) = qfx1_z(i,j,0);
635  } else {
636  zflux(i,j,k) = zero;
637  }
638  } else {
639  zflux(i,j,k) = -rhoAlpha * GradCz;
640  }
641 
642  if (qty_index == RhoTheta_comp) {
643  if (!SurfLayer_on_zlo) {
644  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
645  }
646  } else if (qty_index == RhoQ1_comp) {
647  if (!SurfLayer_on_zlo) {
648  qfx1_z(i,j,k) = zflux(i,j,k) * explicit_fac;
649  }
650  } else if (qty_index == RhoQ2_comp) {
651  qfx2_z(i,j,k) = zflux(i,j,k);
652  }
653  });
654  }
655 
656  // NOTE: With terrain, we implicitly treat the leading order vertical gradient (no metric terms)
657  // This allows us to do semi-implicit discretization of the vertical diffusive terms
658  if (qty_index == RhoTheta_comp ||
659  qty_index == RhoQ1_comp) {
660  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
661  {
662  zflux(i,j,k) *= explicit_fac;
663  });
664  }
665 
666  //-----------------------------------------------------------------------------------
667  //
668  // Modify fluxes by terrain and use fluxes to compute RHS
669  //
670  // Note that we combine all of these operations in order to keep this section
671  // of the loop tiling-safe.
672  //-----------------------------------------------------------------------------------
673  ParallelFor(bx,[=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
674  {
675  Real xfluxbar_lo, yfluxbar_lo;
676  Real met_h_xi_lo = Compute_h_xi_AtKface (i,j,k ,cellSizeInv,z_nd);
677  Real met_h_eta_lo = Compute_h_eta_AtKface(i,j,k ,cellSizeInv,z_nd);
678  if (k == dom_lo.z) {
679  Real xfluxlo = myhalf * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
680  Real xfluxhi = myhalf * ( xflux(i,j,k+1) + xflux(i+1,j,k+1) );
681  xfluxbar_lo = Real(1.5)*xfluxlo - myhalf*xfluxhi;
682 
683  Real yfluxlo = myhalf * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
684  Real yfluxhi = myhalf * ( yflux(i,j,k+1) + yflux(i,j+1,k+1) );
685  yfluxbar_lo = Real(1.5)*yfluxlo - myhalf*yfluxhi;
686  } else {
687  xfluxbar_lo = fourth * ( xflux(i,j,k ) + xflux(i+1,j ,k )
688  + xflux(i,j,k-1) + xflux(i+1,j ,k-1) );
689  yfluxbar_lo = fourth * ( yflux(i,j,k ) + yflux(i ,j+1,k )
690  + yflux(i,j,k-1) + yflux(i ,j+1,k-1) );
691  }
692 
693  Real xfluxbar_hi, yfluxbar_hi;
694  Real met_h_xi_hi = Compute_h_xi_AtKface (i,j,k+1,cellSizeInv,z_nd);
695  Real met_h_eta_hi = Compute_h_eta_AtKface(i,j,k+1,cellSizeInv,z_nd);
696  if (k == dom_hi.z) {
697  Real xfluxlo = myhalf * ( xflux(i,j,k-1) + xflux(i+1,j,k-1) );
698  Real xfluxhi = myhalf * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
699  xfluxbar_hi = Real(1.5)*xfluxhi - myhalf*xfluxlo;
700 
701  Real yfluxlo = myhalf * ( yflux(i,j,k-1) + yflux(i,j+1,k-1) );
702  Real yfluxhi = myhalf * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
703  yfluxbar_hi = Real(1.5)*yfluxhi - myhalf*yfluxlo;
704  } else {
705  xfluxbar_hi = fourth * ( xflux(i,j,k+1) + xflux(i+1,j ,k+1)
706  + xflux(i,j,k ) + xflux(i+1,j ,k ) );
707  yfluxbar_hi = fourth * ( yflux(i,j,k+1) + yflux(i ,j+1,k+1)
708  + yflux(i,j,k ) + yflux(i ,j+1,k ) );
709  }
710 
711  // Allow semi-implicit discretization of the vertical diffusive terms
712  Real zflux_lo;
713  if ( use_SurfLayer &&
714  k == dom_lo.z &&
715  qty_index != RhoTheta_comp &&
716  qty_index != RhoQ1_comp ) {
717  zflux_lo = zero;
718  } else {
719  zflux_lo = zflux(i,j,k )
720  - met_h_xi_lo * mf_mx(i,j,0) * xfluxbar_lo
721  - met_h_eta_lo * mf_my(i,j,0) * yfluxbar_lo;
722  }
723  Real zflux_hi = zflux(i,j,k+1)
724  - met_h_xi_hi * mf_mx(i,j,0) * xfluxbar_hi
725  - met_h_eta_hi * mf_my(i,j,0) * yfluxbar_hi;
726 
727  Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
728  Real stateContrib = ( xflux(i+1,j ,k ) * ax(i+1,j,k) / mf_uy(i+1,j,0)
729  -xflux(i ,j ,k ) * ax(i ,j,k) / mf_uy(i ,j,0) ) * dx_inv * mfsq // Diffusive flux in x-dir
730  +( yflux(i ,j+1,k ) * ay(i,j+1,k) / mf_vx(i,j+1,0)
731  -yflux(i ,j ,k ) * ay(i,j ,k) / mf_vx(i,j ,0) ) * dy_inv * mfsq // Diffusive flux in y-dir
732  +( zflux_hi - zflux_lo) * dz_inv; // Diffusive flux in z-dir
733 
734  stateContrib /= detJ(i,j,k);
735 
736  cell_rhs(i,j,k,qty_index) -= stateContrib;
737  });
738  } // n
739 
740 #include "ERF_AddTKESources.H"
741 #include "ERF_AddQKESources.H"
742 }
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real fourth
Definition: ERF_Constants.H:12
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
void DiffusionSrcForState_T(const Box &bx, const Box &domain, int start_comp, int num_comp, const bool &rotate, 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 Array4< const Real > &z_nd, const Array4< const Real > &z_cc, const Array4< const Real > &ax, const Array4< const Real > &ay, const Array4< const Real > &, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &SmnSmn_a, const Array4< const Real > &mf_mx, const Array4< const Real > &mf_ux, const Array4< const Real > &mf_vx, const Array4< const Real > &mf_my, const Array4< const Real > &mf_uy, const Array4< const Real > &mf_vy, Array4< Real > &hfx_x, Array4< Real > &hfx_y, Array4< Real > &hfx_z, Array4< Real > &qfx1_x, Array4< Real > &qfx1_y, 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_SurfLayer, const Real implicit_fac)
Definition: ERF_DiffusionSrcForState_T.cpp:44
#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 NSCALARS
Definition: ERF_IndexDefines.H:16
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
#define PrimScalar_comp
Definition: ERF_IndexDefines.H:57
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
int * d_eddy_diff_idx
Definition: ERF_SetupDiff.H:45
const Box zbx
Definition: ERF_SetupDiff.H:9
const Real dx_inv
Definition: ERF_SetupDiff.H:4
const Real dy_inv
Definition: ERF_SetupDiff.H:5
int * d_eddy_diff_idy
Definition: ERF_SetupDiff.H:46
const Box xbx
Definition: ERF_SetupDiff.H:7
const Box ybx
Definition: ERF_SetupDiff.H:8
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
bool l_turb
Definition: ERF_SetupVertDiff.H:9
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
bool l_consA
Definition: ERF_SetupVertDiff.H:8
int * d_eddy_diff_idz
Definition: ERF_SetupVertDiff.H:99
Real * d_alpha_eff
Definition: ERF_SetupVertDiff.H:98
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:117
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:184
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:376
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:198
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:170
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:211
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:90
@ ext_dir
Definition: ERF_IndexDefines.H:243
@ ext_dir_prim
Definition: ERF_IndexDefines.H:246
real(c_double), parameter c2
Definition: ERF_module_model_constants.F90:35
real(c_double), private c1
Definition: ERF_module_mp_morr_two_moment.F90:212
Here is the call graph for this function: