Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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
14namespace amrex {
15
16extern "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
28using 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{
37public:
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
51protected:
54 bool m_run_on_gpu = false;
55};
56
62template <class F>
64{
65public:
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
89protected:
91};
92
94{
97 int, int, amrex::GeometryData const&, amrex::Real,
98 const amrex::BCRec*, int, int) const
99 {}
100};
101
104{
105public:
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
115protected:
117};
118
120{
121public:
122 void operator() (MultiFab& /*mf*/, int /*dcomp*/, int /*ncomp*/, IntVect const& /*nghost*/,
123 Real /*time*/, int /*bccomp*/) {}
124};
125
127{
128public:
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
177template <class F>
179{
180public:
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 = amrex::convert(m_geom.Domain(), mf.boxArray().ixType());
210 Box gdomain = domain;
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
249private:
250 Geometry m_geom;
251 Vector<BCRec> m_bcr;
252 F m_f;
253};
254
255template <class F>
256void
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 (bx.ixType().cellCentered()) {
264 ccfcdoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp, FilccCell{});
265 } else if (bx.ixType().nodeCentered()) {
266 nddoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp);
267 } else {
268 ccfcdoit(bx,dest,dcomp,numcomp,geom,time,bcr,bcomp,orig_comp, FilfcFace{});
269 }
270}
271
272template <class F>
273void
275 int dcomp, int numcomp,
276 Geometry const& geom, Real time,
277 const Vector<BCRec>& bcr, int bcomp,
278 int orig_comp)
279{
280 const IntVect& len = bx.length();
281
282 Box const& domain = amrex::convert(geom.Domain(),IntVect::TheNodeVector());
283 Box gdomain = domain;
284 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
285 if (geom.isPeriodic(idim)) {
286 gdomain.grow(idim,len[idim]);
287 }
288 }
289
290 if (gdomain.contains(bx)) { return; }
291
292 Array4<Real> const& fab = dest.array();
293 const auto geomdata = geom.data();
294
295 AsyncArray<BCRec> bcr_aa(bcr.data()+bcomp, numcomp);
296 BCRec* bcr_p = bcr_aa.data();
297
298 const auto f_user = m_user_f;
299
300 // xlo
301 if (!geom.isPeriodic(0) && bx.smallEnd(0) < domain.smallEnd(0)) {
302 Box bndry = bx;
303 int dxlo = domain.smallEnd(0);
304 bndry.setBig(0,dxlo-1);
305 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
306 {
307 for (int n = dcomp; n < dcomp+numcomp; ++n) {
308 fab(i,j,k,n) = fab(dxlo,j,k,n);
309 }
310 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
311 bcr_p, 0, orig_comp);
312 });
313 }
314
315 // xhi
316 if (!geom.isPeriodic(0) && bx.bigEnd(0) > domain.bigEnd(0)) {
317 Box bndry = bx;
318 int dxhi = domain.bigEnd(0);
319 bndry.setSmall(0,dxhi+1);
320 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
321 {
322 for (int n = dcomp; n < dcomp+numcomp; ++n) {
323 fab(i,j,k,n) = fab(dxhi,j,k,n);
324 }
325 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
326 bcr_p, 0, orig_comp);
327 });
328 }
329
330#if (AMREX_SPACEDIM >= 2)
331 // ylo
332 if (!geom.isPeriodic(1) && bx.smallEnd(1) < domain.smallEnd(1)) {
333 Box bndry = bx;
334 int dylo = domain.smallEnd(1);
335 bndry.setBig(1,dylo-1);
336 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
337 {
338 for (int n = dcomp; n < dcomp+numcomp; ++n) {
339 fab(i,j,k,n) = fab(i,dylo,k,n);
340 }
341 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
342 bcr_p, 0, orig_comp);
343 });
344 }
345
346 // yhi
347 if (!geom.isPeriodic(1) && bx.bigEnd(1) > domain.bigEnd(1)) {
348 Box bndry = bx;
349 int dyhi = domain.bigEnd(1);
350 bndry.setSmall(1,dyhi+1);
351 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
352 {
353 for (int n = dcomp; n < dcomp+numcomp; ++n) {
354 fab(i,j,k,n) = fab(i,dyhi,k,n);
355 }
356 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
357 bcr_p, 0, orig_comp);
358 });
359 }
360#endif
361
362#if (AMREX_SPACEDIM == 3)
363 // zlo
364 if (!geom.isPeriodic(2) && bx.smallEnd(2) < domain.smallEnd(2)) {
365 Box bndry = bx;
366 int dzlo = domain.smallEnd(2);
367 bndry.setBig(2,dzlo-1);
368 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
369 {
370 for (int n = dcomp; n < dcomp+numcomp; ++n) {
371 fab(i,j,k,n) = fab(i,j,dzlo,n);
372 }
373 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
374 bcr_p, 0, orig_comp);
375 });
376 }
377
378 // zhi
379 if (!geom.isPeriodic(2) && bx.bigEnd(2) > domain.bigEnd(2)) {
380 Box bndry = bx;
381 int dzhi = domain.bigEnd(2);
382 bndry.setSmall(2,dzhi+1);
383 amrex::ParallelFor(bndry, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
384 {
385 for (int n = dcomp; n < dcomp+numcomp; ++n) {
386 fab(i,j,k,n) = fab(i,j,dzhi,n);
387 }
388 f_user(IntVect(AMREX_D_DECL(i,j,k)), fab, dcomp, numcomp, geomdata, time,
389 bcr_p, 0, orig_comp);
390 });
391 }
392#endif
393}
394
395template <class F>
396template <class FF>
397void
399 int dcomp, int numcomp,
400 Geometry const& geom, Real time,
401 const Vector<BCRec>& bcr, int bcomp,
402 int orig_comp, FF const& fillfunc)
403{
404 const IntVect& len = bx.length();
405
406 IndexType idxType = bx.ixType();
407 Box const& domain = amrex::convert(geom.Domain(),idxType);
408 Box gdomain = domain;
409 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
410 if (geom.isPeriodic(idim)) {
411 gdomain.grow(idim,len[idim]);
412 }
413 }
414
415 if (gdomain.contains(bx)) { return; }
416
417 Array4<Real> const& fab = dest.array();
418 const auto geomdata = geom.data();
419
420#ifdef AMREX_USE_GPU
421 AsyncArray<BCRec> bcr_aa(bcr.data()+bcomp, numcomp);
422 BCRec* bcr_p = bcr_aa.data();
423
424 const auto f_user = m_user_f;
425
426 // fill on the faces first
427 {
428 Array<Box,2*AMREX_SPACEDIM> dom_face_boxes
429 = { AMREX_D_DECL(amrex::convert(amrex::adjCellLo(gdomain, 0, len[0]),idxType),
430 amrex::convert(amrex::adjCellLo(gdomain, 1, len[1]),idxType),
431 amrex::convert(amrex::adjCellLo(gdomain, 2, len[2]),idxType)),
432 AMREX_D_DECL(amrex::convert(amrex::adjCellHi(gdomain, 0, len[0]),idxType),
433 amrex::convert(amrex::adjCellHi(gdomain, 1, len[1]),idxType),
434 amrex::convert(amrex::adjCellHi(gdomain, 2, len[2]),idxType)) };
435
436 Vector<Box> face_boxes;
437 for (const Box& b : dom_face_boxes) {
438 Box tmp = b & bx;
439 if (tmp.ok()) { face_boxes.push_back(tmp); }
440 }
441 const int n_face_boxes = face_boxes.size();
442 if (n_face_boxes == 1) {
443 amrex::ParallelFor(face_boxes[0],
444 [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
445 {
447 IntVect const idx(AMREX_D_DECL(i,j,k));
448 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
449 f_user(idx, fab, dcomp, numcomp, geomdata, time,
450 bcr_p, 0, orig_comp);
451 });
452 } else if (n_face_boxes > 1) {
453 AsyncArray<Box> face_boxes_aa(face_boxes.data(), n_face_boxes);
454 Box* boxes_p = face_boxes_aa.data();
455 Long ncounts = 0;
456 for (const auto& b : face_boxes) {
457 ncounts += b.numPts();
458 }
459 amrex::ParallelFor(ncounts,
460 [=] AMREX_GPU_DEVICE (Long icount) noexcept
461 {
462 const auto& idx = getCell(boxes_p, n_face_boxes, icount);
463 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
464 f_user(idx, fab, dcomp, numcomp, geomdata, time,
465 bcr_p, 0, orig_comp);
466 });
467 }
468 }
469
470#if (AMREX_SPACEDIM >= 2)
471 // fill on the box edges
472 {
473#if (AMREX_SPACEDIM == 2)
474 Array<Box,4> dom_edge_boxes
475 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
476 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
477 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
478 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType) }};
479#else
480 Array<Box,12> dom_edge_boxes
481 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
482 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
483 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
484 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
485 //
486 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
487 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
488 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
489 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
490 //
491 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
492 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType),
493 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
494 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType) }};
495#endif
496 Vector<Box> edge_boxes;
497 for (const Box& b : dom_edge_boxes) {
498 Box tmp = b & bx;
499 if (tmp.ok()) { edge_boxes.push_back(tmp); }
500 }
501 const int n_edge_boxes = edge_boxes.size();
502 if (n_edge_boxes == 1) {
503 amrex::ParallelFor(edge_boxes[0],
504 [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
505 {
507 IntVect const idx(AMREX_D_DECL(i,j,k));
508 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
509 f_user(idx, fab, dcomp, numcomp, geomdata, time,
510 bcr_p, 0, orig_comp);
511 });
512 } else if (n_edge_boxes > 1) {
513 AsyncArray<Box> edge_boxes_aa(edge_boxes.data(), n_edge_boxes);
514 Box* boxes_p = edge_boxes_aa.data();
515 Long ncounts = 0;
516 for (const auto& b : edge_boxes) {
517 ncounts += b.numPts();
518 }
519 amrex::ParallelFor(ncounts,
520 [=] AMREX_GPU_DEVICE (Long icount) noexcept
521 {
522 const auto& idx = getCell(boxes_p, n_edge_boxes, icount);
523 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
524 f_user(idx, fab, dcomp, numcomp, geomdata, time,
525 bcr_p, 0, orig_comp);
526 });
527 }
528 }
529#endif
530
531#if (AMREX_SPACEDIM == 3)
532 // fill on corners
533 {
534 Array<Box,8> dom_corner_boxes
535 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
536 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
537 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
538 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
539 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
540 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
541 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
542 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType) }};
543
544 Vector<Box> corner_boxes;
545 for (const Box& b : dom_corner_boxes) {
546 Box tmp = b & bx;
547 if (tmp.ok()) { corner_boxes.push_back(tmp); }
548 }
549 const int n_corner_boxes = corner_boxes.size();
550 if (n_corner_boxes == 1) {
551 amrex::ParallelFor(corner_boxes[0],
552 [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
553 {
554 IntVect const idx(AMREX_D_DECL(i,j,k));
555 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
556 f_user(idx, fab, dcomp, numcomp, geomdata, time,
557 bcr_p, 0, orig_comp);
558 });
559 } else if (n_corner_boxes > 1) {
560 AsyncArray<Box> corner_boxes_aa(corner_boxes.data(), n_corner_boxes);
561 Box* boxes_p = corner_boxes_aa.data();
562 Long ncounts = 0;
563 for (const auto& b : corner_boxes) {
564 ncounts += b.numPts();
565 }
566 amrex::ParallelFor(ncounts,
567 [=] AMREX_GPU_DEVICE (Long icount) noexcept
568 {
569 const auto& idx = getCell(boxes_p, n_corner_boxes, icount);
570 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
571 f_user(idx, fab, dcomp, numcomp, geomdata, time,
572 bcr_p, 0, orig_comp);
573 });
574 }
575 }
576#endif
577
578#else
579 BCRec const* bcr_p = bcr.data()+bcomp;
580
581 const auto& f_user = m_user_f;
582
583 // fill on the box faces first
584 {
585 Array<Box,2*AMREX_SPACEDIM> dom_face_boxes
586 = {{ AMREX_D_DECL(amrex::convert(amrex::adjCellLo(gdomain, 0, len[0]),idxType),
587 amrex::convert(amrex::adjCellLo(gdomain, 1, len[1]),idxType),
588 amrex::convert(amrex::adjCellLo(gdomain, 2, len[2]),idxType)),
589 AMREX_D_DECL(amrex::convert(amrex::adjCellHi(gdomain, 0, len[0]),idxType),
590 amrex::convert(amrex::adjCellHi(gdomain, 1, len[1]),idxType),
591 amrex::convert(amrex::adjCellHi(gdomain, 2, len[2]),idxType)) }};
592
593 for (const Box& b : dom_face_boxes) {
594 Box tmp = b & bx;
595 amrex::For(tmp, [=] (int i, int j, int k) noexcept
596 {
598 IntVect const idx(AMREX_D_DECL(i,j,k));
599 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
600 f_user(idx, fab, dcomp, numcomp, geomdata, time,
601 bcr_p, 0, orig_comp);
602 });
603 }
604 }
605
606#if (AMREX_SPACEDIM >= 2)
607 // fill on the box edges
608 {
609#if (AMREX_SPACEDIM == 2)
610 Array<Box,4> dom_edge_boxes
611 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
612 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
613 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
614 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType) }};
615#else
616 Array<Box,12> dom_edge_boxes
617 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
618 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
619 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
620 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
621 //
622 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
623 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
624 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
625 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
626 //
627 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
628 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType),
629 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
630 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType) }};
631#endif
632
633 for (const Box& b : dom_edge_boxes) {
634 Box tmp = b & bx;
635 amrex::For(tmp, [=] (int i, int j, int k) noexcept
636 {
638 IntVect const idx(AMREX_D_DECL(i,j,k));
639 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
640 f_user(idx, fab, dcomp, numcomp, geomdata, time,
641 bcr_p, 0, orig_comp);
642 });
643 }
644 }
645#endif
646
647#if (AMREX_SPACEDIM == 3)
648 // fill on box corners
649 {
650 Array<Box,8> dom_corner_boxes
651 = {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
652 amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
653 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
654 amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
655 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
656 amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
657 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
658 amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType) }};
659
660 for (const Box& b : dom_corner_boxes) {
661 Box tmp = b & bx;
662 amrex::For(tmp, [=] (int i, int j, int k) noexcept
663 {
664 IntVect const idx(AMREX_D_DECL(i,j,k));
665 fillfunc(idx, fab, dcomp, numcomp, domain, bcr_p, 0);
666 f_user(idx, fab, dcomp, numcomp, geomdata, time,
667 bcr_p, 0, orig_comp);
668 });
669 }
670 }
671#endif
672
673#endif
674}
675
676}
677#endif
#define AMREX_ARLIM_P(x)
Definition AMReX_ArrayLim.H:30
#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:173
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Boundary Condition Records. Necessary information and functions for computing boundary conditions.
Definition AMReX_BCRec.H:17
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:382
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:854
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:123
__host__ __device__ BoxND & setBig(const IntVectND< dim > &bg) noexcept
Redefine the big end of the BoxND.
Definition AMReX_Box.H:487
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:212
__host__ __device__ IndexTypeND< dim > ixType() const noexcept
Return the indexing type.
Definition AMReX_Box.H:135
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:208
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
__host__ __device__ BoxND & setSmall(const IntVectND< dim > &sm) noexcept
Redefine the small end of the BoxND.
Definition AMReX_Box.H:479
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:231
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:80
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:95
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
GeometryData data() const noexcept
Returns non-static copy of geometry's stored data.
Definition AMReX_Geometry.H:125
bool isAllPeriodic() const noexcept
Is domain periodic in all directions?
Definition AMReX_Geometry.H:344
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:337
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:274
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:398
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:34
T const * data() const noexcept
Definition AMReX_GpuAsyncArray.H:74
__host__ __device__ constexpr CellIndex ixType(int dir) const noexcept
Returns the CellIndex in direction dir.
Definition AMReX_IndexType.H:119
__host__ __device__ constexpr int min() const noexcept
minimum (no absolute values) value
Definition AMReX_IntVect.H:324
__host__ __device__ constexpr 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:542
__host__ static __device__ constexpr 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:801
__host__ __device__ constexpr 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:492
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
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
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
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:28
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:50
__host__ __device__ 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:1735
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
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:193
void(*)(Real *data, const int &, const int &, const int &, const int &, const int &, const int &, 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
__host__ __device__ BoxND< dim > adjCellLo(const BoxND< dim > &b, int dir, int len=1) noexcept
Return the cell centered BoxND of length len adjacent to b on the low end along the coordinate direct...
Definition AMReX_Box.H:1714
IntVect nGrowVect(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2857
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
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
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
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
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
__host__ __device__ IntVectND< dim > getCell(BoxND< dim > const *boxes, int nboxes, Long icell) noexcept
Definition AMReX_Box.H:2103
AMREX_ATTRIBUTE_FLATTEN_FOR void For(T n, L const &f) noexcept
Definition AMReX_GpuLaunchFunctsC.H:136
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_PhysBCFunct.H:94
__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_Geometry.H:26