Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_FillPatchUtil_I.H
Go to the documentation of this file.
1#ifndef AMREX_FillPatchUtil_I_H_
2#define AMREX_FillPatchUtil_I_H_
3#include <AMReX_Config.H>
4
5namespace amrex {
6
8namespace detail {
9
10template <typename F, typename MF>
11auto call_interp_hook (F const& f, MF& mf, int icomp, int ncomp)
12 -> decltype(f(mf[0],Box(),icomp,ncomp))
13{
14#ifdef AMREX_USE_OMP
15#pragma omp parallel if (Gpu::notInLaunchRegion())
16#endif
17 for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
18 auto& dfab = mf[mfi];
19 const Box& dbx = dfab.box();
20 f(dfab, dbx, icomp, ncomp);
21 }
22}
23
24template <typename F, typename MF>
25auto call_interp_hook (F const& f, MF& mf, int icomp, int ncomp)
26 -> decltype(f(mf,icomp,ncomp))
27{
28 f(mf, icomp, ncomp);
29}
30
31}
33
34template <typename Interp>
35bool ProperlyNested (const IntVect& ratio, const IntVect& blocking_factor, int ngrow,
36 const IndexType& boxType, Interp* mapper)
37{
38 int ratio_max = ratio[0];
39#if (AMREX_SPACEDIM > 1)
40 ratio_max = std::max(ratio_max, ratio[1]);
41#endif
42#if (AMREX_SPACEDIM == 3)
43 ratio_max = std::max(ratio_max, ratio[2]);
44#endif
45 // There are at least this many coarse cells outside fine grids
46 // (except at physical boundaries).
47 const IntVect& nbuf = blocking_factor / ratio_max;
48
49 Box crse_box(IntVect(AMREX_D_DECL(0 ,0 ,0 )), IntVect(AMREX_D_DECL(4*nbuf[0]-1,
50 4*nbuf[1]-1,
51 4*nbuf[2]-1)));
52 crse_box.convert(boxType);
53 Box fine_box(nbuf, IntVect(AMREX_D_DECL(3*nbuf[0]-1,3*nbuf[1]-1,3*nbuf[2]-1)));
54 fine_box.convert(boxType);
55 fine_box.refine(ratio_max);
56 fine_box.grow(ngrow);
57 const Box& fine_box_coarsened = mapper->CoarseBox(fine_box, ratio_max);
58 return crse_box.contains(fine_box_coarsened);
59}
60
61template <typename MF, typename BC>
62std::enable_if_t<IsFabArray<MF>::value>
64 const Vector<MF*>& smf, const Vector<Real>& stime,
65 int scomp, int dcomp, int ncomp,
66 const Geometry& geom,
67 BC& physbcf, int bcfcomp)
68{
69 FillPatchSingleLevel(mf, mf.nGrowVect(), time, smf, stime, scomp, dcomp, ncomp,
70 geom, physbcf, bcfcomp);
71}
72
73template <typename MF, typename BC>
74std::enable_if_t<IsFabArray<MF>::value>
75FillPatchSingleLevel (MF& mf, IntVect const& nghost, Real time,
76 const Vector<MF*>& smf, const Vector<Real>& stime,
77 int scomp, int dcomp, int ncomp,
78 const Geometry& geom,
79 BC& physbcf, int bcfcomp)
80{
81 BL_PROFILE("FillPatchSingleLevel");
82
83 AMREX_ASSERT((!smf.empty()) &&
84 (scomp+ncomp <= smf[0]->nComp()) &&
85 (dcomp+ncomp <= mf.nComp()) &&
86 (smf.size() == stime.size()) &&
87 nghost.allLE(mf.nGrowVect()));
88
89 IntVect src_ghost(0);
90 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
91 src_ghost = physbcf.fp1_src_ghost;
92 }
93
94 if (smf.size() == 1)
95 {
96 if (&mf == smf[0] && scomp == dcomp) {
97 mf.FillBoundary(dcomp, ncomp, nghost, geom.periodicity());
98 } else {
99 mf.ParallelCopy(*smf[0], scomp, dcomp, ncomp, src_ghost, nghost, geom.periodicity());
100 }
101 }
102 else if (smf.size() == 2)
103 {
104 BL_ASSERT(smf[0]->boxArray() == smf[1]->boxArray());
105 MF raii;
106 MF * dmf;
107 int destcomp;
108 bool sameba;
109 if (mf.boxArray() == smf[0]->boxArray() &&
110 mf.DistributionMap() == smf[0]->DistributionMap())
111 {
112 dmf = &mf;
113 destcomp = dcomp;
114 sameba = true;
115 } else {
116 raii.define(smf[0]->boxArray(), smf[0]->DistributionMap(), ncomp, src_ghost,
117 MFInfo(), smf[0]->Factory());
118
119 dmf = &raii;
120 destcomp = 0;
121 sameba = false;
122 }
123
124 if ((dmf != smf[0] && dmf != smf[1]) || scomp != dcomp)
125 {
126 IntVect interp_ghost(0);
127 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
128 interp_ghost = physbcf.fp1_src_ghost;
129 if (sameba) {
130 interp_ghost.min(nghost);
131 }
132 }
133#ifdef AMREX_USE_OMP
134#pragma omp parallel if (Gpu::notInLaunchRegion())
135#endif
136 for (MFIter mfi(*dmf,TilingIfNotGPU()); mfi.isValid(); ++mfi)
137 {
138 const Box& bx = mfi.growntilebox(interp_ghost);
139 const Real t0 = stime[0];
140 const Real t1 = stime[1];
141 auto const sfab0 = smf[0]->array(mfi);
142 auto const sfab1 = smf[1]->array(mfi);
143 auto dfab = dmf->array(mfi);
144
145 if (time == t0)
146 {
147 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
148 {
149 dfab(i,j,k,n+destcomp) = sfab0(i,j,k,n+scomp);
150 });
151 }
152 else if (time == t1)
153 {
154 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
155 {
156 dfab(i,j,k,n+destcomp) = sfab1(i,j,k,n+scomp);
157 });
158 }
159 else if (! amrex::almostEqual(t0,t1))
160 {
161 Real alpha = (t1-time)/(t1-t0);
162 Real beta = (time-t0)/(t1-t0);
163 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
164 {
165 dfab(i,j,k,n+destcomp) = alpha*sfab0(i,j,k,n+scomp)
166 + beta*sfab1(i,j,k,n+scomp);
167 });
168 }
169 else
170 {
171 AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( bx, ncomp, i, j, k, n,
172 {
173 dfab(i,j,k,n+destcomp) = sfab0(i,j,k,n+scomp);
174 });
175 }
176 }
177 }
178
179 if (sameba)
180 {
181 // Note that when sameba is true mf's BoxArray is nonoverlapping.
182 // So FillBoundary is safe.
183 mf.FillBoundary(dcomp, ncomp, nghost, geom.periodicity());
184 }
185 else
186 {
187 mf.ParallelCopy(*dmf, 0, dcomp, ncomp, src_ghost, nghost, geom.periodicity());
188 }
189 }
190 else {
191 amrex::Abort("FillPatchSingleLevel: high-order interpolation in time not implemented yet");
192 }
193
194 physbcf(mf, dcomp, ncomp, nghost, time, bcfcomp);
195}
196
214void FillPatchInterp (MultiFab& mf_fine_patch, int fcomp, MultiFab const& mf_crse_patch, int ccomp,
215 int ncomp, IntVect const& ng, const Geometry& cgeom, const Geometry& fgeom,
216 Box const& dest_domain, const IntVect& ratio,
217 MFInterpolater* mapper, const Vector<BCRec>& bcs, int bcscomp);
218
219template <typename MF, typename Interp>
220std::enable_if_t<IsFabArray<MF>::value && !std::is_same_v<Interp,MFInterpolater>>
221FillPatchInterp (MF& mf_fine_patch, int fcomp, MF const& mf_crse_patch, int ccomp,
222 int ncomp, IntVect const& ng, const Geometry& cgeom, const Geometry& fgeom,
223 Box const& dest_domain, const IntVect& ratio,
224 Interp* mapper, const Vector<BCRec>& bcs, int bcscomp)
225{
226 BL_PROFILE("FillPatchInterp(Fab)");
227
228 Box const& cdomain = amrex::convert(cgeom.Domain(), mf_fine_patch.ixType());
229 int idummy=0;
230#ifdef AMREX_USE_OMP
231#pragma omp parallel if (Gpu::notInLaunchRegion())
232#endif
233 {
234 Vector<BCRec> bcr(ncomp);
235 for (MFIter mfi(mf_fine_patch); mfi.isValid(); ++mfi)
236 {
237 auto& sfab = mf_crse_patch[mfi];
238 const Box& sbx = sfab.box();
239
240 auto& dfab = mf_fine_patch[mfi];
241 Box const& dbx = amrex::grow(mfi.validbox(),ng) & dest_domain;
242
243 amrex::setBC(sbx,cdomain,bcscomp,0,ncomp,bcs,bcr);
244 mapper->interp(sfab, ccomp, dfab, fcomp, ncomp, dbx, ratio,
245 cgeom, fgeom, bcr, idummy, idummy, RunOn::Gpu);
246 }
247 }
248}
249
250template <typename MF>
251std::enable_if_t<IsFabArray<MF>::value>
252FillPatchInterp (MF& mf_fine_patch, int fcomp, MF const& mf_crse_patch, int ccomp,
253 int ncomp, IntVect const& ng, const Geometry& cgeom, const Geometry& fgeom,
254 Box const& dest_domain, const IntVect& ratio,
255 InterpBase* mapper, const Vector<BCRec>& bcs, int bcscomp)
256{
257 if (dynamic_cast<MFInterpolater*>(mapper)) {
258 FillPatchInterp(mf_fine_patch, fcomp, mf_crse_patch, ccomp,
259 ncomp, ng, cgeom, fgeom, dest_domain, ratio,
260 static_cast<MFInterpolater*>(mapper), bcs, bcscomp);
261 } else if (dynamic_cast<Interpolater*>(mapper)) {
262 FillPatchInterp(mf_fine_patch, fcomp, mf_crse_patch, ccomp,
263 ncomp, ng, cgeom, fgeom, dest_domain, ratio,
264 static_cast<Interpolater*>(mapper), bcs, bcscomp);
265 } else {
266 amrex::Abort("FillPatchInterp: unknown InterpBase");
267 }
268}
269
270template <typename MF, typename iMF, typename Interp>
271std::enable_if_t<IsFabArray<MF>::value && !std::is_same_v<Interp,MFInterpolater>>
272InterpFace (Interp *interp,
273 MF const& mf_crse_patch, int crse_comp,
274 MF& mf_refined_patch, int fine_comp,
275 int ncomp, const IntVect& ratio,
276 const iMF& solve_mask, const Geometry& crse_geom, const Geometry& fine_geom,
277 int bcscomp, RunOn gpu_or_cpu,
278 const Vector<BCRec>& bcs)
279{
280 Vector<BCRec> bcr(ncomp);
281 Box const& cdomain = amrex::convert(crse_geom.Domain(), mf_crse_patch.ixType());
282 for (MFIter mfi(mf_refined_patch);mfi.isValid(); ++mfi)
283 {
284 auto& sfab = mf_crse_patch[mfi];
285 const Box& sbx = sfab.box();
286 auto& dfab = mf_refined_patch[mfi];
287 Box const& dbx = dfab.box();
288 auto& ifab = solve_mask[mfi];
289 amrex::setBC(sbx,cdomain,bcscomp,0,ncomp,bcs,bcr);
290 interp->interp_face(sfab,crse_comp,dfab,fine_comp,ncomp,
291 dbx, ratio, ifab, crse_geom, fine_geom,
292 bcr, bcscomp, gpu_or_cpu);
293 }
294}
295
296template <typename MF, typename iMF>
297std::enable_if_t<IsFabArray<MF>::value>
299 MF const& mf_crse_patch, int crse_comp,
300 MF& mf_refined_patch, int fine_comp,
301 int ncomp, const IntVect& ratio,
302 const iMF& solve_mask, const Geometry& crse_geom, const Geometry& fine_geom,
303 int bccomp, RunOn gpu_or_cpu,
304 const Vector<BCRec>& bcs)
305{
306 if (dynamic_cast<MFInterpolater*>(interp)){
307 InterpFace(static_cast<MFInterpolater*>(interp),
308 mf_crse_patch, crse_comp,mf_refined_patch, fine_comp,
309 ncomp, ratio, solve_mask, crse_geom, fine_geom, bccomp,
310 gpu_or_cpu, bcs);
311 }
312 else if (dynamic_cast<Interpolater*>(interp)){
313 InterpFace(static_cast<Interpolater*>(interp),
314 mf_crse_patch, crse_comp,mf_refined_patch, fine_comp,
315 ncomp, ratio, solve_mask, crse_geom, fine_geom, bccomp,
316 gpu_or_cpu, bcs);
317 }
318 else {
319 amrex::Abort("InterpFace: unknown InterpBase");
320 }
321}
322
323
325namespace detail {
326
327// ======== FArrayBox
328
329 template <typename MF,
330 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
331 FArrayBox>,
332 int> = 0>
333 MF make_mf_crse_patch (FabArrayBase::FPinfo const& fpc, int ncomp)
334 {
335 MF mf_crse_patch(fpc.ba_crse_patch, fpc.dm_patch, ncomp, 0, MFInfo(),
336 *fpc.fact_crse_patch);
337 return mf_crse_patch;
338 }
339
340 template <typename MF,
341 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
342 FArrayBox>,
343 int> = 0>
344 MF make_mf_crse_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type)
345 {
346 MF mf_crse_patch(amrex::convert(fpc.ba_crse_patch, idx_type), fpc.dm_patch,
347 ncomp, 0, MFInfo(), *fpc.fact_crse_patch);
348 return mf_crse_patch;
349 }
350
351 template <typename MF,
352 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
353 FArrayBox>,
354 int> = 0>
355 MF make_mf_fine_patch (FabArrayBase::FPinfo const& fpc, int ncomp)
356 {
357 MF mf_fine_patch(fpc.ba_fine_patch, fpc.dm_patch, ncomp, 0, MFInfo(),
358 *fpc.fact_fine_patch);
359 return mf_fine_patch;
360 }
361
362 template <typename MF,
363 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
364 FArrayBox>,
365 int> = 0>
366 MF make_mf_fine_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type)
367 {
368 MF mf_fine_patch(amrex::convert(fpc.ba_fine_patch, idx_type), fpc.dm_patch,
369 ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
370 return mf_fine_patch;
371 }
372
373 template <typename MF,
374 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
375 FArrayBox>,
376 int> = 0>
377 MF make_mf_refined_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type, IntVect ratio)
378 {
379 MF mf_refined_patch(amrex::convert( amrex::refine( amrex::coarsen(fpc.ba_fine_patch, ratio), ratio), idx_type),
380 fpc.dm_patch, ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
381 return mf_refined_patch;
382 }
383
384 template <typename MF,
385 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
386 FArrayBox>,
387 int> = 0>
388 MF make_mf_crse_mask (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type, IntVect ratio)
389 {
390 MF mf_crse_mask(amrex::convert(amrex::coarsen(fpc.ba_fine_patch, ratio), idx_type),
391 fpc.dm_patch, ncomp, 0, MFInfo(), *fpc.fact_fine_patch);
392 return mf_crse_mask;
393 }
394
395 template <typename MF,
396 std::enable_if_t<std::is_same_v<typename MF::FABType::value_type,
397 FArrayBox>,
398 int> = 0>
399 void mf_set_domain_bndry (MF &mf, Geometry const & geom)
400 {
401 mf.setDomainBndry(std::numeric_limits<Real>::quiet_NaN(), geom);
402 }
403
404
405// ======== Not FArrayBox
406
407 template <typename MF,
408 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
409 FArrayBox>,
410 int> = 0>
411 MF make_mf_crse_patch (FabArrayBase::FPinfo const& fpc, int ncomp)
412 {
413 return MF(fpc.ba_crse_patch, fpc.dm_patch, ncomp, 0);
414 }
415
416 template <typename MF,
417 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
418 FArrayBox>,
419 int> = 0>
420 MF make_mf_crse_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type)
421 {
422 return MF(amrex::convert(fpc.ba_crse_patch, idx_type), fpc.dm_patch, ncomp, 0);
423 }
424
425 template <typename MF,
426 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
427 FArrayBox>,
428 int> = 0>
429 MF make_mf_fine_patch (FabArrayBase::FPinfo const& fpc, int ncomp)
430 {
431 return MF(fpc.ba_fine_patch, fpc.dm_patch, ncomp, 0);
432 }
433
434 template <typename MF,
435 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
436 FArrayBox>,
437 int> = 0>
438 MF make_mf_fine_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type)
439 {
440 return MF(amrex::convert(fpc.ba_fine_patch, idx_type), fpc.dm_patch, ncomp, 0);
441 }
442
443 template <typename MF,
444 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
445 FArrayBox>,
446 int> = 0>
447 MF make_mf_refined_patch (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type, IntVect ratio)
448 {
449 return MF(amrex::convert( amrex::refine( amrex::coarsen(fpc.ba_fine_patch, ratio), ratio), idx_type), fpc.dm_patch, ncomp, 0);
450 }
451
452 template <typename MF,
453 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
454 FArrayBox>,
455 int> = 0>
456 MF make_mf_crse_mask (FabArrayBase::FPinfo const& fpc, int ncomp, IndexType idx_type, IntVect ratio)
457 {
458 return MF(amrex::convert(amrex::coarsen(fpc.ba_fine_patch, ratio), idx_type), fpc.dm_patch, ncomp, 0);
459 }
460
461 template <typename MF,
462 std::enable_if_t<!std::is_same_v<typename MF::FABType::value_type,
463 FArrayBox>,
464 int> = 0>
465 void mf_set_domain_bndry (MF &/*mf*/, Geometry const & /*geom*/)
466 {
467 // nothing
468 }
469
470 template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
471 std::enable_if_t<IsFabArray<MF>::value,int>
472 FillPatchTwoLevels_doit (MF& mf, IntVect const& nghost, Real time,
473 const Vector<MF*>& cmf, const Vector<Real>& ct,
474 const Vector<MF*>& fmf, const Vector<Real>& ft,
475 int scomp, int dcomp, int ncomp,
476 const Geometry& cgeom, const Geometry& fgeom,
477 BC& cbc, int cbccomp,
478 BC& fbc, int fbccomp,
479 const IntVect& ratio,
480 Interp* mapper,
481 const Vector<BCRec>& bcs, int bcscomp,
482 const PreInterpHook& pre_interp,
483 const PostInterpHook& post_interp,
484 EB2::IndexSpace const* index_space,
485 bool return_error_code = false)
486 {
487 BL_PROFILE("FillPatchTwoLevels");
488
489 int success_code = return_error_code ? 0 : -1;
490 int failure_code = 1;
491
492 if (nghost.max() > 0 || mf.getBDKey() != fmf[0]->getBDKey())
493 {
494 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
495
496 // Test for Face-centered data
497 if ( AMREX_D_TERM( mf.ixType().nodeCentered(0),
498 + mf.ixType().nodeCentered(1),
499 + mf.ixType().nodeCentered(2) ) == 1 )
500 {
501 if ( !dynamic_cast<Interpolater*>(mapper) ){
502 amrex::Abort("This interpolater has not yet implemented a version for face-based data");
503 }
504
505 // Convert to cell-centered MF meta-data for FPInfo.
506 MF mf_cc_dummy( amrex::convert(mf.boxArray(), IntVect::TheZeroVector()),
507 mf.DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(false) );
508 MF fmf_cc_dummy( amrex::convert(fmf[0]->boxArray(), IntVect::TheZeroVector()),
509 fmf[0]->DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(false) );
510
511 const FabArrayBase::FPinfo& fpc = FabArrayBase::TheFPinfo(fmf_cc_dummy, mf_cc_dummy,
512 nghost,
513 coarsener,
514 fgeom,
515 cgeom,
516 index_space);
517
518 if ( ! fpc.ba_crse_patch.empty())
519 {
520 if (return_error_code) {
521 BoxArray const& cba = amrex::convert(cmf[0]->boxArray(), IntVect(0));
522 if (!cba.contains(fpc.ba_crse_patch,cgeom.periodicity())) {
523 return failure_code;
524 }
525 }
526
527 MF mf_crse_patch = make_mf_crse_patch<MF> (fpc, ncomp, mf.boxArray().ixType());
528 // Must make sure fine exists under needed coarse faces.
529 // It stores values for the final (interior) interpolation,
530 // which is done from this fine MF that's been partially filled
531 // (with only faces overlying coarse having valid data).
532 MF mf_refined_patch = make_mf_refined_patch<MF> (fpc, ncomp, mf.boxArray().ixType(), ratio);
533 auto solve_mask = make_mf_crse_mask<iMultiFab>(fpc, ncomp, mf.boxArray().ixType(), ratio);
534
535 mf_set_domain_bndry(mf_crse_patch, cgeom);
536 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
537 cbc.fp1_src_ghost = cbc.cghost;
538 }
539 FillPatchSingleLevel(mf_crse_patch, time, cmf, ct, scomp, 0, ncomp,
540 cgeom, cbc, cbccomp);
541
542 mf_set_domain_bndry(mf_refined_patch, fgeom);
543 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
544 fbc.fp1_src_ghost = IntVect(0);
545 }
546 FillPatchSingleLevel(mf_refined_patch, time, fmf, ft, scomp, 0, ncomp,
547 fgeom, fbc, fbccomp);
548
549 // Aliased MFs, used to allow CPC caching.
550 MF mf_known( amrex::coarsen(fmf[0]->boxArray(), ratio), fmf[0]->DistributionMap(),
551 ncomp, nghost, MFInfo().SetAlloc(false) );
552 MF mf_solution( amrex::coarsen(mf_refined_patch.boxArray(), ratio), mf_refined_patch.DistributionMap(),
553 ncomp, 0, MFInfo().SetAlloc(false) );
554
555 const FabArrayBase::CPC mask_cpc( mf_solution, IntVect::TheZeroVector(),
556 mf_known, IntVect::TheZeroVector(),
557 cgeom.periodicity());
558
559 solve_mask.setVal(1); // Values to solve.
560 solve_mask.setVal(0, mask_cpc, 0, 1); // Known values.
561
562 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
563
564 InterpFace(mapper, mf_crse_patch, 0, mf_refined_patch, 0, ncomp,
565 ratio, solve_mask, cgeom, fgeom, bcscomp, RunOn::Gpu, bcs);
566
567 detail::call_interp_hook(post_interp, mf_refined_patch, 0, ncomp);
568
569 bool aliasing = false;
570 for (auto const& fmf_a : fmf) {
571 aliasing = aliasing || (&mf == fmf_a);
572 }
573 if (aliasing) {
574 mf.ParallelCopyToGhost(mf_refined_patch, 0, dcomp, ncomp,
575 IntVect{0}, nghost);
576 } else {
577 mf.ParallelCopy(mf_refined_patch, 0, dcomp, ncomp,
578 IntVect{0}, nghost);
579 }
580 }
581 }
582 else
583 {
584 const FabArrayBase::FPinfo& fpc = FabArrayBase::TheFPinfo(*fmf[0], mf,
585 nghost,
586 coarsener,
587 fgeom,
588 cgeom,
589 index_space);
590
591 if ( ! fpc.ba_crse_patch.empty())
592 {
593 if (return_error_code) {
594 BoxArray const& cba = cmf[0]->boxArray();
595 if (!cba.contains(fpc.ba_crse_patch,cgeom.periodicity())) {
596 return failure_code;
597 }
598 }
599
600 MF mf_crse_patch = make_mf_crse_patch<MF>(fpc, ncomp);
601 mf_set_domain_bndry (mf_crse_patch, cgeom);
602
603 if constexpr (std::is_same_v<BC,PhysBCFunctUseCoarseGhost>) {
604 cbc.fp1_src_ghost = cbc.cghost;
605 }
606 FillPatchSingleLevel(mf_crse_patch, time, cmf, ct, scomp, 0, ncomp, cgeom, cbc, cbccomp);
607
608 MF mf_fine_patch = make_mf_fine_patch<MF>(fpc, ncomp);
609
610 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
611
612 Box fdomain_g( amrex::convert(fgeom.Domain(),mf.ixType()) );
613 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
614 if (fgeom.isPeriodic(i)) {
615 fdomain_g.grow(i, nghost[i]);
616 } else {
617 if constexpr (std::is_same_v
618 <BC, PhysBCFunctUseCoarseGhost>) {
619 fdomain_g.grow(i, fbc.nghost_outside_domain[i]);
620 }
621 }
622 }
623 FillPatchInterp(mf_fine_patch, 0, mf_crse_patch, 0,
624 ncomp, IntVect(0), cgeom, fgeom,
625 fdomain_g, ratio, mapper, bcs, bcscomp);
626
627 detail::call_interp_hook(post_interp, mf_fine_patch, 0, ncomp);
628
629 mf.ParallelCopy(mf_fine_patch, 0, dcomp, ncomp, IntVect{0}, nghost);
630 }
631 }
632 }
633
634 if constexpr(std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
635 fbc.fp1_src_ghost = IntVect(0);
636 }
637 FillPatchSingleLevel(mf, nghost, time, fmf, ft, scomp, dcomp, ncomp,
638 fgeom, fbc, fbccomp);
639
640 return success_code;
641 }
642
643 template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
644 std::enable_if_t<IsFabArray<MF>::value>
645 FillPatchTwoLevels_doit (Array<MF*, AMREX_SPACEDIM> const& mf, IntVect const& nghost, Real time,
646 const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
647 const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
648 int scomp, int dcomp, int ncomp,
649 const Geometry& cgeom, const Geometry& fgeom,
650 Array<BC, AMREX_SPACEDIM>& cbc, const Array<int, AMREX_SPACEDIM>& cbccomp,
651 Array<BC, AMREX_SPACEDIM>& fbc, const Array<int, AMREX_SPACEDIM>& fbccomp,
652 const IntVect& ratio,
653 Interp* mapper,
654 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, const Array<int, AMREX_SPACEDIM>& bcscomp,
655 const PreInterpHook& pre_interp,
656 const PostInterpHook& post_interp,
657 EB2::IndexSpace const* index_space)
658 {
659 BL_PROFILE("FillPatchTwoLevels (Array<MF*>)");
660
661 using FAB = typename MF::FABType::value_type;
662 using iFAB = typename iMultiFab::FABType::value_type;
663
664 AMREX_ASSERT(AMREX_D_TERM(mf[0]->ixType().nodeCentered(0),
665 && mf[1]->ixType().nodeCentered(1),
666 && mf[2]->ixType().nodeCentered(2)));
667
668 // These need to be true: (ba[0] == ba[1] == ba[2]) & (dm[0] == dm[1] == dm[2]).
669 // Debatable whether these are required, or will be enforced elsewhere prior to this func.
671 && BoxArray::SameRefs(mf[0]->boxArray(), mf[1]->boxArray()),
672 && BoxArray::SameRefs(mf[0]->boxArray(), mf[2]->boxArray())));
673/*
674 AMREX_ASSERT(AMREX_D_TERM(true,
675 && DistributionMapping::SameRefs(mf[0]->DistributionMap(), mf[1]->DistributionMap()),
676 && DistributionMapping::SameRefs(mf[0]->DistributionMap(), mf[2]->DistributionMap())));
677*/
678
679
680 // Test all of them?
681 if (nghost.max() > 0 || mf[0]->getBDKey() != fmf[0][0]->getBDKey())
682 {
683 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
684
685 // Convert to cell-centered MF meta-data for FPInfo.
686 MF mf_cc_dummy( amrex::convert(mf[0]->boxArray(), IntVect::TheZeroVector()),
687 mf[0]->DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(false) );
688 MF fmf_cc_dummy( amrex::convert(fmf[0][0]->boxArray(), IntVect::TheZeroVector()),
689 fmf[0][0]->DistributionMap(), ncomp, nghost, MFInfo().SetAlloc(false) );
690
691 const FabArrayBase::FPinfo& fpc = FabArrayBase::TheFPinfo(fmf_cc_dummy, mf_cc_dummy,
692 nghost,
693 coarsener,
694 fgeom,
695 cgeom,
696 index_space);
697
698 if ( !fpc.ba_crse_patch.empty() )
699 {
700 Array<MF, AMREX_SPACEDIM> mf_crse_patch;
701 Array<MF, AMREX_SPACEDIM> mf_refined_patch;
702 Array<iMultiFab, AMREX_SPACEDIM> solve_mask;
703
704 for (int d=0; d<AMREX_SPACEDIM; ++d)
705 {
706 mf_crse_patch[d] = make_mf_crse_patch<MF> (fpc, ncomp, mf[d]->boxArray().ixType());
707 mf_refined_patch[d] = make_mf_refined_patch<MF> (fpc, ncomp, mf[d]->boxArray().ixType(), ratio);
708 solve_mask[d] = make_mf_crse_mask<iMultiFab>(fpc, ncomp, mf[d]->boxArray().ixType(), ratio);
709
710 mf_set_domain_bndry(mf_crse_patch[d], cgeom);
711 Vector<MF*> cmf_time;
712 for (const auto & mfab : cmf)
713 { cmf_time.push_back(mfab[d]); }
714
715 FillPatchSingleLevel(mf_crse_patch[d], time, cmf_time, ct, scomp, 0, ncomp,
716 cgeom, cbc[d], cbccomp[d]);
717
718 mf_set_domain_bndry(mf_refined_patch[d], fgeom);
719 Vector<MF*> fmf_time;
720 for (const auto & mfab : fmf)
721 { fmf_time.push_back(mfab[d]); }
722
723 FillPatchSingleLevel(mf_refined_patch[d], time, fmf_time, ft, scomp, 0, ncomp,
724 fgeom, fbc[d], fbccomp[d]);
725
726
727 // Aliased MFs, used to allow CPC caching.
728 MF mf_known( amrex::coarsen(fmf[0][d]->boxArray(), ratio), fmf[0][d]->DistributionMap(),
729 ncomp, nghost, MFInfo().SetAlloc(false) );
730 MF mf_solution( amrex::coarsen(mf_refined_patch[d].boxArray(), ratio), mf_refined_patch[d].DistributionMap(),
731 ncomp, 0, MFInfo().SetAlloc(false) );
732
733 const FabArrayBase::CPC mask_cpc( mf_solution, IntVect::TheZeroVector(),
734 mf_known, IntVect::TheZeroVector(),
735 cgeom.periodicity() );
736
737 solve_mask[d].setVal(1); // Values to solve.
738 solve_mask[d].setVal(0, mask_cpc, 0, 1); // Known values.
739 }
740
741 int idummy=0;
742#ifdef AMREX_USE_OMP
743// bool cc = fpc.ba_crse_patch.ixType().cellCentered();
744 bool cc = false; // can anything be done to allow threading, or can the OpenMP just be removed?
745#pragma omp parallel if (cc && Gpu::notInLaunchRegion() )
746#endif
747 {
748 Vector<Array<BCRec, AMREX_SPACEDIM> > bcr(ncomp);
749 for (MFIter mfi(mf_refined_patch[0]); mfi.isValid(); ++mfi)
750 {
751 Array<FAB*, AMREX_SPACEDIM> sfab{ AMREX_D_DECL( &(mf_crse_patch[0][mfi]),
752 &(mf_crse_patch[1][mfi]),
753 &(mf_crse_patch[2][mfi]) )};
754 Array<FAB*, AMREX_SPACEDIM> dfab{ AMREX_D_DECL( &(mf_refined_patch[0][mfi]),
755 &(mf_refined_patch[1][mfi]),
756 &(mf_refined_patch[2][mfi]) )};
757 Array<iFAB*, AMREX_SPACEDIM> mfab{ AMREX_D_DECL( &(solve_mask[0][mfi]),
758 &(solve_mask[1][mfi]),
759 &(solve_mask[2][mfi]) )};
760
761 const Box& sbx_cc = amrex::convert(sfab[0]->box(), IntVect::TheZeroVector());
762 const Box& dbx_cc = amrex::convert(dfab[0]->box(), IntVect::TheZeroVector());
763
764 for (int d=0; d<AMREX_SPACEDIM; ++d)
765 {
766 Vector<BCRec> bcr_d(ncomp);
767
768 amrex::setBC(sfab[d]->box(), amrex::convert(cgeom.Domain(), mf[d]->ixType()),
769 bcscomp[d],0,ncomp,bcs[d],bcr_d);
770
771 for (int n=0; n<ncomp; ++n)
772 { bcr[n][d] = bcr_d[n]; }
773 }
774
775 pre_interp(sfab, sbx_cc, 0, ncomp);
776
777 mapper->interp_arr(sfab, 0, dfab, 0, ncomp, dbx_cc, ratio, mfab,
778 cgeom, fgeom, bcr, idummy, idummy, RunOn::Gpu);
779
780 post_interp(dfab, dbx_cc, 0, ncomp);
781 }
782 }
783
784 for (int d=0; d<AMREX_SPACEDIM; ++d)
785 {
786 bool aliasing = false;
787 for (auto const& fmf_a : fmf) {
788 aliasing = aliasing || (mf[d] == fmf_a[d]);
789 }
790 if (aliasing) {
791 mf[d]->ParallelCopyToGhost(mf_refined_patch[d], 0, dcomp, ncomp,
792 IntVect{0}, nghost);
793 } else {
794 mf[d]->ParallelCopy(mf_refined_patch[d], 0, dcomp, ncomp,
795 IntVect{0}, nghost);
796 }
797 }
798 }
799 }
800
801 for (int d=0; d<AMREX_SPACEDIM; ++d)
802 {
803 Vector<MF*> fmf_time;
804 for (auto const& ffab : fmf)
805 { fmf_time.push_back(ffab[d]); }
806
807 FillPatchSingleLevel(*mf[d], nghost, time, fmf_time, ft, scomp, dcomp, ncomp,
808 fgeom, fbc[d], fbccomp[d]);
809 }
810 }
811
812} // namespace detail
814
815template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
816std::enable_if_t<IsFabArray<MF>::value>
817FillPatchTwoLevels (MF& mf, IntVect const& nghost, Real time,
818 const Vector<MF*>& cmf, const Vector<Real>& ct,
819 const Vector<MF*>& fmf, const Vector<Real>& ft,
820 int scomp, int dcomp, int ncomp,
821 const Geometry& cgeom, const Geometry& fgeom,
822 BC& cbc, int cbccomp,
823 BC& fbc, int fbccomp,
824 const IntVect& ratio,
825 Interp* mapper,
826 const Vector<BCRec>& bcs, int bcscomp,
827 const PreInterpHook& pre_interp,
828 const PostInterpHook& post_interp)
829{
830#ifdef AMREX_USE_EB
831 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
832#else
833 EB2::IndexSpace const* index_space = nullptr;
834#endif
835 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
836 scomp,dcomp,ncomp,cgeom,fgeom,
837 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
838 pre_interp,post_interp,index_space);
839}
840
841template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
842std::enable_if_t<IsFabArray<MF>::value>
844 const Vector<MF*>& cmf, const Vector<Real>& ct,
845 const Vector<MF*>& fmf, const Vector<Real>& ft,
846 int scomp, int dcomp, int ncomp,
847 const Geometry& cgeom, const Geometry& fgeom,
848 BC& cbc, int cbccomp,
849 BC& fbc, int fbccomp,
850 const IntVect& ratio,
851 Interp* mapper,
852 const Vector<BCRec>& bcs, int bcscomp,
853 const PreInterpHook& pre_interp,
854 const PostInterpHook& post_interp)
855{
856#ifdef AMREX_USE_EB
857 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
858#else
859 EB2::IndexSpace const* index_space = nullptr;
860#endif
861
862 detail::FillPatchTwoLevels_doit(mf,mf.nGrowVect(),time,cmf,ct,fmf,ft,
863 scomp,dcomp,ncomp,cgeom,fgeom,
864 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
865 pre_interp,post_interp,index_space);
866}
867
868template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
869std::enable_if_t<IsFabArray<MF>::value>
871 const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
872 const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
873 int scomp, int dcomp, int ncomp,
874 const Geometry& cgeom, const Geometry& fgeom,
877 const IntVect& ratio,
878 Interp* mapper,
879 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, const Array<int, AMREX_SPACEDIM>& bcscomp,
880 const PreInterpHook& pre_interp,
881 const PostInterpHook& post_interp)
882{
883#ifdef AMREX_USE_EB
884 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
885#else
886 EB2::IndexSpace const* index_space = nullptr;
887#endif
888
889 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
890 scomp,dcomp,ncomp,cgeom,fgeom,
891 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
892 pre_interp,post_interp,index_space);
893}
894
895template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
896std::enable_if_t<IsFabArray<MF>::value>
898 const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
899 const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
900 int scomp, int dcomp, int ncomp,
901 const Geometry& cgeom, const Geometry& fgeom,
902 Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
903 Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
904 const IntVect& ratio,
905 Interp* mapper,
906 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
907 const PreInterpHook& pre_interp,
908 const PostInterpHook& post_interp)
909{
910#ifdef AMREX_USE_EB
911 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
912#else
913 EB2::IndexSpace const* index_space = nullptr;
914#endif
915
916 Array<int, AMREX_SPACEDIM> cbccomp_arr = {AMREX_D_DECL(cbccomp,cbccomp,cbccomp)};
917 Array<int, AMREX_SPACEDIM> fbccomp_arr = {AMREX_D_DECL(fbccomp,fbccomp,fbccomp)};
918 Array<int, AMREX_SPACEDIM> bcscomp_arr = {AMREX_D_DECL(bcscomp,bcscomp,bcscomp)};
919
920 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
921 scomp,dcomp,ncomp,cgeom,fgeom,
922 cbc,cbccomp_arr,fbc,fbccomp_arr,ratio,mapper,bcs,bcscomp_arr,
923 pre_interp,post_interp,index_space);
924}
925
926template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
927std::enable_if_t<IsFabArray<MF>::value>
929 const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
930 const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
931 int scomp, int dcomp, int ncomp,
932 const Geometry& cgeom, const Geometry& fgeom,
933 Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
934 Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
935 const IntVect& ratio,
936 Interp* mapper,
937 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
938 const PreInterpHook& pre_interp,
939 const PostInterpHook& post_interp)
940{
941#ifdef AMREX_USE_EB
942 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
943#else
944 EB2::IndexSpace const* index_space = nullptr;
945#endif
946
947 Array<int, AMREX_SPACEDIM> cbccomp_arr = {AMREX_D_DECL(cbccomp,cbccomp,cbccomp)};
948 Array<int, AMREX_SPACEDIM> fbccomp_arr = {AMREX_D_DECL(fbccomp,fbccomp,fbccomp)};
949 Array<int, AMREX_SPACEDIM> bcscomp_arr = {AMREX_D_DECL(bcscomp,bcscomp,bcscomp)};
950
951 detail::FillPatchTwoLevels_doit(mf,mf[0]->nGrowVect(),time,cmf,ct,fmf,ft,
952 scomp,dcomp,ncomp,cgeom,fgeom,
953 cbc,cbccomp_arr,fbc,fbccomp_arr,ratio,mapper,bcs,bcscomp_arr,
954 pre_interp,post_interp,index_space);
955}
956
957#ifdef AMREX_USE_EB
958template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
959std::enable_if_t<IsFabArray<MF>::value>
960FillPatchTwoLevels (MF& mf, IntVect const& nghost, Real time,
961 const EB2::IndexSpace& index_space,
962 const Vector<MF*>& cmf, const Vector<Real>& ct,
963 const Vector<MF*>& fmf, const Vector<Real>& ft,
964 int scomp, int dcomp, int ncomp,
965 const Geometry& cgeom, const Geometry& fgeom,
966 BC& cbc, int cbccomp,
967 BC& fbc, int fbccomp,
968 const IntVect& ratio,
969 Interp* mapper,
970 const Vector<BCRec>& bcs, int bcscomp,
971 const PreInterpHook& pre_interp,
972 const PostInterpHook& post_interp)
973{
974 detail::FillPatchTwoLevels_doit(mf,nghost,time,cmf,ct,fmf,ft,
975 scomp,dcomp,ncomp,cgeom,fgeom,
976 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
977 pre_interp,post_interp,&index_space);
978}
979
980template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
981std::enable_if_t<IsFabArray<MF>::value>
983 const EB2::IndexSpace& index_space,
984 const Vector<MF*>& cmf, const Vector<Real>& ct,
985 const Vector<MF*>& fmf, const Vector<Real>& ft,
986 int scomp, int dcomp, int ncomp,
987 const Geometry& cgeom, const Geometry& fgeom,
988 BC& cbc, int cbccomp,
989 BC& fbc, int fbccomp,
990 const IntVect& ratio,
991 Interp* mapper,
992 const Vector<BCRec>& bcs, int bcscomp,
993 const PreInterpHook& pre_interp,
994 const PostInterpHook& post_interp)
995{
996 detail::FillPatchTwoLevels_doit(mf,mf.nGrowVect(),time,cmf,ct,fmf,ft,
997 scomp,dcomp,ncomp,cgeom,fgeom,
998 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
999 pre_interp,post_interp,&index_space);
1000}
1001#endif
1002
1003template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
1004std::enable_if_t<IsFabArray<MF>::value>
1006 const MF& cmf, int scomp, int dcomp, int ncomp,
1007 const Geometry& cgeom, const Geometry& fgeom,
1008 BC& cbc, int cbccomp,
1009 BC& fbc, int fbccomp,
1010 const IntVect& ratio,
1011 Interp* mapper,
1012 const Vector<BCRec>& bcs, int bcscomp,
1013 const PreInterpHook& pre_interp,
1014 const PostInterpHook& post_interp)
1015{
1016#ifdef AMREX_USE_EB
1017 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
1018#else
1019 EB2::IndexSpace const* index_space = nullptr;
1020#endif
1021
1022 InterpFromCoarseLevel(mf,mf.nGrowVect(),time,index_space,cmf,scomp,dcomp,ncomp,cgeom,fgeom,
1023 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
1024 pre_interp,post_interp);
1025}
1026
1027template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
1028std::enable_if_t<IsFabArray<MF>::value>
1030 const Array<MF*, AMREX_SPACEDIM>& cmf, int scomp, int dcomp, int ncomp,
1031 const Geometry& cgeom, const Geometry& fgeom,
1032 Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
1033 Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
1034 const IntVect& ratio,
1035 Interp* mapper,
1036 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
1037 const PreInterpHook& pre_interp,
1038 const PostInterpHook& post_interp)
1039{
1040 InterpFromCoarseLevel(mf,mf[0]->nGrowVect(),time,cmf,scomp,dcomp,ncomp,cgeom,fgeom,
1041 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
1042 pre_interp,post_interp);
1043}
1044
1045template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
1046std::enable_if_t<IsFabArray<MF>::value>
1047InterpFromCoarseLevel (MF& mf, IntVect const& nghost, Real time,
1048 const MF& cmf, int scomp, int dcomp, int ncomp,
1049 const Geometry& cgeom, const Geometry& fgeom,
1050 BC& cbc, int cbccomp,
1051 BC& fbc, int fbccomp,
1052 const IntVect& ratio,
1053 Interp* mapper,
1054 const Vector<BCRec>& bcs, int bcscomp,
1055 const PreInterpHook& pre_interp,
1056 const PostInterpHook& post_interp)
1057{
1058#ifdef AMREX_USE_EB
1059 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
1060#else
1061 EB2::IndexSpace const* index_space = nullptr;
1062#endif
1063
1064 InterpFromCoarseLevel(mf,nghost,time,index_space,cmf,scomp,dcomp,ncomp,cgeom,fgeom,
1065 cbc,cbccomp,fbc,fbccomp,ratio,mapper,bcs,bcscomp,
1066 pre_interp,post_interp);
1067}
1068
1069template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
1070std::enable_if_t<IsFabArray<MF>::value>
1071InterpFromCoarseLevel (MF& mf, IntVect const& nghost, Real time,
1072 const EB2::IndexSpace* index_space,
1073 const MF& cmf, int scomp, int dcomp, int ncomp,
1074 const Geometry& cgeom, const Geometry& fgeom,
1075 BC& cbc, int cbccomp,
1076 BC& fbc, int fbccomp,
1077 const IntVect& ratio,
1078 Interp* mapper,
1079 const Vector<BCRec>& bcs, int bcscomp,
1080 const PreInterpHook& pre_interp,
1081 const PostInterpHook& post_interp)
1082{
1083 BL_PROFILE("InterpFromCoarseLevel");
1084
1085 using FAB = typename MF::FABType::value_type;
1086
1087 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
1088
1089 const BoxArray& ba = mf.boxArray();
1090 const DistributionMapping& dm = mf.DistributionMap();
1091
1092 const IndexType& typ = ba.ixType();
1093
1094 BL_ASSERT(typ == cmf.boxArray().ixType());
1095
1096 Box fdomain_g( amrex::convert(fgeom.Domain(),mf.ixType()) );
1097 for (int i = 0; i < AMREX_SPACEDIM; ++i) {
1098 if (fgeom.isPeriodic(i)) {
1099 fdomain_g.grow(i, nghost[i]);
1100 } else {
1101 if constexpr (std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
1102 fdomain_g.grow(i, fbc.nghost_outside_domain[i]);
1103 }
1104 }
1105 }
1106
1107 MF mf_crse_patch;
1108 IntVect send_ghost(0), recv_ghost(0);
1109 if constexpr (std::is_same_v<BC, PhysBCFunctUseCoarseGhost>) {
1110 mf_crse_patch.define(amrex::coarsen(ba,ratio), dm, ncomp, fbc.src_ghost);
1111 send_ghost = fbc.cghost;
1112 recv_ghost = fbc.src_ghost;
1113 } else {
1114 BoxArray ba_crse_patch(ba.size());
1115 { // TODO: later we might want to cache this
1116 for (int i = 0, N = ba.size(); i < N; ++i)
1117 {
1118 Box bx = amrex::convert(amrex::grow(ba[i],nghost), typ);
1119 bx &= fdomain_g;
1120 ba_crse_patch.set(i, coarsener.doit(bx));
1121 }
1122 }
1123
1124#ifndef AMREX_USE_EB
1125 amrex::ignore_unused(index_space);
1126#else
1127 if (index_space) {
1128 auto factory = makeEBFabFactory(index_space, cgeom, ba_crse_patch, dm,
1129 {0,0,0}, EBSupport::basic);
1130 mf_crse_patch.define(ba_crse_patch, dm, ncomp, 0, MFInfo(), *factory);
1131 } else
1132#endif
1133 {
1134 mf_crse_patch.define(ba_crse_patch, dm, ncomp, 0);
1135 }
1136 detail::mf_set_domain_bndry (mf_crse_patch, cgeom);
1137 }
1138
1139 mf_crse_patch.ParallelCopy(cmf, scomp, 0, ncomp, send_ghost, recv_ghost,
1140 cgeom.periodicity());
1141
1142 cbc(mf_crse_patch, 0, ncomp, mf_crse_patch.nGrowVect(), time, cbccomp);
1143
1144 detail::call_interp_hook(pre_interp, mf_crse_patch, 0, ncomp);
1145
1146 FillPatchInterp(mf, dcomp, mf_crse_patch, 0, ncomp, nghost, cgeom, fgeom, fdomain_g,
1147 ratio, mapper, bcs, bcscomp);
1148
1149#ifdef AMREX_USE_OMP
1150#pragma omp parallel if (Gpu::notInLaunchRegion())
1151#endif
1152 for (MFIter mfi(mf); mfi.isValid(); ++mfi)
1153 {
1154 FAB& dfab = mf[mfi];
1155 Box dfab_bx = dfab.box();
1156 dfab_bx.grow(nghost-mf.nGrowVect());
1157 const Box& dbx = dfab_bx & fdomain_g;
1158
1159 post_interp(dfab, dbx, dcomp, ncomp);
1160 }
1161
1162 fbc(mf, dcomp, ncomp, nghost, time, fbccomp);
1163}
1164
1165template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
1166std::enable_if_t<IsFabArray<MF>::value>
1168 const Array<MF*, AMREX_SPACEDIM>& cmf, int scomp, int dcomp, int ncomp,
1169 const Geometry& cgeom, const Geometry& fgeom,
1170 Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
1171 Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
1172 const IntVect& ratio,
1173 Interp* mapper,
1174 const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
1175 const PreInterpHook& pre_interp,
1176 const PostInterpHook& post_interp)
1177{
1178 BL_PROFILE("InterpFromCoarseLevel(array)");
1179
1180 using FAB = typename MF::FABType::value_type;
1181 using iFAB = typename iMultiFab::FABType::value_type;
1182
1183 const InterpolaterBoxCoarsener& coarsener = mapper->BoxCoarsener(ratio);
1184 const BoxArray& ba = mf[0]->boxArray();
1185 const DistributionMapping& dm = mf[0]->DistributionMap();
1186
1187 AMREX_ASSERT(AMREX_D_TERM(mf[0]->ixType().nodeCentered(0),
1188 && mf[1]->ixType().nodeCentered(1),
1189 && mf[2]->ixType().nodeCentered(2)));
1190
1191 // These need to be true: (ba[0] == ba[1] == ba[2]) & (dm[0] == dm[1] == dm[2]).
1192 // Debatable whether these are required, or will be enforced elsewhere prior to this func.
1194 && BoxArray::SameRefs(mf[0]->boxArray(), mf[1]->boxArray()),
1195 && BoxArray::SameRefs(mf[0]->boxArray(), mf[2]->boxArray())));
1196/*
1197 AMREX_ASSERT(AMREX_D_TERM(true,
1198 && DistributionMapping::SameRefs(mf[0]->DistributionMap(), mf[1]->DistributionMap()),
1199 && DistributionMapping::SameRefs(mf[0]->DistributionMap(), mf[2]->DistributionMap())));
1200*/
1201
1202 // If needed, adjust to fully overlap the coarse cells.
1203 IntVect nghost_adj = nghost;
1204 for (int d=0; d<AMREX_SPACEDIM; ++d) {
1205 if (nghost[d] % ratio[d] != 0) {
1206 nghost_adj[d] += ratio[d] - (nghost[d] % ratio[d]);
1207 }
1208 }
1209
1210 Array<MF*, AMREX_SPACEDIM> mf_local = mf;
1211 int dcomp_adj = dcomp;
1212 Array<std::unique_ptr<MF>, AMREX_SPACEDIM> mf_temp;
1213 if (! nghost.allGE(nghost_adj)) {
1214 for (int d=0; d<AMREX_SPACEDIM; ++d) {
1215 mf_temp[d] = std::make_unique<MF>(mf[d]->boxArray(),
1216 mf[d]->DistributionMap(), ncomp, nghost_adj);
1217 mf_local[d] = mf_temp[d].get();
1218 }
1219 dcomp_adj = 0;
1220 }
1221
1222 // Create a cell-centered boxArray of the region to interp.
1223 // Convert this boxArray and domain as needed.
1224 Box fdomain = amrex::convert(fgeom.Domain(), IntVect::TheZeroVector());
1225 Box fdomain_g(fdomain);
1226 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1227 if (fgeom.isPeriodic(d)) {
1228 fdomain_g.grow(d,nghost_adj[d]);
1229 }
1230 }
1231
1232 // Build patches, using domain to account for periodic bcs.
1233 BoxArray ba_crse_patch(ba.size());
1234 { // TODO: later we might want to cache this
1235 for (int i = 0, N = ba.size(); i < N; ++i)
1236 {
1237 Box bx = amrex::convert(amrex::grow(ba[i], nghost_adj), IntVect::TheZeroVector());
1238 bx &= fdomain_g;
1239 ba_crse_patch.set(i, coarsener.doit(bx));
1240 }
1241 }
1242
1243 Array<MF, AMREX_SPACEDIM> mf_crse_patch;
1244 for (int d = 0; d<AMREX_SPACEDIM; ++d)
1245 {
1246 IndexType typ = mf[d]->boxArray().ixType();
1247 BoxArray ba_crse_idxed = amrex::convert(ba_crse_patch, typ);
1248
1249#ifdef AMREX_USE_EB
1250 auto crse_factory = makeEBFabFactory(cgeom, ba_crse_idxed, dm, {0,0,0}, EBSupport::basic);
1251 mf_crse_patch[d].define(ba_crse_idxed, dm, ncomp, 0, MFInfo(), *crse_factory);
1252#else
1253 mf_crse_patch[d].define(ba_crse_idxed, dm, ncomp, 0);
1254#endif
1255 detail::mf_set_domain_bndry(mf_crse_patch[d], cgeom);
1256
1257 mf_crse_patch[d].ParallelCopy(*(cmf[d]), scomp, 0, ncomp, cgeom.periodicity());
1258 cbc[d](mf_crse_patch[d], 0, ncomp, mf_crse_patch[d].nGrowVect(), time, cbccomp);
1259 }
1260
1261 int idummy1=0, idummy2=0;
1262#ifdef AMREX_USE_OMP
1263#pragma omp parallel if (Gpu::notInLaunchRegion())
1264#endif
1265 {
1267
1268 // Empty containers describing that all points must be solved (no mask).
1269 Array<iFAB*, AMREX_SPACEDIM> mfab{ AMREX_D_DECL( nullptr, nullptr, nullptr ) };
1270
1271 for (MFIter mfi(mf_crse_patch[0]); mfi.isValid(); ++mfi)
1272 {
1273 Array<FAB*, AMREX_SPACEDIM> sfab{ AMREX_D_DECL( &(mf_crse_patch[0][mfi]),
1274 &(mf_crse_patch[1][mfi]),
1275 &(mf_crse_patch[2][mfi]) )};
1276 Array<FAB*, AMREX_SPACEDIM> dfab{ AMREX_D_DECL( &(*mf_local[0])[mfi],
1277 &(*mf_local[1])[mfi],
1278 &(*mf_local[2])[mfi] )};
1279
1280 const Box& sbx_cc = amrex::convert(sfab[0]->box(), IntVect::TheZeroVector());
1281 Box dfab_cc = amrex::convert(dfab[0]->box(), IntVect::TheZeroVector());
1282 const Box& dbx_cc = dfab_cc & fdomain_g;
1283
1284 for (int d=0; d<AMREX_SPACEDIM; ++d)
1285 {
1286 Vector<BCRec> bcr_d(ncomp);
1287
1288 amrex::setBC(sfab[d]->box(),
1289 amrex::convert(cgeom.Domain(), sfab[d]->box().ixType()),
1290 bcscomp,0,ncomp,bcs[d],bcr_d);
1291
1292 for (int n=0; n<ncomp; ++n)
1293 { bcr[n][d] = bcr_d[n]; }
1294 }
1295
1296 pre_interp(sfab, sbx_cc, 0, ncomp);
1297
1298 mapper->interp_arr(sfab, 0, dfab, 0, ncomp, dbx_cc, ratio, mfab,
1299 cgeom, fgeom, bcr, idummy1, idummy2, RunOn::Gpu);
1300
1301 post_interp(dfab, dbx_cc, 0, ncomp);
1302 }
1303 }
1304
1305 for (int d=0; d<AMREX_SPACEDIM; ++d)
1306 {
1307 if (mf[d] != mf_local[d]) {
1308 amrex::Copy(*mf[d], *mf_local[d], 0, dcomp_adj, ncomp, nghost);
1309 }
1310
1311 fbc[d](*mf[d], dcomp, ncomp, nghost, time, fbccomp);
1312 }
1313}
1314
1315template <typename MF, typename Interp>
1316std::enable_if_t<IsFabArray<MF>::value>
1317InterpFromCoarseLevel (MF& mf, IntVect const& nghost,
1318 IntVect const& nghost_outside_domain,
1319 const MF& cmf, int scomp, int dcomp, int ncomp,
1320 const Geometry& cgeom, const Geometry& fgeom,
1321 const IntVect& ratio, Interp* mapper,
1322 const Vector<BCRec>& bcs, int bcscomp)
1323{
1324 PhysBCFunctUseCoarseGhost erfbc(cmf,nghost,nghost_outside_domain,ratio,mapper);
1325 InterpFromCoarseLevel(mf, nghost, Real(0.0), cmf, scomp, dcomp, ncomp,
1326 cgeom, fgeom, erfbc, 0, erfbc, 0, ratio, mapper,
1327 bcs, bcscomp);
1328}
1329
1330template <typename MF>
1331std::enable_if_t<IsFabArray<MF>::value>
1332FillPatchSingleLevel (MF& mf, IntVect const& nghost, Real time,
1333 const Vector<MF*>& smf, IntVect const& snghost,
1334 const Vector<Real>& stime, int scomp, int dcomp, int ncomp,
1335 const Geometry& geom)
1336{
1337 PhysBCFunctUseCoarseGhost erfbc(snghost);
1338 FillPatchSingleLevel(mf, nghost, time, smf, stime, scomp, dcomp, ncomp, geom,
1339 erfbc, 0);
1340}
1341
1342template <typename MF, typename Interp>
1343std::enable_if_t<IsFabArray<MF>::value>
1344FillPatchTwoLevels (MF& mf, IntVect const& nghost,
1345 IntVect const& nghost_outside_domain, Real time,
1346 const Vector<MF*>& cmf, const Vector<Real>& ct,
1347 const Vector<MF*>& fmf, const Vector<Real>& ft,
1348 int scomp, int dcomp, int ncomp,
1349 const Geometry& cgeom, const Geometry& fgeom,
1350 const IntVect& ratio, Interp* mapper,
1351 const Vector<BCRec>& bcs, int bcscomp)
1352{
1353 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(nghost_outside_domain == 0, "TODO");
1354 PhysBCFunctUseCoarseGhost erfbc(*cmf[0], nghost, nghost_outside_domain, ratio,
1355 mapper);
1356 FillPatchTwoLevels(mf, nghost, time, cmf, ct, fmf, ft, scomp, dcomp, ncomp,
1357 cgeom, fgeom, erfbc, 0, erfbc, 0, ratio, mapper,
1358 bcs, bcscomp);
1359}
1360
1361template <typename MF, typename BC, typename Interp>
1362std::enable_if_t<IsFabArray<MF>::value>
1363FillPatchNLevels (MF& mf, int level, const IntVect& nghost, Real time,
1364 const Vector<Vector<MF*>>& smf, const Vector<Vector<Real>>& st,
1365 int scomp, int dcomp, int ncomp,
1366 const Vector<Geometry>& geom,
1367 Vector<BC>& bc, int bccomp,
1368 const Vector<IntVect>& ratio,
1369 Interp* mapper,
1370 const Vector<BCRec>& bcr, int bcrcomp)
1371{
1372 BL_PROFILE("FillPatchNLevels");
1373
1374 // FillPatchTwolevels relies on that mf's valid region is inside the
1375 // domain at periodic boundaries. But when we create coarsen boxarray
1376 // using mapper->CoarseBox, the resulting boxarray might violate the
1377 // requirement. If that happens, we need to create a second version of
1378 // the boxarray that is safe for FillPatchTwolevels.
1379
1380 auto get_clayout = [&] () -> std::tuple<BoxArray,BoxArray,DistributionMapping>
1381 {
1382 if (level == 0) {
1383 return std::make_tuple(BoxArray(),BoxArray(),DistributionMapping());
1384 } else {
1385 BoxArray const& ba = mf.boxArray();
1386 auto const& typ = ba.ixType();
1387 std::map<int,Vector<Box>> extra_boxes_map;
1388 BoxList cbl(typ);
1389 cbl.reserve(ba.size());
1390 for (int i = 0, N = int(ba.size()); i < N; ++i) {
1391 Box const& cbox = mapper->CoarseBox(amrex::grow(ba[i],nghost),ratio[level-1]);
1392 cbl.push_back(cbox);
1393 Box gdomain = geom[level-1].growNonPeriodicDomain(cbox.length());
1394 gdomain.convert(typ);
1395 if (!gdomain.contains(cbox)) {
1396 auto& extra_boxes = extra_boxes_map[i];
1397 auto const& pshift = geom[level-1].periodicity().shiftIntVect();
1398 for (auto const& piv : pshift) {
1399 auto const& ibox = amrex::shift(cbox,piv) & gdomain;
1400 if (ibox.ok()) {
1401 extra_boxes.push_back(ibox);
1402 }
1403 }
1404 }
1405 }
1406
1407 BoxArray cba2;
1409 if (!extra_boxes_map.empty()) {
1410 BoxList cbl2 = cbl;
1411 auto& lbox = cbl2.data();
1412 DistributionMapping const& dm = mf.DistributionMap();
1413 Vector<int> procmap2 = dm.ProcessorMap();
1414 for (auto const& [i, vb] : extra_boxes_map) {
1415 lbox[i] = vb[0];
1416 for (int j = 1, nj = int(vb.size()); j < nj; ++j) {
1417 lbox.push_back(vb[j]);
1418 procmap2.push_back(dm[i]);
1419 }
1420 }
1421 cba2 = BoxArray(std::move(cbl2));
1422 dm2 = DistributionMapping(std::move(procmap2));
1423 }
1424
1425 return std::make_tuple(BoxArray(std::move(cbl)), cba2, dm2);
1426 }
1427 };
1428
1429#ifdef AMREX_USE_EB
1430 EB2::IndexSpace const* index_space = EB2::TopIndexSpaceIfPresent();
1431#else
1432 EB2::IndexSpace const* index_space = nullptr;
1433#endif
1434
1435 AMREX_ALWAYS_ASSERT(level < int(geom.size()) &&
1436 level < int(bc.size()) &&
1437 level < int(ratio.size()+1));
1438 if (level == 0) {
1439 FillPatchSingleLevel(mf, nghost, time, smf[0], st[0], scomp, dcomp, ncomp, geom[0],
1440 bc[0], bccomp);
1441 } else if (level >= int(smf.size()))
1442 {
1443 auto const& [ba1, ba2, dm2] = get_clayout();
1444 MF cmf1, cmf2;
1445#ifdef AMREX_USE_EB
1446 if (index_space) {
1447 auto factory = makeEBFabFactory(index_space, geom[level-1], ba1,
1448 mf.DistributionMap(), {0,0,0},
1450 cmf1.define(ba1, mf.DistributionMap(), ncomp, 0, MFInfo(), *factory);
1451 if (!ba2.empty()) {
1452 auto factory2 = makeEBFabFactory(index_space, geom[level-1], ba2,
1453 dm2, {0,0,0},
1455 cmf2.define(ba2, dm2, ncomp, 0, MFInfo(), *factory2);
1456 }
1457 } else
1458#endif
1459 {
1460 cmf1.define(ba1, mf.DistributionMap(), ncomp, 0);
1461 if (!ba2.empty()) {
1462 cmf2.define(ba2, dm2, ncomp, 0);
1463 }
1464 }
1465
1466 MF* p_mf_inside = (ba2.empty()) ? &cmf1 : &cmf2;
1467 FillPatchNLevels(*p_mf_inside, level-1, IntVect(0), time, smf, st, scomp, 0, ncomp,
1468 geom, bc, bccomp, ratio, mapper, bcr, bcrcomp);
1469 if (&cmf1 != p_mf_inside) {
1470 cmf1.ParallelCopy(*p_mf_inside, geom[level-1].periodicity());
1471 }
1472 Box domain_g = geom[level].growPeriodicDomain(nghost);
1473 domain_g.convert(mf.ixType());
1474 FillPatchInterp(mf, dcomp, cmf1, 0, ncomp, nghost, geom[level-1], geom[level],
1475 domain_g, ratio[level-1], mapper, bcr, bcrcomp);
1476 } else {
1478 int error_code = detail::FillPatchTwoLevels_doit(mf, nghost, time,
1479 smf[level-1], st[level-1],
1480 smf[level ], st[level ],
1481 scomp, dcomp, ncomp,
1482 geom[level-1], geom[level],
1483 bc[level-1], bccomp,
1484 bc[level ], bccomp,
1485 ratio[level-1], mapper, bcr, bcrcomp,
1486 hook, hook, index_space, true);
1487 if (error_code == 0) { return; }
1488
1489 auto const& [ba1, ba2, dm2] = get_clayout();
1490 MF cmf_tmp;
1491#ifdef AMREX_USE_EB
1492 if (index_space) {
1493 if (ba2.empty()) {
1494 auto factory = makeEBFabFactory(index_space, geom[level-1], ba1,
1495 mf.DistributionMap(), {0,0,0},
1497 cmf_tmp.define(ba1, mf.DistributionMap(), ncomp, 0, MFInfo(), *factory);
1498 } else {
1499 auto factory = makeEBFabFactory(index_space, geom[level-1], ba2,
1500 dm2, {0,0,0},
1502 cmf_tmp.define(ba2, dm2, ncomp, 0, MFInfo(), *factory);
1503 }
1504 } else
1505#endif
1506 {
1507 if (ba2.empty()) {
1508 cmf_tmp.define(ba1, mf.DistributionMap(), ncomp, 0);
1509 } else {
1510 cmf_tmp.define(ba2, dm2, ncomp, 0);
1511 }
1512 }
1513
1514 FillPatchNLevels(cmf_tmp, level-1, IntVect(0), time, smf, st, scomp, 0, ncomp,
1515 geom, bc, bccomp, ratio, mapper, bcr, bcrcomp);
1516
1517 Vector<MF*> cmf{&cmf_tmp};
1518 Vector<MF*> fmf = smf[level];
1519 Vector<MF> fmf_raii;
1520 if (scomp != 0) {
1521 for (auto const* p : fmf) {
1522 fmf_raii.emplace_back(*p, amrex::make_alias, scomp, ncomp);
1523 }
1524 }
1525
1526 detail::FillPatchTwoLevels_doit(mf, nghost, time,
1527 cmf, {time},
1528 fmf, st[level],
1529 0, dcomp, ncomp,
1530 geom[level-1], geom[level],
1531 bc[level-1], bccomp,
1532 bc[level ], bccomp,
1533 ratio[level-1], mapper, bcr, bccomp,
1534 hook, hook, index_space);
1535 }
1536}
1537
1538}
1539
1540#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
IndexType ixType() const noexcept
Return index type of this BoxArray.
Definition AMReX_BoxArray.H:858
static bool SameRefs(const BoxArray &lhs, const BoxArray &rhs)
whether two BoxArrays share the same data
Definition AMReX_BoxArray.H:841
Long size() const noexcept
Return the number of boxes in the BoxArray.
Definition AMReX_BoxArray.H:615
void set(int i, const Box &ibox)
Set element i in this BoxArray to Box ibox.
Definition AMReX_BoxArray.cpp:878
A class for managing a List of Boxes that share a common IndexType. This class implements operations ...
Definition AMReX_BoxList.H:52
Vector< Box > & data() noexcept
Returns a reference to the Vector<Box>.
Definition AMReX_BoxList.H:215
void reserve(std::size_t n)
Definition AMReX_BoxList.H:90
void push_back(const Box &bn)
Append a Box to this BoxList.
Definition AMReX_BoxList.H:93
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ IntVectND< dim > length() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:154
__host__ __device__ BoxND & convert(IndexTypeND< dim > typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:974
__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__ BoxND & refine(int ref_ratio) noexcept
Refine BoxND by given (positive) refinement ratio. NOTE: if type(dir) = CELL centered: lo <- lo*ratio...
Definition AMReX_Box.H:698
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
const Vector< int > & ProcessorMap() const noexcept
Returns a constant reference to the mapping of boxes in the underlying BoxArray to the CPU that holds...
Definition AMReX_DistributionMapping.cpp:47
Definition AMReX_EB2.H:28
static const FPinfo & TheFPinfo(const FabArrayBase &srcfa, const FabArrayBase &dstfa, const IntVect &dstng, const BoxConverter &coarsener, const Geometry &fgeom, const Geometry &cgeom, const EB2::IndexSpace *)
Definition AMReX_FabArrayBase.cpp:2065
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:332
__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:232
__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:450
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:679
__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:400
Definition AMReX_InterpBase.H:34
Definition AMReX_InterpBase.H:20
Box doit(const Box &fine) const override
Apply the coarse-box logic to the supplied fine box.
Definition AMReX_InterpBase.cpp:10
Virtual base class for interpolaters.
Definition AMReX_Interpolater.H:27
Definition AMReX_MFInterpolater.H:20
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
Definition AMReX_PhysBCFunct.H:127
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
std::unique_ptr< EBFArrayBoxFactory > makeEBFabFactory(const Geometry &a_geom, const BoxArray &a_ba, const DistributionMapping &a_dm, const Vector< int > &a_ngrow, EBSupport a_support)
Definition AMReX_EBFabFactory.cpp:215
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__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
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
std::array< T, N > Array
Definition AMReX_Array.H:26
const IndexSpace * TopIndexSpaceIfPresent() noexcept
Definition AMReX_EB2.cpp:93
Definition AMReX_Amr.cpp:49
@ make_alias
Definition AMReX_MakeType.H:7
__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
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2866
IntVect nGrowVect(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2856
std::enable_if_t< IsFabArray< MF >::value > FillPatchNLevels(MF &mf, int level, const IntVect &nghost, Real time, const Vector< Vector< MF * > > &smf, const Vector< Vector< Real > > &st, int scomp, int dcomp, int ncomp, const Vector< Geometry > &geom, Vector< BC > &bc, int bccomp, const Vector< IntVect > &ratio, Interp *mapper, const Vector< BCRec > &bcr, int bcrcomp)
FillPatch with data from AMR levels.
Definition AMReX_FillPatchUtil_I.H:1363
RunOn
Definition AMReX_GpuControl.H:69
bool ProperlyNested(const IntVect &ratio, const IntVect &blocking_factor, int ngrow, const IndexType &boxType, Interp *mapper)
Test if AMR grids are properly nested.
Definition AMReX_FillPatchUtil_I.H:35
std::enable_if_t< IsFabArray< MF >::value > FillPatchSingleLevel(MF &mf, IntVect const &nghost, Real time, const Vector< MF * > &smf, const Vector< Real > &stime, int scomp, int dcomp, int ncomp, const Geometry &geom, BC &physbcf, int bcfcomp)
FillPatch with data from the current level.
Definition AMReX_FillPatchUtil_I.H:75
void Copy(FabArray< DFAB > &dst, FabArray< SFAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:180
std::enable_if_t< IsFabArray< MF >::value > FillPatchTwoLevels(MF &mf, IntVect const &nghost, Real time, const Vector< MF * > &cmf, const Vector< Real > &ct, const Vector< MF * > &fmf, const Vector< Real > &ft, int scomp, int dcomp, int ncomp, const Geometry &cgeom, const Geometry &fgeom, BC &cbc, int cbccomp, BC &fbc, int fbccomp, const IntVect &ratio, Interp *mapper, const Vector< BCRec > &bcs, int bcscomp, const PreInterpHook &pre_interp={}, const PostInterpHook &post_interp={})
FillPatch with data from the current level and the level below.
Definition AMReX_FillPatchUtil_I.H:817
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
std::enable_if_t< IsFabArray< MF >::value > InterpFromCoarseLevel(MF &mf, Real time, const MF &cmf, int scomp, int dcomp, int ncomp, const Geometry &cgeom, const Geometry &fgeom, BC &cbc, int cbccomp, BC &fbc, int fbccomp, const IntVect &ratio, Interp *mapper, const Vector< BCRec > &bcs, int bcscomp, const PreInterpHook &pre_interp={}, const PostInterpHook &post_interp={})
Fill with interpolation of coarse level data.
Definition AMReX_FillPatchUtil_I.H:1005
__host__ __device__ BoxND< dim > shift(const BoxND< dim > &b, int dir, int nzones) noexcept
Return a BoxND with indices shifted by nzones in dir direction.
Definition AMReX_Box.H:1495
@ basic
EBCellFlag.
IndexTypeND< 3 > IndexType
IndexType is an alias for amrex::IndexTypeND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:36
std::enable_if_t< IsFabArray< MF >::value &&!std::is_same_v< Interp, MFInterpolater > > InterpFace(Interp *interp, MF const &mf_crse_patch, int crse_comp, MF &mf_refined_patch, int fine_comp, int ncomp, const IntVect &ratio, const iMF &solve_mask, const Geometry &crse_geom, const Geometry &fine_geom, int bcscomp, RunOn gpu_or_cpu, const Vector< BCRec > &bcs)
Definition AMReX_FillPatchUtil_I.H:272
__host__ __device__ std::enable_if_t< std::is_floating_point_v< T >, bool > almostEqual(T x, T y, int ulp=2)
Definition AMReX_Algorithm.H:122
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
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
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2861
void FillPatchInterp(MultiFab &mf_fine_patch, int fcomp, MultiFab const &mf_crse_patch, int ccomp, int ncomp, IntVect const &ng, const Geometry &cgeom, const Geometry &fgeom, Box const &dest_domain, const IntVect &ratio, MFInterpolater *mapper, const Vector< BCRec > &bcs, int bcscomp)
Helper that applies a MFInterpolater to fill a fine patch from a coarse patch.
Definition AMReX_FillPatchUtil.cpp:140
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
Definition AMReX_FillPatchUtil.H:39