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 = 1.0 - 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 = 0.5 * ( 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 += 0.5 * ( 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;
123 
124  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
125 
126  Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
127  Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
128  Real GradCz = 0.5 * ( 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 = 0.5 * ( 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 += 0.5 * ( 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;
155  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
156 
157  Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
158  Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
159  Real GradCz = 0.5 * ( 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 = 0.5 * ( 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 += 0.5 * ( 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;
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].lo(5) == ERFBCType::ext_dir) ||
189  (bc_ptr[bc_comp].lo(5) == 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 = 1.0 / dz0;
199  Real f = (dz1 / dz0) + 2.0;
200  Real f2 = f*f;
201  Real c3 = 2.0 / (f - f2);
202  Real c2 = -f2*c3;
203  Real c1 = -(1.0-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 = 1.0 / dz0;
213  Real f = (dz1 / dz0) + 2.0;
214  Real f2 = f*f;
215  Real c3 = 2.0 / (f - f2);
216  Real c2 = -f2*c3;
217  Real c1 = -(1.0-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 && (qty_index == RhoTheta_comp)) {
227  zflux(i,j,k) = hfx_z(i,j,0);
228  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
229  zflux(i,j,k) = qfx1_z(i,j,0);
230  } else {
231  zflux(i,j,k) = -rhoAlpha * GradCz;
232  }
233 
234  if (qty_index == RhoTheta_comp) {
235  if (!SurfLayer_on_zlo) {
236  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
237  }
238  } else if (qty_index == RhoQ1_comp) {
239  if (!SurfLayer_on_zlo) {
240  qfx1_z(i,j,k) = zflux(i,j,k);
241  }
242  } else if (qty_index == RhoQ2_comp) {
243  qfx2_z(i,j,k) = zflux(i,j,k);
244  }
245  });
246  // Constant rho*alpha & Turb model
247  } else if (l_turb) {
248  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
249  {
250  const int prim_index = qty_index - 1;
251 
252  Real rhoAlpha = d_alpha_eff[prim_index];
253  rhoAlpha += 0.5 * ( mu_turb(i , j, k, d_eddy_diff_idx[prim_index])
254  + mu_turb(i-1, j, k, d_eddy_diff_idx[prim_index]) );
255 
256  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
257 
258  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
259  BCVars::RhoScalar_bc_comp : qty_index;
261  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
262 
263  Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
264  Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
265  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
266  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
267  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
268 
269  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
270  xflux(i,j,k) = hfx_x(i,j,0);
271  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
272  xflux(i,j,k) = qfx1_x(i,j,0);
273  } else {
274  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
275  }
276  });
277  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
278  {
279  const int prim_index = qty_index - 1;
280 
281  Real rhoAlpha = d_alpha_eff[prim_index];
282  rhoAlpha += 0.5 * ( mu_turb(i, j , k, d_eddy_diff_idy[prim_index])
283  + mu_turb(i, j-1, k, d_eddy_diff_idy[prim_index]) );
284 
285  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
286 
287  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
288  BCVars::RhoScalar_bc_comp : qty_index;
290  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
291 
292  Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
293  Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
294  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
295  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
296  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
297 
298  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
299  yflux(i,j,k) = hfx_y(i,j,0);
300  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
301  yflux(i,j,k) = qfx1_y(i,j,0);
302  } else {
303  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
304  }
305  });
306  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
307  {
308  const int prim_index = qty_index - 1;
309 
310  Real rhoAlpha = d_alpha_eff[prim_index];
311  rhoAlpha += 0.5 * ( mu_turb(i, j, k , d_eddy_diff_idz[prim_index])
312  + mu_turb(i, j, k-1, d_eddy_diff_idz[prim_index]) );
313 
314  Real GradCz;
315  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
316  BCVars::RhoScalar_bc_comp : qty_index;
318  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
319  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
320  && k == dom_lo.z);
321  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir) ||
322  (bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir_prim))
323  && k == dom_hi.z+1);
324  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
325 
326  if (ext_dir_on_zlo) {
327  // Third order stencil with variable dz
328  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
329  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
330  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
331  Real idz0 = 1.0 / dz0;
332  Real f = (dz1 / dz0) + 2.0;
333  Real f2 = f*f;
334  Real c3 = 2.0 / (f - f2);
335  Real c2 = -f2*c3;
336  Real c1 = -(1.0-f2)*c3;
337  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
338  + c2 * cell_prim(i, j, k , prim_index)
339  + c3 * cell_prim(i, j, k+1, prim_index) );
340  } else if (ext_dir_on_zhi) {
341  // Third order stencil with variable dz
342  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
343  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
344  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
345  Real idz0 = 1.0 / dz0;
346  Real f = (dz1 / dz0) + 2.0;
347  Real f2 = f*f;
348  Real c3 = 2.0 / (f - f2);
349  Real c2 = -f2*c3;
350  Real c1 = -(1.0-f2)*c3;
351  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
352  + c2 * cell_prim(i, j, k-1, prim_index)
353  + c3 * cell_prim(i, j, k-2, prim_index) ) );
354  } else {
355  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
356  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
357  }
358 
359  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
360  zflux(i,j,k) = hfx_z(i,j,0);
361  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
362  zflux(i,j,k) = qfx1_z(i,j,0);
363  } else {
364  zflux(i,j,k) = -rhoAlpha * GradCz;
365  }
366 
367  if (qty_index == RhoTheta_comp) {
368  if (!SurfLayer_on_zlo) {
369  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
370  }
371  } else if (qty_index == RhoQ1_comp) {
372  if (!SurfLayer_on_zlo) {
373  qfx1_z(i,j,k) = zflux(i,j,k);
374  }
375  } else if (qty_index == RhoQ2_comp) {
376  qfx2_z(i,j,k) = zflux(i,j,k);
377  }
378 
379  });
380  // Constant alpha & no LES/PBL model
381  } else if(l_consA) {
382  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
383  {
384  const int prim_index = qty_index - 1;
385 
386  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i-1, j, k, Rho_comp) );
387  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
388 
389  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
390 
391  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
392  BCVars::RhoScalar_bc_comp : qty_index;
394  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
395 
396  Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
397  Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
398  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
399  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
400  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
401 
402  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
403  xflux(i,j,k) = hfx_x(i,j,0);
404  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
405  xflux(i,j,k) = qfx1_x(i,j,0);
406  } else {
407  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
408  }
409  });
410  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
411  {
412  const int prim_index = qty_index - 1;
413 
414  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
415  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
416 
417  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
418 
419  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
420  BCVars::RhoScalar_bc_comp : qty_index;
422  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
423 
424  Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
425  Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
426  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
427  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
428  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
429 
430  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
431  yflux(i,j,k) = hfx_y(i,j,0);
432  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
433  yflux(i,j,k) = qfx1_y(i,j,0);
434  } else {
435  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
436  }
437  });
438  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
439  {
440  const int prim_index = qty_index - 1;
441 
442  Real rhoFace = 0.5 * ( cell_data(i, j, k, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
443  Real rhoAlpha = rhoFace * d_alpha_eff[prim_index];
444 
445  Real GradCz;
446  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
447  BCVars::RhoScalar_bc_comp : qty_index;
449  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
450  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
451  && k == dom_lo.z);
452  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir) ||
453  (bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir_prim))
454  && k == dom_hi.z+1);
455  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
456 
457  if (ext_dir_on_zlo) {
458  // Third order stencil with variable dz
459  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
460  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
461  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
462  Real idz0 = 1.0 / dz0;
463  Real f = (dz1 / dz0) + 2.0;
464  Real f2 = f*f;
465  Real c3 = 2.0 / (f - f2);
466  Real c2 = -f2*c3;
467  Real c1 = -(1.0-f2)*c3;
468  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
469  + c2 * cell_prim(i, j, k , prim_index)
470  + c3 * cell_prim(i, j, k+1, prim_index) );
471  } else if (ext_dir_on_zhi) {
472  // Third order stencil with variable dz
473  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
474  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
475  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
476  Real idz0 = 1.0 / dz0;
477  Real f = (dz1 / dz0) + 2.0;
478  Real f2 = f*f;
479  Real c3 = 2.0 / (f - f2);
480  Real c2 = -f2*c3;
481  Real c1 = -(1.0-f2)*c3;
482  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
483  + c2 * cell_prim(i, j, k-1, prim_index)
484  + c3 * cell_prim(i, j, k-2, prim_index) ) );
485  } else {
486  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
487  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
488  }
489 
490  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
491  zflux(i,j,k) = hfx_z(i,j,0);
492  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
493  zflux(i,j,k) = qfx1_z(i,j,0);
494  } else {
495  zflux(i,j,k) = -rhoAlpha * GradCz;
496  }
497 
498  if (qty_index == RhoTheta_comp) {
499  if (!SurfLayer_on_zlo) {
500  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
501  }
502  } else if (qty_index == RhoQ1_comp) {
503  if (!SurfLayer_on_zlo) {
504  qfx1_z(i,j,k) = zflux(i,j,k);
505  }
506  } else if (qty_index == RhoQ2_comp) {
507  qfx2_z(i,j,k) = zflux(i,j,k);
508  }
509  });
510  // Constant rho*alpha & no LES/PBL model
511  } else {
512  ParallelFor(xbx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
513  {
514  const int prim_index = qty_index - 1;
515 
516  Real rhoAlpha = d_alpha_eff[prim_index];
517 
518  Real met_h_xi = Compute_h_xi_AtIface (i,j,k,cellSizeInv,z_nd);
519 
520  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
521  BCVars::RhoScalar_bc_comp : qty_index;
523  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
524 
525  Real idz_hi = 1.0 / (z_cc(i ,j,k+1) - z_cc(i ,j,k-1));
526  Real idz_lo = 1.0 / (z_cc(i-1,j,k+1) - z_cc(i-1,j,k-1));
527  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i-1, j, k+1, prim_index)*idz_lo
528  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i-1, j, k-1, prim_index)*idz_lo );
529  Real GradCx = dx_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i-1, j, k , prim_index) );
530 
531  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
532  xflux(i,j,k) = hfx_x(i,j,0);
533  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
534  xflux(i,j,k) = qfx1_x(i,j,0);
535  } else {
536  xflux(i,j,k) = -rhoAlpha * mf_ux(i,j,0) * ( GradCx - met_h_xi*GradCz );
537  }
538  });
539  ParallelFor(ybx_g1, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
540  {
541  const int prim_index = qty_index - 1;
542 
543  Real rhoAlpha = d_alpha_eff[prim_index];
544 
545  Real met_h_eta = Compute_h_eta_AtJface (i,j,k,cellSizeInv,z_nd);
546 
547  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
548  BCVars::RhoScalar_bc_comp : qty_index;
550  bool SurfLayer_on_zlo = ( use_SurfLayer && rotate && k == dom_lo.z);
551 
552  Real idz_hi = 1.0 / (z_cc(i,j ,k+1) - z_cc(i,j ,k-1));
553  Real idz_lo = 1.0 / (z_cc(i,j-1,k+1) - z_cc(i,j-1,k-1));
554  Real GradCz = 0.5 * ( cell_prim(i, j, k+1, prim_index)*idz_hi + cell_prim(i, j-1, k+1, prim_index)*idz_lo
555  - cell_prim(i, j, k-1, prim_index)*idz_hi - cell_prim(i, j-1, k-1, prim_index)*idz_lo );
556  Real GradCy = dy_inv * ( cell_prim(i, j, k , prim_index) - cell_prim(i, j-1, k , prim_index) );
557 
558  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
559  yflux(i,j,k) = hfx_y(i,j,0);
560  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
561  yflux(i,j,k) = qfx1_y(i,j,0);
562  } else {
563  yflux(i,j,k) = -rhoAlpha * mf_vy(i,j,0) * ( GradCy - met_h_eta*GradCz );
564  }
565  });
566  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
567  {
568  const int prim_index = qty_index - 1;
569 
570  Real rhoAlpha = d_alpha_eff[prim_index];
571 
572 
573  Real GradCz;
574  int bc_comp = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ?
575  BCVars::RhoScalar_bc_comp : qty_index;
577  bool ext_dir_on_zlo = ( ((bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir) ||
578  (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim))
579  && k == dom_lo.z);
580  bool ext_dir_on_zhi = ( ((bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir) ||
581  (bc_ptr[bc_comp].lo(5) == ERFBCType::ext_dir_prim))
582  && k == dom_hi.z+1);
583  bool SurfLayer_on_zlo = ( use_SurfLayer && k == dom_lo.z);
584 
585  if (ext_dir_on_zlo) {
586  // Third order stencil with variable dz
587  Real zm = Compute_Z_AtWFace(i,j,k+1,z_nd);
588  Real dz0 = zm - Compute_Z_AtWFace(i,j,k,z_nd);
589  Real dz1 = Compute_Z_AtWFace(i,j,k+2,z_nd) - zm;
590  Real idz0 = 1.0 / dz0;
591  Real f = (dz1 / dz0) + 2.0;
592  Real f2 = f*f;
593  Real c3 = 2.0 / (f - f2);
594  Real c2 = -f2*c3;
595  Real c1 = -(1.0-f2)*c3;
596  GradCz = idz0 * ( c1 * cell_prim(i, j, k-1, prim_index)
597  + c2 * cell_prim(i, j, k , prim_index)
598  + c3 * cell_prim(i, j, k+1, prim_index) );
599  } else if (ext_dir_on_zhi) {
600  // Third order stencil with variable dz
601  Real zm = Compute_Z_AtWFace(i,j,k-1,z_nd);
602  Real dz0 = Compute_Z_AtWFace(i,j,k,z_nd) - zm;
603  Real dz1 = zm - Compute_Z_AtWFace(i,j,k-2,z_nd);
604  Real idz0 = 1.0 / dz0;
605  Real f = (dz1 / dz0) + 2.0;
606  Real f2 = f*f;
607  Real c3 = 2.0 / (f - f2);
608  Real c2 = -f2*c3;
609  Real c1 = -(1.0-f2)*c3;
610  GradCz = idz0 * ( -( c1 * cell_prim(i, j, k , prim_index)
611  + c2 * cell_prim(i, j, k-1, prim_index)
612  + c3 * cell_prim(i, j, k-2, prim_index) ) );
613  } else {
614  Real met_h_zeta = Compute_h_zeta_AtKface(i,j,k,cellSizeInv,z_nd);
615  GradCz = (dz_inv/met_h_zeta) * ( cell_prim(i, j, k, prim_index) - cell_prim(i, j, k-1, prim_index) );
616  }
617 
618  if (SurfLayer_on_zlo && (qty_index == RhoTheta_comp)) {
619  zflux(i,j,k) = hfx_z(i,j,0);
620  } else if (SurfLayer_on_zlo && (qty_index == RhoQ1_comp)) {
621  zflux(i,j,k) = qfx1_z(i,j,0);
622  } else {
623  zflux(i,j,k) = -rhoAlpha * GradCz;
624  }
625 
626  if (qty_index == RhoTheta_comp) {
627  if (!SurfLayer_on_zlo) {
628  hfx_z(i,j,k) = zflux(i,j,k) * explicit_fac;
629  }
630  } else if (qty_index == RhoQ1_comp) {
631  if (!SurfLayer_on_zlo) {
632  qfx1_z(i,j,k) = zflux(i,j,k);
633  }
634  } else if (qty_index == RhoQ2_comp) {
635  qfx2_z(i,j,k) = zflux(i,j,k);
636  }
637  });
638  }
639 
640  //-----------------------------------------------------------------------------------
641  //
642  // Modify fluxes by terrain and use fluxes to compute RHS
643  //
644  // Note that we combine all of these operations in order to keep this section
645  // of the loop tiling-safe.
646  //-----------------------------------------------------------------------------------
647  ParallelFor(bx,[=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
648  {
649  Real met_h_xi_lo = Compute_h_xi_AtKface (i,j,k ,cellSizeInv,z_nd);
650  Real met_h_xi_hi = Compute_h_xi_AtKface (i,j,k+1,cellSizeInv,z_nd);
651 
652  Real met_h_eta_lo = Compute_h_eta_AtKface(i,j,k ,cellSizeInv,z_nd);
653  Real met_h_eta_hi = Compute_h_eta_AtKface(i,j,k+1,cellSizeInv,z_nd);
654 
655  Real xfluxbar_lo, yfluxbar_lo;
656  if (k == dom_lo.z) {
657  Real xfluxlo = 0.5 * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
658  Real xfluxhi = 0.5 * ( xflux(i,j,k+1) + xflux(i+1,j,k+1) );
659  xfluxbar_lo = 1.5*xfluxlo - 0.5*xfluxhi;
660 
661  Real yfluxlo = 0.5 * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
662  Real yfluxhi = 0.5 * ( yflux(i,j,k+1) + yflux(i,j+1,k+1) );
663  yfluxbar_lo = 1.5*yfluxlo - 0.5*yfluxhi;
664  } else {
665  xfluxbar_lo = 0.25 * ( xflux(i,j,k ) + xflux(i+1,j ,k )
666  + xflux(i,j,k-1) + xflux(i+1,j ,k-1) );
667  yfluxbar_lo = 0.25 * ( yflux(i,j,k ) + yflux(i ,j+1,k )
668  + yflux(i,j,k-1) + yflux(i ,j+1,k-1) );
669  }
670 
671  Real xfluxbar_hi, yfluxbar_hi;
672  if (k == dom_hi.z) {
673  Real xfluxlo = 0.5 * ( xflux(i,j,k-1) + xflux(i+1,j,k-1) );
674  Real xfluxhi = 0.5 * ( xflux(i,j,k ) + xflux(i+1,j,k ) );
675  xfluxbar_hi = 1.5*xfluxhi - 0.5*xfluxlo;
676 
677  Real yfluxlo = 0.5 * ( yflux(i,j,k-1) + yflux(i,j+1,k-1) );
678  Real yfluxhi = 0.5 * ( yflux(i,j,k ) + yflux(i,j+1,k ) );
679  yfluxbar_hi = 1.5*yfluxhi - 0.5*yfluxlo;
680  } else {
681  xfluxbar_hi = 0.25 * ( xflux(i,j,k+1) + xflux(i+1,j ,k+1)
682  + xflux(i,j,k ) + xflux(i+1,j ,k ) );
683  yfluxbar_hi = 0.25 * ( yflux(i,j,k+1) + yflux(i ,j+1,k+1)
684  + yflux(i,j,k ) + yflux(i ,j+1,k ) );
685  }
686 
687  // Allow semi-implicit discretization of the vertical diffusive terms
688  Real zflux_lo = explicit_fac * zflux(i,j,k )
689  - met_h_xi_lo * mf_mx(i,j,0) * xfluxbar_lo
690  - met_h_eta_lo * mf_my(i,j,0) * yfluxbar_lo;
691  Real zflux_hi = explicit_fac * zflux(i,j,k+1)
692  - met_h_xi_hi * mf_mx(i,j,0) * xfluxbar_hi
693  - met_h_eta_hi * mf_my(i,j,0) * yfluxbar_hi;
694 
695  Real mfsq = mf_mx(i,j,0) * mf_my(i,j,0);
696  Real stateContrib = ( xflux(i+1,j ,k ) * ax(i+1,j,k) / mf_uy(i+1,j,0)
697  -xflux(i ,j ,k ) * ax(i ,j,k) / mf_uy(i ,j,0) ) * dx_inv * mfsq // Diffusive flux in x-dir
698  +( yflux(i ,j+1,k ) * ay(i,j+1,k) / mf_vx(i,j+1,0)
699  -yflux(i ,j ,k ) * ay(i,j ,k) / mf_vx(i,j ,0) ) * dy_inv * mfsq // Diffusive flux in y-dir
700  +( zflux_hi - zflux_lo) * dz_inv; // Diffusive flux in z-dir
701 
702  stateContrib /= detJ(i,j,k);
703 
704  cell_rhs(i,j,k,qty_index) -= stateContrib;
705  });
706  } // n
707 
708 #include "ERF_AddTKESources.H"
709 #include "ERF_AddQKESources.H"
710 }
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
const int bc_comp
Definition: ERF_Implicit.H:8
#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:52
int * d_eddy_diff_idx
Definition: ERF_SetupDiff.H:47
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:48
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:8
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
bool l_consA
Definition: ERF_SetupVertDiff.H:7
int * d_eddy_diff_idz
Definition: ERF_SetupVertDiff.H:107
Real * d_alpha_eff
Definition: ERF_SetupVertDiff.H:104
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:115
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:182
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:374
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:196
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:168
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:209
@ RhoScalar_bc_comp
Definition: ERF_IndexDefines.H:80
@ ext_dir
Definition: ERF_IndexDefines.H:209
@ ext_dir_prim
Definition: ERF_IndexDefines.H:211
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: