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

Macros

#define INSTANTIATE_IMPLICIT_DIFF_FOR_MOM_LU(STAGDIR)
 

Functions

void ImplicitDiffForStateLU_S (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 Gpu::DeviceVector< Real > &stretched_dz_d, 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_S (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 Gpu::DeviceVector< Real > &stretched_dz_d, 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_S<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 Gpu::DeviceVector<Real>&, \
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_S()

template<int stagdir>
void ImplicitDiffForMomLU_S ( 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 Gpu::DeviceVector< Real > &  stretched_dz_d,
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, with a vertically stretched grid over flat terrain.

This function (explicitly instantiated below) handles staggering in x, y, or z through the template parameter, stagdir.

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]stretched_dz_darray over z of dz[k]
[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
228 {
229  BL_PROFILE_VAR("ImplicitDiffForMom_S()",ImplicitDiffForMom_S);
230 
231  // setup quantities for getRhoAlphaAtFaces()
232  DiffChoice dc = solverChoice.diffChoice;
233  TurbChoice tc = solverChoice.turbChoice[level];
235  bool l_turb = tc.use_kturb;
236  Real mu_eff = (l_consA) ? two * dc.dynamic_viscosity / dc.rho0_trans
237  : two * dc.dynamic_viscosity;
238 
239  // g(S*) coefficient
240  // stagdir==0: tau_corr = myhalf * du/dz * mu_tot
241  // stagdir==1: tau_corr = myhalf * dv/dz * mu_tot
242  // stagdir==2: tau_corr = dw/dz * mu_tot
243  constexpr Real gfac = (stagdir == 2) ? two/three : one;
244 
245  // offsets used to average to faces
246  constexpr int ioff = (stagdir == 0) ? 1 : 0;
247  constexpr int joff = (stagdir == 1) ? 1 : 0;
248 
249  // Box bounds
250  int ilo = bx.smallEnd(0);
251  int ihi = bx.bigEnd(0);
252  int jlo = bx.smallEnd(1);
253  int jhi = bx.bigEnd(1);
254  int klo = bx.smallEnd(2);
255  int khi = bx.bigEnd(2);
256  amrex::ignore_unused(ilo, ihi, jlo, jhi);
257 
258  // Temporary FABs for tridiagonal solve (allocated on column)
259  // A[k] * x[k-1] + B[k] * x[k] + C[k+1] = RHS[k]
260  amrex::FArrayBox RHS_fab, soln_fab, coeffG_fab;
261  RHS_fab.resize(bx,1, amrex::The_Async_Arena());
262  soln_fab.resize(bx,1, amrex::The_Async_Arena());
263  coeffG_fab.resize(bx,1, amrex::The_Async_Arena());
264  auto const& RHS_a = RHS_fab.array();
265  auto const& soln_a = soln_fab.array();
266  auto const& coeffG_a = coeffG_fab.array();
267 
268  auto dz_ptr = stretched_dz_d.data();
269 
270  int bc_comp = BCVars::xvel_bc + stagdir;
271  bool ext_dir_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir ||
272  bc_ptr[bc_comp].lo(2) == ERFBCType::ext_dir_prim);
273  bool ext_dir_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir ||
274  bc_ptr[bc_comp].hi(2) == ERFBCType::ext_dir_prim);
275  bool foextrap_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::foextrap);
276  amrex::ignore_unused(foextrap_on_zhi);
277 
278  AMREX_ASSERT_WITH_MESSAGE(ext_dir_on_zlo || use_SurfLayer,
279  "Unexpected lower BC used with implicit vertical diffusion");
280  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zhi || ext_dir_on_zhi,
281  "Unexpected upper BC used with implicit vertical diffusion");
282  if (stagdir < 2 && (ext_dir_on_zlo || ext_dir_on_zhi)) {
283  amrex::Warning("No-slip walls have not been fully tested");
284  }
285 
286  Real Fact = implicit_fac * dt;
287 
288 #ifdef AMREX_USE_GPU
289  ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int)
290  {
291 #else
292  for (int j(jlo); j<=jhi; ++j) {
293  for (int i(ilo); i<=ihi; ++i) {
294 #endif
295  // Notes:
296  //
297  // - In DiffusionSrcForMom (e.g., for x-mom)
298  //
299  // Real diffContrib = ...
300  // + (tau13(i,j,k+1) - tau13(i,j,k)) / dzinv
301  // rho_u_rhs(i,j,k) -= diffContrib; // note the negative sign
302  //
303  // - We need to scale the explicit _part_ of `tau13` (for x-mom) by (1 - implicit_fac)
304  // The part that needs to be scaled is stored in `tau_corr`.
305  // E.g., tau13 = myhalf * (du/dz + dw/dx)
306  // tau13_corr = myhalf * du/dz
307  //
308  // - The momentum (`face_data`) was set to `S_old + S_rhs * dt`
309  // prior to including "ERF_Implicit.H". Recall that S_rhs includes
310  // sources from advection and other forcings, not just diffusion.
311  //
312  // - To correct momentum, we need to subtract `implicit_fac * diffContrib_corr`
313  // from S_rhs to recover `(1 - implicit_fac) * diffContrib_corr`,
314  // where `diffContrib_corr = -d(tau_corr)/dz`. The negative sign
315  // comes from our convention for the RHS diffusion source.
316  //
317  // Subtracting a negative gives the += below; multiply by dt to
318  // get the intermediate momentum on the RHS of the tridiagonal
319  // system.
320  //
321  // - With a surface_layer BC, tau13/23 holds the vertical flux -d_z(k*u_i)
322  // directly. We must use tau at klo (not tau_corr) with SL BCs.
323  //
324  // - Finally, the terms ~ RHS += (tau_corr_hi - tau_corr_lo) / dz (below)
325  // essentially undo the explicit diffusion update that will be
326  // handled here implicitly.
327 
328  // Bottom boundary coefficients and RHS for L decomp
329  //===================================================
330  Real rhoface, rhoAlpha_lo, rhoAlpha_hi;
331  Real dz_inv, dz_inv_lo, dz_inv_hi;
332  Real a_tmp, b_tmp, c_tmp, inv_b2_tmp;
333  {
334  rhoface = myhalf * (cell_data(i,j,klo,Rho_comp) + cell_data(i-ioff,j-joff,klo,Rho_comp));
335  getRhoAlphaForFaces(i, j, klo, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
336  cell_data, mu_turb, mu_eff,
337  l_consA, l_turb);
338 
339  dz_inv = one / dz_ptr[klo];
340  dz_inv_lo = dz_inv;
341  dz_inv_hi = two / (dz_ptr[klo] + dz_ptr[klo+1]);
342 
343  a_tmp = zero;
344  c_tmp = -Fact * gfac * rhoAlpha_hi * dz_inv_hi * dz_inv;
345 
346  RHS_a(i,j,klo) = face_data(i,j,klo); // NOTE: this is momenta; solution is velocity
347 
348  // BCs: Dirichlet (u_i = val), slip wall (w = 0), or surface layer (w = 0)
349  if (ext_dir_on_zlo) {
350  RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau_corr(i,j,klo)) * dz_inv;
351  if (stagdir==2) {
352  c_tmp = 0.;
353  RHS_a(i,j,klo) = 0.;
354  } else {
355  a_tmp = -two * Fact * rhoAlpha_lo * dz_inv_lo * dz_inv;
356  RHS_a(i,j,klo) += two * rhoAlpha_lo * face_data(i,j,klo-1) * dz_inv_lo * dz_inv;
357  }
358  } else if (use_SurfLayer) {
359  // NOTE: tau = -mu*d_z(u_i) w/ SL
360  RHS_a(i,j,klo) += Fact * gfac * (tau_corr(i,j,klo+1) - tau(i,j,klo)) * dz_inv;
361  RHS_a(i,j,klo) += Fact * dz_inv * tau(i,j,klo);
362  }
363 
364  b_tmp = rhoface - a_tmp - c_tmp;
365  inv_b2_tmp = one;
366 
367  RHS_a(i,j,klo) /= b_tmp; // NOTE: this is now "rho"
368  coeffG_a(i,j,klo) = c_tmp / b_tmp; // NOTE: this is now "gamma"
369  }
370 
371  // Build the coefficients and RHS for L decomp
372  //===================================================
373  for (int k(klo+1); k < khi; k++) {
374  rhoface = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i-ioff,j-joff,k,Rho_comp));
375  getRhoAlphaForFaces(i, j, k, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
376  cell_data, mu_turb, mu_eff,
377  l_consA, l_turb);
378 
379  dz_inv = one / dz_ptr[k];
380  dz_inv_lo = two / (dz_ptr[k] + dz_ptr[k-1]);
381  dz_inv_hi = two / (dz_ptr[k] + dz_ptr[k+1]);
382 
383  a_tmp = -Fact * rhoAlpha_lo * dz_inv_lo * dz_inv;
384  c_tmp = -Fact * rhoAlpha_hi * dz_inv_hi * dz_inv;
385  b_tmp = rhoface - a_tmp - c_tmp;
386  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,k-1));
387 
388  RHS_a(i,j,k) = face_data(i,j,k); // NOTE: this is momenta; solution is velocity
389  RHS_a(i,j,k) += Fact * gfac * (tau_corr(i,j,k+1) - tau_corr(i,j,k)) * dz_inv;
390 
391  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"
392  coeffG_a(i,j,k) = c_tmp * inv_b2_tmp; // NOTE: this is now "gamma"
393  } // k
394 
395  // Top boundary coefficients and RHS for L decomp
396  //===================================================
397  {
398  rhoface = myhalf * (cell_data(i,j,khi,Rho_comp) + cell_data(i-ioff,j-joff,khi,Rho_comp));
399  getRhoAlphaForFaces(i, j, khi, ioff, joff, rhoAlpha_lo, rhoAlpha_hi,
400  cell_data, mu_turb, mu_eff,
401  l_consA, l_turb);
402 
403  dz_inv = one / dz_ptr[khi];
404  dz_inv_lo = two / (dz_ptr[khi] + dz_ptr[khi-1]);
405  dz_inv_hi = dz_inv;
406 
407  a_tmp = -Fact * gfac * rhoAlpha_lo * dz_inv_lo * dz_inv;
408  c_tmp = zero;
409 
410  RHS_a(i,j,khi) = face_data(i,j,khi); // NOTE: this is momenta; solution is velocity
411  RHS_a(i,j,khi) += Fact * gfac * (tau_corr(i,j,khi+1) - tau_corr(i,j,khi)) * dz_inv;
412 
413  // BCs: Dirichlet (u_i = val), slip wall (w = 0)
414  if (ext_dir_on_zhi) {
415  if (stagdir==2) {
416  a_tmp = zero;
417  RHS_a(i,j,khi) = zero;
418  } else {
419  c_tmp = -two * Fact * rhoAlpha_hi * dz_inv_hi * dz_inv;
420  RHS_a(i,j,khi) += two * rhoAlpha_hi * face_data(i,j,khi+1) * dz_inv_hi * dz_inv;
421  }
422  }
423 
424  b_tmp = rhoface - a_tmp - c_tmp;
425  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,khi-1));
426 
427  // First solve
428  soln_a(i,j,khi) = (RHS_a(i,j,khi) - a_tmp * RHS_a(i,j,khi-1)) * inv_b2_tmp;
429  }
430 
431  // Back sweep the U decomp solution
432  //===================================================
433  for (int k(khi-1); k>=klo; --k) {
434  soln_a(i,j,k) = RHS_a(i,j,k) - coeffG_a(i,j,k) * soln_a(i,j,k+1);
435  }
436 
437  // Convert back to momenta
438  //===================================================
439  for (int k(klo); k<=khi; ++k) {
440  rhoface = myhalf * (cell_data(i,j,k,Rho_comp) + cell_data(i-ioff,j-joff,k,Rho_comp));
441  face_data(i,j,k) = rhoface * soln_a(i,j,k);
442  }
443 
444 #ifdef AMREX_USE_GPU
445  });
446 #else
447  } // i
448  } // j
449 #endif
450 }
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_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_S()

