ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_ImplicitDiff_T.cpp File Reference
Include dependency graph for ERF_ImplicitDiff_T.cpp:

Macros

#define INSTANTIATE_IMPLICIT_DIFF_FOR_MOM_LU(STAGDIR)
 

Functions

void ImplicitDiffForStateLU_T (const Box &bx, const Box &domain, const int level, const int n, const Real dt, const GpuArray< Real, AMREX_SPACEDIM *2 > &bc_neumann_vals, const Array4< Real > &cell_data, const Array4< const Real > &z_nd, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &scalar_zflux, const Array4< const Real > &mu_turb, const SolverChoice &solverChoice, const BCRec *bc_ptr, const bool use_SurfLayer, const Real implicit_fac)
 
template<int stagdir>
void ImplicitDiffForMomLU_T (const Box &bx, const Box &, const int level, const Real dt, const Array4< const Real > &cell_data, const Array4< Real > &face_data, const Array4< const Real > &tau, const Array4< const Real > &tau_corr, const Array4< const Real > &z_nd, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mu_turb, const SolverChoice &solverChoice, const BCRec *bc_ptr, const bool use_SurfLayer, const Real implicit_fac)
 

Macro Definition Documentation

◆ INSTANTIATE_IMPLICIT_DIFF_FOR_MOM_LU

#define INSTANTIATE_IMPLICIT_DIFF_FOR_MOM_LU (   STAGDIR)
Value:
template void ImplicitDiffForMomLU_T<STAGDIR> ( \
const Box&, \
const Box&, \
const int, \
const Real, \
const Array4<const Real>&, \
const Array4< Real>&, \
const Array4<const Real>&, \
const Array4<const Real>&, \
const Array4<const Real>&, \
const Array4<const Real>&, \
const GpuArray<Real, AMREX_SPACEDIM>&, \
const Array4<const Real>&, \
const SolverChoice&, \
const BCRec*, \
const bool, \
const Real);
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_DataStruct.H:130

Function Documentation

◆ ImplicitDiffForMomLU_T()

template<int stagdir>
void ImplicitDiffForMomLU_T ( const Box &  bx,
const Box &  ,
const int  level,
const Real  dt,
const Array4< const Real > &  cell_data,
const Array4< Real > &  face_data,
const Array4< const Real > &  tau,
const Array4< const Real > &  tau_corr,
const Array4< const Real > &  z_nd,
const Array4< const Real > &  detJ,
const GpuArray< Real, AMREX_SPACEDIM > &  cellSizeInv,
const Array4< const Real > &  mu_turb,
const SolverChoice solverChoice,
const BCRec *  bc_ptr,
const bool  use_SurfLayer,
const Real  implicit_fac 
)

Function for computing the implicit contribution to the vertical diffusion of momentum, over terrain.

This function (explicitly instantiated below) handles staggering in x, y, or z through the template parameter, stagdir. NOTE: implicit diffusion of w has remains an experimental feature and has not been tested yet with terrain.

