chombo-discharge
CD_ParticleContainer.H
Go to the documentation of this file.
1 /* chombo-discharge
2  * Copyright © 2021 SINTEF Energy Research.
3  * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4  */
5 
12 #ifndef CD_ParticleContainer_H
13 #define CD_ParticleContainer_H
14 
15 // Chombo includes
16 #include <Particle.H>
17 #include <ParticleData.H>
18 #include <ParticleValidRegion.H>
19 #include <BinItem.H>
20 #include <BinFab.H>
21 #include <ProblemDomain.H>
22 #include <DisjointBoxLayout.H>
23 #include <BaseEBCellFAB.H>
24 
25 // Our includes
26 #include <CD_OpenMP.H>
27 #include <CD_LevelTiles.H>
28 #include <CD_NamespaceHeader.H>
29 
33 template <class P>
34 using AMRParticles = Vector<RefCountedPtr<ParticleData<P>>>;
35 
39 template <class P>
40 using AMRCellParticles = Vector<RefCountedPtr<LayoutData<BinFab<P>>>>;
41 
48 template <class P>
50 {
51 public:
52  using ValidMask = RefCountedPtr<LevelData<BaseFab<bool>>>;
53 
58 
72  ParticleContainer(const Vector<DisjointBoxLayout>& a_grids,
73  const Vector<ProblemDomain>& a_domains,
74  const Vector<Real>& a_dx,
75  const Vector<int>& a_refRat,
76  const Vector<ValidMask>& a_validMask,
77  const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
78  const RealVect& a_probLo,
79  const int a_blockingFactor,
80  const int a_finestLevel,
81  const std::string a_realm);
82 
86  virtual ~ParticleContainer();
87 
101  void
102  define(const Vector<DisjointBoxLayout>& a_grids,
103  const Vector<ProblemDomain>& a_domains,
104  const Vector<Real>& a_dx,
105  const Vector<int>& a_refRat,
106  const Vector<ValidMask>& a_validMask,
107  const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
108  const RealVect& a_probLo,
109  const int a_blockingFactor,
110  const int a_finestLevel,
111  const std::string a_realm);
112 
124  void
125  regrid(const Vector<DisjointBoxLayout>& a_grids,
126  const Vector<ProblemDomain>& a_domains,
127  const Vector<Real>& a_dx,
128  const Vector<int>& a_refRat,
129  const Vector<ValidMask>& a_validMask,
130  const Vector<RefCountedPtr<LevelTiles>>& a_levelTiles,
131  const int a_base,
132  const int a_newFinestLevel);
133 
138  void
139  preRegrid(const int a_base);
140 
146  void
147  copyMaskParticles(const Vector<RefCountedPtr<LevelData<BaseFab<bool>>>>& a_mask) const;
148 
154  void
155  copyMaskParticles(const int a_level, const LevelData<BaseFab<bool>>& a_mask) const;
156 
162  void
163  transferMaskParticles(const Vector<RefCountedPtr<LevelData<BaseFab<bool>>>>& a_mask);
164 
170  void
171  transferMaskParticles(const int a_level, const LevelData<BaseFab<bool>>& a_mask);
172 
176  void
177  sortParticles() noexcept;
178 
182  void
183  clearParticles();
184 
188  void
189  clearBufferParticles() const;
190 
194  void
195  clearMaskParticles() const;
196 
200  void
201  clearOutcast() noexcept;
202 
207  void
208  clear(AMRParticles<P>& a_particles) const;
209 
213  bool
214  isOrganizedByCell() const;
215 
219  int
220  getFinestLevel() const;
221 
225  const std::string
226  getRealm() const;
227 
233  getParticles();
234 
239  const AMRParticles<P>&
240  getParticles() const;
241 
248 
253  const AMRParticles<P>&
254  getBufferParticles() const;
255 
262 
267  const AMRParticles<P>&
268  getMaskParticles() const;
269 
273  const Vector<DisjointBoxLayout>&
274  getGrids() const;
275 
280  const RealVect
281  getProbLo() const noexcept;
282 
287  const Vector<RealVect>
288  getDx() const noexcept;
289 
294  ParticleData<P>&
295  operator[](const int a_level);
296 
301  const ParticleData<P>&
302  operator[](const int a_level) const;
303 
310 
315  const AMRCellParticles<P>&
316  getCellParticles() const;
317 
323  LayoutData<BinFab<P>>&
324  getCellParticles(const int a_level);
325 
331  const LayoutData<BinFab<P>>&
332  getCellParticles(const int a_level) const;
333 
340  BinFab<P>&
341  getCellParticles(const int a_level, const DataIndex a_dit);
342 
349  const BinFab<P>&
350  getCellParticles(const int a_level, const DataIndex a_dit) const;
351 
359  void
360  getCellParticles(BinFab<P>& a_cellParticles, const int a_lvl, const DataIndex a_dit) const;
361 
369  void
370  getCellParticlesDestructive(BinFab<P>& a_cellParticles, const int a_lvl, const DataIndex a_dit);
371 
376  void
378 
383  void
385 
390  void
391  addParticles(const List<P>& a_particles);
392 
397  void
398  addParticlesDestructive(List<P>& a_particles);
399 
406  void
407  addParticles(const BinFab<P>& a_particles, const int a_lvl, const DataIndex a_dit);
408 
415  void
416  addParticlesDestructive(BinFab<P>& a_particles, const int a_lvl, const DataIndex a_dit);
417 
422  void
423  addParticles(const ParticleContainer<P>& a_otherContainer);
424 
429  void
431 
435  void
436  sanityCheck() const noexcept;
437 
441  void
442  remap();
443 
449  void
450  transferParticles(ParticleContainer<P>& a_otherContainer);
451 
457  void
458  transferParticles(AMRParticles<P>& a_otherContainer);
459 
466  void
467  admitValidParticles(List<P>& a_evictedParticles, ParticleData<P>& a_particles, const int a_coarseLevel);
468 
476  template <Real& (P::*particleScalarField)()>
477  inline void
478  setValue(const Real a_value);
479 
487  template <RealVect& (P::*particleVectorField)()>
488  inline void
489  setValue(const RealVect a_value);
490 
494  unsigned long long
496 
500  unsigned long long
502 
506  unsigned long long
508 
512  unsigned long long
514 
518  unsigned long long
520 
524  unsigned long long
526 
527 protected:
528  // Reductions for performing catenation and join operations in a thread-safe manner.
529 #pragma omp declare reduction(catenate : List <P> : ThreadSafeCatenation <P>(omp_out, omp_in))
530 #pragma omp declare reduction(join : List <P> : ThreadSafeJoin <P>(omp_out, omp_in))
531 
535  std::string m_realm;
536 
541 
546 
550  RealVect m_probLo;
551 
556 
561 
565  bool m_profile;
566 
570  bool m_debug;
571 
575  bool m_verbose;
576 
580  Vector<RefCountedPtr<LevelTiles>> m_levelTiles;
581 
585  Vector<DisjointBoxLayout> m_grids;
586 
590  Vector<DisjointBoxLayout> m_oldGrids;
591 
595  Vector<ProblemDomain> m_domains;
596 
600  Vector<RefCountedPtr<LevelData<BaseFab<bool>>>> m_validRegion;
601 
605  Vector<BoxLayout> m_grownGrids;
606 
610  Vector<RealVect> m_dx;
611 
615  Vector<int> m_refRat;
616 
621 
630 
636 
642 
647 
653  void
654  setupGrownGrids(const int a_base, const int a_finestLevel);
655 
661  void
662  setupParticleData(const int a_base, const int a_finestLevel);
663 
670  inline void
671  transferParticlesToSingleList(List<P>& a_list, AMRParticles<P>& a_particles) const noexcept;
672 
679  inline void
680  copyParticlesToSingleList(List<P>& a_list, const AMRParticles<P>& a_particles) const noexcept;
681 
687  inline void
688  mapParticlesToAMRGrid(std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_mappedParticles,
689  List<P>& a_unmappedParticles) const noexcept;
690 
696  inline void
698  std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_globalParticles,
699  std::vector<std::map<std::pair<unsigned int, unsigned int>, List<P>>>& a_localParticles) const noexcept;
700 
706  inline void
707  assignLocalParticles(std::map<std::pair<unsigned int, unsigned int>, List<P>>& a_mappedParticles,
708  AMRParticles<P>& a_particleData) const noexcept;
709 };
710 
711 #include <CD_NamespaceFooter.H>
712 
714 
715 #endif
Declaration of LevelTiles.
Declaration of various useful OpenMP-related utilities.
Implementation of CD_ParticleContainer.H.
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
Templated class for holding particles on an AMR hierarchy with particle remapping.
Definition: CD_ParticleContainer.H:50
Vector< ProblemDomain > m_domains
Problem domains.
Definition: CD_ParticleContainer.H:595
AMRCellParticles< P > & getCellParticles()
Get all cell particles.
Definition: CD_ParticleContainerImplem.H:358
Vector< BoxLayout > m_grownGrids
Grown layouts.
Definition: CD_ParticleContainer.H:605
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
bool m_verbose
Verbose or not.
Definition: CD_ParticleContainer.H:575
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
Vector< RefCountedPtr< LevelData< BaseFab< bool > > > > m_validRegion
Map of valid cells; evaluates to false if a cell is covered by a finer grid.
Definition: CD_ParticleContainer.H:600
AMRCellParticles< P > m_cellSortedParticles
Cell particles.
Definition: CD_ParticleContainer.H:646
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
AMRParticles< P > m_particles
The actual particles that this ParticleContainer represents.
Definition: CD_ParticleContainer.H:620
const RealVect getProbLo() const noexcept
Get lower-left corner.
Definition: CD_ParticleContainerImplem.H:256
void admitValidParticles(List< P > &a_evictedParticles, ParticleData< P > &a_particles, const int a_coarseLevel)
Evict particles if they move out of the valid region.
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
bool m_debug
Debug or not.
Definition: CD_ParticleContainer.H:570
bool m_isOrganizedByCell
Check if particle container is "cell sorted".
Definition: CD_ParticleContainer.H:560
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
Vector< RefCountedPtr< LevelTiles > > m_levelTiles
Tiled AMR space.
Definition: CD_ParticleContainer.H:580
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
int m_blockingFactor
Blocking factor, aka grid size.
Definition: CD_ParticleContainer.H:540
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
Vector< int > m_refRat
Refinement ratios. Entry at index 'l' is the refinement between level 'l' and level 'l+1'.
Definition: CD_ParticleContainer.H:615
bool m_profile
Profile or not.
Definition: CD_ParticleContainer.H:565
void clearParticles()
Clear all particles.
Definition: CD_ParticleContainerImplem.H:1670
int m_finestLevel
Finest grid level.
Definition: CD_ParticleContainer.H:545
bool m_isDefined
Check if particle container is defined.
Definition: CD_ParticleContainer.H:555
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
AMRParticles< P > m_bufferParticles
Special data holder that holds particles on a grown grid layout.
Definition: CD_ParticleContainer.H:629
int getFinestLevel() const
Get finest AMR level.
Definition: CD_ParticleContainerImplem.H:231
std::string m_realm
Realm on which the ParticleContainer is defined.
Definition: CD_ParticleContainer.H:535
void clearBufferParticles() const
Clear the buffer particles.
Definition: CD_ParticleContainerImplem.H:1679
Vector< DisjointBoxLayout > m_grids
AMR grids.
Definition: CD_ParticleContainer.H:585
const std::string getRealm() const
Get the realm where this ParticleContainer lives.
Definition: CD_ParticleContainerImplem.H:240
RealVect m_probLo
Lower left corner of the physical domain.
Definition: CD_ParticleContainer.H:550
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
Vector< DisjointBoxLayout > m_oldGrids
Old AMR grids. Used during regrids.
Definition: CD_ParticleContainer.H:590
AMRParticles< P > m_maskParticles
Special data holder that holds copies of particles on the coarse level that are within a specified ra...
Definition: CD_ParticleContainer.H:635
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
Vector< RealVect > m_dx
Resolutions on each grid level.
Definition: CD_ParticleContainer.H:610
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
AMRParticles< P > m_cacheParticles
Particles stored on the old grid.
Definition: CD_ParticleContainer.H:641
void transferParticles(ParticleContainer< P > &a_otherContainer)
Move particles into this container.
Definition: CD_ParticleContainerImplem.H:914