1#ifndef AMREX_NONLOCAL_BC_H_
5#ifndef AMREX_NONLOCAL_BC_IMPL_H_
6#define AMREX_NONLOCAL_BC_IMPL_H_
7#include <AMReX_Config.H>
11struct Rotate90ClockWise {
13 IntVect operator() (IntVect
const& iv)
const noexcept {
18 Dim3 operator() (Dim3
const& a)
const noexcept {
19 return Dim3{.x = a.y, .y = -1-a.x, .z = a.z};
22 Box operator() (Box
const& box)
const noexcept {
32struct Rotate90CounterClockWise {
34 IntVect operator() (IntVect
const& iv)
const noexcept {
39 Dim3 operator() (Dim3
const& a)
const noexcept {
40 return Dim3{.x = -1-a.y, .y = a.x, .z = a.z};
43 Box operator() (Box
const& box)
const noexcept {
53struct Rotate90DstToSrc
56 Dim3 operator() (Dim3
const& a)
const noexcept {
58 return Rotate90ClockWise()(a);
60 return Rotate90CounterClockWise()(a);
69 IntVect operator() (IntVect
const& iv)
const noexcept {
74 Dim3 operator() (Dim3
const& a)
const noexcept {
75 return Dim3{.x = -1-a.x, .y = Ly-1-a.y, .z = a.z};
78 Box operator() (Box
const& box)
const noexcept {
92 int i_index (
int i)
const noexcept {
93 return (i < Lx/2) ? -1-i : 2*Lx-1-i;
97 int j_index (
int j)
const noexcept {
98 return (j < Ly/2) ? j+Ly/2 : j-Ly/2;
102 IntVect operator() (IntVect
const& iv)
const noexcept {
107 Dim3 operator() (Dim3
const& a)
const noexcept {
108 return Dim3{.x = i_index(a.x), .y = j_index(a.y), .z = a.z};
111 [[nodiscard]]
Box operator() (Box
const& box)
const noexcept {
125 int i_index (
int i)
const noexcept {
126 return (i < Lx/2) ? -1-i : 2*Lx-1-i;
130 int j_index (
int j)
const noexcept {
133 }
else if (j >= Ly) {
135 }
else if (j < Ly/2) {
143 IntVect operator() (IntVect
const& iv)
const noexcept {
148 Dim3 operator() (Dim3
const& a)
const noexcept {
149 return Dim3{.x = i_index(a.x), .y = j_index(a.y), .z = a.z};
152 [[nodiscard]]
Box operator() (Box
const& box)
const noexcept {
185template <BaseFabType FAB,
class DTOS,
class Proj>
186requires (IsCallableR<Dim3, DTOS, Dim3>::value && IsFabProjection<Proj, FAB>::value)
188local_copy_cpu (FabArray<FAB>& dest,
const FabArray<FAB>& src,
int dcomp,
int scomp,
int ncomp,
191 const auto N_locs =
static_cast<int>(local_tags.size());
192 if (N_locs == 0) {
return; }
194#pragma omp parallel for
196 for (
int itag = 0; itag < N_locs; ++itag) {
197 const auto& tag = local_tags[itag];
198 auto const& sfab = src.const_array(tag.srcIndex);
199 auto const& dfab = dest.array (tag.dstIndex);
202 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
203 dfab(i,j,k,dcomp+n) = proj(sfab,si,scomp+n);
213template <BaseFabType FAB,
class DTOS,
class Proj>
214requires (IsCallableR<Dim3, DTOS, Dim3>::value && IsFabProjection<Proj, FAB>::value)
216unpack_recv_buffer_cpu (FabArray<FAB>& mf,
int dcomp,
int ncomp, Vector<char*>
const& recv_data,
217 Vector<std::size_t>
const& recv_size,
218 Vector<FabArrayBase::CopyComTagsContainer const*>
const& recv_cctc,
219 DTOS
const& dtos, Proj
const& proj)
noexcept
223 const auto N_rcvs =
static_cast<int>(recv_cctc.size());
224 if (N_rcvs == 0) {
return; }
226 using T =
typename FAB::value_type;
228#pragma omp parallel for
230 for (
int ircv = 0; ircv < N_rcvs; ++ircv) {
231 const char* dptr = recv_data[ircv];
232 auto const& cctc = *recv_cctc[ircv];
233 for (
auto const& tag : cctc) {
234 auto const& dfab = mf.array(tag.dstIndex);
237 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
238 dfab(i, j, k, dcomp + n) = proj(sfab, si, n);
240 dptr += tag.sbox.numPts() * ncomp *
sizeof(T);
241 AMREX_ASSERT(dptr <= recv_data[ircv] + recv_size[ircv]);
248struct Array4Array4Box {
250 Array4<T const> sfab;
254 Box const& box () const noexcept {
return dbox; }
257template <BaseFabType FAB,
class DTOS,
class Proj>
258requires (IsCallableR<Dim3, DTOS, Dim3>::value && IsFabProjection<Proj, FAB>::value)
260local_copy_gpu (FabArray<FAB>& dest,
const FabArray<FAB>& src,
int dcomp,
int scomp,
int ncomp,
261 FabArrayBase::CopyComTagsContainer
const& local_tags, DTOS
const& dtos, Proj
const& proj)
noexcept
263 int N_locs = local_tags.
size();
264 if (N_locs == 0) {
return; }
266 using T =
typename FAB::value_type;
267 Vector<Array4Array4Box<T> > loc_copy_tags;
268 loc_copy_tags.reserve(N_locs);
269 for (
auto const& tag : local_tags) {
270 loc_copy_tags.push_back(Array4Array4Box<T>{.dfab = dest.array(tag.dstIndex),
271 .sfab = src.const_array(tag.srcIndex),
276 Array4Array4Box<T>
const& tag)
noexcept
278 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
279 tag.dfab(i,j,k,dcomp+n) = proj(tag.sfab, si, scomp+n);
283template <BaseFabType FAB,
class DTOS,
class Proj>
284requires (IsCallableR<Dim3, DTOS, Dim3>::value && IsFabProjection<Proj, FAB>::value)
286unpack_recv_buffer_gpu (FabArray<FAB>& mf,
int scomp,
int ncomp,
287 Vector<char*>
const& recv_data,
288 Vector<std::size_t>
const& recv_size,
289 Vector<FabArrayBase::CopyComTagsContainer const*>
const& recv_cctc,
290 DTOS
const& dtos, Proj
const& proj)
294 const int N_rcvs = recv_cctc.size();
295 if (N_rcvs == 0) {
return; }
297 char* pbuffer = recv_data[0];
299 std::size_t szbuffer = 0;
302 if (not ParallelDescriptor::UseGpuAwareMpi()) {
304 szbuffer = (recv_data[N_rcvs-1]-recv_data[0]) + recv_size[N_rcvs-1];
306 Gpu::copyAsync(Gpu::hostToDevice,recv_data[0],recv_data[0]+szbuffer,pbuffer);
307 Gpu::streamSynchronize();
311 using T =
typename FAB::value_type;
312 using TagType = Array4Array4Box<T>;
313 Vector<TagType> tags;
314 tags.reserve(N_rcvs);
316 for (
int k = 0; k < N_rcvs; ++k)
318 std::size_t
offset = recv_data[k]-recv_data[0];
319 const char* dptr = pbuffer +
offset;
320 auto const& cctc = *recv_cctc[k];
321 for (
auto const& tag : cctc)
323 tags.emplace_back(TagType{mf.array(tag.dstIndex),
326 dptr += tag.dbox.numPts() * ncomp *
sizeof(T);
332 Array4Array4Box<T>
const& tag)
noexcept
334 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
335 tag.dfab(i,j,k,scomp+n) = proj(tag.sfab, si ,n);
340 if (pbuffer != recv_data[0]) {
346template <
typename DTOS>
347requires (IsIndexMapping<DTOS>::value)
348MultiBlockCommMetaData::MultiBlockCommMetaData (
const FabArrayBase& dst,
const Box& dstbox,
const FabArrayBase& src,
349 const IntVect& ngrow, DTOS
const& dtos)
350 : MultiBlockCommMetaData(dst.boxArray(), dst.DistributionMap(), dstbox, src.boxArray(),
351 src.DistributionMap(), ngrow, dtos) {}
353template <
typename DTOS>
354requires (IsIndexMapping<DTOS>::value)
355MultiBlockCommMetaData::MultiBlockCommMetaData (
const BoxArray& dstba,
const DistributionMapping& dstdm,
356 const Box& dstbox,
const BoxArray& srcba,
357 const DistributionMapping& srcdm,
const IntVect& ngrow, DTOS
const& dtos)
359 define(dstba, dstdm, dstbox, srcba, srcdm, ngrow, dtos);
362template <
typename DTOS>
363requires (IsIndexMapping<DTOS>::value)
365MultiBlockCommMetaData::define (
const BoxArray& dstba,
const DistributionMapping& dstdm,
const Box& dstbox,
366 const BoxArray& srcba,
const DistributionMapping& srcdm,
const IntVect& ngrow,
369 m_LocTags = std::make_unique<FabArrayBase::CopyComTagsContainer>();
370 m_SndTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
371 m_RcvTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
372 const int myproc = ParallelDescriptor::MyProc();
373 for (
int i = 0, N =
static_cast<int>(dstba.size()); i < N; ++i) {
374 const int dest_owner = dstdm[i];
375 const Box partial_dstbox =
grow(dstba[i], ngrow) & dstbox;
376 if (partial_dstbox.isEmpty()) {
379 const Box partial_dstbox_mapped_in_src =
Image(dtos, partial_dstbox).
setType(srcba.ixType());
380 enum { not_first_only = 0, first_only = 1 };
381 std::vector<std::pair<int, Box>> boxes_from_src =
382 srcba.intersections(partial_dstbox_mapped_in_src, not_first_only, ngrow);
383 for (std::pair<int, Box> counted_box : boxes_from_src) {
384 const int k = counted_box.first;
385 const Box src_box = counted_box.second;
387 const int src_owner = srcdm[k];
388 if (dest_owner == myproc || src_owner == myproc) {
389 if (src_owner == dest_owner) {
390 const BoxList tilelist(src_box, FabArrayBase::comm_tile_size);
391 for (
const Box& tilebox : tilelist) {
393 if ((inverse_image & partial_dstbox).ok()) {
394 m_LocTags->emplace_back(inverse_image, tilebox, i, k);
399 if ((inverse_image & partial_dstbox).ok()) {
400 FabArrayBase::CopyComTagsContainer& copy_tags =
401 (src_owner == myproc) ? (*m_SndTags)[dest_owner]
402 : (*m_RcvTags)[src_owner];
403 copy_tags.emplace_back(inverse_image, src_box, i, k);
411template <BaseFabType FAB,
class DTOS,
class Proj>
416Comm_nowait (FabArray<FAB>& mf,
int scomp,
int ncomp, FabArrayBase::CommMetaData
const& cmd,
417 DTOS
const& dtos, Proj
const& proj)
420 if (ParallelContext::NProcsSub() == 1)
423 if (cmd.m_LocTags->empty()) {
return CommHandler{}; }
425 if (Gpu::inLaunchRegion()) {
426 local_copy_gpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
430 local_copy_cpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
432 return CommHandler{};
440 int SeqNum = ParallelDescriptor::SeqNum();
442 const auto N_locs = cmd.m_LocTags->size();
443 const auto N_rcvs = cmd.m_RcvTags->size();
444 const auto N_snds = cmd.m_SndTags->size();
446 if (N_locs == 0 && N_rcvs == 0 && N_snds == 0) {
448 return CommHandler{};
451 CommHandler handler{};
455 handler.recv.the_data = FabArray<FAB>::PostRcvs(*cmd.m_RcvTags, handler.recv.data, handler.recv.size,
456 handler.recv.rank, handler.recv.request, ncomp, SeqNum);
460 handler.send.the_data =
461 FabArray<FAB>::PrepareSendBuffers(*cmd.m_SndTags, handler.send.data, handler.send.size,
462 handler.send.rank, handler.send.request, handler.send.cctc, ncomp);
465 if (Gpu::inLaunchRegion()) {
466 FabArray<FAB>::pack_send_buffer_gpu(mf, scomp, ncomp, handler.send.data,
467 handler.send.size, handler.send.cctc, handler.send.id);
471 FabArray<FAB>::pack_send_buffer_cpu(mf, scomp, ncomp, handler.send.data,
472 handler.send.size, handler.send.cctc);
475 FabArray<FAB>::PostSnds(handler.send.data, handler.send.size, handler.send.rank, handler.send.request, SeqNum);
481 if (Gpu::inLaunchRegion()) {
482 local_copy_gpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
486 local_copy_cpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
495template <BaseFabType FAB,
class DTOS,
class Proj>
497Comm_finish (FabArray<FAB>& mf,
int scomp,
int ncomp, FabArrayBase::CommMetaData
const& cmd,
498 CommHandler handler, DTOS
const& dtos, Proj
const& proj)
500 if (ParallelContext::NProcsSub() == 1) {
return; }
502 const auto N_rcvs =
static_cast<int>(cmd.m_RcvTags->size());
505 handler.recv.cctc.resize(N_rcvs,
nullptr);
506 for (
int k = 0; k < N_rcvs; ++k) {
507 auto const& cctc = cmd.m_RcvTags->at(handler.recv.rank[k]);
508 handler.recv.cctc[k] = &cctc;
510 handler.recv.stats.resize(handler.recv.request.size());
511 ParallelDescriptor::Waitall(handler.recv.request, handler.recv.stats);
513 if (!CheckRcvStats(handler.recv.stats, handler.recv.size, handler.mpi_tag)) {
514 amrex::Abort(
"NonLocalBC::Comm_finish failed with wrong message size");
519 if (Gpu::inLaunchRegion())
522 handler.recv.size, handler.recv.cctc, dtos, proj);
527 handler.recv.size, handler.recv.cctc, dtos, proj);
531 if ( ! cmd.m_SndTags->empty() ) {
532 handler.send.stats.resize(handler.send.request.size());
533 ParallelDescriptor::Waitall(handler.send.request, handler.send.stats);
538template <BaseFabType FAB>
540Rotate90 (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
548 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
549 AMREX_ASSERT(nghost.allLE(mf.nGrowVect()) && nghost[0] == nghost[1]);
551 if (nghost[0] <= 0) {
return; }
553 const FabArrayBase::RB90& TheRB90 = mf.getRB90(nghost, domain);
555 auto handler = Comm_nowait(mf, scomp, ncomp, TheRB90,Rotate90DstToSrc{},
558 Box corner(-nghost, IntVect{
AMREX_D_DECL(-1,-1,domain.bigEnd(2)+nghost[2])});
560#pragma omp parallel if (Gpu::notInLaunchRegion())
562 for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
563 Box const& bx = corner & mfi.fabbox();
565 auto const& fab = mf.array(mfi);
568 fab(i,j,k,n) = fab(-i-1,-j-1,k,n);
574 Comm_finish(mf, scomp, ncomp, TheRB90, std::move(handler), Rotate90DstToSrc{},
581template <BaseFabType FAB>
583Rotate90 (FabArray<FAB>& mf, Box
const& domain)
585 Rotate90(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
588template <BaseFabType FAB>
590Rotate180 (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
598 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
601 if (nghost[0] <= 0) {
return; }
603 const FabArrayBase::RB180& TheRB180 = mf.getRB180(nghost, domain);
605 auto handler = Comm_nowait(mf, scomp, ncomp, TheRB180,
606 Rotate180Fn{domain.length(1)}, Identity{});
609 Comm_finish(mf, scomp, ncomp, TheRB180, std::move(handler),
610 Rotate180Fn{domain.length(1)}, Identity{});
616template <BaseFabType FAB>
618Rotate180 (FabArray<FAB>& mf, Box
const& domain)
620 Rotate180(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
623template <BaseFabType FAB>
625FillPolar (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
633 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
636 if (nghost[0] <= 0) {
return; }
638 const FabArrayBase::PolarB& ThePolarB = mf.getPolarB(nghost, domain);
640 auto handler = Comm_nowait(mf, scomp, ncomp, ThePolarB,
641 PolarFn{.Lx = domain.length(0), .Ly = domain.length(1)},
645 Comm_finish(mf, scomp, ncomp, ThePolarB, std::move(handler),
646 PolarFn{.Lx = domain.length(0), .Ly = domain.length(1)}, Identity{});
652template <BaseFabType FAB>
654FillPolar (FabArray<FAB>& mf, Box
const& domain)
656 FillPolar(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
659template <BaseFabType FAB,
typename DTOS,
typename Proj>
660requires (IsCallableR<Dim3,DTOS,Dim3>::value && IsFabProjection<Proj,FAB>::value)
662FillBoundary_nowait (FabArray<FAB>& mf,
const FabArrayBase::CommMetaData& cmd,
663 int scomp,
int ncomp, DTOS
const& dtos, Proj
const& proj)
666 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
667 return Comm_nowait(mf, scomp, ncomp, cmd, dtos, proj);
670template <BaseFabType FAB,
typename DTOS,
typename Proj>
671requires (IsCallableR<Dim3,DTOS,Dim3>::value && IsFabProjection<Proj,FAB>::value)
673FillBoundary_finish (CommHandler handler,
674 FabArray<FAB>& mf,
const FabArrayBase::CommMetaData& cmd,
675 int scomp,
int ncomp, DTOS
const& dtos, Proj
const& proj)
679 Comm_finish(mf, scomp, ncomp, cmd, std::move(handler), dtos, proj);
685template <
typename DTOS>
686Vector<std::pair<Box,Box>>
687get_src_dst_boxes (DTOS
const& dtos, Box
const& dstbox, Box
const& domain)
689 Vector<std::pair<Box,Box>> r;
692 if (!domain.contains(mapped_smallend) || !domain.contains(mapped_bigend)) {
698 auto dtype = dstbox.type();
702 Array<Array<std::pair<int,int>,2>,AMREX_SPACEDIM> ends;
703 Array<Array<std::pair<int,int>,2>,AMREX_SPACEDIM> dst_ends;
704 Array<int,AMREX_SPACEDIM> nboxes;
705 for (
int ddim = 0; ddim < AMREX_SPACEDIM; ++ddim) {
706 int sdim = perm[ddim];
707 auto mm = std::minmax(mapped_smallend[sdim],mapped_bigend[sdim]);
708 if (((sign[ddim] > 0) && (mapped_smallend[sdim] <= mapped_bigend[sdim])) ||
709 ((sign[ddim] < 0) && (mapped_bigend[sdim] <= mapped_smallend[sdim])))
713 dst_ends[ddim][0] = std::make_pair(dstbox.smallEnd(ddim),
714 dstbox.bigEnd(ddim));
717 ends[sdim][0].first = domain.smallEnd(sdim);
718 ends[sdim][0].second = mm.first;
719 ends[sdim][1].first = mm.second;
720 ends[sdim][1].second = domain.bigEnd(sdim);
721 int n0 = ends[sdim][0].second - ends[sdim][0].first;
722 int n1 = ends[sdim][1].second - ends[sdim][1].first;
723 if (mm.first == mapped_smallend[sdim]) {
724 dst_ends[ddim][0] = std::make_pair(dstbox.smallEnd(ddim),
725 dstbox.smallEnd(ddim)+n0);
726 dst_ends[ddim][1] = std::make_pair(dstbox.bigEnd(ddim)-n1,
727 dstbox.bigEnd(ddim));
729 dst_ends[ddim][0] = std::make_pair(dstbox.bigEnd(ddim)-n0,
730 dstbox.bigEnd(ddim));
731 dst_ends[ddim][1] = std::make_pair(dstbox.smallEnd(ddim),
732 dstbox.smallEnd(ddim)+n1);
737 r.reserve(
AMREX_D_TERM(nboxes[0],*nboxes[1],*nboxes[2]));
739#if (AMREX_SPACEDIM == 3)
740 for (
int kbox = 0; kbox < nboxes[2]; ++kbox) {
742#if (AMREX_SPACEDIM >=2 )
743 for (
int jbox = 0; jbox < nboxes[1]; ++jbox) {
745 for (
int ibox = 0; ibox < nboxes[0]; ++ibox)
751 ends[2][kbox].first)),
753 ends[1][jbox].second,
754 ends[2][kbox].second)),
757 dst_ends[1][div[1]].first,
758 dst_ends[2][div[2]].first)),
760 dst_ends[1][div[1]].second,
761 dst_ends[2][div[2]].second)),
768template <
typename DTOS>
769Box get_dst_subbox (DTOS
const& dtos, std::pair<Box,Box>
const& sdboxes,
770 Box
const& srcsubbox)
772 Box const& srcbox = sdboxes.first;
773 Box const& dstbox = sdboxes.second;
774 if (srcbox == srcsubbox) {
779 Box dstsubbox = dstbox;
780 for (
int ddim = 0; ddim < AMREX_SPACEDIM; ++ddim) {
781 int sdim = perm[ddim];
782 if (sign[ddim] > 0) {
783 dstsubbox.
growLo(ddim, srcbox.smallEnd(sdim)-srcsubbox.smallEnd(sdim));
784 dstsubbox.
growHi(ddim, srcsubbox.bigEnd(sdim)-srcbox.bigEnd(sdim));
786 dstsubbox.
growLo(ddim, srcsubbox.bigEnd(sdim)-srcbox.bigEnd(sdim));
787 dstsubbox.
growHi(ddim, srcbox.smallEnd(sdim)-srcsubbox.smallEnd(sdim));
796 void split_boxes (BoxList& bl, Box
const& domain);
800template <BaseFabType FAB,
typename DTOS>
801requires (IsCallableR<Dim3,DTOS,Dim3>::value)
802FabArrayBase::CommMetaData
803makeFillBoundaryMetaData (FabArray<FAB>& mf, IntVect
const& nghost,
804 Geometry
const& geom, DTOS
const& dtos)
806 FabArrayBase::CommMetaData cmd;
807 cmd.m_LocTags = std::make_unique<FabArrayBase::CopyComTagsContainer>();
808 cmd.m_SndTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
809 cmd.m_RcvTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
812 mf.define_fb_metadata(cmd, nghost,
false, geom.periodicity(),
false);
814 BoxArray
const& ba = mf.boxArray();
815 DistributionMapping
const& dm = mf.DistributionMap();
819 const int myproc = ParallelDescriptor::MyProc();
820 const auto nboxes =
static_cast<int>(ba.size());
821 std::vector<std::pair<int,Box> > isects;
823 for (
int i = 0; i < nboxes; ++i) {
826 if (bl.isEmpty()) {
continue; }
828 detail::split_boxes(bl, dombox);
830 const int dst_owner = dm[i];
831 for (
auto const& dst_box : bl) {
832 auto const& src_dst_boxes = get_src_dst_boxes(dtos, dst_box, dombox);
833 for (
auto const& sd_box_pair : src_dst_boxes) {
834 ba.intersections(sd_box_pair.first, isects);
835 for (
auto const& is : isects) {
836 int const k = is.first;
837 Box const src_b = is.second;
838 int const src_owner = dm[k];
839 if (dst_owner == myproc || src_owner == myproc) {
840 Box const& dst_b = get_dst_subbox(dtos, sd_box_pair, src_b);
841 if (src_owner == dst_owner) {
842 cmd.m_LocTags->emplace_back(dst_b, src_b, i, k);
844 auto& tags = (dst_owner == myproc) ?
845 (*cmd.m_RcvTags)[src_owner] :
846 (*cmd.m_SndTags)[dst_owner];
847 tags.emplace_back(dst_b, src_b, i, k);
858struct SphThetaPhiRIndexMapping
860 SphThetaPhiRIndexMapping (Box
const& a_domain)
869 Dim3 operator() (Dim3
const& ijk)
const noexcept
876 bool imd = i >= 0 && i < nx;
879 bool jmd = j >= 0 && j < ny;
881 bool kmd = k >= 0 && k < nz;
884 if (ilo && jmd && kmd)
886 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
888 else if (ihi && jmd && kmd)
890 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
892 else if (imd && jlo && kmd)
894 return Dim3{.x = i, .y = j+ny, .z = k};
896 else if (imd && jhi && kmd)
898 return Dim3{.x = i, .y = j-ny, .z = k};
900 else if (imd && jmd && klo)
902 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
904 else if (ilo && jlo && kmd)
906 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
908 else if (ihi && jlo && kmd)
910 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
912 else if (ilo && jhi && kmd)
914 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
916 else if (ihi && jhi && kmd)
918 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
920 else if (imd && jlo && klo)
922 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
924 else if (imd && jhi && klo)
926 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
934 [[nodiscard]]
IntVect sign (Dim3
const& ijk)
const noexcept
938 }
else if (ijk.z >=0 && ijk.z < nz &&
939 (ijk.x < 0 || ijk.x >= nx)) {
946 [[nodiscard]]
IntVect permutation (Dim3
const& )
const noexcept
955struct SphThetaPhiRComponentMapping
957 SphThetaPhiRComponentMapping (Box
const& a_domain,
int a_start_index)
961 scomp(a_start_index) {}
963 template <
typename T>
965 T operator()(Array4<const T>
const& a, Dim3
const& ijk,
int n)
const noexcept
972 if ((i >= 0 && i < nx) &&
973 (j < 0 || j >= ny) &&
974 (k >= 0 && k < nz)) {
981 }
else if (n == scomp+2) {
996extern template MultiBlockCommMetaData
ParallelCopy(FabArray<FArrayBox>& dest,
const Box& destbox,
997 const FabArray<FArrayBox>& src,
int destcomp,
998 int srccomp,
int numcomp,
const IntVect& ngrow,
999 MultiBlockIndexMapping
const&, Identity
const&);
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_NODISCARD
Definition AMReX_Extension.H:252
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
virtual void free(void *pt)=0
A pure virtual function for deleting the arena pointed to by pt.
virtual void * alloc(std::size_t sz)=0
__host__ __device__ IntVectND< dim > size() const noexcept
Return the length of the BoxND.
Definition AMReX_Box.H:148
__host__ __device__ BoxND & growLo(int idir, int n_cell=1) noexcept
Grow the BoxND on the low end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the Bo...
Definition AMReX_Box.H:670
__host__ __device__ BoxND & growHi(int idir, int n_cell=1) noexcept
Grow the BoxND on the high end by n_cell cells in direction idir. NOTE: n_cell negative shrinks the B...
Definition AMReX_Box.H:681
__host__ __device__ BoxND & setType(const IndexTypeND< dim > &t) noexcept
Set indexing type.
Definition AMReX_Box.H:513
CopyComTag::CopyComTagsContainer CopyComTagsContainer
Definition AMReX_FabArrayBase.H:220
__host__ __device__ Dim3 ubound(Array4< T > const &a) noexcept
Return the inclusive upper bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1359
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
Return the spatial extents of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1373
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Return the inclusive lower bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1345
__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
Arena * The_Arena()
Definition AMReX_Arena.cpp:820
Definition AMReX_NonLocalBC.cpp:39
void Rotate90(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
void unpack_recv_buffer_cpu(FabArray< FAB > &mf, int dcomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< FabArrayBase::CopyComTagsContainer const * > const &recv_cctc, DTOS const &dtos=DTOS{}, Proj const &proj=Proj{}) noexcept
Box InverseImage(DTOS const &dtos, const Box &box)
Applies the inverse Dim3 to Dim3 mapping onto Boxes without changing the index type.
Definition AMReX_NonLocalBC.H:189
void Rotate180(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
void local_copy_cpu(FabArray< FAB > &dest, const FabArray< FAB > &src, int dcomp, int scomp, int ncomp, FabArrayBase::CopyComTagsContainer const &local_tags, DTOS const &dtos=DTOS{}, Proj const &proj=Proj{}) noexcept
void FillPolar(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
void local_copy_gpu(FabArray< FAB > &dest, const FabArray< FAB > &src, int dcomp, int scomp, int ncomp, FabArrayBase::CopyComTagsContainer const &local_tags, DTOS const &dtos=DTOS{}, Proj const &proj=Proj{}) noexcept
Box Image(DTOS const &dtos, const Box &box)
Applies the Dim3 to Dim3 mapping onto Boxes but does not change the index type.
Definition AMReX_NonLocalBC.H:123
void unpack_recv_buffer_gpu(FabArray< FAB > &mf, int scomp, int ncomp, Vector< char * > const &recv_data, Vector< std::size_t > const &recv_size, Vector< FabArrayBase::CopyComTagsContainer const * > const &recv_cctc, DTOS const &dtos=DTOS{}, Proj const &proj=Proj{})
int SeqNum() noexcept
Returns sequential message sequence numbers, usually used as tags for send/recv.
Definition AMReX_ParallelDescriptor.H:678
__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
__host__ __device__ Array4< T > makeArray4(T *p, Box const &bx, int ncomp) noexcept
Definition AMReX_BaseFab.H:92
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
BoxList boxDiff(const Box &b1in, const Box &b2)
Returns BoxList defining the compliment of b2 in b1in.
Definition AMReX_BoxList.cpp:599
void ParallelCopy(MF &dst, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &ng_src=IntVect(0), IntVect const &ng_dst=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic())
dst = src w/ MPI communication
Definition AMReX_FabArrayUtility.H:1951
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void LoopConcurrentOnCpu(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:388
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241