ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_DiffusionSrcForState_T.cpp File Reference
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)
 

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 
)

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