Parameters
[in]bxcell-centered box to loop over
[in]domainbox of the whole domain
[in]dttime step
[in]cell_dataconserved cell-centered rho
[in,out]face_dataconserved momentum
[in]tau_corrstress contribution to momentum that will be corrected by the implicit solve
[in]z_ndnodal array of z
[in]detJJacobian determinant
[in]cellSizeInvinverse cell size array
[in]mu_turbturbulent viscosity
[in]solverChoicecontainer of parameters
[in]bc_ptrcontainer with boundary conditions
[in]use_SurfLayerwhether we have turned on subgrid diffusion
[in]implicit_facif 1 then fully implicit; if 0 then fully explicit
238 {
239  BL_PROFILE_VAR("ImplicitDiffForMom_T()",ImplicitDiffForMom_T);
240 
241  // setup quantities for getRhoAlphaAtFaces()
242  DiffChoice dc = solverChoice.diffChoice;
243  TurbChoice tc = solverChoice.turbChoice[level];
245  bool l_turb = tc.use_kturb;
246  Real mu_eff = (l_consA) ? two * dc.dynamic_viscosity / dc.rho0_trans
247  : two * dc.dynamic_viscosity;
248 
249  // g(S*) coefficient
250  // stagdir==0: tau_corr = myhalf * du/dz * mu_tot
251  // stagdir==1: tau_corr = myhalf * dv/dz * mu_tot
252  // stagdir==2: tau_corr = dw/dz * mu_tot
253  constexpr Real gfac = (stagdir == 2) ? two/three : one;
254 
255  // offsets used to average to faces
256  constexpr int ioff = (stagdir == 0) ? 1 : 0;
257  constexpr int joff = (stagdir == 1) ? 1 : 0;
258 
259  // Box bounds
260  int ilo = bx.smallEnd(0);
261  int ihi = bx.bigEnd(0);
262  int jlo = bx.smallEnd(1);
263  int jhi = bx.bigEnd(1);
264  int klo = bx.smallEnd(2);
265  int khi = bx.bigEnd(2);
266  amrex::ignore_unused(ilo, ihi, jlo, jhi);
267 
268  // Temporary FABs for tridiagonal solve (allocated on column)
269  // A[k] * x[k-1] + B[k] * x[k] + C[k+1] = RHS[k]
270  amrex::FArrayBox RHS_fab, soln_fab, coeffG_fab;
271  RHS_fab.resize(bx,1, amrex::The_Async_Arena());
272  soln_fab.resize(bx,1, amrex::The_Async_Arena());
273  coeffG_fab.resize(bx,1, amrex::The_Async_Arena());
274  auto const& RHS_a = RHS_fab.array();
275  auto const& soln_a = soln_fab.array();
276  auto const& coeffG_a = coeffG_fab.array();
277 
278  Real dz_inv = cellSizeInv[2];
279 
280  int bc_comp = BCVars::xvel_bc + stagdir;
281  bool ext_dir_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir ||
282  bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim);
283  bool ext_dir_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir ||
284  bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim);
285  bool foextrap_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::foextrap);
286  amrex::ignore_unused(foextrap_on_zhi);
287 
288  AMREX_ASSERT_WITH_MESSAGE(ext_dir_on_zlo || use_SurfLayer,
289  "Unexpected lower BC used with implicit vertical diffusion");
290  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zhi || ext_dir_on_zhi,
291  "Unexpected upper BC used with implicit vertical diffusion");
292  if (stagdir < 2 && (ext_dir_on_zlo || ext_dir_on_zhi)) {
293  amrex::Warning("No-slip walls have not been fully tested");
294  }
295 
296  Real Fact = implicit_fac * dt * dz_inv;
297 
298 #ifdef AMREX_USE_GPU
299  ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int)
300  {
301 #else
302  for (int j(jlo); j<=jhi; ++j) {
303  for (int i(ilo); i<=ihi; ++i) {
304 #endif
305  // Notes:
306  //
307  // - In DiffusionSrcForMom (e.g., for x-mom)
308  //
309  // Real diffContrib = ...
310  // + (tau13(i,j,k+1) - tau13(i,j,k)) / dzinv
311  // rho_u_rhs(i,j,k) -= diffContrib; // note the negative sign
312  //
313  // - We need to scale the explicit _part_ of `tau13` (for x-mom) by (1 - implicit_fac)
314  // The part that needs to be scaled is stored in `tau_corr`.
315  // E.g., tau13 = myhalf * (du/dz + dw/dx)
316  // tau13_corr = myhalf * du/dz
317  //
318  // - The momentum (`face_data`) was set to `S_old + S_rhs * dt`
319  // prior to including "ERF_Implicit.H". Recall that S_rhs includes
320  // sources from advection and other forcings, not just diffusion.
321  //
322  // - To correct momentum, we need to subtract `implicit_fac * diffContrib_corr`
323  // from S_rhs to recover `(1 - implicit_fac) * diffContrib_corr`,
324  // where `diffContrib_corr = -d(tau_corr)/dz`. The negative sign
325  // comes from our convention for the RHS diffusion source.
326  //
327  // Subtracting a negative gives the += below; multiply by dt to
328  // get the intermediate momentum on the RHS of the tridiagonal
329  // system.
330  //
331  // - With a surface_layer BC, tau13/23 holds the vertical flux -d_z(k*u_i)
332  // directly. We must use tau at klo (not tau_corr) with SL BCs.
333  //
334  // - The detJ for divergence was multiplied through.
335  // Therefore, it doesn't show up in the A/B denominator,
336  // but it does modify B and the RHS.
337  //
338  // - Finally, the terms ~ RHS += (tau_corr_hi - tau_corr_lo) / dz (below)
339  // essentially undo the explicit diffusion update that will be
340  // handled here implicitly.
341 
342  // Bottom boundary coefficients and RHS for L decomp
343  //===================================================
344  Real rhoface, rhoAlpha_lo, rhoAlpha_hi;
345  Real detJface, met_h_zeta_lo, met_h_zeta_hi;
346  Real a_tmp, b_tmp, c_tmp, inv_b2_tmp;
347  {
348  detJface = myhalf * (detJ(i,j,klo) + detJ(i-ioff,j-joff,klo));
349  rhoface = myhalf * (cell_data(i,j,klo,Rho_comp) + cell_data(i-ioff,j-joff,klo,Rho_comp));
350  getRhoAlphaForFaces(i, j, klo, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
351  cell_data, mu_turb, mu_eff,
352  l_consA, l_turb);
353 
354  met_h_zeta_lo = myhalf * ( Compute_h_zeta_AtKface(i ,j ,klo ,cellSizeInv,z_nd)
355  + Compute_h_zeta_AtKface(i-ioff,j-joff,klo ,cellSizeInv,z_nd) );
356  met_h_zeta_hi = myhalf * ( Compute_h_zeta_AtKface(i ,j ,klo+1,cellSizeInv,z_nd)
357  + Compute_h_zeta_AtKface(i-ioff,j-joff,klo+1,cellSizeInv,z_nd) );
358 
359  a_tmp = zero;
360  c_tmp = -Fact * gfac * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
361 
362  RHS_a(i,j,klo) = detJface * face_data(i,j,klo); // NOTE: this is momenta; solution is velocity
363 
364  // BCs: Dirichlet (u_i = val), slip wall (w = 0), or surface layer (w = 0)
365  if (ext_dir_on_zlo) {
366  RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau_corr(i,j,klo));
367  if (stagdir==2) {
368  c_tmp = zero;
369  RHS_a(i,j,klo) = zero;
370  } else {
371  a_tmp = -two * Fact * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
372  RHS_a(i,j,klo) += two * rhoAlpha_lo * face_data(i,j,klo-1) * dz_inv * dz_inv / met_h_zeta_lo;
373  }
374  } else if (use_SurfLayer) {
375  // NOTE: tau = -mu*d_z(u_i) w/ SL
376  RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau(i,j,klo));
377  RHS_a(i,j,klo) += Fact * tau(i,j,klo);
378  }
379 
380  b_tmp = detJface * rhoface - a_tmp - c_tmp;
381  inv_b2_tmp = one;
382 
383  RHS_a(i,j,klo) /= b_tmp; // NOTE: this is now "rho"
384  coeffG_a(i,j,klo) = c_tmp / b_tmp; // NOTE: this is now "gamma"
385  }
386 
387  // Build the coefficients and RHS for L decomp
388  //===================================================
389  for (int k(klo+1); k < khi; k++) {
390  detJface = myhalf * (detJ(i,j,k) + detJ(i-ioff,j-joff,k));
391  rhoface = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i-ioff,j-joff,k,Rho_comp));
392  getRhoAlphaForFaces(i, j, k, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
393  cell_data, mu_turb, mu_eff,
394  l_consA, l_turb);
395 
396  met_h_zeta_lo = myhalf * ( Compute_h_zeta_AtKface(i ,j ,k ,cellSizeInv,z_nd)
397  + Compute_h_zeta_AtKface(i-ioff,j-joff,k ,cellSizeInv,z_nd) );
398  met_h_zeta_hi = myhalf * ( Compute_h_zeta_AtKface(i ,j ,k+1,cellSizeInv,z_nd)
399  + Compute_h_zeta_AtKface(i-ioff,j-joff,k+1,cellSizeInv,z_nd) );
400 
401  a_tmp = -Fact * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
402  c_tmp = -Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
403  b_tmp = detJface * rhoface - a_tmp - c_tmp;
404  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,k-1));
405 
406  RHS_a(i,j,k) = detJface * face_data(i,j,k); // NOTE: this is momenta; solution is velocity
407  RHS_a(i,j,k) += Fact * gfac * (tau_corr(i,j,k+1) - tau_corr(i,j,k));
408 
409  RHS_a(i,j,k) = (RHS_a(i,j,k) - a_tmp * RHS_a(i,j,k-1)) * inv_b2_tmp; // NOTE: This is now "rho"
410  coeffG_a(i,j,k) = c_tmp * inv_b2_tmp; // NOTE: this is now "gamma"
411  } // k
412 
413  // Top boundary coefficients and RHS for L decomp
414  //===================================================
415  {
416  detJface = myhalf * (detJ(i,j,khi) + detJ(i-ioff,j-joff,khi));
417  rhoface = myhalf * (cell_data(i,j,khi,Rho_comp) + cell_data(i-ioff,j-joff,khi,Rho_comp));
418  getRhoAlphaForFaces(i, j, khi, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
419  cell_data, mu_turb, mu_eff,
420  l_consA, l_turb);
421 
422  met_h_zeta_lo = myhalf * ( Compute_h_zeta_AtKface(i ,j ,khi ,cellSizeInv,z_nd)
423  + Compute_h_zeta_AtKface(i-ioff,j-joff,khi ,cellSizeInv,z_nd) );
424  met_h_zeta_hi = myhalf * ( Compute_h_zeta_AtKface(i ,j ,khi+1,cellSizeInv,z_nd)
425  + Compute_h_zeta_AtKface(i-ioff,j-joff,khi+1,cellSizeInv,z_nd) );
426 
427  a_tmp = -Fact * gfac * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
428  c_tmp = zero;
429 
430  RHS_a(i,j,khi) = detJface * face_data(i,j,khi); // NOTE: this is momenta; solution is velocity
431  RHS_a(i,j,khi) += Fact * gfac * (tau_corr(i,j,khi+1) - tau_corr(i,j,khi));
432 
433  // BCs: Dirichlet (u_i = val), slip wall (w = 0)
434  if (ext_dir_on_zhi) {
435  if (stagdir==2) {
436  a_tmp = zero;
437  RHS_a(i,j,khi) = zero;
438  } else {
439  c_tmp = -two * Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
440  RHS_a(i,j,khi) += two * rhoAlpha_hi * face_data(i,j,khi+1) * dz_inv * dz_inv / met_h_zeta_hi;
441  }
442  }
443 
444  b_tmp = detJface * rhoface - a_tmp - c_tmp;
445  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,khi-1));
446 
447  // First solve
448  soln_a(i,j,khi) = (RHS_a(i,j,khi) - a_tmp * RHS_a(i,j,khi-1)) * inv_b2_tmp;
449  }
450 
451  // Back sweep the U decomp solution
452  //===================================================
453  for (int k(khi-1); k>=klo; --k) {
454  soln_a(i,j,k) = RHS_a(i,j,k) - coeffG_a(i,j,k) * soln_a(i,j,k+1);
455  }
456 
457  // Convert back to momenta
458  //===================================================
459  for (int k(klo); k<=khi; ++k) {
460  rhoface = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i-ioff,j-joff,k,Rho_comp));
461  face_data(i,j,k) = rhoface * soln_a(i,j,k);
462  }
463 
464 #ifdef AMREX_USE_GPU
465  });
466 #else
467  } // i
468  } // j
469 #endif
470 }
constexpr amrex::Real three
Definition: ERF_Constants.H:9
constexpr amrex::Real two
Definition: ERF_Constants.H:8
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
constexpr amrex::Real myhalf
Definition: ERF_Constants.H:11
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void getRhoAlphaForFaces(int i, int j, int k, int ioff, int joff, amrex::Real &rhoAlpha_lo, amrex::Real &rhoAlpha_hi, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &mu_turb, const amrex::Real mu_eff, bool l_consA, bool l_turb)
Definition: ERF_GetRhoAlphaForFaces.H:5
#define Rho_comp
Definition: ERF_IndexDefines.H:36
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+myhalf) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=one) { Real dT=T_pert *(std::cos(PI *L)+one)/two;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
bool l_turb
Definition: ERF_SetupVertDiff.H:9
bool l_consA
Definition: ERF_SetupVertDiff.H:8
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_ASSERT_WITH_MESSAGE(wbar_cutoff_min > wbar_cutoff_max, "ERROR: wbar_cutoff_min < wbar_cutoff_max")
@ xvel_bc
Definition: ERF_IndexDefines.H:87
@ foextrap
Definition: ERF_IndexDefines.H:226
@ ext_dir
Definition: ERF_IndexDefines.H:227
@ ext_dir_prim
Definition: ERF_IndexDefines.H:229
Definition: ERF_DiffStruct.H:19
amrex::Real rho0_trans
Definition: ERF_DiffStruct.H:91
MolecDiffType molec_diff_type
Definition: ERF_DiffStruct.H:84
amrex::Real dynamic_viscosity
Definition: ERF_DiffStruct.H:96
DiffChoice diffChoice
Definition: ERF_DataStruct.H:1088
amrex::Vector< TurbChoice > turbChoice
Definition: ERF_DataStruct.H:1091
Definition: ERF_TurbStruct.H:42
bool use_kturb
Definition: ERF_TurbStruct.H:424
Here is the call graph for this function:

◆ ImplicitDiffForStateLU_T()

void ImplicitDiffForStateLU_T ( const Box &  bx,
const Box &  domain,
const int  level,
const int  n,
const Real  dt,
const GpuArray< Real, AMREX_SPACEDIM *2 > &  bc_neumann_vals,
const Array4< Real > &  cell_data,
const Array4< const Real > &  z_nd,
const Array4< const Real > &  detJ,
const GpuArray< Real, AMREX_SPACEDIM > &  cellSizeInv,
const Array4< const Real > &  scalar_zflux,
const Array4< const Real > &  mu_turb,
const SolverChoice solverChoice,
const BCRec *  bc_ptr,
const bool  use_SurfLayer,
const Real  implicit_fac 
)

Function for computing the implicit contribution to the vertical diffusion of theta, with terrain.

Parameters
[in]bxcell-centered box to loop over
[in]domainbox of the whole domain
[in]dttime step
[in]bc_neumann_valsvalues of derivatives if bc_type == Neumann
[in,out]cell_dataconserved cell-centered rho, rho theta
[in]z_ndnodal array of z
[in]detJJacobian determinant
[in]cellSizeInvinverse cell size array
[in,out]hfx_zheat flux in z-dir
[in]mu_turbturbulent viscosity
[in]solverChoicecontainer of parameters
[in]bc_ptrcontainer with boundary conditions
[in]use_SurfLayerwhether we have turned on subgrid diffusion
[in]implicit_facif 1 then fully implicit; if 0 then fully explicit
45 {
46  BL_PROFILE_VAR("ImplicitDiffForState_T()",ImplicitDiffForState_T);
47 
48  // setup quantities for getRhoAlpha()
49 #include "ERF_SetupVertDiff.H"
50  const int qty_index = n;
51  const int prim_index = qty_index - 1;
52  const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
53 
54  // Box bounds
55  int ilo = bx.smallEnd(0);
56  int ihi = bx.bigEnd(0);
57  int jlo = bx.smallEnd(1);
58  int jhi = bx.bigEnd(1);
59  int klo = bx.smallEnd(2);
60  int khi = bx.bigEnd(2);
61  amrex::ignore_unused(ilo, ihi, jlo, jhi);
62 
63  // Temporary FABs for tridiagonal solve (allocated on column)
64  // A[k] * x[k-1] + B[k] * x[k] + C[k+1] = RHS[k]
65 
66  // With LU decomposition, M * x = r is written as L * U * x = r with U * x = rho
67  // We then first have L * rho = r and U * x = rho
68  amrex::FArrayBox RHS_fab, soln_fab, coeffG_fab;
69  RHS_fab.resize(bx,1, amrex::The_Async_Arena());
70  soln_fab.resize(bx,1, amrex::The_Async_Arena());
71  coeffG_fab.resize(bx,1, amrex::The_Async_Arena());
72  auto const& RHS_a = RHS_fab.array();
73  auto const& soln_a = soln_fab.array();
74  auto const& coeffG_a = coeffG_fab.array();
75 
76  Real dz_inv = cellSizeInv[2];
77 
78  int bc_comp = qty_index;
79  bool foextrap_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::foextrap);
80  bool foextrap_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::foextrap);
81  bool neumann_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::neumann);
82  bool neumann_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::neumann);
83  amrex::ignore_unused(foextrap_on_zlo, foextrap_on_zhi);
84 
85  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zlo || neumann_on_zlo || use_SurfLayer,
86  "Unexpected lower BC used with implicit vertical diffusion");
87  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zhi || neumann_on_zhi,
88  "Unexpected upper BC used with implicit vertical diffusion");
89 
90  Real Fact = implicit_fac * dt * dz_inv;
91 
92 #ifdef AMREX_USE_GPU
93  ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int)
94  {
95 #else
96  for (int j(jlo); j<=jhi; ++j) {
97  for (int i(ilo); i<=ihi; ++i) {
98 #endif
99  // Notes: The detJ for divergence was multiplied through.
100  // Therefore, it doesn't show up in the A/B denominator,
101  // but it does modify B and the RHS.
102 
103  // Bottom boundary coefficients and RHS for L decomp
104  //===================================================
105  Real rhoAlpha_lo, rhoAlpha_hi;
106  Real met_h_zeta_lo, met_h_zeta_hi;
107  Real a_tmp, b_tmp, c_tmp, inv_b2_tmp;
108  {
109  getRhoAlpha(i, j, klo, rhoAlpha_lo, rhoAlpha_hi,
110  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
111  prim_index, prim_scal_index, l_consA, l_turb);
112 
113  met_h_zeta_lo = Compute_h_zeta_AtKface(i,j,klo ,cellSizeInv,z_nd);
114  met_h_zeta_hi = Compute_h_zeta_AtKface(i,j,klo+1,cellSizeInv,z_nd);
115 
116  a_tmp = zero;
117  c_tmp = -Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
118  b_tmp = detJ(i,j,klo) * cell_data(i,j,klo,Rho_comp) - a_tmp - c_tmp;
119  inv_b2_tmp = one;
120 
121  RHS_a(i,j,klo) = detJ(i,j,klo) * cell_data(i,j,klo,n); // NOTE: this is rho*phi; solution is phi
122  if (use_SurfLayer && scalar_zflux) {
123  RHS_a(i,j,klo) += Fact * scalar_zflux(i,j,klo); // NOTE: scalar_zflux = -K*d_z(\phi)
124  } else if (neumann_on_zlo) {
125  RHS_a(i,j,klo) += -Fact * rhoAlpha_lo * bc_neumann_vals[2]; // NOTE: N_val = d_z(\phi)
126  }
127 
128  RHS_a(i,j,klo) /= b_tmp; // NOTE: this is now "rho"
129  coeffG_a(i,j,klo) = c_tmp / b_tmp; // NOTE: this is now "gamma"
130  }
131 
132  // Build the coefficients and RHS for L decomp
133  //===================================================
134  for (int k(klo+1); k < khi; k++) {
135  getRhoAlpha(i, j, k, rhoAlpha_lo, rhoAlpha_hi,
136  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
137  prim_index, prim_scal_index, l_consA, l_turb);
138 
139  met_h_zeta_lo = Compute_h_zeta_AtKface(i,j,k ,cellSizeInv,z_nd);
140  met_h_zeta_hi = Compute_h_zeta_AtKface(i,j,k+1,cellSizeInv,z_nd);
141 
142  a_tmp = -Fact * rhoAlpha_lo * dz_inv / met_h_zeta_lo;
143  c_tmp = -Fact * rhoAlpha_hi * dz_inv / met_h_zeta_hi;
144  b_tmp = detJ(i,j,k) * cell_data(i,j,k,Rho_comp) - a_tmp - c_tmp;
145  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,k-1));
146 
147  RHS_a(i,j,k) = detJ(i,j,k) * cell_data(i,j,k,n); // NOTE: this is rho*phi; solution is phi
148 
149  RHS_a(i,j,k) = (RHS_a(i,j,k) - a_tmp * RHS_a(i,j,k-1)) * inv_b2_tmp; // NOTE: This is now "rho"
150  coeffG_a(i,j,k) = c_tmp * inv_b2_tmp; // NOTE: this is now "gamma"
151  } // k
152 
153  // Top boundary coefficients and RHS for L decomp
154  //===================================================
155  {
156  getRhoAlpha(i, j, khi, rhoAlpha_lo, rhoAlpha_hi,
157  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
158  prim_index, prim_scal_index, l_consA, l_turb);
159 
160  met_h_zeta_lo = Compute_h_zeta_AtKface(i,j,khi ,cellSizeInv,z_nd);
161  met_h_zeta_hi = Compute_h_zeta_AtKface(i,j,khi+1,cellSizeInv,z_nd);
162 
163  a_tmp = -Fact * rhoAlpha_lo * dz_inv / met_h_zeta_hi;
164  c_tmp = zero;
165  b_tmp = detJ(i,j,khi) * cell_data(i,j,khi,Rho_comp) - a_tmp - c_tmp;
166  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,khi-1));
167 
168  RHS_a(i,j,khi) = detJ(i,j,khi) * cell_data(i,j,khi,n); // NOTE: this is rho*phi; solution is phi
169  if (neumann_on_zhi) {
170  RHS_a(i,j,khi) -= -Fact * rhoAlpha_hi * bc_neumann_vals[5]; // NOTE: N_val = d_z(\phi)
171  }
172 
173  // First solve
174  soln_a(i,j,khi) = (RHS_a(i,j,khi) - a_tmp * RHS_a(i,j,khi-1)) * inv_b2_tmp;
175  }
176 
177  // Back sweep the U decomp solution
178  //===================================================
179  for (int k(khi-1); k>=klo; --k) {
180  soln_a(i,j,k) = RHS_a(i,j,k) - coeffG_a(i,j,k) * soln_a(i,j,k+1);
181  }
182 
183  // Convert back to rho*theta
184  //===================================================
185  for (int k(klo); k<=khi; ++k) {
186  cell_data(i,j,k,n) = cell_data(i,j,k,Rho_comp) * soln_a(i,j,k);
187  }
188 
189 #ifdef AMREX_USE_GPU
190  });
191 #else
192  } // i
193  } // j
194 #endif
195 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void getRhoAlpha(int i, int j, int k, amrex::Real &rhoAlpha_lo, amrex::Real &rhoAlpha_hi, const amrex::Array4< const amrex::Real > &cell_data, const amrex::Array4< const amrex::Real > &mu_turb, const amrex::Real *d_alpha_eff, const int *d_eddy_diff_idz, int prim_index, int prim_scal_index, bool l_consA, bool l_turb)
Definition: ERF_GetRhoAlpha.H:5
#define RhoScalar_comp
Definition: ERF_IndexDefines.H:40
#define NSCALARS
Definition: ERF_IndexDefines.H:16
#define PrimScalar_comp
Definition: ERF_IndexDefines.H:52
int * d_eddy_diff_idz
Definition: ERF_SetupVertDiff.H:108
Real * d_alpha_eff
Definition: ERF_SetupVertDiff.H:105
@ neumann
Definition: ERF_IndexDefines.H:231
Here is the call graph for this function: