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_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 (MFIter const &mfi, int ncomp, const Geometry &geom, MultiFab &result, MultiFab &result_tmp, MultiFab const &state, EBFArrayBoxFactory 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,
EBFArrayBoxFactory const &  ebfact,
BCRec const *  bc,
Real const  local_dt 
)
19 {
20  // ************************************************************************
21  // Redistribute result_tmp and pass out result
22  // ************************************************************************
23  AMREX_ASSERT(result.nComp() == state.nComp());
24 
25  result_tmp.FillBoundary(geom.periodicity());
26 
27 #ifdef _OPENMP
28 #pragma omp parallel if (Gpu::notInLaunchRegion())
29 #endif
30  for (MFIter mfi(state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
31  {
32  redistribute_term(mfi, ncomp, geom, result, result_tmp, state, ebfact, bc, local_dt);
33  }
34 }
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)
Definition: ERF_EBRedistribute.cpp:11

Referenced by erf_slow_rhs_post().

Here is the caller graph for this function:

◆ redistribute_term() [2/2]

void redistribute_term ( MFIter const &  mfi,
int  ncomp,
const Geometry &  geom,
MultiFab &  result,
MultiFab &  result_tmp,
MultiFab const &  state,
EBFArrayBoxFactory const &  ebfact,
BCRec const *  bc,
Real const  local_dt 
)
46 {
47  AMREX_ASSERT(result.nComp() == state.nComp());
48 
49  Box const& bx = mfi.tilebox();
50 
51  // EBFArrayBoxFactory const& ebfact = EBFactory(lev);
52  EBCellFlagFab const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
53  Array4<EBCellFlag const> const& flag = flagfab.const_array();
54 
55  bool regular = (flagfab.getType(amrex::grow(bx,4)) == FabType::regular);
56  bool covered = (flagfab.getType(bx) == FabType::covered);
57 
58  Array4<Real> out = result.array(mfi);
59  Array4<Real> in = result_tmp.array(mfi);
60 
61  if (!regular && !covered)
62  {
63  auto const& vfrac = ebfact.getVolFrac().const_array(mfi);
64  auto const& ccc = ebfact.getCentroid().const_array(mfi);
65 
66  auto const& apx = ebfact.getAreaFrac()[0]->const_array(mfi);
67  auto const& apy = ebfact.getAreaFrac()[1]->const_array(mfi);
68  auto const& apz = ebfact.getAreaFrac()[2]->const_array(mfi);
69 
70  auto const& fcx = ebfact.getFaceCent()[0]->const_array(mfi);
71  auto const& fcy = ebfact.getFaceCent()[1]->const_array(mfi);
72  auto const& fcz = ebfact.getFaceCent()[2]->const_array(mfi);
73 
74  Box gbx = bx; gbx.grow(3);
75 
76  FArrayBox scratch_fab(gbx,ncomp);
77  Array4<Real> scratch = scratch_fab.array();
78  Elixir eli_scratch = scratch_fab.elixir();
79 
80  // This is scratch space if calling StateRedistribute
81  // ParallelFor(Box(scratch), [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
82  // {
83  // scratch(i,j,k) = 1.;
84  // });
85 
86  std::string redistribution_type = "StateRedist";
87 
88  // State redist acts on the state.
89  Array4<Real const> state_arr = state.const_array(mfi);
90  ApplyRedistribution( bx, ncomp, out, in, state_arr,
91  scratch, flag,
92  apx, apy, apz, vfrac,
93  fcx, fcy, fcz, ccc,
94  bc, geom, local_dt, redistribution_type);
95  }
96  else
97  {
98  ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
99  {
100  out(i,j,k,n) = in(i,j,k,n);
101  });
102  }
103 }