ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_EBRedistribute.cpp File Reference
#include <AMReX_Config.H>
#include <AMReX_Geometry.H>
#include <ERF.H>
#include <ERF_EB.H>
#include <ERF_EBRedistribute.H>
#include "AMReX_EB_Redistribution.H"
Include dependency graph for ERF_EBRedistribute.cpp:

Functions

void redistribute_term (int ncomp, const Geometry &geom, MultiFab &result, MultiFab &result_tmp, MultiFab const &state, EBFArrayBoxFactory const &ebfact, BCRec const *bc, Real const local_dt)
 
void redistribute_term (int ncomp, const Geometry &geom, MultiFab &result, MultiFab &result_tmp, MultiFab const &state, eb_aux_ const &ebfact, BCRec const *bc, Real const local_dt)
 

Function Documentation

◆ redistribute_term() [1/2]

void redistribute_term ( int  ncomp,
const Geometry &  geom,
MultiFab &  result,
MultiFab &  result_tmp,
MultiFab const &  state,
eb_aux_ const &  ebfact,
BCRec const *  bc,
Real const  local_dt 
)
100 {
101  // ************************************************************************
102  // Redistribute result_tmp and pass out result
103  // ************************************************************************
104  AMREX_ASSERT(result.nComp() == state.nComp());
105 
106  result_tmp.FillBoundary(geom.periodicity());
107 
108 #ifdef _OPENMP
109 #pragma omp parallel if (Gpu::notInLaunchRegion())
110 #endif
111  for (MFIter mfi(state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
112  {
113  Box const& bx = mfi.tilebox();
114 
115  EBCellFlagFab const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
116  Array4<EBCellFlag const> const& flag = flagfab.const_array();
117 
118  bool regular = (flagfab.getType(amrex::grow(bx,4)) == FabType::regular);
119  bool covered = (flagfab.getType(bx) == FabType::covered);
120 
121  Array4<Real> out = result.array(mfi);
122  Array4<Real> in = result_tmp.array(mfi);
123 
124  if (!regular && !covered)
125  {
126  auto const& vfrac = ebfact.getVolFrac().const_array(mfi);
127  auto const& ccc = ebfact.getCentroid().const_array(mfi);
128 
129  auto const& apx = ebfact.getAreaFrac()[0]->const_array(mfi);
130  auto const& apy = ebfact.getAreaFrac()[1]->const_array(mfi);
131  auto const& apz = ebfact.getAreaFrac()[2]->const_array(mfi);
132 
133  auto const& fcx = ebfact.getFaceCent()[0]->const_array(mfi);
134  auto const& fcy = ebfact.getFaceCent()[1]->const_array(mfi);
135  auto const& fcz = ebfact.getFaceCent()[2]->const_array(mfi);
136 
137  Box bx_cc = bx;
138  bx_cc = bx_cc.convert(IntVect::TheZeroVector());
139  // bx_cc = bx_cc.setBig(0, bx_cc.bigEnd(0) + 1);
140  Box gbx = bx_cc; gbx.grow(3);
141 
142  FArrayBox scratch_fab(gbx,ncomp);
143  Array4<Real> scratch = scratch_fab.array();
144  Elixir eli_scratch = scratch_fab.elixir();
145 
146  // This is scratch space if calling StateRedistribute
147 
148  std::string redistribution_type = "StateRedist";
149 
150  // State redist acts on the state.
151  Array4<Real const> state_arr = state.const_array(mfi);
152  ApplyRedistribution( bx_cc, ncomp, out, in, state_arr,
153  scratch, flag,
154  apx, apy, apz, vfrac,
155  fcx, fcy, fcz, ccc,
156  bc, geom, local_dt, redistribution_type,
157  false, 2, 0.25_rt, {});
158  }
159  else
160  {
161  ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
162  {
163  out(i,j,k,n) = in(i,j,k,n);
164  });
165  }
166  } // MFIter
167 }
real(c_double), private bc
Definition: ERF_module_mp_morr_two_moment.F90:181
Here is the call graph for this function:

◆ redistribute_term() [2/2]

void redistribute_term ( int  ncomp,
const Geometry &  geom,
MultiFab &  result,
MultiFab &  result_tmp,
MultiFab const &  state,
EBFArrayBoxFactory const &  ebfact,
BCRec const *  bc,
Real const  local_dt 
)
21 {
22  // ************************************************************************
23  // Redistribute result_tmp and pass out result
24  // ************************************************************************
25  AMREX_ASSERT(result.nComp() == state.nComp());
26 
27  result_tmp.FillBoundary(geom.periodicity());
28 
29 #ifdef _OPENMP
30 #pragma omp parallel if (Gpu::notInLaunchRegion())
31 #endif
32  for (MFIter mfi(state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
33  {
34  Box const& bx = mfi.tilebox();
35 
36  EBCellFlagFab const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
37  Array4<EBCellFlag const> const& flag = flagfab.const_array();
38 
39  bool regular = (flagfab.getType(amrex::grow(bx,4)) == FabType::regular);
40  bool covered = (flagfab.getType(bx) == FabType::covered);
41 
42  Array4<Real> out = result.array(mfi);
43  Array4<Real> in = result_tmp.array(mfi);
44 
45  if (!regular && !covered)
46  {
47  auto const& vfrac = ebfact.getVolFrac().const_array(mfi);
48  auto const& ccc = ebfact.getCentroid().const_array(mfi);
49 
50  auto const& apx = ebfact.getAreaFrac()[0]->const_array(mfi);
51  auto const& apy = ebfact.getAreaFrac()[1]->const_array(mfi);
52  auto const& apz = ebfact.getAreaFrac()[2]->const_array(mfi);
53 
54  auto const& fcx = ebfact.getFaceCent()[0]->const_array(mfi);
55  auto const& fcy = ebfact.getFaceCent()[1]->const_array(mfi);
56  auto const& fcz = ebfact.getFaceCent()[2]->const_array(mfi);
57 
58  Box gbx = bx; gbx.grow(3);
59 
60  FArrayBox scratch_fab(gbx,ncomp);
61  Array4<Real> scratch = scratch_fab.array();
62  Elixir eli_scratch = scratch_fab.elixir();
63 
64  // This is scratch space if calling StateRedistribute
65  // ParallelFor(Box(scratch), [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
66  // {
67  // scratch(i,j,k) = 1.;
68  // });
69 
70  std::string redistribution_type = "StateRedist";
71 
72  // State redist acts on the state.
73  Array4<Real const> state_arr = state.const_array(mfi);
74  ApplyRedistribution( bx, ncomp, out, in, state_arr,
75  scratch, flag,
76  apx, apy, apz, vfrac,
77  fcx, fcy, fcz, ccc,
78  bc, geom, local_dt, redistribution_type,
79  false, 2, 0.25_rt, {});
80  }
81  else
82  {
83  ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
84  {
85  out(i,j,k,n) = in(i,j,k,n);
86  });
87  }
88  } // MFIter
89 }

Referenced by erf_slow_rhs_post().

Here is the caller graph for this function: