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