Block-Structured AMR Software Framework
AMReX_PhysBCFunct.H
Go to the documentation of this file.
1 #ifndef AMREX_PhysBCFunct_H_
2 #define AMREX_PhysBCFunct_H_
3 #include <AMReX_Config.H>
4 
5 #include <AMReX_BCRec.H>
6 #include <AMReX_Geometry.H>
7 #include <AMReX_MultiFab.H>
8 #include <AMReX_ArrayLim.H>
9 #include <AMReX_FilCC_C.H>
10 #include <AMReX_FilND_C.H>
11 #include <AMReX_FilFC_C.H>
12 #include <AMReX_TypeTraits.H>
13 
14 namespace amrex {
15 
16 extern "C"
17 {
18  using BndryFuncDefault = void (*)(Real* data, AMREX_ARLIM_P(lo), AMREX_ARLIM_P(hi),
19  const int* dom_lo, const int* dom_hi,
20  const Real* dx, const Real* grd_lo,
21  const Real* time, const int* bc);
22  using BndryFunc3DDefault = void (*)(Real* data, const int* lo, const int* hi,
23  const int* dom_lo, const int* dom_hi,
24  const Real* dx, const Real* grd_lo,
25  const Real* time, const int* bc);
26 }
27 
28 using UserFillBox = void (*)(Box const& bx, Array4<Real> const& dest,
29  int dcomp, int numcomp,
30  GeometryData const& geom, Real time,
31  const BCRec* bcr, int bcomp,
32  int orig_comp);
33 
36 {
37 public:
38  BndryFuncArray () noexcept = default;
39  BndryFuncArray (BndryFuncDefault inFunc) noexcept : m_func(inFunc) {}
40  BndryFuncArray (BndryFunc3DDefault inFunc) noexcept : m_func3D(inFunc) {}
41 
42  void operator() (Box const& bx, FArrayBox& dest,
43  int dcomp, int numcomp,
44  Geometry const& geom, Real time,
45  const Vector<BCRec>& bcr, int bcomp,
46  int orig_comp);
47 
48  [[nodiscard]] bool RunOnGPU () const noexcept { return m_run_on_gpu; }
49 // void setRunOnGPU (bool b) noexcept { m_run_on_gpu = b; }
50 
51 protected:
54  bool m_run_on_gpu = false;
55 };
56 
62 template <class F>
64 {
65 public:
66  GpuBndryFuncFab () = default;
67  GpuBndryFuncFab (F const& a_f) : m_user_f(a_f) {}
68  GpuBndryFuncFab (F&& a_f) : m_user_f(std::move(a_f)) {}
69 
70  void operator() (Box const& bx, FArrayBox& dest,
71  int dcomp, int numcomp,
72  Geometry const& geom, Real time,
73  const Vector<BCRec>& bcr, int bcomp,
74  int orig_comp);
75 
76  template <class FF>
77  void ccfcdoit (Box const& bx, FArrayBox& dest,
78  int dcomp, int numcomp,
79  Geometry const& geom, Real time,
80  const Vector<BCRec>& bcr, int bcomp,
81  int orig_comp, FF const& fillfunc);
82 
83  void nddoit (Box const& bx, FArrayBox& dest,
84  int dcomp, int numcomp,
85  Geometry const& geom, Real time,
86  const Vector<BCRec>& bcr, int bcomp,
87  int orig_comp);
88 
89 protected:
91 };
92 
94 {
97  int, int, amrex::GeometryData const&, amrex::Real,
98  const amrex::BCRec*, int, int) const
99  {}
100 };
101 
104 {
105 public:
106  CpuBndryFuncFab () = default;
108 
109  void operator() (Box const& bx, FArrayBox& dest,
110  int dcomp, int numcomp,
111  Geometry const& geom, Real time,
112  const Vector<BCRec>& bcr, int bcomp,
113  int orig_comp);
114 
115 protected:
116  UserFillBox f_user = nullptr;
117 };
118 
120 {
121 public:
122  void operator() (MultiFab& /*mf*/, int /*dcomp*/, int /*ncomp*/, IntVect const& /*nghost*/,
123  Real /*time*/, int /*bccomp*/) {}
124 };
125 
127 {
128 public:
129 
130  PhysBCFunctUseCoarseGhost (IntVect const& a_fp1_src_ghost)
131  : fp1_src_ghost(a_fp1_src_ghost) {}
132 
133  template <typename MF, typename Interp>
134  PhysBCFunctUseCoarseGhost (MF const& cmf, IntVect const& a_nghost,
135  IntVect const& a_nghost_outside_domain,
136  IntVect const& ratio, Interp* mapper)
137  : nghost(a_nghost),
138  nghost_outside_domain(a_nghost_outside_domain),
139  cghost(cmf.nGrowVect())
140  {
141  IndexType typ = cmf.ixType();
142 
143  const auto& coarsener = mapper->BoxCoarsener(ratio);
144  Box tmp(-nghost, IntVect(32), typ);
145  Box tmp2 = coarsener.doit(tmp);
146  src_ghost = -tmp2.smallEnd();
147 
148  tmp = Box(-nghost_outside_domain, IntVect(32), typ);
149  tmp2 = coarsener.doit(tmp);
151 
154  }
155 
156  void operator() (MultiFab& /*mf*/, int /*dcomp*/, int /*ncomp*/, IntVect const& /*nghost*/,
157  Real /*time*/, int /*bccomp*/) {}
158 
159  IntVect nghost; // # of fine ghosts inside domain to be filled
160 
161  IntVect nghost_outside_domain; // # of fine ghosts outside domain to be filled
162 
163  // This is the number of coarse ghost cells needed to interpolate fine
164  // ghost cells inside the domain
166 
167  // This is the number of coarse ghost cells needed to interpolate
168  // nghost_outside_domain fine ghost cells outside the domain
170 
171  // This is the minimum number of ghost cells needed in coarse MF
173 
174  IntVect fp1_src_ghost; // Used to pass information into FillPatchSingleLevel
175 };
176 
177 template <class F>
179 {
180 public:
181  PhysBCFunct () = default;
182 
183  PhysBCFunct (const Geometry& geom, const Vector<BCRec>& bcr, F const& f)
184  : m_geom(geom), m_bcr(bcr), m_f(f)
185  {}
186 
187  PhysBCFunct (const Geometry& geom, const Vector<BCRec>& bcr, F&& f)
188  : m_geom(geom), m_bcr(bcr), m_f(std::move(f))
189  {}
190 
191  void define (const Geometry& geom, const Vector<BCRec>& bcr, F const& f) {
192  m_geom = geom; m_bcr = bcr; m_f = f;
193  }
194 
195  void define (const Geometry& geom, const Vector<BCRec>& bcr, F&& f) {
196  m_geom = geom; m_bcr = bcr; m_f = std::move(f);
197  }
198 
199  void operator() (MultiFab& mf, int icomp, int ncomp, IntVect const& nghost,
200  Real time, int bccomp)
201  {
202  if (m_geom.isAllPeriodic()) { return; }
203 
204  BL_PROFILE("PhysBCFunct::()");
205 
206  AMREX_ASSERT(nghost.allLE(mf.nGrowVect()));
207 
209  const Box& domain = m_geom.Domain();
210  Box gdomain = amrex::convert(domain, mf.boxArray().ixType());
211  for (int i = 0; i < AMREX_SPACEDIM; ++i) {
212  if (m_geom.isPeriodic(i)) {
213  gdomain.grow(i, nghost[i]);
214  }
215  }
216 
217 #ifdef AMREX_USE_OMP
218 #pragma omp parallel if (Gpu::notInLaunchRegion())
219 #endif
220  {
221  Vector<BCRec> bcrs(ncomp);
222  for (MFIter mfi(mf); mfi.isValid(); ++mfi)
223  {
224  FArrayBox& dest = mf[mfi];
225  const Box& bx = grow(mfi.validbox(),nghost);
226 
230  if (!gdomain.contains(bx))
231  {
233  amrex::setBC(bx, domain, bccomp, 0, ncomp, m_bcr, bcrs);
234 
236  m_f(bx, dest, icomp, ncomp, m_geom, time, bcrs, 0, bccomp);
237  }
238  }
239  }
240  }
241 
242  // For backward compatibility
243  void FillBoundary (MultiFab& mf, int icomp, int ncomp, IntVect const& nghost,
244  Real time, int bccomp) {
245  // NOLINTNEXTLINE(readability-suspicious-call-argument)
246  this->operator()(mf,icomp,ncomp,nghost,time,bccomp);
247  }
248 
249 private:
253 };
254 
255 template <class F>
256 void
258  int dcomp, int numcomp,
259  Geometry const& geom, Real time,
260  const Vector<BCRec>& bcr, int bcomp,
261  int orig_comp)
262 {
263 #if defined(__CUDACC__) && (__CUDACC_VER_MAJOR__ == 11) && (__CUDACC_VER_MINOR__ == 6)
264  // The compiler is allowed to always reject static_assert(false,..) even
265  // without instantiating this function. The following is to work around
266  // the issue. Now the compiler is not allowed to reject the code unless
267  // GpuBndryFuncFab::operator() is instantiated.
268  static_assert(std::is_same<F,int>::value,
269  "CUDA 11.6 bug: https://github.com/AMReX-Codes/amrex/issues/2607");
270 #endif
271  if (bx.ixType().cellCentered()) {
272  ccfcdoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp, FilccCell{});
273  } else if (bx.ixType().nodeCentered()) {
274  nddoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp);
275  } else {
276  ccfcdoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp, FilfcFace{});
277  }
278 }
279 
280 template <class F>
281 void
283  int dcomp, int numcomp,
284  Geometry const& geom, Real time,
285  const Vector<BCRec>& bcr, int bcomp,
286  int orig_comp)
287 {
288  const IntVect& len = bx.length();
289 
290  Box const& domain = amrex::convert(geom.Domain(),IntVect::TheNodeVector());
291  Box gdomain = domain;
292  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
293  if (geom.isPeriodic(idim)) {
294  gdomain.grow(idim,len[idim]);
295  }
296  }
297 
298  if (gdomain.contains(bx)) { return; }
299 
300  Array4<Real> const& fab = dest.array();
301  const auto geomdata = geom.data();
302 
303  AsyncArray<BCRec> bcr_aa(bcr.data()+bcomp, numcomp);
304  BCRec* bcr_p = bcr_aa.data();
305 
306  const auto f_user = m_user_f;
307 
308  // xlo
309  if (!geom.isPeriodic(0) && bx.smallEnd(0) < domain.smallEnd(0)) {
310  Box bndry = bx;
311  int dxlo = domain.smallEnd(0);
312  bndry.setBig(0,dxlo-1);
313  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
314  {
315  for (int n = dcomp; n < dcomp+numcomp; ++n) {
316  fab(i,j,k,n) = fab(dxlo,j,k,n);
317  }
318  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
319  bcr_p, 0, orig_comp);
320  });
321  }
322 
323  // xhi
324  if (!geom.isPeriodic(0) && bx.bigEnd(0) > domain.bigEnd(0)) {
325  Box bndry = bx;
326  int dxhi = domain.bigEnd(0);
327  bndry.setSmall(0,dxhi+1);
328  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
329  {
330  for (int n = dcomp; n < dcomp+numcomp; ++n) {
331  fab(i,j,k,n) = fab(dxhi,j,k,n);
332  }
333  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
334  bcr_p, 0, orig_comp);
335  });
336  }
337 
338 #if (AMREX_SPACEDIM >= 2)
339  // ylo
340  if (!geom.isPeriodic(1) && bx.smallEnd(1) < domain.smallEnd(1)) {
341  Box bndry = bx;
342  int dylo = domain.smallEnd(1);
343  bndry.setBig(1,dylo-1);
344  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
345  {
346  for (int n = dcomp; n < dcomp+numcomp; ++n) {
347  fab(i,j,k,n) = fab(i,dylo,k,n);
348  }
349  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
350  bcr_p, 0, orig_comp);
351  });
352  }
353 
354  // yhi
355  if (!geom.isPeriodic(1) && bx.bigEnd(1) > domain.bigEnd(1)) {
356  Box bndry = bx;
357  int dyhi = domain.bigEnd(1);
358  bndry.setSmall(1,dyhi+1);
359  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
360  {
361  for (int n = dcomp; n < dcomp+numcomp; ++n) {
362  fab(i,j,k,n) = fab(i,dyhi,k,n);
363  }
364  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
365  bcr_p, 0, orig_comp);
366  });
367  }
368 #endif
369 
370 #if (AMREX_SPACEDIM == 3)
371  // zlo
372  if (!geom.isPeriodic(2) && bx.smallEnd(2) < domain.smallEnd(2)) {
373  Box bndry = bx;
374  int dzlo = domain.smallEnd(2);
375  bndry.setBig(2,dzlo-1);
376  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
377  {
378  for (int n = dcomp; n < dcomp+numcomp; ++n) {
379  fab(i,j,k,n) = fab(i,j,dzlo,n);
380  }
381  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
382  bcr_p, 0, orig_comp);
383  });
384  }
385 
386  // zhi
387  if (!geom.isPeriodic(2) && bx.bigEnd(2) > domain.bigEnd(2)) {
388  Box bndry = bx;
389  int dzhi = domain.bigEnd(2);
390  bndry.setSmall(2,dzhi+1);
391  amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
392  {
393  for (int n = dcomp; n < dcomp+numcomp; ++n) {
394  fab(i,j,k,n) = fab(i,j,dzhi,n);
395  }
396  f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
397  bcr_p, 0, orig_comp);
398  });
399  }
400 #endif
401 }
402 
403 template <class F>
404 template <class FF>
405 void
407  int dcomp, int numcomp,
408  Geometry const& geom, Real time,
409  const Vector<BCRec>& bcr, int bcomp,
410  int orig_comp, FF const& fillfunc)
411 {
412  const IntVect& len = bx.length();
413 
414  IndexType idxType = bx.ixType();
415  Box const& domain = amrex::convert(geom.Domain(),idxType);
416  Box gdomain = domain;
417  for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
418  if (geom.isPeriodic(idim)) {
419  gdomain.grow(idim,len[idim]);
420  }
421  }
422 
423  if (gdomain.contains(bx)) { return; }
424 
425  Array4<Real> const& fab = dest.array();
426  const auto geomdata = geom.data();
427 
428 #ifdef AMREX_USE_GPU
429  AsyncArray<BCRec> bcr_aa(bcr.data()+bcomp, numcomp);
430  BCRec* bcr_p = bcr_aa.data();
431 
432  const auto f_user = m_user_f;
433 
434  // fill on the faces first
435  {
436  Array<Box,2*AMREX_SPACEDIM> dom_face_boxes
437  = { AMREX_D_DECL(amrex::convert(amrex::adjCellLo(gdomain, 0, len[0]),idxType),
438  amrex::convert(amrex::adjCellLo(gdomain, 1, len[1]),idxType),
439  amrex::convert(amrex::adjCellLo(gdomain, 2, len[2]),idxType)),
440  AMREX_D_DECL(amrex::convert(amrex::adjCellHi(gdomain, 0, len[0]),idxType),
441  amrex::convert(amrex::adjCellHi(gdomain, 1, len[1]),idxType),
442  amrex::convert(amrex::adjCellHi(gdomain, 2, len[2]),idxType)) };
443 
444  Vector<Box> face_boxes;
445  for (const Box& b : dom_face_boxes) {
446  Box tmp = b & bx;
447  if (tmp.ok()) { face_boxes.push_back(tmp); }
448  }
449  const int n_face_boxes = face_boxes.size();
450  if (n_face_boxes == 1) {
451  amrex::ParallelFor(face_boxes[0],
452  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
453  {
455  IntVect const idx(AMREX_D_DECL(i,j,k));
456  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
457  f_user(idx, fab, dcomp, numcomp, geomdata, time,
458  bcr_p, 0, orig_comp);
459  });
460  } else if (n_face_boxes > 1) {
461  AsyncArray<Box> face_boxes_aa(face_boxes.data(), n_face_boxes);
462  Box* boxes_p = face_boxes_aa.data();
463  Long ncounts = 0;
464  for (const auto& b : face_boxes) {
465  ncounts += b.numPts();
466  }
467  amrex::ParallelFor(ncounts,
468  [=] AMREX_GPU_DEVICE (Long icount) noexcept
469  {
470  const auto& idx = getCell(boxes_p, n_face_boxes, icount);
471  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
472  f_user(idx, fab, dcomp, numcomp, geomdata, time,
473  bcr_p, 0, orig_comp);
474  });
475  }
476  }
477 
478 #if (AMREX_SPACEDIM >= 2)
479  // fill on the box edges
480  {
481 #if (AMREX_SPACEDIM == 2)
482  Array<Box,4> dom_edge_boxes
483  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
484  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
485  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
486  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType) }};
487 #else
488  Array<Box,12> dom_edge_boxes
489  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
490  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
491  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
492  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
493  //
494  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
495  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
496  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
497  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
498  //
499  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
500  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType),
501  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
502  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType) }};
503 #endif
504  Vector<Box> edge_boxes;
505  for (const Box& b : dom_edge_boxes) {
506  Box tmp = b & bx;
507  if (tmp.ok()) { edge_boxes.push_back(tmp); }
508  }
509  const int n_edge_boxes = edge_boxes.size();
510  if (n_edge_boxes == 1) {
511  amrex::ParallelFor(edge_boxes[0],
512  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
513  {
515  IntVect const idx(AMREX_D_DECL(i,j,k));
516  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
517  f_user(idx, fab, dcomp, numcomp, geomdata, time,
518  bcr_p, 0, orig_comp);
519  });
520  } else if (n_edge_boxes > 1) {
521  AsyncArray<Box> edge_boxes_aa(edge_boxes.data(), n_edge_boxes);
522  Box* boxes_p = edge_boxes_aa.data();
523  Long ncounts = 0;
524  for (const auto& b : edge_boxes) {
525  ncounts += b.numPts();
526  }
527  amrex::ParallelFor(ncounts,
528  [=] AMREX_GPU_DEVICE (Long icount) noexcept
529  {
530  const auto& idx = getCell(boxes_p, n_edge_boxes, icount);
531  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
532  f_user(idx, fab, dcomp, numcomp, geomdata, time,
533  bcr_p, 0, orig_comp);
534  });
535  }
536  }
537 #endif
538 
539 #if (AMREX_SPACEDIM == 3)
540  // fill on corners
541  {
542  Array<Box,8> dom_corner_boxes
543  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
544  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
545  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
546  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
547  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
548  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
549  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
550  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType) }};
551 
552  Vector<Box> corner_boxes;
553  for (const Box& b : dom_corner_boxes) {
554  Box tmp = b & bx;
555  if (tmp.ok()) { corner_boxes.push_back(tmp); }
556  }
557  const int n_corner_boxes = corner_boxes.size();
558  if (n_corner_boxes == 1) {
559  amrex::ParallelFor(corner_boxes[0],
560  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
561  {
562  IntVect const idx(AMREX_D_DECL(i,j,k));
563  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
564  f_user(idx, fab, dcomp, numcomp, geomdata, time,
565  bcr_p, 0, orig_comp);
566  });
567  } else if (n_corner_boxes > 1) {
568  AsyncArray<Box> corner_boxes_aa(corner_boxes.data(), n_corner_boxes);
569  Box* boxes_p = corner_boxes_aa.data();
570  Long ncounts = 0;
571  for (const auto& b : corner_boxes) {
572  ncounts += b.numPts();
573  }
574  amrex::ParallelFor(ncounts,
575  [=] AMREX_GPU_DEVICE (Long icount) noexcept
576  {
577  const auto& idx = getCell(boxes_p, n_corner_boxes, icount);
578  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
579  f_user(idx, fab, dcomp, numcomp, geomdata, time,
580  bcr_p, 0, orig_comp);
581  });
582  }
583  }
584 #endif
585 
586 #else
587  BCRec const* bcr_p = bcr.data()+bcomp;
588 
589  const auto& f_user = m_user_f;
590 
591  // fill on the box faces first
592  {
593  Array<Box,2*AMREX_SPACEDIM> dom_face_boxes
594  = {{ AMREX_D_DECL(amrex::convert(amrex::adjCellLo(gdomain, 0, len[0]),idxType),
595  amrex::convert(amrex::adjCellLo(gdomain, 1, len[1]),idxType),
596  amrex::convert(amrex::adjCellLo(gdomain, 2, len[2]),idxType)),
597  AMREX_D_DECL(amrex::convert(amrex::adjCellHi(gdomain, 0, len[0]),idxType),
598  amrex::convert(amrex::adjCellHi(gdomain, 1, len[1]),idxType),
599  amrex::convert(amrex::adjCellHi(gdomain, 2, len[2]),idxType)) }};
600 
601  for (const Box& b : dom_face_boxes) {
602  Box tmp = b & bx;
603  amrex::For(tmp, [=] (int i, int j, int k) noexcept
604  {
606  IntVect const idx(AMREX_D_DECL(i,j,k));
607  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
608  f_user(idx, fab, dcomp, numcomp, geomdata, time,
609  bcr_p, 0, orig_comp);
610  });
611  }
612  }
613 
614 #if (AMREX_SPACEDIM >= 2)
615  // fill on the box edges
616  {
617 #if (AMREX_SPACEDIM == 2)
618  Array<Box,4> dom_edge_boxes
619  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
620  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
621  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
622  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType) }};
623 #else
624  Array<Box,12> dom_edge_boxes
625  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
626  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
627  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
628  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
629  //
630  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
631  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
632  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
633  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
634  //
635  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
636  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType),
637  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
638  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType) }};
639 #endif
640 
641  for (const Box& b : dom_edge_boxes) {
642  Box tmp = b & bx;
643  amrex::For(tmp, [=] (int i, int j, int k) noexcept
644  {
646  IntVect const idx(AMREX_D_DECL(i,j,k));
647  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
648  f_user(idx, fab, dcomp, numcomp, geomdata, time,
649  bcr_p, 0, orig_comp);
650  });
651  }
652  }
653 #endif
654 
655 #if (AMREX_SPACEDIM == 3)
656  // fill on box corners
657  {
658  Array<Box,8> dom_corner_boxes
659  = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
660  amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
661  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
662  amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
663  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
664  amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
665  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
666  amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType) }};
667 
668  for (const Box& b : dom_corner_boxes) {
669  Box tmp = b & bx;
670  amrex::For(tmp, [=] (int i, int j, int k) noexcept
671  {
672  IntVect const idx(AMREX_D_DECL(i,j,k));
673  fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
674  f_user(idx, fab, dcomp, numcomp, geomdata, time,
675  bcr_p, 0, orig_comp);
676  });
677  }
678  }
679 #endif
680 
681 #endif
682 }
683 
684 }
685 #endif
#define BL_PROFILE(a)
Definition: AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition: AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition: AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition: AMReX_GpuQualifiers.H:18
#define AMREX_D_PICK(a, b, c)
Definition: AMReX_SPACE.H:151
#define AMREX_D_DECL(a, b, c)
Definition: AMReX_SPACE.H:104
Boundary Condition Records. Necessary information and functions for computing boundary conditions.
Definition: AMReX_BCRec.H:17
AMREX_FORCE_INLINE Array4< T const > array() const noexcept
Definition: AMReX_BaseFab.H:379
This version calls function working on array.
Definition: AMReX_PhysBCFunct.H:36
BndryFuncArray() noexcept=default
bool m_run_on_gpu
Definition: AMReX_PhysBCFunct.H:54
void operator()(Box const &bx, FArrayBox &dest, int dcomp, int numcomp, Geometry const &geom, Real time, const Vector< BCRec > &bcr, int bcomp, int orig_comp)
Definition: AMReX_PhysBCFunct.cpp:7
BndryFuncArray(BndryFunc3DDefault inFunc) noexcept
Definition: AMReX_PhysBCFunct.H:40
BndryFunc3DDefault m_func3D
Definition: AMReX_PhysBCFunct.H:53
BndryFuncDefault m_func
Definition: AMReX_PhysBCFunct.H:52
bool RunOnGPU() const noexcept
Definition: AMReX_PhysBCFunct.H:48
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition: AMReX_BoxArray.H:836
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & smallEnd() const &noexcept
Get the smallend of the BoxND.
Definition: AMReX_Box.H:105
AMREX_GPU_HOST_DEVICE BoxND & grow(int i) noexcept
Definition: AMReX_Box.H:627
AMREX_GPU_HOST_DEVICE IndexTypeND< dim > ixType() const noexcept
Returns the indexing type.
Definition: AMReX_Box.H:127
AMREX_GPU_HOST_DEVICE const IntVectND< dim > & bigEnd() const &noexcept
Get the bigend.
Definition: AMReX_Box.H:116
AMREX_GPU_HOST_DEVICE IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition: AMReX_Box.H:146
AMREX_GPU_HOST_DEVICE BoxND & setSmall(const IntVectND< dim > &sm) noexcept
Redefine the small end of the BoxND.
Definition: AMReX_Box.H:466
AMREX_GPU_HOST_DEVICE bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition: AMReX_Box.H:200
AMREX_GPU_HOST_DEVICE BoxND & setBig(const IntVectND< dim > &bg) noexcept
Redefine the big end of the BoxND.
Definition: AMReX_Box.H:474
AMREX_GPU_HOST_DEVICE bool contains(const IntVectND< dim > &p) const noexcept
Returns true if argument is contained within BoxND.
Definition: AMReX_Box.H:204
This cpu version calls function working on FArrayBox.
Definition: AMReX_PhysBCFunct.H:104
UserFillBox f_user
Definition: AMReX_PhysBCFunct.H:116
CpuBndryFuncFab(UserFillBox a_f)
Definition: AMReX_PhysBCFunct.H:107
void operator()(Box const &bx, FArrayBox &dest, int dcomp, int numcomp, Geometry const &geom, Real time, const Vector< BCRec > &bcr, int bcomp, int orig_comp)
Definition: AMReX_PhysBCFunct.cpp:48
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
IntVect nGrowVect() const noexcept
Definition: AMReX_FabArrayBase.H:79
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition: AMReX_FabArrayBase.H:94
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition: AMReX_Geometry.H:210
GeometryData data() const noexcept
Returns non-static copy of geometry's stored data.
Definition: AMReX_Geometry.H:123
bool isAllPeriodic() const noexcept
Is domain periodic in all directions?
Definition: AMReX_Geometry.H:338
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition: AMReX_Geometry.H:331
Definition: AMReX_PhysBCFunct.H:64
GpuBndryFuncFab(F const &a_f)
Definition: AMReX_PhysBCFunct.H:67
void nddoit(Box const &bx, FArrayBox &dest, int dcomp, int numcomp, Geometry const &geom, Real time, const Vector< BCRec > &bcr, int bcomp, int orig_comp)
Definition: AMReX_PhysBCFunct.H:282
GpuBndryFuncFab(F &&a_f)
Definition: AMReX_PhysBCFunct.H:68
void ccfcdoit(Box const &bx, FArrayBox &dest, int dcomp, int numcomp, Geometry const &geom, Real time, const Vector< BCRec > &bcr, int bcomp, int orig_comp, FF const &fillfunc)
Definition: AMReX_PhysBCFunct.H:406
void operator()(Box const &bx, FArrayBox &dest, int dcomp, int numcomp, Geometry const &geom, Real time, const Vector< BCRec > &bcr, int bcomp, int orig_comp)
Definition: AMReX_PhysBCFunct.H:257
F m_user_f
Definition: AMReX_PhysBCFunct.H:90
Definition: AMReX_GpuAsyncArray.H:29
T const * data() const noexcept
Definition: AMReX_GpuAsyncArray.H:69
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE CellIndex ixType(int dir) const noexcept
Returns the CellIndex in direction dir.
Definition: AMReX_IndexType.H:116
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool allGE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is greater than or equal to argument for all components. NOTE: This is NOT a str...
Definition: AMReX_IntVect.H:443
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE bool allLE(const IntVectND< dim > &rhs) const noexcept
Returns true if this is less than or equal to argument for all components. NOTE: This is NOT a strict...
Definition: AMReX_IntVect.H:393
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVectND< dim > TheNodeVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition: AMReX_IntVect.H:702
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int min() const noexcept
minimum (no absolute values) value
Definition: AMReX_IntVect.H:225
Definition: AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition: AMReX_MFIter.H:141
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Definition: AMReX_PhysBCFunct.H:120
void operator()(MultiFab &, int, int, IntVect const &, Real, int)
Definition: AMReX_PhysBCFunct.H:122
Definition: AMReX_PhysBCFunct.H:127
IntVect src_ghost_outside_domain
Definition: AMReX_PhysBCFunct.H:169
IntVect fp1_src_ghost
Definition: AMReX_PhysBCFunct.H:174
PhysBCFunctUseCoarseGhost(IntVect const &a_fp1_src_ghost)
Definition: AMReX_PhysBCFunct.H:130
IntVect cghost
Definition: AMReX_PhysBCFunct.H:172
void operator()(MultiFab &, int, int, IntVect const &, Real, int)
Definition: AMReX_PhysBCFunct.H:156
IntVect src_ghost
Definition: AMReX_PhysBCFunct.H:165
IntVect nghost
Definition: AMReX_PhysBCFunct.H:159
PhysBCFunctUseCoarseGhost(MF const &cmf, IntVect const &a_nghost, IntVect const &a_nghost_outside_domain, IntVect const &ratio, Interp *mapper)
Definition: AMReX_PhysBCFunct.H:134
IntVect nghost_outside_domain
Definition: AMReX_PhysBCFunct.H:161
Definition: AMReX_PhysBCFunct.H:179
void define(const Geometry &geom, const Vector< BCRec > &bcr, F const &f)
Definition: AMReX_PhysBCFunct.H:191
void define(const Geometry &geom, const Vector< BCRec > &bcr, F &&f)
Definition: AMReX_PhysBCFunct.H:195
PhysBCFunct(const Geometry &geom, const Vector< BCRec > &bcr, F const &f)
Definition: AMReX_PhysBCFunct.H:183
Vector< BCRec > m_bcr
Definition: AMReX_PhysBCFunct.H:251
void operator()(MultiFab &mf, int icomp, int ncomp, IntVect const &nghost, Real time, int bccomp)
Definition: AMReX_PhysBCFunct.H:199
PhysBCFunct(const Geometry &geom, const Vector< BCRec > &bcr, F &&f)
Definition: AMReX_PhysBCFunct.H:187
PhysBCFunct()=default
Geometry m_geom
Definition: AMReX_PhysBCFunct.H:250
F m_f
Definition: AMReX_PhysBCFunct.H:252
void FillBoundary(MultiFab &mf, int icomp, int ncomp, IntVect const &nghost, Real time, int bccomp)
Definition: AMReX_PhysBCFunct.H:243
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition: AMReX_Vector.H:27
Long size() const noexcept
Definition: AMReX_Vector.H:50
static int f(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition: AMReX_SundialsIntegrator.H:44
Definition: AMReX_Amr.cpp:49
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition: AMReX_CTOParallelForImpl.H:200
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Returns the cell centered BoxND of length len adjacent to b on the low end along the coordinate direc...
Definition: AMReX_Box.H:1591
BoxND< AMREX_SPACEDIM > Box
Definition: AMReX_BaseFwd.H:27
IntVect nGrowVect(FabArrayBase const &fa)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > adjCellHi(const BoxND< dim > &b, int dir, int len=1) noexcept
Similar to adjCellLo but builds an adjacent BoxND on the high end.
Definition: AMReX_Box.H:1612
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition: AMReX_Box.H:1211
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Returns a BoxND with different type.
Definition: AMReX_Box.H:1435
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE IntVectND< dim > getCell(BoxND< dim > const *boxes, int nboxes, Long icell) noexcept
Definition: AMReX_Box.H:1978
void(*)(Real *data, const int *lo, const int *hi, const int *dom_lo, const int *dom_hi, const Real *dx, const Real *grd_lo, const Real *time, const int *bc) BndryFunc3DDefault
Definition: AMReX_PhysBCFunct.H:25
IntVectND< AMREX_SPACEDIM > IntVect
Definition: AMReX_BaseFwd.H:30
void(*)(Box const &bx, Array4< Real > const &dest, int dcomp, int numcomp, GeometryData const &geom, Real time, const BCRec *bcr, int bcomp, int orig_comp) UserFillBox
Definition: AMReX_PhysBCFunct.H:32
void(*)(Real *data, AMREX_ARLIM_P(lo), AMREX_ARLIM_P(hi), const int *dom_lo, const int *dom_hi, const Real *dx, const Real *grd_lo, const Real *time, const int *bc) BndryFuncDefault
Definition: AMReX_PhysBCFunct.H:21
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition: AMReX.H:111
void setBC(const Box &bx, const Box &domain, int src_comp, int dest_comp, int ncomp, const Vector< BCRec > &bc_dom, Vector< BCRec > &bcr) noexcept
Function for setting array of BCs.
Definition: AMReX_BCRec.cpp:8
AMREX_ATTRIBUTE_FLATTEN_FOR void For(T n, L const &f) noexcept
Definition: AMReX_GpuLaunchFunctsC.H:134
std::array< T, N > Array
Definition: AMReX_Array.H:24
Definition: AMReX_Array4.H:61
Definition: AMReX_PhysBCFunct.H:94
AMREX_GPU_DEVICE void operator()(const amrex::IntVect &, amrex::Array4< amrex::Real > const &, int, int, amrex::GeometryData const &, amrex::Real, const amrex::BCRec *, int, int) const
Definition: AMReX_PhysBCFunct.H:96
Definition: AMReX_FilCC_1D_C.H:12
Definition: AMReX_FilFC_1D_C.H:12
Definition: AMReX_Geometry.H:25