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 {
186template <
class FAB,
class DTOS,
class Proj>
187std::enable_if_t<IsBaseFab<FAB>() && IsCallableR<Dim3, DTOS, Dim3>() && IsFabProjection<Proj, FAB>()>
188local_copy_cpu (FabArray<FAB>& dest,
const FabArray<FAB>& src,
int dcomp,
int scomp,
int ncomp,
190 const auto N_locs =
static_cast<int>(local_tags.size());
191 if (N_locs == 0) {
return; }
193#pragma omp parallel for
195 for (
int itag = 0; itag < N_locs; ++itag) {
196 const auto& tag = local_tags[itag];
197 auto const& sfab = src.const_array(tag.srcIndex);
198 auto const& dfab = dest.array (tag.dstIndex);
201 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
202 dfab(i,j,k,dcomp+n) = proj(sfab,si,scomp+n);
212template <
class FAB,
class DTOS,
class Proj>
213std::enable_if_t<IsBaseFab<FAB>() && IsCallableR<Dim3, DTOS, Dim3>() && IsFabProjection<Proj, FAB>()>
215 Vector<std::size_t>
const& recv_size,
216 Vector<FabArrayBase::CopyComTagsContainer const*>
const& recv_cctc,
217 DTOS
const& dtos, Proj
const& proj)
noexcept {
220 const auto N_rcvs =
static_cast<int>(recv_cctc.size());
221 if (N_rcvs == 0) {
return; }
223 using T =
typename FAB::value_type;
225#pragma omp parallel for
227 for (
int ircv = 0; ircv < N_rcvs; ++ircv) {
228 const char* dptr = recv_data[ircv];
229 auto const& cctc = *recv_cctc[ircv];
230 for (
auto const& tag : cctc) {
231 auto const& dfab = mf.array(tag.dstIndex);
234 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
235 dfab(i, j, k, dcomp + n) = proj(sfab, si, n);
237 dptr += tag.sbox.numPts() * ncomp *
sizeof(T);
238 AMREX_ASSERT(dptr <= recv_data[ircv] + recv_size[ircv]);
245struct Array4Array4Box {
247 Array4<T const> sfab;
251 Box const& box () const noexcept {
return dbox; }
254template <
class FAB,
class DTOS,
class Proj>
255std::enable_if_t<IsBaseFab<FAB>() && IsCallableR<Dim3, DTOS, Dim3>() && IsFabProjection<Proj, FAB>()>
256local_copy_gpu (FabArray<FAB>& dest,
const FabArray<FAB>& src,
int dcomp,
int scomp,
int ncomp,
257 FabArrayBase::CopyComTagsContainer
const& local_tags, DTOS
const& dtos, Proj
const& proj)
noexcept {
258 int N_locs = local_tags.size();
259 if (N_locs == 0) {
return; }
261 using T =
typename FAB::value_type;
262 Vector<Array4Array4Box<T> > loc_copy_tags;
263 loc_copy_tags.reserve(N_locs);
264 for (
auto const& tag : local_tags) {
265 loc_copy_tags.push_back(Array4Array4Box<T>{.dfab = dest.array(tag.dstIndex),
266 .sfab = src.const_array(tag.srcIndex),
271 Array4Array4Box<T>
const& tag)
noexcept
273 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
274 tag.dfab(i,j,k,dcomp+n) = proj(tag.sfab, si, scomp+n);
278template <
class FAB,
class DTOS,
class Proj>
279std::enable_if_t<IsBaseFab<FAB>() && IsCallableR<Dim3, DTOS, Dim3>() && IsFabProjection<Proj, FAB>()>
281 Vector<char*>
const& recv_data,
282 Vector<std::size_t>
const& recv_size,
283 Vector<FabArrayBase::CopyComTagsContainer const*>
const& recv_cctc,
284 DTOS
const& dtos, Proj
const& proj)
288 const int N_rcvs = recv_cctc.size();
289 if (N_rcvs == 0) {
return; }
291 char* pbuffer = recv_data[0];
293 std::size_t szbuffer = 0;
296 if (not ParallelDescriptor::UseGpuAwareMpi()) {
298 szbuffer = (recv_data[N_rcvs-1]-recv_data[0]) + recv_size[N_rcvs-1];
300 Gpu::copyAsync(Gpu::hostToDevice,recv_data[0],recv_data[0]+szbuffer,pbuffer);
301 Gpu::streamSynchronize();
305 using T =
typename FAB::value_type;
306 using TagType = Array4Array4Box<T>;
307 Vector<TagType> tags;
308 tags.reserve(N_rcvs);
310 for (
int k = 0; k < N_rcvs; ++k)
312 std::size_t
offset = recv_data[k]-recv_data[0];
313 const char* dptr = pbuffer +
offset;
314 auto const& cctc = *recv_cctc[k];
315 for (
auto const& tag : cctc)
317 tags.emplace_back(TagType{mf.array(tag.dstIndex),
320 dptr += tag.dbox.numPts() * ncomp *
sizeof(T);
326 Array4Array4Box<T>
const& tag)
noexcept
328 auto const si = dtos(Dim3{.x = i, .y = j, .z = k});
329 tag.dfab(i,j,k,scomp+n) = proj(tag.sfab, si ,n);
334 if (pbuffer != recv_data[0]) {
340template <
typename DTOS,
typename>
341MultiBlockCommMetaData::MultiBlockCommMetaData (
const FabArrayBase& dst,
const Box& dstbox,
const FabArrayBase& src,
342 const IntVect& ngrow, DTOS
const& dtos)
346template <
typename DTOS,
typename>
347MultiBlockCommMetaData::MultiBlockCommMetaData (
const BoxArray& dstba,
const DistributionMapping& dstdm,
348 const Box& dstbox,
const BoxArray& srcba,
349 const DistributionMapping& srcdm,
const IntVect& ngrow, DTOS
const& dtos) {
350 define(dstba, dstdm, dstbox, srcba, srcdm, ngrow, dtos);
353template <
typename DTOS>
354std::enable_if_t<IsIndexMapping<DTOS>::value>
355MultiBlockCommMetaData::define (
const BoxArray& dstba,
const DistributionMapping& dstdm,
const Box& dstbox,
356 const BoxArray& srcba,
const DistributionMapping& srcdm,
const IntVect& ngrow,
358 m_LocTags = std::make_unique<FabArrayBase::CopyComTagsContainer>();
359 m_SndTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
360 m_RcvTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
361 const int myproc = ParallelDescriptor::MyProc();
362 for (
int i = 0, N =
static_cast<int>(dstba.size()); i < N; ++i) {
363 const int dest_owner = dstdm[i];
364 const Box partial_dstbox =
grow(dstba[i], ngrow) & dstbox;
365 if (partial_dstbox.isEmpty()) {
368 const Box partial_dstbox_mapped_in_src =
Image(dtos, partial_dstbox).setType(srcba.ixType());
369 enum { not_first_only = 0, first_only = 1 };
370 std::vector<std::pair<int, Box>> boxes_from_src =
371 srcba.intersections(partial_dstbox_mapped_in_src, not_first_only, ngrow);
372 for (std::pair<int, Box> counted_box : boxes_from_src) {
373 const int k = counted_box.first;
374 const Box src_box = counted_box.second;
376 const int src_owner = srcdm[k];
377 if (dest_owner == myproc || src_owner == myproc) {
378 if (src_owner == dest_owner) {
379 const BoxList tilelist(src_box, FabArrayBase::comm_tile_size);
380 for (
const Box& tilebox : tilelist) {
381 const Box inverse_image =
InverseImage(dtos, tilebox).setType(dstba.ixType());
382 if ((inverse_image & partial_dstbox).ok()) {
383 m_LocTags->emplace_back(inverse_image, tilebox, i, k);
387 const Box inverse_image =
InverseImage(dtos, src_box).setType(dstba.ixType());
388 if ((inverse_image & partial_dstbox).ok()) {
389 FabArrayBase::CopyComTagsContainer& copy_tags =
390 (src_owner == myproc) ? (*m_SndTags)[dest_owner]
391 : (*m_RcvTags)[src_owner];
392 copy_tags.emplace_back(inverse_image, src_box, i, k);
400template <
class FAB,
class DTOS,
class Proj>
405Comm_nowait (FabArray<FAB>& mf,
int scomp,
int ncomp, FabArrayBase::CommMetaData
const& cmd,
406 DTOS
const& dtos, Proj
const& proj)
409 if (ParallelContext::NProcsSub() == 1)
412 if (cmd.m_LocTags->empty()) {
return CommHandler{}; }
414 if (Gpu::inLaunchRegion()) {
415 local_copy_gpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
419 local_copy_cpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
421 return CommHandler{};
429 int SeqNum = ParallelDescriptor::SeqNum();
431 const auto N_locs = cmd.m_LocTags->size();
432 const auto N_rcvs = cmd.m_RcvTags->size();
433 const auto N_snds = cmd.m_SndTags->size();
435 if (N_locs == 0 && N_rcvs == 0 && N_snds == 0) {
437 return CommHandler{};
440 CommHandler handler{};
444 handler.recv.the_data = FabArray<FAB>::PostRcvs(*cmd.m_RcvTags, handler.recv.data, handler.recv.size,
445 handler.recv.rank, handler.recv.request, ncomp, SeqNum);
449 handler.send.the_data =
450 FabArray<FAB>::PrepareSendBuffers(*cmd.m_SndTags, handler.send.data, handler.send.size,
451 handler.send.rank, handler.send.request, handler.send.cctc, ncomp);
454 if (Gpu::inLaunchRegion()) {
455 FabArray<FAB>::pack_send_buffer_gpu(mf, scomp, ncomp, handler.send.data,
456 handler.send.size, handler.send.cctc, handler.send.id);
460 FabArray<FAB>::pack_send_buffer_cpu(mf, scomp, ncomp, handler.send.data,
461 handler.send.size, handler.send.cctc);
464 FabArray<FAB>::PostSnds(handler.send.data, handler.send.size, handler.send.rank, handler.send.request, SeqNum);
470 if (Gpu::inLaunchRegion()) {
471 local_copy_gpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
475 local_copy_cpu(mf, mf, scomp, scomp, ncomp, *cmd.m_LocTags, dtos, proj);
484template <
class FAB,
class DTOS,
class Proj>
486Comm_finish (FabArray<FAB>& mf,
int scomp,
int ncomp, FabArrayBase::CommMetaData
const& cmd,
487 CommHandler handler, DTOS
const& dtos, Proj
const& proj)
489 if (ParallelContext::NProcsSub() == 1) {
return; }
491 const auto N_rcvs =
static_cast<int>(cmd.m_RcvTags->size());
494 handler.recv.cctc.resize(N_rcvs,
nullptr);
495 for (
int k = 0; k < N_rcvs; ++k) {
496 auto const& cctc = cmd.m_RcvTags->at(handler.recv.rank[k]);
497 handler.recv.cctc[k] = &cctc;
499 handler.recv.stats.resize(handler.recv.request.size());
500 ParallelDescriptor::Waitall(handler.recv.request, handler.recv.stats);
502 if (!CheckRcvStats(handler.recv.stats, handler.recv.size, handler.mpi_tag)) {
503 amrex::Abort(
"NonLocalBC::Comm_finish failed with wrong message size");
508 if (Gpu::inLaunchRegion())
511 handler.recv.size, handler.recv.cctc, dtos, proj);
516 handler.recv.size, handler.recv.cctc, dtos, proj);
520 if ( ! cmd.m_SndTags->empty() ) {
521 handler.send.stats.resize(handler.send.request.size());
522 ParallelDescriptor::Waitall(handler.send.request, handler.send.stats);
528std::enable_if_t<IsBaseFab<FAB>::value>
529Rotate90 (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
537 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
538 AMREX_ASSERT(nghost.allLE(mf.nGrowVect()) && nghost[0] == nghost[1]);
540 if (nghost[0] <= 0) {
return; }
542 const FabArrayBase::RB90& TheRB90 = mf.getRB90(nghost, domain);
544 auto handler = Comm_nowait(mf, scomp, ncomp, TheRB90,Rotate90DstToSrc{},
547 Box corner(-nghost, IntVect{
AMREX_D_DECL(-1,-1,domain.bigEnd(2)+nghost[2])});
549#pragma omp parallel if (Gpu::notInLaunchRegion())
551 for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
552 Box const& bx = corner & mfi.fabbox();
554 auto const& fab = mf.array(mfi);
557 fab(i,j,k,n) = fab(-i-1,-j-1,k,n);
563 Comm_finish(mf, scomp, ncomp, TheRB90, std::move(handler), Rotate90DstToSrc{},
571std::enable_if_t<IsBaseFab<FAB>::value>
572Rotate90 (FabArray<FAB>& mf, Box
const& domain)
574 Rotate90(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
578std::enable_if_t<IsBaseFab<FAB>::value>
579Rotate180 (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
587 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
590 if (nghost[0] <= 0) {
return; }
592 const FabArrayBase::RB180& TheRB180 = mf.getRB180(nghost, domain);
594 auto handler = Comm_nowait(mf, scomp, ncomp, TheRB180,
595 Rotate180Fn{domain.length(1)}, Identity{});
598 Comm_finish(mf, scomp, ncomp, TheRB180, std::move(handler),
599 Rotate180Fn{domain.length(1)}, Identity{});
606std::enable_if_t<IsBaseFab<FAB>::value>
607Rotate180 (FabArray<FAB>& mf, Box
const& domain)
609 Rotate180(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
613std::enable_if_t<IsBaseFab<FAB>::value>
614FillPolar (FabArray<FAB>& mf,
int scomp,
int ncomp, IntVect
const& nghost, Box
const& domain)
622 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
625 if (nghost[0] <= 0) {
return; }
627 const FabArrayBase::PolarB& ThePolarB = mf.getPolarB(nghost, domain);
629 auto handler = Comm_nowait(mf, scomp, ncomp, ThePolarB,
630 PolarFn{.Lx = domain.length(0), .Ly = domain.length(1)},
634 Comm_finish(mf, scomp, ncomp, ThePolarB, std::move(handler),
635 PolarFn{.Lx = domain.length(0), .Ly = domain.length(1)}, Identity{});
642std::enable_if_t<IsBaseFab<FAB>::value>
643FillPolar (FabArray<FAB>& mf, Box
const& domain)
645 FillPolar(mf, 0, mf.nComp(), mf.nGrowVect(), domain);
648template <
typename FAB,
typename DTOS,
typename Proj>
649std::enable_if_t<IsBaseFab<FAB>() &&
650 IsCallableR<Dim3,DTOS,Dim3>() &&
651 IsFabProjection<Proj,FAB>(),
654 int scomp,
int ncomp, DTOS
const& dtos, Proj
const& proj)
657 AMREX_ASSERT(scomp < mf.nComp() && scomp+ncomp <= mf.nComp());
658 return Comm_nowait(mf, scomp, ncomp, cmd, dtos, proj);
661template <
typename FAB,
typename DTOS,
typename Proj>
662std::enable_if_t<IsBaseFab<FAB>() &&
663 IsCallableR<Dim3,DTOS,Dim3>() &&
664 IsFabProjection<Proj,FAB>()>
666 FabArray<FAB>& mf,
const FabArrayBase::CommMetaData& cmd,
667 int scomp,
int ncomp, DTOS
const& dtos, Proj
const& proj)
671 Comm_finish(mf, scomp, ncomp, cmd, std::move(handler), dtos, proj);
677template <
typename DTOS>
678Vector<std::pair<Box,Box>>
679get_src_dst_boxes (DTOS
const& dtos, Box
const& dstbox, Box
const& domain)
681 Vector<std::pair<Box,Box>> r;
684 if (!domain.contains(mapped_smallend) || !domain.contains(mapped_bigend)) {
690 auto dtype = dstbox.type();
694 Array<Array<std::pair<int,int>,2>,AMREX_SPACEDIM> ends;
695 Array<Array<std::pair<int,int>,2>,AMREX_SPACEDIM> dst_ends;
696 Array<int,AMREX_SPACEDIM> nboxes;
697 for (
int ddim = 0; ddim < AMREX_SPACEDIM; ++ddim) {
698 int sdim = perm[ddim];
699 auto mm = std::minmax(mapped_smallend[sdim],mapped_bigend[sdim]);
700 if (((sign[ddim] > 0) && (mapped_smallend[sdim] <= mapped_bigend[sdim])) ||
701 ((sign[ddim] < 0) && (mapped_bigend[sdim] <= mapped_smallend[sdim])))
705 dst_ends[ddim][0] = std::make_pair(dstbox.smallEnd(ddim),
706 dstbox.bigEnd(ddim));
709 ends[sdim][0].first = domain.smallEnd(sdim);
710 ends[sdim][0].second = mm.first;
711 ends[sdim][1].first = mm.second;
712 ends[sdim][1].second = domain.bigEnd(sdim);
713 int n0 = ends[sdim][0].second - ends[sdim][0].first;
714 int n1 = ends[sdim][1].second - ends[sdim][1].first;
715 if (mm.first == mapped_smallend[sdim]) {
716 dst_ends[ddim][0] = std::make_pair(dstbox.smallEnd(ddim),
717 dstbox.smallEnd(ddim)+n0);
718 dst_ends[ddim][1] = std::make_pair(dstbox.bigEnd(ddim)-n1,
719 dstbox.bigEnd(ddim));
721 dst_ends[ddim][0] = std::make_pair(dstbox.bigEnd(ddim)-n0,
722 dstbox.bigEnd(ddim));
723 dst_ends[ddim][1] = std::make_pair(dstbox.smallEnd(ddim),
724 dstbox.smallEnd(ddim)+n1);
729 r.reserve(
AMREX_D_TERM(nboxes[0],*nboxes[1],*nboxes[2]));
731#if (AMREX_SPACEDIM == 3)
732 for (
int kbox = 0; kbox < nboxes[2]; ++kbox) {
734#if (AMREX_SPACEDIM >=2 )
735 for (
int jbox = 0; jbox < nboxes[1]; ++jbox) {
737 for (
int ibox = 0; ibox < nboxes[0]; ++ibox)
743 ends[2][kbox].first)),
745 ends[1][jbox].second,
746 ends[2][kbox].second)),
749 dst_ends[1][div[1]].first,
750 dst_ends[2][div[2]].first)),
752 dst_ends[1][div[1]].second,
753 dst_ends[2][div[2]].second)),
760template <
typename DTOS>
761Box get_dst_subbox (DTOS
const& dtos, std::pair<Box,Box>
const& sdboxes,
762 Box
const& srcsubbox)
764 Box const& srcbox = sdboxes.first;
765 Box const& dstbox = sdboxes.second;
766 if (srcbox == srcsubbox) {
771 Box dstsubbox = dstbox;
772 for (
int ddim = 0; ddim < AMREX_SPACEDIM; ++ddim) {
773 int sdim = perm[ddim];
774 if (sign[ddim] > 0) {
775 dstsubbox.
growLo(ddim, srcbox.smallEnd(sdim)-srcsubbox.smallEnd(sdim));
776 dstsubbox.
growHi(ddim, srcsubbox.bigEnd(sdim)-srcbox.bigEnd(sdim));
778 dstsubbox.
growLo(ddim, srcsubbox.bigEnd(sdim)-srcbox.bigEnd(sdim));
779 dstsubbox.
growHi(ddim, srcbox.smallEnd(sdim)-srcsubbox.smallEnd(sdim));
788 void split_boxes (BoxList& bl, Box
const& domain);
792template <
typename FAB,
typename DTOS>
793std::enable_if_t<IsBaseFab<FAB>() && IsCallableR<Dim3,DTOS,Dim3>(),
794 FabArrayBase::CommMetaData>
796 Geometry
const& geom, DTOS
const& dtos)
798 FabArrayBase::CommMetaData cmd;
799 cmd.m_LocTags = std::make_unique<FabArrayBase::CopyComTagsContainer>();
800 cmd.m_SndTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
801 cmd.m_RcvTags = std::make_unique<FabArrayBase::MapOfCopyComTagContainers>();
804 mf.define_fb_metadata(cmd, nghost,
false, geom.periodicity(),
false);
806 BoxArray
const& ba = mf.boxArray();
807 DistributionMapping
const& dm = mf.DistributionMap();
811 const int myproc = ParallelDescriptor::MyProc();
812 const auto nboxes =
static_cast<int>(ba.size());
813 std::vector<std::pair<int,Box> > isects;
815 for (
int i = 0; i < nboxes; ++i) {
818 if (bl.isEmpty()) {
continue; }
820 detail::split_boxes(bl, dombox);
822 const int dst_owner = dm[i];
823 for (
auto const& dst_box : bl) {
824 auto const& src_dst_boxes = get_src_dst_boxes(dtos, dst_box, dombox);
825 for (
auto const& sd_box_pair : src_dst_boxes) {
826 ba.intersections(sd_box_pair.first, isects);
827 for (
auto const& is : isects) {
828 int const k = is.first;
829 Box const src_b = is.second;
830 int const src_owner = dm[k];
831 if (dst_owner == myproc || src_owner == myproc) {
832 Box const& dst_b = get_dst_subbox(dtos, sd_box_pair, src_b);
833 if (src_owner == dst_owner) {
834 cmd.m_LocTags->emplace_back(dst_b, src_b, i, k);
836 auto& tags = (dst_owner == myproc) ?
837 (*cmd.m_RcvTags)[src_owner] :
838 (*cmd.m_SndTags)[dst_owner];
839 tags.emplace_back(dst_b, src_b, i, k);
850struct SphThetaPhiRIndexMapping
852 SphThetaPhiRIndexMapping (Box
const& a_domain)
861 Dim3 operator() (Dim3
const& ijk)
const noexcept
868 bool imd = i >= 0 && i < nx;
871 bool jmd = j >= 0 && j < ny;
873 bool kmd = k >= 0 && k < nz;
876 if (ilo && jmd && kmd)
878 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
880 else if (ihi && jmd && kmd)
882 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
884 else if (imd && jlo && kmd)
886 return Dim3{.x = i, .y = j+ny, .z = k};
888 else if (imd && jhi && kmd)
890 return Dim3{.x = i, .y = j-ny, .z = k};
892 else if (imd && jmd && klo)
894 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
896 else if (ilo && jlo && kmd)
898 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
900 else if (ihi && jlo && kmd)
902 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
904 else if (ilo && jhi && kmd)
906 return Dim3{.x = -1-i, .y = (j+ny/2)%ny, .z = k};
908 else if (ihi && jhi && kmd)
910 return Dim3{.x = 2*nx-1-i, .y = (j+ny/2)%ny, .z = k};
912 else if (imd && jlo && klo)
914 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
916 else if (imd && jhi && klo)
918 return Dim3{.x = nx-1-i, .y = (j+ny/2)%ny, .z = -1-k};
926 [[nodiscard]]
IntVect sign (Dim3
const& ijk)
const noexcept
930 }
else if (ijk.z >=0 && ijk.z < nz &&
931 (ijk.x < 0 || ijk.x >= nx)) {
938 [[nodiscard]]
IntVect permutation (Dim3
const& )
const noexcept
947struct SphThetaPhiRComponentMapping
949 SphThetaPhiRComponentMapping (Box
const& a_domain,
int a_start_index)
953 scomp(a_start_index) {}
955 template <
typename T>
957 T operator()(Array4<const T>
const& a, Dim3
const& ijk,
int n)
const noexcept
964 if ((i >= 0 && i < nx) &&
965 (j < 0 || j >= ny) &&
966 (k >= 0 && k < nz)) {
973 }
else if (n == scomp+2) {
988extern template MultiBlockCommMetaData
ParallelCopy(FabArray<FArrayBox>& dest,
const Box& destbox,
989 const FabArray<FArrayBox>& src,
int destcomp,
990 int srccomp,
int numcomp,
const IntVect& ngrow,
991 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__ 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:662
__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:673
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:1331
__host__ __device__ Dim3 length(Array4< T > const &a) noexcept
Return the spatial extents of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1345
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Return the inclusive lower bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1317
__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
Arena * The_Arena()
Definition AMReX_Arena.cpp:820
Definition AMReX_NonLocalBC.cpp:39
std::enable_if_t< IsBaseFab< FAB >() &&IsCallableR< Dim3, DTOS, Dim3 >() &&IsFabProjection< Proj, FAB >()> 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
std::enable_if_t< IsBaseFab< FAB >::value > Rotate90(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
std::enable_if_t< IsBaseFab< FAB >::value > Rotate180(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
std::enable_if_t< IsCallableR< Dim3, DTOS, Dim3 >::value &&!IsCallableR< IndexType, DTOS, IndexType >::value, 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:120
std::enable_if_t< IsBaseFab< FAB >() &&IsCallableR< Dim3, DTOS, Dim3 >(), FabArrayBase::CommMetaData > makeFillBoundaryMetaData(FabArray< FAB > &mf, IntVect const &nghost, Geometry const &geom, DTOS const &dtos)
Make metadata for FillBoundary.
std::enable_if_t< IsBaseFab< FAB >() &&IsCallableR< Dim3, DTOS, Dim3 >() &&IsFabProjection< Proj, FAB >()> 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
std::enable_if_t< IsBaseFab< FAB >::value > FillPolar(FabArray< FAB > &mf, int scomp, int ncomp, IntVect const &nghost, Box const &domain)
std::enable_if_t< HasInverseMemFn< DTOS >::value &&!IsCallableR< IndexType, DTOS, IndexType >::value, 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:183
std::enable_if_t< IsBaseFab< FAB >() &&IsCallableR< Dim3, DTOS, Dim3 >() &&IsFabProjection< Proj, FAB >()> 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
std::enable_if_t< IsBaseFab< FAB >() &&IsCallableR< Dim3, DTOS, Dim3 >() &&IsFabProjection< Proj, FAB >()> 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:696
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:193
__host__ __device__ Array4< T > makeArray4(T *p, Box const &bx, int ncomp) noexcept
Definition AMReX_BaseFab.H:92
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2867
std::enable_if_t< IsFabArray< MF >::value > FillBoundary_finish(Vector< MF * > const &mf)
Wait for outstanding FillBoundary_nowait operations launched with the vector helper to complete.
Definition AMReX_FabArrayCommI.H:1113
std::enable_if_t< IsFabArray< MF >::value > FillBoundary_nowait(Vector< MF * > const &mf, Vector< int > const &scomp, Vector< int > const &ncomp, Vector< IntVect > const &nghost, Vector< Periodicity > const &period, Vector< int > const &cross={})
Launch FillBoundary_nowait across a vector of FabArrays.
Definition AMReX_FabArrayCommI.H:1067
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
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
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:2019
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2862