ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_FillZeroAreaFaceFluxes.cpp File Reference
#include "ERF.H"
#include "ERF_Utils.H"
Include dependency graph for ERF_FillZeroAreaFaceFluxes.cpp:

Functions

void FillZeroAreaFaceFluxes (MultiFab &phi, Array< MultiFab, AMREX_SPACEDIM > &fluxes, const Geometry &geom, EBFArrayBoxFactory const &ebfact, eb_aux_ const &ebfact_u, eb_aux_ const &ebfact_v, eb_aux_ const &ebfact_w)
 

Function Documentation

◆ FillZeroAreaFaceFluxes()

void FillZeroAreaFaceFluxes ( MultiFab &  phi,
Array< MultiFab, AMREX_SPACEDIM > &  fluxes,
const Geometry &  geom,
EBFArrayBoxFactory const &  ebfact,
eb_aux_ const &  ebfact_u,
eb_aux_ const &  ebfact_v,
eb_aux_ const &  ebfact_w 
)

Compute phi gradients where the area of the face is zero

13 {
14  BL_PROFILE("ERF::FillZeroAreaFaceFluxes()");
15 
16  const GpuArray<Real, AMREX_SPACEDIM> dxInv = geom.InvCellSizeArray();
17 
18  for (MFIter mfi(phi,TileNoZ()); mfi.isValid(); ++mfi)
19  {
20  const Box& tbx = mfi.tilebox();
21  const Box& xbx = mfi.nodaltilebox(0);
22  const Box& ybx = mfi.nodaltilebox(1);
23  const Box& zbx = mfi.nodaltilebox(2);
24 
25  // EBCellFlagFab const& cflag_fab = (ebfact.get_const_factory())->getMultiEBCellFlagFab()[mfi];
26  EBCellFlagFab const& cflag_fab = (ebfact.getMultiEBCellFlagFab())[mfi];
27  Array4<const EBCellFlag> cflag = cflag_fab.const_array();
28 
29  if (cflag_fab.getType(tbx) == FabType::singlevalued)
30  {
31  Array4<const Real> apx = ebfact.getAreaFrac()[0]->const_array(mfi);
32  Array4<const Real> apy = ebfact.getAreaFrac()[1]->const_array(mfi);
33  Array4<const Real> apz = ebfact.getAreaFrac()[2]->const_array(mfi);
34 
35  Array4<const EBCellFlag> u_cflag = ebfact_u.getMultiEBCellFlagFab()[mfi].const_array();
36  Array4<const EBCellFlag> v_cflag = ebfact_v.getMultiEBCellFlagFab()[mfi].const_array();
37  Array4<const EBCellFlag> w_cflag = ebfact_w.getMultiEBCellFlagFab()[mfi].const_array();
38 
39  Array4<Real const> const& p_arr = phi.const_array(mfi);
40  Array4<Real> const& fx = fluxes[0].array(mfi);
41  Array4<Real> const& fy = fluxes[1].array(mfi);
42  Array4<Real> const& fz = fluxes[2].array(mfi);
43 
44  ParallelFor(xbx, ybx, zbx,
45  // x-face
46  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
47  {
48  if (apx(i,j,k) == Real(0.0)) {
49  if (!u_cflag(i,j,k).isCovered()) {
50  if (cflag(i,j,k).isCovered() && !cflag(i-1,j,k).isCovered()) {
51  fx(i,j,k) = dxInv[0] * (p_arr(i-3,j,k) - 3.*p_arr(i-2,j,k) + 2.*p_arr(i-1,j,k));
52  } else if (cflag(i-1,j,k).isCovered() && !cflag(i,j,k).isCovered()) {
53  fx(i,j,k) = dxInv[0] * (3.*p_arr(i+1,j,k) - p_arr(i+2,j,k) - 2.*p_arr(i,j,k));
54  }
55  }
56  }
57  },
58  // y-face
59  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
60  {
61  if (apy(i,j,k) == Real(0.0)) {
62  if (!v_cflag(i,j,k).isCovered()) {
63  if (cflag(i,j,k).isCovered() && !cflag(i,j-1,k).isCovered()) {
64  fy(i,j,k) = dxInv[1] * (p_arr(i,j-3,k) - 3.*p_arr(i,j-2,k) + 2.*p_arr(i,j-1,k));
65  } else if (cflag(i,j-1,k).isCovered() && !cflag(i,j,k).isCovered()) {
66  fy(i,j,k) = dxInv[1] * (3.*p_arr(i,j+1,k) - p_arr(i,j+2,k) - 2.*p_arr(i,j,k));
67  }
68  }
69  }
70  },
71  // z-face
72  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
73  {
74  if (apz(i,j,k) == Real(0.0)) {
75  if (!w_cflag(i,j,k).isCovered()) {
76  if (cflag(i,j,k).isCovered() && !cflag(i,j,k-1).isCovered()) {
77  fz(i,j,k) = dxInv[2] * (p_arr(i,j,k-3) - 3.*p_arr(i,j,k-2) + 2.*p_arr(i,j,k-1));
78  } else if (cflag(i,j,k-1).isCovered() && !cflag(i,j,k).isCovered()) {
79  fz(i,j,k) = dxInv[2] * (3.*p_arr(i,j,k+1) - p_arr(i,j,k+2) - 2.*p_arr(i,j,k));
80  }
81  }
82  }
83  });
84  } // single-valued
85  } // mfi
86 }
const Box zbx
Definition: ERF_SetupDiff.H:9
const Box xbx
Definition: ERF_SetupDiff.H:7
const Box ybx
Definition: ERF_SetupDiff.H:8
amrex::Real Real
Definition: ERF_ShocInterface.H:19
AMREX_FORCE_INLINE amrex::IntVect TileNoZ()
Definition: ERF_TileNoZ.H:11

Referenced by solve_with_EB_mlmg().

Here is the call graph for this function:
Here is the caller graph for this function: