ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
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, int const igrid)
 

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,
int const  igrid 
)
96 {
97  BL_PROFILE_VAR("redistribute_term2", redistribute_term2);
98  // ************************************************************************
99  // Redistribute result_tmp and pass out result
100  // ************************************************************************
101  AMREX_ASSERT(result.nComp() == state.nComp());
102 
103  result_tmp.FillBoundary(geom.periodicity());
104 
105 #ifdef _OPENMP
106 #pragma omp parallel if (Gpu::notInLaunchRegion())
107 #endif
108  for (MFIter mfi(state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
109  {
110  Box const& bx = mfi.tilebox();
111 
112  EBCellFlagFab const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
113  Array4<EBCellFlag const> const& flag = flagfab.const_array();
114 
115  bool regular = (flagfab.getType(amrex::grow(bx,4)) == FabType::regular);
116  bool covered = (flagfab.getType(bx) == FabType::covered);
117 
118  Array4<Real> out = result.array(mfi);
119  Array4<Real> in = result_tmp.array(mfi);
120 
121  if (!regular && !covered)
122  {
123  auto const& vfrac = ebfact.getVolFrac().const_array(mfi);
124  auto const& ccc = ebfact.getCentroid().const_array(mfi);
125 
126  auto const& apx = ebfact.getAreaFrac()[0]->const_array(mfi);
127  auto const& apy = ebfact.getAreaFrac()[1]->const_array(mfi);
128  auto const& apz = ebfact.getAreaFrac()[2]->const_array(mfi);
129 
130  auto const& fcx = ebfact.getFaceCent()[0]->const_array(mfi);
131  auto const& fcy = ebfact.getFaceCent()[1]->const_array(mfi);
132  auto const& fcz = ebfact.getFaceCent()[2]->const_array(mfi);
133 
134  Box bx_cc = bx;
135  bx_cc = bx_cc.convert(IntVect::TheZeroVector());
136 
137  // Extend box for staggered grids
138  if (igrid == IntVars::xmom) bx_cc = bx_cc.setBig(0, bx_cc.bigEnd(0) + 1);
139  if (igrid == IntVars::ymom) bx_cc = bx_cc.setBig(1, bx_cc.bigEnd(1) + 1);
140  if (igrid == IntVars::zmom) bx_cc = bx_cc.setBig(2, bx_cc.bigEnd(2) + 1);
141 
142  Box gbx = bx_cc; gbx.grow(3);
143 
144  // Extended geometry domain
145  Box domain_grown = geom.Domain();
146  domain_grown.grow(igrid-1, 1); // Extend geometry domain by 1 in the staggering direction
147  Geometry geom_new(domain_grown, geom.ProbDomain(), geom.Coord(), geom.isPeriodic());
148 
149  FArrayBox scratch_fab(gbx,ncomp);
150  Array4<Real> scratch = scratch_fab.array();
151  Elixir eli_scratch = scratch_fab.elixir();
152 
153  // This is scratch space if calling StateRedistribute
154 
155  std::string redistribution_type = "StateRedist";
156 
157  // State redist acts on the state.
158  Array4<Real const> state_arr = state.const_array(mfi);
159  ApplyRedistribution( bx_cc, ncomp, out, in, state_arr,
160  scratch, flag,
161  apx, apy, apz, vfrac,
162  fcx, fcy, fcz, ccc,
163  bc, geom_new, local_dt, redistribution_type,
164  false, 2, 0.5_rt, {});
165  }
166  else
167  {
168  ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
169  {
170  out(i,j,k,n) = in(i,j,k,n);
171  });
172  }
173  } // MFIter
174 }
struct @19 out
struct @19 in
@ ymom
Definition: ERF_IndexDefines.H:160
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
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  BL_PROFILE_VAR("redistribute_term1", redistribute_term1);
23  // ************************************************************************
24  // Redistribute result_tmp and pass out result
25  // ************************************************************************
26  AMREX_ASSERT(result.nComp() == state.nComp());
27 
28  result_tmp.FillBoundary(geom.periodicity());
29 
30 #ifdef _OPENMP
31 #pragma omp parallel if (Gpu::notInLaunchRegion())
32 #endif
33  for (MFIter mfi(state,TilingIfNotGPU()); mfi.isValid(); ++mfi)
34  {
35  Box const& bx = mfi.tilebox();
36 
37  EBCellFlagFab const& flagfab = ebfact.getMultiEBCellFlagFab()[mfi];
38  Array4<EBCellFlag const> const& flag = flagfab.const_array();
39 
40  bool regular = (flagfab.getType(amrex::grow(bx,4)) == FabType::regular);
41  bool covered = (flagfab.getType(bx) == FabType::covered);
42 
43  Array4<Real> out = result.array(mfi);
44  Array4<Real> in = result_tmp.array(mfi);
45 
46  if (!regular && !covered)
47  {
48  auto const& vfrac = ebfact.getVolFrac().const_array(mfi);
49  auto const& ccc = ebfact.getCentroid().const_array(mfi);
50 
51  auto const& apx = ebfact.getAreaFrac()[0]->const_array(mfi);
52  auto const& apy = ebfact.getAreaFrac()[1]->const_array(mfi);
53  auto const& apz = ebfact.getAreaFrac()[2]->const_array(mfi);
54 
55  auto const& fcx = ebfact.getFaceCent()[0]->const_array(mfi);
56  auto const& fcy = ebfact.getFaceCent()[1]->const_array(mfi);
57  auto const& fcz = ebfact.getFaceCent()[2]->const_array(mfi);
58 
59  Box gbx = bx; gbx.grow(3);
60 
61  FArrayBox scratch_fab(gbx,ncomp);
62  Array4<Real> scratch = scratch_fab.array();
63  Elixir eli_scratch = scratch_fab.elixir();
64 
65  std::string redistribution_type = "StateRedist";
66 
67  // State redist acts on the state.
68  Array4<Real const> state_arr = state.const_array(mfi);
69  ApplyRedistribution( bx, ncomp, out, in, state_arr,
70  scratch, flag,
71  apx, apy, apz, vfrac,
72  fcx, fcy, fcz, ccc,
73  bc, geom, local_dt, redistribution_type,
74  false, 2, 0.5_rt, {});
75  }
76  else
77  {
78  ParallelFor(bx, ncomp, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
79  {
80  out(i,j,k,n) = in(i,j,k,n);
81  });
82  }
83  } // MFIter
84 }