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

Referenced by erf_slow_rhs_post().

Here is the caller graph for this function: