12 #ifndef CD_ParticleContainerImplem_H
13 #define CD_ParticleContainerImplem_H
19 #include <ParmParse.H>
22 #include <ParticleDataI.H>
31 #include <CD_NamespaceHeader.H>
37 m_isOrganizedByCell =
false;
45 const Vector<ProblemDomain>& a_domains,
46 const Vector<Real>& a_dx,
47 const Vector<int>& a_refRat,
48 const Vector<ValidMask>& a_validMask,
49 const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
50 const RealVect& a_probLo,
51 const int a_blockingFactor,
52 const int a_finestLevel,
53 const std::string a_realm)
55 CH_TIME(
"ParticleContainer::ParticleContainer");
72 CH_TIME(
"ParticleContainer::~ParticleContainer");
78 const Vector<ProblemDomain>& a_domains,
79 const Vector<Real>& a_dx,
80 const Vector<int>& a_refRat,
81 const Vector<ValidMask>& a_validMask,
82 const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
83 const RealVect& a_probLo,
84 const int a_blockingFactor,
85 const int a_finestLevel,
86 const std::string a_realm)
88 CH_TIME(
"ParticleContainer::~ParticleContainer");
91 m_domains = a_domains;
93 m_validRegion = a_validMask;
95 m_finestLevel = a_finestLevel;
97 m_blockingFactor = a_blockingFactor;
98 m_levelTiles = a_levelTiles;
100 m_dx.resize(1 + m_finestLevel);
101 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
102 m_dx[lvl] = a_dx[lvl] * RealVect::Unit;
105 constexpr
int base = 0;
108 this->setupGrownGrids(base, m_finestLevel);
109 this->setupParticleData(base, m_finestLevel);
112 m_isOrganizedByCell =
false;
115 ParmParse pp(
"ParticleContainer");
116 pp.query(
"profile", m_profile);
117 pp.query(
"debug", m_debug);
118 pp.query(
"verbose", m_verbose);
125 CH_TIME(
"ParticleContainer::setupGrownGrids");
127 pout() <<
"ParticleContainer::setupGrownGrids" << endl;
133 m_grownGrids.resize(1 + a_finestLevel);
135 for (
int lvl = 0; lvl <= a_finestLevel; lvl++) {
136 const DisjointBoxLayout& dbl = m_grids[lvl];
137 const ProblemDomain& domain = m_domains[lvl];
140 Vector<Box> boxes = dbl.boxArray();
142 for (
auto& box : boxes.stdVector()) {
143 box.grow(m_refRat[lvl - 1]);
148 m_grownGrids[lvl] = BoxLayout(boxes, dbl.procIDs());
156 CH_TIME(
"ParticleContainer::setupParticleData");
158 pout() <<
"ParticleContainer::setupParticleData" << endl;
176 m_particles.resize(1 + a_finestLevel);
177 m_bufferParticles.resize(1 + a_finestLevel);
178 m_maskParticles.resize(1 + a_finestLevel);
179 m_cellSortedParticles.resize(1 + a_finestLevel);
181 for (
int lvl = a_base; lvl <= a_finestLevel; lvl++) {
182 m_particles[lvl] = RefCountedPtr<ParticleData<P>>(
183 new ParticleData<P>(m_grids[lvl], m_domains[lvl], m_blockingFactor, m_dx[lvl], m_probLo));
185 m_bufferParticles[lvl] = RefCountedPtr<ParticleData<P>>(
186 new ParticleData<P>(m_grownGrids[lvl], m_domains[lvl], m_blockingFactor, m_dx[lvl], m_probLo));
188 m_maskParticles[lvl] = RefCountedPtr<ParticleData<P>>(
189 new ParticleData<P>(m_grids[lvl], m_domains[lvl], m_blockingFactor, m_dx[lvl], m_probLo));
191 m_cellSortedParticles[lvl] = RefCountedPtr<LayoutData<BinFab<P>>>(
new LayoutData<BinFab<P>>(m_grids[lvl]));
199 CH_TIME(
"ParticleContainer::sortParticles");
201 pout() <<
"ParticleContainer::sortParticles" << endl;
204 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
206 const DisjointBoxLayout& dbl = m_grids[lvl];
207 const DataIterator& dit = dbl.dataIterator();
209 const int nbox = dit.size();
211 #pragma omp parallel for schedule(runtime)
212 for (
int mybox = 0; mybox < nbox; mybox++) {
213 const DataIndex& din = dit[mybox];
215 List<P>& particles = (*m_particles[lvl])[din].listItems();
226 return m_isOrganizedByCell;
233 CH_assert(m_isDefined);
235 return m_finestLevel;
242 CH_assert(m_isDefined);
248 const Vector<DisjointBoxLayout>&
262 const Vector<RealVect>
272 CH_assert(m_isDefined);
274 if (m_isOrganizedByCell) {
275 MayDay::Error(
"ParticleContainer::getParticles - particles are sorted by cell!");
285 CH_assert(m_isDefined);
287 if (m_isOrganizedByCell) {
288 MayDay::Abort(
"ParticleContainer::getParticles - particles are sorted by cell!");
298 CH_assert(m_isDefined);
300 return m_bufferParticles;
307 CH_assert(m_isDefined);
309 return m_bufferParticles;
316 CH_assert(m_isDefined);
318 return m_maskParticles;
325 CH_assert(m_isDefined);
327 return m_maskParticles;
334 CH_assert(m_isDefined);
336 if (m_isOrganizedByCell) {
337 MayDay::Error(
"ParticleContainer::operator[](a_lvl) - particles are sorted by cell!");
340 return *m_particles[a_lvl];
344 const ParticleData<P>&
347 CH_assert(m_isDefined);
349 if (m_isOrganizedByCell) {
350 MayDay::Error(
"ParticleContainer::operator[](a_lvl) - particles are sorted by cell!");
353 return *m_particles[a_level];
360 CH_assert(m_isDefined);
362 if (!m_isOrganizedByCell) {
363 MayDay::Error(
"ParticleContainer::getCellParticles()- particles are not sorted by cell!");
366 return m_cellSortedParticles;
373 CH_assert(m_isDefined);
375 if (!m_isOrganizedByCell) {
376 MayDay::Error(
"ParticleContainer::getCellParticles()- particles are not sorted by cell!");
379 return m_cellSortedParticles;
383 LayoutData<BinFab<P>>&
386 CH_assert(m_isDefined);
388 if (!m_isOrganizedByCell) {
389 MayDay::Error(
"ParticleContainer::getCellParticles(level)- particles are not sorted by cell!");
392 return *m_cellSortedParticles[a_level];
396 const LayoutData<BinFab<P>>&
399 CH_TIME(
"ParticleContainer::getCellParticles(int)");
401 CH_assert(m_isDefined);
403 if (!m_isOrganizedByCell) {
404 MayDay::Error(
"ParticleContainer::getCellParticles(level)- particles are not sorted by cell!");
407 return *m_cellSortedParticles[a_level];
414 CH_TIME(
"ParticleContainer::getCellParticles(BinFab)");
416 CH_assert(m_isDefined);
418 if (!m_isOrganizedByCell) {
419 MayDay::Error(
"ParticleContainer::getCellParticles - particles are not sorted by cell!");
422 cellParticles.define(m_grids[a_lvl][a_dit], m_dx[a_lvl], m_probLo);
423 cellParticles.addItems((*m_particles[a_lvl])[a_dit].listItems());
430 CH_TIME(
"ParticleContainer::getCellParticlesDestructive");
432 CH_assert(m_isDefined);
434 if (!m_isOrganizedByCell) {
435 MayDay::Error(
"ParticleContainer::getCellParticlesDestructive - particles are not sorted by cell!");
438 cellParticles.define(m_grids[a_lvl].get(a_dit), m_dx[a_lvl], m_probLo);
439 cellParticles.addItemsDestructive((*m_particles[a_lvl])[a_dit].listItems());
446 CH_TIME(
"ParticleContainer::getCellParticles(int, DataIndex)");
448 CH_assert(m_isDefined);
450 if (!m_isOrganizedByCell) {
451 MayDay::Error(
"ParticleContainer::getCellParticles(int, dit) - particles are not sorted by cell!");
454 return (*m_cellSortedParticles[a_level])[a_dit];
461 CH_TIME(
"ParticleContainer::getCellParticles(int, DataIndex)");
463 CH_assert(m_isDefined);
465 if (!m_isOrganizedByCell) {
466 MayDay::Error(
"ParticleContainer::getCellParticles(int, dit) - particles are not sorted by cell!");
469 return (*m_cellSortedParticles[a_level])[a_dit];
476 CH_TIME(
"ParticleContainer::organizeParticlesByCell");
478 pout() <<
"ParticleContainer::organizeParticlesByCell" << endl;
481 CH_assert(m_isDefined);
483 if (!m_isOrganizedByCell) {
485 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
486 const DisjointBoxLayout& dbl = m_grids[lvl];
487 const DataIterator& dit = dbl.dataIterator();
489 const int nbox = dit.size();
491 #pragma omp parallel for schedule(runtime)
492 for (
int mybox = 0; mybox < nbox; mybox++) {
493 const DataIndex& din = dit[mybox];
495 BinFab<P>& cellParticles = (*m_cellSortedParticles[lvl])[din];
497 cellParticles.define(dbl[din], m_dx[lvl], m_probLo);
498 cellParticles.addItemsDestructive((*m_particles[lvl])[din].listItems());
502 m_isOrganizedByCell =
true;
510 CH_TIME(
"ParticleContainer::organizeParticlesByPatch");
512 pout() <<
"ParticleContainer::organizeParticlesByPatch" << endl;
515 CH_assert(m_isDefined);
517 if (m_isOrganizedByCell) {
519 constexpr
int comp = 0;
521 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
522 const DisjointBoxLayout& dbl = m_grids[lvl];
523 const DataIterator& dit = dbl.dataIterator();
525 const int nbox = dit.size();
527 #pragma omp parallel for schedule(runtime)
528 for (
int mybox = 0; mybox < nbox; mybox++) {
529 const DataIndex& din = dit[mybox];
531 ListBox<P>& patchParticles = (*m_particles[lvl])[din];
532 BinFab<P>& cellParticles = (*m_cellSortedParticles[lvl])[din];
535 auto kernel = [&](
const IntVect& iv) ->
void {
536 patchParticles.addItemsDestructive(cellParticles(iv, comp));
541 cellParticles.clear();
545 m_isOrganizedByCell =
false;
553 CH_TIME(
"ParticleContainer::addParticles");
555 pout() <<
"ParticleContainer::addParticles" << endl;
558 CH_assert(m_isDefined);
560 if (m_isOrganizedByCell) {
561 MayDay::Error(
"ParticleContainer::addParticles(List<P>) - particles are sorted by cell!");
567 const unsigned int numRanks = numProc();
568 const unsigned int myRank = procID();
570 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
572 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
580 outcasts.join(a_particles);
583 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
586 this->transferParticlesToSingleList(outcasts, m_particles);
589 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
593 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
597 this->assignLocalParticles(particlesToSend[myRank], m_particles);
601 std::map<LevelAndIndex, List<P>> receivedParticles;
603 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
607 this->assignLocalParticles(receivedParticles, m_particles);
619 CH_TIME(
"ParticleContainer::addParticlesDestructive");
621 pout() <<
"ParticleContainer::addParticlesDestructive" << endl;
624 CH_assert(m_isDefined);
626 if (m_isOrganizedByCell) {
627 MayDay::Error(
"ParticleContainer::addParticlesDestructive(List<P>) - particles are sorted by cell!");
633 const unsigned int numRanks = numProc();
634 const unsigned int myRank = procID();
636 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
638 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
646 outcasts.catenate(a_particles);
649 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
652 this->transferParticlesToSingleList(outcasts, m_particles);
655 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
659 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
663 this->assignLocalParticles(particlesToSend[myRank], m_particles);
667 std::map<LevelAndIndex, List<P>> receivedParticles;
669 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
673 this->assignLocalParticles(receivedParticles, m_particles);
685 CH_TIME(
"ParticleContainer::addParticles(BinFab)");
687 pout() <<
"ParticleContainer::addParticles(BinFab)" << endl;
690 CH_assert(m_isDefined);
691 CH_assert(m_grids[a_lvl].get(a_dit) == a_particles.getRegion());
693 if (!m_isOrganizedByCell) {
694 MayDay::Abort(
"ParticleContainer::addParticles(BinFab<P>) - particles are not sorted by cell!");
697 constexpr
int comp = 0;
699 BinFab<P>& boxParticles = (*m_cellSortedParticles[a_lvl])[a_dit];
701 auto kernel = [&](
const IntVect& iv) ->
void {
702 List<P>& myParticles = boxParticles(iv, comp);
703 const List<P>& inputParticles = a_particles(iv, comp);
705 myParticles.join(inputParticles);
715 CH_TIME(
"ParticleContainer::addParticlesDestructive(BinFab)");
717 pout() <<
"ParticleContainer::addParticlesDestructive(BinFab)" << endl;
720 CH_assert(m_isDefined);
721 CH_assert(m_grids[a_lvl].get(a_dit) == a_particles.getRegion());
723 if (!m_isOrganizedByCell) {
724 MayDay::Error(
"ParticleContainer::addParticles(BinFab<P>) - particles are not sorted by cell!");
727 constexpr
int comp = 0;
729 BinFab<P>& boxParticles = (*m_cellSortedParticles[a_lvl])[a_dit];
731 auto kernel = [&](
const IntVect& iv) ->
void {
732 List<P>& myParticles = boxParticles(iv, comp);
733 const List<P>& inputParticles = a_particles(iv, comp);
735 myParticles.catenate(inputParticles);
745 CH_TIME(
"ParticleContainer::addParticles(ParticleContainer)");
747 pout() <<
"ParticleContainer::addParticles(ParticleContainer)" << endl;
750 CH_assert(m_isDefined);
752 if (m_isOrganizedByCell) {
753 MayDay::Error(
"ParticleContainer::addParticles(ParticleContainer<P>) - particles are sorted by cell!");
759 if (m_realm == a_otherContainer.
getRealm()) {
760 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
761 const DisjointBoxLayout& dbl = m_grids[lvl];
762 const DataIterator& dit = dbl.dataIterator();
764 const int nbox = dit.size();
766 #pragma omp parallel for schedule(runtime)
767 for (
int mybox = 0; mybox < nbox; mybox++) {
768 const DataIndex& din = dit[mybox];
770 ListBox<P>& myParticles = (*m_particles[lvl])[din];
771 const ListBox<P>& otherParticles = a_otherContainer[lvl][din];
773 myParticles.addItems(otherParticles.listItems());
781 const unsigned int numRanks = numProc();
782 const unsigned int myRank = procID();
784 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
786 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
792 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
798 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
802 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
806 this->assignLocalParticles(particlesToSend[myRank], m_particles);
810 std::map<LevelAndIndex, List<P>> receivedParticles;
812 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
816 this->assignLocalParticles(receivedParticles, m_particles);
829 CH_TIME(
"ParticleContainer::addParticlesDestructive(ParticleContainer)");
831 pout() <<
"ParticleContainer::addParticlesDestructive(ParticleContainer)" << endl;
837 CH_assert(m_isDefined);
839 if (m_isOrganizedByCell) {
840 MayDay::Error(
"ParticleContainer::addParticles(ParticleContainer<P>) - particles are sorted by cell!");
846 if (m_realm == a_otherContainer.
getRealm()) {
847 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
848 const DisjointBoxLayout& dbl = m_grids[lvl];
849 const DataIterator& dit = dbl.dataIterator();
851 const int nbox = dit.size();
853 #pragma omp parallel for schedule(runtime)
854 for (
int mybox = 0; mybox < nbox; mybox++) {
855 const DataIndex& din = dit[mybox];
857 ListBox<P>& myParticles = (*m_particles[lvl])[din];
858 const ListBox<P>& otherParticles = a_otherContainer[lvl][din];
860 myParticles.addItems(otherParticles.listItems());
868 const unsigned int numRanks = numProc();
869 const unsigned int myRank = procID();
871 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
873 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
879 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
885 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
889 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
893 this->assignLocalParticles(particlesToSend[myRank], m_particles);
897 std::map<LevelAndIndex, List<P>> receivedParticles;
899 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
903 this->assignLocalParticles(receivedParticles, m_particles);
916 CH_TIME(
"ParticleContainer::transferParticles");
918 CH_assert(m_isDefined);
919 CH_assert(this->getRealm() == a_otherContainer.
getRealm());
921 if (m_isOrganizedByCell) {
922 MayDay::Error(
"ParticleContainer::transferParticles(ParticleContainer<P>) - particles are sorted by cell!");
926 MayDay::Error(
"ParticleContainer::transferParticles(ParticleContainer<P>) - other particles are sorted by cell!");
929 if (a_otherContainer.
getRealm() != m_realm) {
931 "ParticleContainer::transferParticles(ParticleContainer<P>) - other container defined over a different realm");
934 this->transferParticles(a_otherContainer.
getParticles());
941 CH_TIME(
"ParticleContainer::transferParticles(AMRParticles)");
943 pout() <<
"ParticleContainer::transferParticle(AMRParticles)" << endl;
946 CH_assert(m_isDefined);
948 if (m_isOrganizedByCell) {
949 MayDay::Error(
"ParticleContainer::transferParticles(ParticleContainer<P>) - particles are sorted by cell!");
953 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
954 const DisjointBoxLayout& dbl = m_grids[lvl];
955 const DataIterator& dit = dbl.dataIterator();
957 const int nbox = dit.size();
959 #pragma omp parallel for schedule(runtime)
960 for (
int mybox = 0; mybox < nbox; mybox++) {
961 const DataIndex& din = dit[mybox];
963 List<P>& myParticles = (*m_particles[lvl])[din].listItems();
964 List<P>& otherParticles = (*a_particles[lvl])[din].listItems();
966 myParticles.catenate(otherParticles);
975 CH_TIME(
"ParticleContainer<P>::remap");
977 pout() <<
"ParticleContainer::remap" << endl;
990 const unsigned int numRanks = numProc();
991 const unsigned int myRank = procID();
993 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
998 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
1000 #pragma omp parallel
1004 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
1007 this->transferParticlesToSingleList(outcasts, m_particles);
1010 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
1013 #pragma omp critical
1014 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
1018 this->assignLocalParticles(particlesToSend[myRank], m_particles);
1022 std::map<LevelAndIndex, List<P>> receivedParticles;
1024 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
1028 this->assignLocalParticles(receivedParticles, m_particles);
1032 this->sanityCheck();
1036 template <
typename P>
1040 CH_TIME(
"ParticleContainer::sanityCheck");
1042 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1043 const DisjointBoxLayout& dbl = m_grids[lvl];
1044 const RealVect dx = m_dx[lvl];
1046 for (DataIterator dit(dbl); dit.ok(); ++dit) {
1047 const List<P>& particles = (*m_particles[lvl])[dit()].listItems();
1048 const Box box = dbl[dit()];
1050 for (ListIterator<P> lit(particles); lit.ok(); ++lit) {
1051 const RealVect particlePosition = lit().position();
1054 if (!(box.contains(particleCell))) {
1055 MayDay::Error(
"ParticleContainer::remap -- error2: particle is not contained in grid box!");
1062 template <
typename P>
1066 CH_TIME(
"ParticleContainer::transferParticlesToSingleList");
1068 for (
int lvl = 0; lvl < a_particles.size(); lvl++) {
1069 const BoxLayout& dbl = a_particles[lvl]->getBoxes();
1070 const DataIterator& dit = dbl.dataIterator();
1072 const int nbox = dit.size();
1074 #pragma omp for schedule(runtime)
1075 for (
int mybox = 0; mybox < nbox; mybox++) {
1076 const DataIndex& din = dit[mybox];
1078 List<P>& particles = (*a_particles[lvl])[din].listItems();
1080 a_list.catenate(particles);
1085 template <
typename P>
1089 CH_TIME(
"ParticleContainer::copyParticlesToSingleList");
1091 for (
int lvl = 0; lvl < a_particles.size(); lvl++) {
1092 const BoxLayout& dbl = a_particles[lvl]->getBoxes();
1093 const DataIterator& dit = dbl.dataIterator();
1095 const int nbox = dit.size();
1097 #pragma omp for schedule(runtime)
1098 for (
int mybox = 0; mybox < nbox; mybox++) {
1099 const DataIndex& din = dit[mybox];
1101 const List<P>& particles = (*a_particles[lvl])[din].listItems();
1103 a_list.join(particles);
1108 template <
typename P>
1111 std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_mappedParticles,
1112 List<P>& a_unmappedParticles)
const noexcept
1114 CH_TIME(
"ParticleContainer::mapParticlesToAMRGrid");
1116 const unsigned int numRanks = numProc();
1117 const unsigned int myRank = procID();
1119 std::vector<RealVect> quasiDx(1 + m_finestLevel);
1121 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1122 quasiDx[lvl] = m_blockingFactor * m_dx[lvl];
1127 for (ListIterator<P> lit(a_unmappedParticles); lit.ok();) {
1128 bool foundTile =
false;
1130 for (
int lvl = m_finestLevel; lvl >= 0 && !foundTile; lvl--) {
1131 const IntVect particleTile = locateBin(lit().
position(), quasiDx[lvl], m_probLo);
1133 const auto& myLevelTiles = m_levelTiles[lvl]->getMyTiles();
1134 const auto& otherLevelTiles = m_levelTiles[lvl]->getOtherTiles();
1136 if (myLevelTiles.find(particleTile) != myLevelTiles.end()) {
1138 const unsigned int gridIndex = myLevelTiles.at(particleTile);
1139 const unsigned int toRank = myRank;
1141 a_mappedParticles[toRank][std::pair<unsigned int, unsigned int>(lvl, gridIndex)].transfer(lit);
1146 else if (otherLevelTiles.find(particleTile) != otherLevelTiles.end()) {
1148 const unsigned int gridIndex = otherLevelTiles.at(particleTile).first;
1149 const unsigned int& toRank = otherLevelTiles.at(particleTile).second;
1151 a_mappedParticles[toRank][std::pair<unsigned int, unsigned int>(lvl, gridIndex)].transfer(lit);
1165 a_unmappedParticles.clear();
1168 template <
typename P>
1171 std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_globalParticles,
1172 std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_localParticles)
const noexcept
1174 CH_TIME(
"ParticleContainer::catenateParticleMaps");
1176 const unsigned int numRanks = numProc();
1178 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
1180 for (
int irank = 0; irank < numRanks; irank++) {
1181 std::map<LevelAndIndex, List<P>>& globalParticles = a_globalParticles[irank];
1183 for (
auto& m : a_localParticles[irank]) {
1184 globalParticles[m.first].catenate(m.second);
1189 template <
typename P>
1194 CH_TIME(
"ParticleContainer::assignLocalParticles");
1196 #pragma omp parallel
1200 for (
auto mapIter = a_mappedParticles.begin(); mapIter != a_mappedParticles.end(); mapIter++) {
1201 #pragma omp task firstprivate(mapIter)
1203 const unsigned int gridLevel = mapIter->first.first;
1204 const unsigned int gridIndex = mapIter->first.second;
1205 const DataIndex din = m_levelTiles[gridLevel]->getMyGrids().at(gridIndex);
1207 List<P>& particles = mapIter->second;
1208 List<P>& myParticles = (*a_particleData[gridLevel])[din].listItems();
1210 myParticles.catenate(particles);
1216 a_mappedParticles.clear();
1223 CH_TIME(
"ParticleContainer::preRegrid");
1225 pout() <<
"ParticleContainer::preRegrid" << endl;
1231 CH_assert(m_isDefined);
1232 if (m_isOrganizedByCell) {
1233 MayDay::Error(
"ParticleContainer::preRegrid - particles are sorted by cell!");
1237 m_cacheParticles.resize(1 + m_finestLevel);
1239 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1240 m_cacheParticles[lvl] = RefCountedPtr<ParticleData<P>>(
1241 new ParticleData<P>(m_grids[lvl], m_domains[lvl], m_blockingFactor, m_dx[lvl], m_probLo));
1243 const DisjointBoxLayout& dbl = m_grids[lvl];
1244 const DataIterator& dit = dbl.dataIterator();
1246 const int nbox = dit.size();
1248 #pragma omp parallel for schedule(runtime)
1249 for (
int mybox = 0; mybox < nbox; mybox++) {
1250 const DataIndex& din = dit[mybox];
1252 List<P>& cacheParticles = (*m_cacheParticles[lvl])[din].listItems();
1253 List<P>& solverParticles = (*m_particles[lvl])[din].listItems();
1255 cacheParticles.catenate(solverParticles);
1263 const Vector<ProblemDomain>& a_domains,
1264 const Vector<Real>& a_dx,
1265 const Vector<int>& a_refRat,
1266 const Vector<ValidMask>& a_validMask,
1267 const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
1269 const int a_newFinestLevel)
1271 CH_TIME(
"ParticleContainer::regrid");
1273 pout() <<
"ParticleContainer::regrid" << endl;
1277 CH_assert(m_isDefined);
1279 if (m_isOrganizedByCell) {
1280 MayDay::Error(
"ParticleContainer::regrid(...) - particles are sorted by cell!");
1285 m_domains = a_domains;
1286 m_refRat = a_refRat;
1287 m_validRegion = a_validMask;
1288 m_finestLevel = a_newFinestLevel;
1289 m_levelTiles = a_levelTiles;
1291 m_dx.resize(1 + m_finestLevel);
1292 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1293 m_dx[lvl] = a_dx[lvl] * RealVect::Unit;
1296 this->setupGrownGrids(a_lmin, m_finestLevel);
1297 this->setupParticleData(a_lmin, m_finestLevel);
1300 const unsigned int numRanks = numProc();
1301 const unsigned int myRank = procID();
1303 using LevelAndIndex = std::pair<unsigned int, unsigned int>;
1308 std::vector<std::map<LevelAndIndex, List<P>>> particlesToSend(numRanks);
1310 #pragma omp parallel
1314 std::vector<std::map<LevelAndIndex, List<P>>> threadLocalParticlesToSend(numRanks);
1317 this->transferParticlesToSingleList(outcasts, m_cacheParticles);
1320 this->mapParticlesToAMRGrid(threadLocalParticlesToSend, outcasts);
1323 #pragma omp critical
1324 this->catenateParticleMaps(particlesToSend, threadLocalParticlesToSend);
1328 this->assignLocalParticles(particlesToSend[myRank], m_particles);
1332 std::map<LevelAndIndex, List<P>> receivedParticles;
1334 ParticleOps::scatterParticles(receivedParticles, particlesToSend);
1338 this->assignLocalParticles(receivedParticles, m_particles);
1342 this->sanityCheck();
1346 m_cacheParticles.resize(0);
1350 template <Real& (P::*particleScalarField)()>
1354 CH_TIME(
"ParticleContainer::setValue");
1356 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1357 const DisjointBoxLayout& dbl = m_grids[lvl];
1358 const DataIterator& dit = dbl.dataIterator();
1360 const int nbox = dit.size();
1362 #pragma omp parallel for schedule(runtime)
1363 for (
int mybox = 0; mybox < nbox; mybox++) {
1364 const DataIndex& din = dit[mybox];
1366 List<P>& patchParticles = (*m_particles[lvl])[din].listItems();
1368 for (ListIterator<P> lit(patchParticles); lit.ok(); ++lit) {
1371 (p.*particleScalarField)() = a_value;
1378 template <RealVect& (P::*particleVectorField)()>
1382 CH_TIME(
"ParticleContainer::setValue");
1384 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1385 const DisjointBoxLayout& dbl = m_grids[lvl];
1386 const DataIterator& dit = dbl.dataIterator();
1388 const int nbox = dit.size();
1390 #pragma omp parallel for schedule(runtime)
1391 for (
int mybox = 0; mybox < nbox; mybox++) {
1392 const DataIndex& din = dit[mybox];
1394 List<P>& patchParticles = (*m_particles[lvl])[din].listItems();
1396 for (ListIterator<P> lit(patchParticles); lit.ok(); ++lit) {
1399 (p.*particleVectorField)() = a_value;
1409 CH_TIME(
"ParticleContainer::getNumberOfValidParticlesLocal");
1411 CH_assert(m_isDefined);
1413 unsigned long long n = 0;
1415 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1416 const DisjointBoxLayout& dbl = m_grids[lvl];
1417 const DataIterator& dit = dbl.dataIterator();
1419 const int nbox = dit.size();
1421 #pragma omp parallel for schedule(runtime) reduction(+ : n)
1422 for (
int mybox = 0; mybox < nbox; mybox++) {
1423 const DataIndex& din = dit[mybox];
1425 const List<P>& particles = (*m_particles[lvl])[din].listItems();
1427 n += (
unsigned long long)particles.length();
1438 CH_TIME(
"ParticleContainer::getNumberOfValidParticlesGlobal");
1440 CH_assert(m_isDefined);
1442 const unsigned long long n = this->getNumberOfValidParticlesLocal();
1451 CH_TIME(
"ParticleContainer::getNumberOfOutcastParticlesLocal");
1453 CH_assert(m_isDefined);
1455 unsigned long long n = 0;
1457 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1458 const DisjointBoxLayout& dbl = m_grids[lvl];
1459 const DataIterator& dit = dbl.dataIterator();
1461 const int nbox = dit.size();
1463 #pragma omp parallel for schedule(runtime) reduction(+ : n)
1464 for (
int mybox = 0; mybox < nbox; mybox++) {
1465 const DataIndex& din = dit[mybox];
1467 const List<P>& outcast = m_particles[lvl]->outcast();
1469 n += (
unsigned long long)outcast.length();
1480 CH_TIME(
"ParticleContainer::getNumberOfOutcastParticlesGlobal");
1482 CH_assert(m_isDefined);
1484 const unsigned long long n = this->getNumberOfOutcastParticlesLocal();
1493 CH_TIME(
"ParticleContainer::getNumberOfMaskParticlesLocal");
1495 CH_assert(m_isDefined);
1497 unsigned long long n = 0;
1499 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1500 const DisjointBoxLayout& dbl = m_grids[lvl];
1501 const DataIterator& dit = dbl.dataIterator();
1503 const int nbox = dit.size();
1505 #pragma omp parallel for schedule(runtime) reduction(+ : n)
1506 for (
int mybox = 0; mybox < nbox; mybox++) {
1507 const DataIndex& din = dit[mybox];
1509 const List<P>& maskParticles = (*m_maskParticles[lvl])[din].listItems();
1511 n += (
unsigned long long)maskParticles.length();
1522 CH_TIME(
"ParticleContainer::getNumberOfMaskParticlesGlobal");
1524 CH_assert(m_isDefined);
1526 const unsigned long long n = this->getNumberOfMaskParticlesLocal();
1535 CH_TIME(
"ParticleContainer::copyMaskParticles(amr)");
1537 CH_assert(m_isDefined);
1539 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1540 if (!a_mask[lvl].isNull()) {
1541 this->copyMaskParticles(lvl, *a_mask[lvl]);
1550 CH_TIME(
"ParticleContainer::copyMaskParticles(level)");
1552 pout() <<
"ParticleContainer::copyMaskParticles(level)" << endl;
1555 CH_assert(m_isDefined);
1557 m_maskParticles[a_level]->clear();
1559 const RealVect dx = m_dx[a_level];
1562 const DisjointBoxLayout& dbl = m_grids[a_level];
1563 const DataIterator& dit = dbl.dataIterator();
1565 const int nbox = dit.size();
1567 #pragma omp parallel for schedule(runtime)
1568 for (
int mybox = 0; mybox < nbox; mybox++) {
1569 const DataIndex& din = dit[mybox];
1571 const BaseFab<bool>& mask = a_mask[din];
1572 const Box gridBox = m_grids[a_level][din];
1573 const Box maskBox = mask.box();
1575 CH_assert(gridBox == maskBox);
1577 List<P>& maskParticles = (*m_maskParticles[a_level])[din].listItems();
1578 const List<P>& particles = (*m_particles[a_level])[din].listItems();
1580 for (ListIterator<P> lit(particles); lit.ok(); ++lit) {
1581 const RealVect x = lit().position();
1582 const IntVect iv = IntVect(D_DECL(std::floor((x[0] - m_probLo[0]) / dx[0]),
1583 std::floor((x[1] - m_probLo[1]) / dx[1]),
1584 std::floor((x[2] - m_probLo[2]) / dx[2])));
1585 if (!(maskBox.contains(iv))) {
1586 std::cout << a_level <<
"\t" << iv <<
"\t" << x <<
"\t" << maskBox << std::endl;
1588 MayDay::Error(
"ParticleContainer::copyMaskParticles -- logic bust. Particle has fallen off grid");
1590 else if (mask(iv)) {
1591 maskParticles.add(lit());
1601 CH_TIME(
"ParticleContainer::transferMaskParticles(amr)");
1603 CH_assert(m_isDefined);
1605 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1606 if (!a_mask[lvl].isNull()) {
1607 this->transferMaskParticles(lvl, *a_mask[lvl]);
1616 CH_TIME(
"ParticleContainer::transferMaskParticles(level)");
1618 pout() <<
"ParticleContainer::transferMaskParticles(level)" << endl;
1621 CH_assert(m_isDefined);
1623 const RealVect dx = m_dx[a_level];
1626 const DisjointBoxLayout& dbl = m_grids[a_level];
1627 const DataIterator& dit = dbl.dataIterator();
1629 const int nbox = dit.size();
1631 #pragma omp parallel for schedule(runtime)
1632 for (
int mybox = 0; mybox < nbox; mybox++) {
1633 const DataIndex& din = dit[mybox];
1635 const BaseFab<bool>& mask = a_mask[din];
1636 const Box gridBox = m_grids[a_level][din];
1637 const Box maskBox = mask.box();
1639 CH_assert(gridBox == maskBox);
1641 List<P>& maskParticles = (*m_maskParticles[a_level])[din].listItems();
1642 const List<P>& particles = (*m_particles[a_level])[din].listItems();
1644 for (ListIterator<P> lit(particles); lit.ok();) {
1645 const RealVect x = lit().position();
1646 const IntVect iv = IntVect(D_DECL(std::floor((x[0] - m_probLo[0]) / dx[0]),
1647 std::floor((x[1] - m_probLo[1]) / dx[1]),
1648 std::floor((x[2] - m_probLo[2]) / dx[2])));
1649 if (!(maskBox.contains(iv))) {
1652 MayDay::Warning(
"ParticleContainer::transferMaskParticles -- logic bust. Particle has fallen off grid");
1658 maskParticles.transfer(lit);
1672 CH_assert(m_isDefined);
1674 this->clear(m_particles);
1681 CH_assert(m_isDefined);
1683 this->clear(m_bufferParticles);
1690 CH_assert(m_isDefined);
1692 this->clear(m_maskParticles);
1699 CH_assert(m_isDefined);
1701 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1702 ParticleData<P>& particleData = *m_particles[lvl];
1704 List<P>& outcast = particleData.outcast();
1715 pout() <<
"ParticleContainer::clear(AMRParticles)" << endl;
1718 CH_assert(m_isDefined);
1720 for (
int lvl = 0; lvl <= m_finestLevel; lvl++) {
1722 ParticleData<P>& levelParticles = *a_particles[lvl];
1725 const BoxLayout& grids = levelParticles.getBoxes();
1726 const DataIterator& dit = grids.dataIterator();
1728 const int nbox = dit.size();
1730 #pragma omp parallel for schedule(runtime)
1731 for (
int mybox = 0; mybox < nbox; mybox++) {
1732 const DataIndex din = dit[mybox];
1734 List<P>& patchParticles = levelParticles[din].listItems();
1736 patchParticles.clear();
1740 List<P>& outcast = levelParticles.outcast();
1745 #include <CD_NamespaceFooter.H>
Declaration of a namespace for proto-typing grid and EB loops.
Declaration of various useful OpenMP-related utilities.
Agglomeration of basic MPI reductions.
Declaration of a class for holding particles on an AMR hierarchy.
Vector< RefCountedPtr< LayoutData< BinFab< P > >> > AMRCellParticles
Alias for cell-sorted particles on each AMR level.
Definition: CD_ParticleContainer.H:40
Vector< RefCountedPtr< ParticleData< P > >> AMRParticles
Alias for ParticleData on each AMR level.
Definition: CD_ParticleContainer.H:34
Declaration of a static class containing some common useful particle routines that would otherwise be...
Implementation of CD_Timer.H.
Templated class for holding particles on an AMR hierarchy with particle remapping.
Definition: CD_ParticleContainer.H:50
AMRCellParticles< P > & getCellParticles()
Get all cell particles.
Definition: CD_ParticleContainerImplem.H:358
void getCellParticlesDestructive(BinFab< P > &a_cellParticles, const int a_lvl, const DataIndex a_dit)
Fill a cell-sorted particle data holder with all the particles in the grid patch. The original partic...
Definition: CD_ParticleContainerImplem.H:428
void organizeParticlesByPatch()
Sort particles by cell.
Definition: CD_ParticleContainerImplem.H:508
void define(const Vector< DisjointBoxLayout > &a_grids, const Vector< ProblemDomain > &a_domains, const Vector< Real > &a_dx, const Vector< int > &a_refRat, const Vector< ValidMask > &a_validMask, const Vector< RefCountedPtr< LevelTiles >> &a_levelTiles, const RealVect &a_probLo, const int a_blockingFactor, const int a_finestLevel, const std::string a_realm)
Define the container. This will do a clear-out of all particles.
Definition: CD_ParticleContainerImplem.H:77
const Vector< RealVect > getDx() const noexcept
Get grid resolutions.
Definition: CD_ParticleContainerImplem.H:263
void clearMaskParticles() const
Clear the "mask" particles.
Definition: CD_ParticleContainerImplem.H:1688
unsigned long long getNumberOfOutcastParticlesGlobal() const
Get global number of particles.
Definition: CD_ParticleContainerImplem.H:1478
void setupParticleData(const int a_base, const int a_finestLevel)
Setup function for the particle data (m_particles and m_maskParticles)
Definition: CD_ParticleContainerImplem.H:154
const RealVect getProbLo() const noexcept
Get lower-left corner.
Definition: CD_ParticleContainerImplem.H:256
unsigned long long getNumberOfValidParticlesLocal() const
Get local number of particles.
Definition: CD_ParticleContainerImplem.H:1407
void catenateParticleMaps(std::vector< std::map< std::pair< unsigned int, unsigned int >, List< P >>> &a_globalParticles, std::vector< std::map< std::pair< unsigned int, unsigned int >, List< P >>> &a_localParticles) const noexcept
Catenate the particles. This is usually called within OpenMP parallel regions.
Definition: CD_ParticleContainerImplem.H:1170
void setValue(const Real a_value)
Set the particle member to the input value.
Definition: CD_ParticleContainerImplem.H:1352
void preRegrid(const int a_base)
Cache particles before calling regrid.
Definition: CD_ParticleContainerImplem.H:1221
virtual ~ParticleContainer()
Destructor ( does nothing)
Definition: CD_ParticleContainerImplem.H:70
void sortParticles() noexcept
Sort particles according to the < operator in the particle.
Definition: CD_ParticleContainerImplem.H:197
const Vector< DisjointBoxLayout > & getGrids() const
Get the AMR grids.
Definition: CD_ParticleContainerImplem.H:249
unsigned long long getNumberOfOutcastParticlesLocal() const
Get local number of particles.
Definition: CD_ParticleContainerImplem.H:1449
ParticleData< P > & operator[](const int a_level)
Get particle data on a level.
Definition: CD_ParticleContainerImplem.H:332
void addParticles(const List< P > &a_particles)
Add particles to container.
Definition: CD_ParticleContainerImplem.H:551
void setupGrownGrids(const int a_base, const int a_finestLevel)
Set up grown grids.
Definition: CD_ParticleContainerImplem.H:123
unsigned long long getNumberOfMaskParticlesGlobal() const
Get the number particles in the halo cells.
Definition: CD_ParticleContainerImplem.H:1520
void transferMaskParticles(const Vector< RefCountedPtr< LevelData< BaseFab< bool >>>> &a_mask)
Copy particles to the mask particle data holder.
Definition: CD_ParticleContainerImplem.H:1599
void mapParticlesToAMRGrid(std::vector< std::map< std::pair< unsigned int, unsigned int >, List< P >>> &a_mappedParticles, List< P > &a_unmappedParticles) const noexcept
Iterate through the unmapped particles and map them to proper level, grid indices,...
Definition: CD_ParticleContainerImplem.H:1110
void clear(AMRParticles< P > &a_particles) const
Clear particles on input data holder.
Definition: CD_ParticleContainerImplem.H:1712
void remap()
Remap over the entire AMR hierarchy.
Definition: CD_ParticleContainerImplem.H:973
void assignLocalParticles(std::map< std::pair< unsigned int, unsigned int >, List< P >> &a_mappedParticles, AMRParticles< P > &a_particleData) const noexcept
Gather particles locally.
Definition: CD_ParticleContainerImplem.H:1191
void clearParticles()
Clear all particles.
Definition: CD_ParticleContainerImplem.H:1670
void clearOutcast() noexcept
Clear outcast particles.
Definition: CD_ParticleContainerImplem.H:1697
unsigned long long getNumberOfMaskParticlesLocal() const
Get the number particles in the halo cells.
Definition: CD_ParticleContainerImplem.H:1491
void organizeParticlesByCell()
Sort particles by cell.
Definition: CD_ParticleContainerImplem.H:474
int getFinestLevel() const
Get finest AMR level.
Definition: CD_ParticleContainerImplem.H:231
void clearBufferParticles() const
Clear the buffer particles.
Definition: CD_ParticleContainerImplem.H:1679
const std::string getRealm() const
Get the realm where this ParticleContainer lives.
Definition: CD_ParticleContainerImplem.H:240
void transferParticlesToSingleList(List< P > &a_list, AMRParticles< P > &a_particles) const noexcept
Gather the particles onto a single list.
Definition: CD_ParticleContainerImplem.H:1064
bool isOrganizedByCell() const
Is cell-sorted or not.
Definition: CD_ParticleContainerImplem.H:224
void regrid(const Vector< DisjointBoxLayout > &a_grids, const Vector< ProblemDomain > &a_domains, const Vector< Real > &a_dx, const Vector< int > &a_refRat, const Vector< ValidMask > &a_validMask, const Vector< RefCountedPtr< LevelTiles >> &a_levelTiles, const int a_base, const int a_newFinestLevel)
Regrid function. a_base is the coarsest grid level which did not change.
Definition: CD_ParticleContainerImplem.H:1262
void copyParticlesToSingleList(List< P > &a_list, const AMRParticles< P > &a_particles) const noexcept
Copy the input particles onto a single list.
Definition: CD_ParticleContainerImplem.H:1087
ParticleContainer()
Default constructor. Leaves object in undefined state.
Definition: CD_ParticleContainerImplem.H:34
void sanityCheck() const noexcept
Run a sanity check and make sure all particles are in their correctly assigned box.
Definition: CD_ParticleContainerImplem.H:1038
AMRParticles< P > & getMaskParticles()
Get the mask particles.
Definition: CD_ParticleContainerImplem.H:314
unsigned long long getNumberOfValidParticlesGlobal() const
Get global number of particles.
Definition: CD_ParticleContainerImplem.H:1436
void addParticlesDestructive(List< P > &a_particles)
Add particles to container. The input particles are destroyed by this routine.
Definition: CD_ParticleContainerImplem.H:617
AMRParticles< P > & getBufferParticles()
Get buffer particles on all levels.
Definition: CD_ParticleContainerImplem.H:296
void copyMaskParticles(const Vector< RefCountedPtr< LevelData< BaseFab< bool >>>> &a_mask) const
Copy particles to mask particle data holder.
Definition: CD_ParticleContainerImplem.H:1533
AMRParticles< P > & getParticles()
Get all particles on all levels.
Definition: CD_ParticleContainerImplem.H:270
void transferParticles(ParticleContainer< P > &a_otherContainer)
Move particles into this container.
Definition: CD_ParticleContainerImplem.H:914
static IntVect getParticleCellIndex(const RealVect &a_particlePosition, const RealVect &a_probLo, const Real &a_dx) noexcept
Get the cell index corresponding to the particle position.
Definition: CD_ParticleOpsImplem.H:26
ALWAYS_INLINE void loop(const Box &a_computeBox, Functor &&kernel, const IntVect &a_stride=IntVect::Unit)
Launch a C++ kernel over a regular grid.
Definition: CD_BoxLoopsImplem.H:20
RealVect position(const Location::Cell a_location, const VolIndex &a_vof, const EBISBox &a_ebisbox, const Real &a_dx)
Compute the position (ignoring the "origin) of a Vof.
Definition: CD_LocationImplem.H:20
Real sum(const Real &a_value) noexcept
Compute the sum across all MPI ranks.
Definition: CD_ParallelOpsImplem.H:353