void ImplicitDiffForStateLU_S ( 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 Gpu::DeviceVector< Real > &  stretched_dz_d,
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 a vertically stretched grid over flat 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]stretched_dz_darray over z of dz[k]
[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
41 {
42  BL_PROFILE_VAR("ImplicitDiffForState_S()",ImplicitDiffForState_S);
43 
44  // setup quantities for getRhoAlpha()
45 #include "ERF_SetupVertDiff.H"
46  const int qty_index = n;
47  const int prim_index = qty_index - 1;
48  const int prim_scal_index = (qty_index >= RhoScalar_comp && qty_index < RhoScalar_comp+NSCALARS) ? PrimScalar_comp : prim_index;
49 
50  // Box bounds
51  int ilo = bx.smallEnd(0);
52  int ihi = bx.bigEnd(0);
53  int jlo = bx.smallEnd(1);
54  int jhi = bx.bigEnd(1);
55  int klo = bx.smallEnd(2);
56  int khi = bx.bigEnd(2);
57  amrex::ignore_unused(ilo, ihi, jlo, jhi);
58 
59  // Temporary FABs for tridiagonal solve (allocated on column)
60  // A[k] * x[k-1] + B[k] * x[k] + C[k+1] = RHS[k]
61 
62  // With LU decomposition, M * x = r is written as L * U * x = r with U * x = rho
63  // We then first have L * rho = r and U * x = rho
64  amrex::FArrayBox RHS_fab, soln_fab, coeffG_fab;
65  RHS_fab.resize(bx,1, amrex::The_Async_Arena());
66  soln_fab.resize(bx,1, amrex::The_Async_Arena());
67  coeffG_fab.resize(bx,1, amrex::The_Async_Arena());
68  auto const& RHS_a = RHS_fab.array();
69  auto const& soln_a = soln_fab.array();
70  auto const& coeffG_a = coeffG_fab.array();
71 
72  auto dz_ptr = stretched_dz_d.data();
73 
74  int bc_comp = qty_index;
75  bool foextrap_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::foextrap);
76  bool foextrap_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::foextrap);
77  bool neumann_on_zlo = (bc_ptr[bc_comp].lo(2) == ERFBCType::neumann);
78  bool neumann_on_zhi = (bc_ptr[bc_comp].hi(2) == ERFBCType::neumann);
79  amrex::ignore_unused(foextrap_on_zlo, foextrap_on_zhi);
80 
81  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zlo || neumann_on_zlo || use_SurfLayer,
82  "Unexpected lower BC used with implicit vertical diffusion");
83  AMREX_ASSERT_WITH_MESSAGE(foextrap_on_zhi || neumann_on_zhi,
84  "Unexpected upper BC used with implicit vertical diffusion");
85 
86  Real Fact = implicit_fac * dt;
87 
88 #ifdef AMREX_USE_GPU
89  ParallelFor(makeSlab(bx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int)
90  {
91 #else
92  for (int j(jlo); j<=jhi; ++j) {
93  for (int i(ilo); i<=ihi; ++i) {
94 #endif
95  // Bottom boundary coefficients and RHS for L decomp
96  //===================================================
97  Real rhoAlpha_lo, rhoAlpha_hi;
98  Real dz_inv, dz_inv_lo, dz_inv_hi;
99  Real a_tmp, b_tmp, c_tmp, inv_b2_tmp;
100  {
101  getRhoAlpha(i, j, klo, rhoAlpha_lo, rhoAlpha_hi,
102  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
103  prim_index, prim_scal_index, l_consA, l_turb);
104 
105  dz_inv = one / dz_ptr[klo];
106  dz_inv_lo = dz_inv;
107  dz_inv_hi = two / (dz_ptr[klo] + dz_ptr[klo+1]);
108 
109  a_tmp = zero;
110  c_tmp = -Fact * rhoAlpha_hi * dz_inv_hi * dz_inv;
111  b_tmp = cell_data(i,j,klo,Rho_comp) - a_tmp - c_tmp;
112  inv_b2_tmp = one;
113 
114  RHS_a(i,j,klo) = cell_data(i,j,klo,n); // NOTE: this is rho*phi; solution is phi
115  if (use_SurfLayer && scalar_zflux) {
116  RHS_a(i,j,klo) += Fact * dz_inv * scalar_zflux(i,j,klo); // NOTE: scalar_zflux = -K*d_z(\phi)
117  } else if (neumann_on_zlo) {
118  RHS_a(i,j,klo) += -Fact * dz_inv * rhoAlpha_lo * bc_neumann_vals[2]; // NOTE: N_val = d_z(\phi)
119  }
120 
121  RHS_a(i,j,klo) /= b_tmp; // NOTE: this is now "rho"
122  coeffG_a(i,j,klo) = c_tmp / b_tmp; // NOTE: this is now "gamma"
123  }
124 
125  // Build the coefficients and RHS for L decomp
126  //===================================================
127  for (int k(klo+1); k < khi; k++) {
128  getRhoAlpha(i, j, k, rhoAlpha_lo, rhoAlpha_hi,
129  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
130  prim_index, prim_scal_index, l_consA, l_turb);
131 
132  dz_inv = one / dz_ptr[k];
133  dz_inv_lo = two / (dz_ptr[k] + dz_ptr[k-1]);
134  dz_inv_hi = two / (dz_ptr[k] + dz_ptr[k+1]);
135 
136  a_tmp = -Fact * rhoAlpha_lo * dz_inv_lo * dz_inv;
137  c_tmp = -Fact * rhoAlpha_hi * dz_inv_hi * dz_inv;
138  b_tmp = cell_data(i,j,k,Rho_comp) - a_tmp - c_tmp;
139  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,k-1));
140 
141  RHS_a(i,j,k) = cell_data(i,j,k,n); // NOTE: this is rho*phi; solution is phi
142 
143  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"
144  coeffG_a(i,j,k) = c_tmp * inv_b2_tmp; // NOTE: this is now "gamma"
145  } // k
146 
147  // Top boundary coefficients and RHS for L decomp
148  //===================================================
149  {
150  getRhoAlpha(i, j, khi, rhoAlpha_lo, rhoAlpha_hi,
151  cell_data, mu_turb, d_alpha_eff, d_eddy_diff_idz,
152  prim_index, prim_scal_index, l_consA, l_turb);
153 
154  dz_inv = one / dz_ptr[khi];
155  dz_inv_lo = two / (dz_ptr[khi] + dz_ptr[khi-1]);
156  dz_inv_hi = dz_inv;
157 
158  a_tmp = -Fact * rhoAlpha_lo * dz_inv_lo * dz_inv;
159  c_tmp = zero;
160  b_tmp = cell_data(i,j,khi,Rho_comp) - a_tmp - c_tmp;
161  inv_b2_tmp = one / (b_tmp - a_tmp * coeffG_a(i,j,khi-1));
162 
163  RHS_a(i,j,khi) = cell_data(i,j,khi,n); // NOTE: this is rho*phi; solution is phi
164  if (neumann_on_zhi) {
165  RHS_a(i,j,khi) -= -Fact * dz_inv * rhoAlpha_hi * bc_neumann_vals[5]; // NOTE: N_val = d_z(\phi)
166  }
167 
168  // First solve
169  soln_a(i,j,khi) = (RHS_a(i,j,khi) - a_tmp * RHS_a(i,j,khi-1)) * inv_b2_tmp;
170  }
171 
172  // Back sweep the U decomp solution
173  //===================================================
174  for (int k(khi-1); k>=klo; --k) {
175  soln_a(i,j,k) = RHS_a(i,j,k) - coeffG_a(i,j,k) * soln_a(i,j,k+1);
176  }
177 
178  // Convert back to rho*theta
179  //===================================================
180  for (int k(klo); k<=khi; ++k) {
181  cell_data(i,j,k,n) = cell_data(i,j,k,Rho_comp) * soln_a(i,j,k);
182  }
183 
184 #ifdef AMREX_USE_GPU
185  });
186 #else
187  } // i
188  } // j
189 #endif
190 }
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: