ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERFFillPatcher Class Reference

#include <ERF_FillPatcher.H>

Collaboration diagram for ERFFillPatcher:

Public Member Functions

 ERFFillPatcher (amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
 
void Define (amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
 
void BuildMask (amrex::BoxArray const &fba, int nghost, int nghost_set)
 
void RegisterCoarseData (amrex::Vector< amrex::MultiFab const * > const &crse_data, amrex::Vector< amrex::Real > const &crse_time)
 
void InterpFace (amrex::MultiFab &fine, amrex::MultiFab const &crse, int mask_val)
 
void InterpCell (amrex::MultiFab &fine, amrex::MultiFab const &crse, amrex::Vector< amrex::BCRec > const &bcr, int mask_val)
 
int GetSetMaskVal ()
 
int GetRelaxMaskVal ()
 
amrex::iMultiFab * GetMask ()
 
template<typename BC >
void FillSet (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
 
template<typename BC >
void FillRelax (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
 
template<typename BC >
void Fill (amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs, int mask_val)
 

Private Attributes

amrex::BoxArray m_fba
 
amrex::BoxArray m_cba
 
amrex::DistributionMapping m_fdm
 
amrex::DistributionMapping m_cdm
 
amrex::Geometry m_fgeom
 
amrex::Geometry m_cgeom
 
int m_nghost
 
int m_nghost_subset
 
int m_ncomp
 
amrex::InterpBase * m_interp
 
amrex::IntVect m_ratio
 
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_old
 
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_new
 
std::unique_ptr< amrex::iMultiFab > m_cf_mask
 
amrex::Vector< amrex::Real > m_crse_times
 
amrex::Real m_dt_crse
 
int m_set_mask {2}
 
int m_relax_mask {1}
 

Constructor & Destructor Documentation

◆ ERFFillPatcher()

ERFFillPatcher::ERFFillPatcher ( amrex::BoxArray const &  fba,
amrex::DistributionMapping const &  fdm,
amrex::Geometry const &  fgeom,
amrex::BoxArray const &  cba,
amrex::DistributionMapping const &  cdm,
amrex::Geometry const &  cgeom,
int  nghost,
int  nghost_set,
int  ncomp,
amrex::InterpBase *  interp 
)
27  : m_fba(fba),
28  m_cba(cba),
29  m_fdm(fdm),
30  m_cdm(cdm),
31  m_fgeom(fgeom),
32  m_cgeom(cgeom)
33 {
34  AMREX_ALWAYS_ASSERT(fba.ixType() == cba.ixType());
35 
36  // Vector to hold times for coarse data
37  m_crse_times.resize(2);
38 
39  // Define the coarse and fine MFs
40  Define(fba, fdm, fgeom, cba, cdm, cgeom,nghost, nghost_set, ncomp, interp);
41 }
amrex::Geometry m_fgeom
Definition: ERF_FillPatcher.H:63
amrex::BoxArray m_fba
Definition: ERF_FillPatcher.H:59
amrex::BoxArray m_cba
Definition: ERF_FillPatcher.H:60
void Define(amrex::BoxArray const &fba, amrex::DistributionMapping const &fdm, amrex::Geometry const &fgeom, amrex::BoxArray const &cba, amrex::DistributionMapping const &cdm, amrex::Geometry const &cgeom, int nghost, int nghost_set, int ncomp, amrex::InterpBase *interp)
Definition: ERF_FillPatcher.cpp:57
amrex::DistributionMapping m_cdm
Definition: ERF_FillPatcher.H:62
amrex::Vector< amrex::Real > m_crse_times
Definition: ERF_FillPatcher.H:73
amrex::Geometry m_cgeom
Definition: ERF_FillPatcher.H:64
amrex::DistributionMapping m_fdm
Definition: ERF_FillPatcher.H:61
Here is the call graph for this function:

Member Function Documentation

◆ BuildMask()

void ERFFillPatcher::BuildMask ( amrex::BoxArray const &  fba,
int  nghost,
int  nghost_set 
)
122 {
123  // Minimal bounding box of fine BA plus a halo cell
124  Box fba_bnd = grow(fba.minimalBox(), IntVect(1,1,1));
125 
126  // BoxList and BoxArray to store complement
127  BoxList com_bl; BoxArray com_ba;
128 
129  // Compute the complement
130  fba.complementIn(com_bl,fba_bnd);
131 
132  // com_bl cannot be null since we grew with halo cells
133  AMREX_ALWAYS_ASSERT(com_bl.size() > 0);
134 
135  IntVect box_grow_vect(-nghost,-nghost,0);
136 
137  // cf_set_width = cf_width = 0 is a special case
138  // In this case we set only the normal velocities
139  // (not any cell-centered quantities) and only
140  // on the coarse-fine boundary itself
141  if (nghost == 0) {
142  if (fba.ixType()[0] == IndexType::NODE) {
143  box_grow_vect = IntVect(1,0,0);
144  } else if (fba.ixType()[1] == IndexType::NODE) {
145  box_grow_vect = IntVect(0,1,0);
146  } else if (fba.ixType()[2] == IndexType::NODE) {
147  box_grow_vect = IntVect(0,0,1);
148  }
149  }
150 
151  // Grow the complement boxes and trim with the bounding box
152  Vector<Box>& com_bl_v = com_bl.data();
153  for (int i(0); i<com_bl.size(); ++i) {
154  Box& bx = com_bl_v[i];
155  bx.grow(box_grow_vect);
156  bx &= fba_bnd;
157  }
158 
159 
160  // Do second complement with the grown boxes
161  com_ba.define(std::move(com_bl));
162  com_ba.complementIn(com_bl, fba_bnd);
163 
164  // Fill mask based upon the com_bl BoxList
165  for (MFIter mfi(*m_cf_mask); mfi.isValid(); ++mfi) {
166  const Box& vbx = mfi.validbox();
167  const Array4<int>& mask_arr = m_cf_mask->array(mfi);
168 
169  for (auto const& b : com_bl) {
170  Box com_bx = vbx & b;
171  ParallelFor(com_bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
172  {
173  mask_arr(i,j,k) = mask_val;
174  });
175  }
176  }
177 }
std::unique_ptr< amrex::iMultiFab > m_cf_mask
Definition: ERF_FillPatcher.H:72

Referenced by Define().

Here is the caller graph for this function:

◆ Define()

void ERFFillPatcher::Define ( amrex::BoxArray const &  fba,
amrex::DistributionMapping const &  fdm,
amrex::Geometry const &  fgeom,
amrex::BoxArray const &  cba,
amrex::DistributionMapping const &  cdm,
amrex::Geometry const &  cgeom,
int  nghost,
int  nghost_set,
int  ncomp,
amrex::InterpBase *  interp 
)
63 {
64  AMREX_ALWAYS_ASSERT(nghost <= 0);
65  AMREX_ALWAYS_ASSERT(nghost_set <= 0);
66  AMREX_ALWAYS_ASSERT(nghost <= nghost_set);
67 
68  // Set data members
69  m_fba = fba; m_cba = cba;
70  m_fdm = fdm; m_cdm = cdm;
71  m_fgeom = fgeom; m_cgeom = cgeom;
72  m_nghost = nghost; m_nghost_subset = nghost_set;
73  m_ncomp = ncomp; m_interp = interp;
74 
75  // Delete old MFs if they exist
78  if (m_cf_mask) m_cf_mask.reset();
79 
80  // Index type for the BL/BA
81  IndexType m_ixt = fba.ixType();
82 
83  // Refinement ratios
84  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
85  m_ratio[idim] = m_fgeom.Domain().length(idim) / m_cgeom.Domain().length(idim);
86  }
87 
88  // Coarse box list
89  // NOTE: if we use face_cons_linear_interp then CoarseBox returns the grown box
90  // so we don't need to manually grow it here
91  BoxList cbl;
92  cbl.set(m_ixt);
93  cbl.reserve(fba.size());
94  for (int i(0); i < fba.size(); ++i) {
95  Box coarse_box(interp->CoarseBox(fba[i], m_ratio));
96  cbl.push_back(coarse_box);
97  }
98 
99  // Box arrays for the coarse data
100  BoxArray cf_cba(std::move(cbl));
101 
102  // Two coarse patches to hold the data to be interpolated
103  m_cf_crse_data_old = std::make_unique<MultiFab> (cf_cba, fdm, m_ncomp, 0);
104  m_cf_crse_data_new = std::make_unique<MultiFab> (cf_cba, fdm, m_ncomp, 0);
105 
106  // Integer masking array
107  m_cf_mask = std::make_unique<iMultiFab> (fba, fdm, 1, 0);
108 
109  // Populate mask array
110  if (nghost_set <= 0) {
111  m_cf_mask->setVal(m_set_mask);
112  BuildMask(fba,nghost_set,m_set_mask-1);
113  } else {
114  m_cf_mask->setVal(m_relax_mask);
115  BuildMask(fba,nghost,m_relax_mask-1);
116  }
117 }
amrex::InterpBase * m_interp
Definition: ERF_FillPatcher.H:68
int m_ncomp
Definition: ERF_FillPatcher.H:67
int m_nghost_subset
Definition: ERF_FillPatcher.H:66
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_new
Definition: ERF_FillPatcher.H:71
int m_set_mask
Definition: ERF_FillPatcher.H:75
amrex::IntVect m_ratio
Definition: ERF_FillPatcher.H:69
int m_relax_mask
Definition: ERF_FillPatcher.H:76
void BuildMask(amrex::BoxArray const &fba, int nghost, int nghost_set)
Definition: ERF_FillPatcher.cpp:119
std::unique_ptr< amrex::MultiFab > m_cf_crse_data_old
Definition: ERF_FillPatcher.H:70
int m_nghost
Definition: ERF_FillPatcher.H:65

Referenced by ERFFillPatcher().

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

◆ Fill()

template<typename BC >
void ERFFillPatcher::Fill ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs,
int  mask_val 
)
124 {
125  constexpr amrex::Real eps = std::numeric_limits<float>::epsilon();
126 
127  AMREX_ALWAYS_ASSERT((time >= m_crse_times[0]-eps) && (time <= m_crse_times[1]+eps));
128 
129  // Time interpolation factors
130  amrex::Real fac_new = (time - m_crse_times[0]) / m_dt_crse;
131  amrex::Real fac_old = 1.0 - fac_new;
132 
133  // Boundary condition operator
134  cbc(*(m_cf_crse_data_old), 0, m_ncomp, amrex::IntVect(0), time, 0);
135 
136  // Coarse MF to hold time interpolated data
137  amrex::MultiFab crse_data_time_interp(m_cf_crse_data_old->boxArray(),
138  m_cf_crse_data_old->DistributionMap(),
139  m_ncomp, amrex::IntVect{0});
140 
141  // Time interpolate the coarse data
142  amrex::MultiFab::LinComb(crse_data_time_interp,
143  fac_old, *(m_cf_crse_data_old), 0,
144  fac_new, *(m_cf_crse_data_new), 0,
145  0, m_ncomp, amrex::IntVect{0});
146 
147  // Call correct spatial interpolation type
148  amrex::IndexType m_ixt = mf.boxArray().ixType();
149  int ixt_sum = m_ixt[0]+m_ixt[1]+m_ixt[2];
150  if (ixt_sum == 0) {
151  InterpCell(mf,crse_data_time_interp,bcs,mask_val);
152  } else if (ixt_sum == 1) {
153  InterpFace(mf,crse_data_time_interp,mask_val);
154  } else {
155  amrex::Abort("ERF_FillPatcher only supports face linear and cell cons linear interp!");
156  }
157 }
void InterpFace(amrex::MultiFab &fine, amrex::MultiFab const &crse, int mask_val)
Definition: ERF_FillPatcher.cpp:210
void InterpCell(amrex::MultiFab &fine, amrex::MultiFab const &crse, amrex::Vector< amrex::BCRec > const &bcr, int mask_val)
Definition: ERF_FillPatcher.cpp:328
amrex::Real m_dt_crse
Definition: ERF_FillPatcher.H:74

Referenced by FillRelax(), and FillSet().

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

◆ FillRelax()

template<typename BC >
void ERFFillPatcher::FillRelax ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs 
)
107 {
108  Fill(mf,time,cbc,bcs,m_relax_mask);
109 }
void Fill(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs, int mask_val)
Definition: ERF_FillPatcher.H:122

Referenced by fine_compute_interior_ghost_rhs().

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

◆ FillSet()

template<typename BC >
void ERFFillPatcher::FillSet ( amrex::MultiFab &  mf,
amrex::Real  time,
BC &  cbc,
amrex::Vector< amrex::BCRec > const &  bcs 
)
91 {
92  Fill(mf,time,cbc,bcs,m_set_mask);
93 }
Here is the call graph for this function:

◆ GetMask()

amrex::iMultiFab* ERFFillPatcher::GetMask ( )
inline
43 { return m_cf_mask.get(); }

Referenced by fine_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ GetRelaxMaskVal()

int ERFFillPatcher::GetRelaxMaskVal ( )
inline
41 { return m_relax_mask; }

Referenced by fine_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ GetSetMaskVal()

int ERFFillPatcher::GetSetMaskVal ( )
inline
39 { return m_set_mask; }

Referenced by fine_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ InterpCell()

void ERFFillPatcher::InterpCell ( amrex::MultiFab &  fine,
amrex::MultiFab const &  crse,
amrex::Vector< amrex::BCRec > const &  bcr,
int  mask_val 
)
332 {
333  int ncomp = m_ncomp;
334  IntVect ratio = m_ratio;
335  IndexType m_ixt = fine.boxArray().ixType();
336  Box const& cdomain = convert(m_cgeom.Domain(), m_ixt);
337 
338  for (MFIter mfi(fine); mfi.isValid(); ++mfi) {
339  Box const& fbx = mfi.validbox();
340 
341  Array4<Real> const& fine_arr = fine.array(mfi);
342  Array4<Real const> const& crse_arr = crse.const_array(mfi);
343  Array4<int const> const& mask_arr = m_cf_mask->const_array(mfi);
344 
345  bool run_on_gpu = Gpu::inLaunchRegion();
346  amrex::ignore_unused(run_on_gpu);
347 
348  amrex::ignore_unused(m_fgeom);
349 
350  const Box& crse_region = m_interp->CoarseBox(fbx,ratio);
351  Box cslope_bx(crse_region);
352  for (int dim = 0; dim < AMREX_SPACEDIM; dim++) {
353  if (ratio[dim] > 1) {
354  cslope_bx.grow(dim,-1);
355  }
356  }
357 
358  FArrayBox ccfab(cslope_bx, ncomp*AMREX_SPACEDIM, The_Async_Arena());
359  Array4<Real> const& tmp = ccfab.array();
360  Array4<Real const> const& ctmp = ccfab.const_array();
361 
362 #ifdef AMREX_USE_GPU
363  AsyncArray<BCRec> async_bcr(bcr.data(), (run_on_gpu) ? ncomp : 0);
364  BCRec const* bcrp = (run_on_gpu) ? async_bcr.data() : bcr.data();
365 #else
366  BCRec const* bcrp = bcr.data();
367 #endif
368 
369  AMREX_HOST_DEVICE_PARALLEL_FOR_4D_FLAG(RunOn::Gpu, cslope_bx, ncomp, i, j, k, n,
370  {
371  mf_cell_cons_lin_interp_mcslope(i,j,k,n, tmp, crse_arr, 0, ncomp,
372  cdomain, ratio, bcrp);
373  });
374 
375  AMREX_HOST_DEVICE_PARALLEL_FOR_4D_FLAG(RunOn::Gpu, fbx, ncomp, i, j, k, n,
376  {
377  if (mask_arr(i,j,k) == mask_val) mf_cell_cons_lin_interp(i,j,k,n, fine_arr, 0, ctmp,
378  crse_arr, 0, ncomp, ratio);
379  });
380  } // MFIter
381 }

Referenced by Fill().

Here is the caller graph for this function:

◆ InterpFace()

void ERFFillPatcher::InterpFace ( amrex::MultiFab &  fine,
amrex::MultiFab const &  crse,
int  mask_val 
)
213 {
214  int ncomp = 1;
215  IntVect ratio = m_ratio;
216 
217  FArrayBox slope;
218 
219  //
220  // This box is only used to make sure we don't look outside
221  // the domain for computing the slopes in the interpolation
222  // We need it to be of the type of the faces being filled
223  //
224  // IndexType ixt = fine.boxArray().ixType();
225  // Box const& domface = convert(m_cgeom.Domain(), ixt);
226 
227  // We don't need to worry about face-based domain because this is only used in the tangential interpolation
228  Box per_grown_domain = m_cgeom.Domain();
229  for (int dim = 0; dim < AMREX_SPACEDIM; dim++) {
230  if (m_cgeom.isPeriodic(dim)) {
231  per_grown_domain.grow(dim,1);
232  }
233  }
234 
235  for (MFIter mfi(fine); mfi.isValid(); ++mfi)
236  {
237  Box const& fbx = mfi.validbox();
238 
239  slope.resize(fbx,ncomp,The_Async_Arena());
240 
241  Array4<Real> const& fine_arr = fine.array(mfi);
242  Array4<Real> const& slope_arr = slope.array();
243  Array4<Real const> const& crse_arr = crse.const_array(mfi);
244  Array4<int const> const& mask_arr = m_cf_mask->const_array(mfi);
245 
246  if (fbx.type(0) == IndexType::NODE) // x-faces
247  {
248  // Here do interpolation in the tangential directions
249  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
250  {
251  if (mask_arr(i,j,k) == mask_val) { // x-faces
252  const int ii = coarsen(i,ratio[0]);
253  if (i-ii*ratio[0] == 0) {
254  interp_face_reg(i,j,k,ratio,fine_arr,0,crse_arr,slope_arr,ncomp,per_grown_domain,0);
255  }
256  }
257  });
258 
259  // Here do interpolation in the normal direction
260  // using the fine values that have already been filled
261  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
262  {
263  if (mask_arr(i,j,k) == mask_val) {
264  const int ii = coarsen(i,ratio[0]);
265  if (i-ii*ratio[0] != 0) {
266  Real const w = static_cast<Real>(i-ii*ratio[0]) * (Real(1.)/Real(ratio[0]));
267  fine_arr(i,j,k,0) = (Real(1.)-w) * fine_arr(ii*ratio[0],j,k,0) + w * fine_arr((ii+1)*ratio[0],j,k,0);
268  }
269  }
270  });
271 
272  }
273  else if (fbx.type(1) == IndexType::NODE) // y-faces
274  {
275  // Here do interpolation in the tangential directions
276  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
277  {
278  if (mask_arr(i,j,k) == mask_val) {
279  const int jj = coarsen(j,ratio[1]);
280  if (j-jj*ratio[1] == 0) {
281  interp_face_reg(i,j,k,ratio,fine_arr,0,crse_arr,slope_arr,ncomp,per_grown_domain,1);
282  }
283  }
284  });
285 
286  // Here do interpolation in the normal direction
287  // using the fine values that have already been filled
288  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
289  {
290  if (mask_arr(i,j,k) == mask_val) {
291  const int jj = coarsen(j,ratio[1]);
292  if (j-jj*ratio[1] != 0) {
293  Real const w = static_cast<Real>(j-jj*ratio[1]) * (Real(1.)/Real(ratio[1]));
294  fine_arr(i,j,k,0) = (Real(1.)-w) * fine_arr(i,jj*ratio[1],k,0) + w * fine_arr(i,(jj+1)*ratio[1],k,0);
295  }
296  }
297  });
298  }
299  else // z-faces
300  {
301  // Here do interpolation in the tangential directions
302  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
303  {
304  if (mask_arr(i,j,k) == mask_val) {
305  const int kk = coarsen(k,ratio[2]);
306  if (k-kk*ratio[2] == 0) {
307  interp_face_reg(i,j,k,ratio,fine_arr,0,crse_arr,slope_arr,1,per_grown_domain,2);
308  }
309  }
310  });
311 
312  // Here do interpolation in the normal direction
313  // using the fine values that have already been filled
314  AMREX_HOST_DEVICE_PARALLEL_FOR_3D_FLAG(RunOn::Gpu,fbx,i,j,k,
315  {
316  if (mask_arr(i,j,k) == mask_val) {
317  const int kk = coarsen(k,ratio[2]);
318  if (k-kk*ratio[2] != 0) {
319  Real const w = static_cast<Real>(k-kk*ratio[2]) * (Real(1.)/Real(ratio[2]));
320  fine_arr(i,j,k,0) = (Real(1.)-w) * fine_arr(i,j,kk*ratio[2],0) + w * fine_arr(i,j,(kk+1)*ratio[2],0);
321  }
322  }
323  });
324  } // IndexType::NODE
325  } // MFiter
326 }

