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

Macros

#define INSTANTIATE_IMPLICIT_DIFF_FOR_MOM_LU(STAGDIR)
 

Functions

void ImplicitDiffForStateLU_N (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 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_N (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 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_N<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 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_N()

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

void ImplicitDiffForStateLU_N ( 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 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 a uniform grid, no terrain, and LU decomposition.

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