Referenced by Fill().

Here is the caller graph for this function:

◆ RegisterCoarseData()

void ERFFillPatcher::RegisterCoarseData ( amrex::Vector< amrex::MultiFab const * > const &  crse_data,
amrex::Vector< amrex::Real > const &  crse_time 
)
188 {
189  AMREX_ALWAYS_ASSERT(crse_data.size() == 2); // old and new
190  AMREX_ALWAYS_ASSERT(crse_time[1] >= crse_time[0]);
191 
192  // NOTE: CoarseBox with CellConsLinear interpolation grows the
193  // box by 1 in all directions. This pushes the domain for
194  // m_cf_crse_data into ghost cells in the z-dir. So we need
195  // to include ghost cells for crse_data when doing the copy
196  IntVect src_ng = crse_data[0]->nGrowVect();
197  IntVect dst_ng = m_cf_crse_data_old->nGrowVect();
198 
199  m_cf_crse_data_old->ParallelCopy(*(crse_data[0]), 0, 0, m_ncomp,
200  src_ng, dst_ng, m_cgeom.periodicity()); // old data
201  m_cf_crse_data_new->ParallelCopy(*(crse_data[1]), 0, 0, m_ncomp,
202  src_ng, dst_ng, m_cgeom.periodicity()); // new data
203 
204  m_crse_times[0] = crse_time[0]; // time of "old" coarse data
205  m_crse_times[1] = crse_time[1]; // time of "new" coarse data
206 
207  m_dt_crse = crse_time[1] - crse_time[0];
208 }

Member Data Documentation

◆ m_cba

amrex::BoxArray ERFFillPatcher::m_cba
private

Referenced by Define().

◆ m_cdm

amrex::DistributionMapping ERFFillPatcher::m_cdm
private

Referenced by Define().

◆ m_cf_crse_data_new

std::unique_ptr<amrex::MultiFab> ERFFillPatcher::m_cf_crse_data_new
private

Referenced by Define(), Fill(), and RegisterCoarseData().

◆ m_cf_crse_data_old

std::unique_ptr<amrex::MultiFab> ERFFillPatcher::m_cf_crse_data_old
private

Referenced by Define(), Fill(), and RegisterCoarseData().

◆ m_cf_mask

std::unique_ptr<amrex::iMultiFab> ERFFillPatcher::m_cf_mask
private

◆ m_cgeom

amrex::Geometry ERFFillPatcher::m_cgeom
private

◆ m_crse_times

amrex::Vector<amrex::Real> ERFFillPatcher::m_crse_times
private

◆ m_dt_crse

amrex::Real ERFFillPatcher::m_dt_crse
private

Referenced by Fill(), and RegisterCoarseData().

◆ m_fba

amrex::BoxArray ERFFillPatcher::m_fba
private

Referenced by Define().

◆ m_fdm

amrex::DistributionMapping ERFFillPatcher::m_fdm
private

Referenced by Define().

◆ m_fgeom

amrex::Geometry ERFFillPatcher::m_fgeom
private

Referenced by Define(), and InterpCell().

◆ m_interp

amrex::InterpBase* ERFFillPatcher::m_interp
private

Referenced by Define(), and InterpCell().

◆ m_ncomp

int ERFFillPatcher::m_ncomp
private

◆ m_nghost

int ERFFillPatcher::m_nghost
private

Referenced by Define().

◆ m_nghost_subset

int ERFFillPatcher::m_nghost_subset
private

Referenced by Define().

◆ m_ratio

amrex::IntVect ERFFillPatcher::m_ratio
private

Referenced by Define(), InterpCell(), and InterpFace().

◆ m_relax_mask

int ERFFillPatcher::m_relax_mask {1}
private

Referenced by Define(), FillRelax(), and GetRelaxMaskVal().

◆ m_set_mask

int ERFFillPatcher::m_set_mask {2}
private

Referenced by Define(), FillSet(), and GetSetMaskVal().


The documentation for this class was generated from the following files: