Block-Structured AMR Software Framework
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator > Class Template Reference

A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured AMR hierarchy. More...

#include <AMReX_ParticleContainer.H>

Inheritance diagram for amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >:
amrex::ParticleContainerBase

Public Types

template<typename T >
using AllocatorType = Allocator< T >
 The memory allocator in use. More...
 
using ParticleType = Particle< NStructReal, NStructInt >
 The type of Particles we hold. More...
 
using SuperParticleType = Particle< NStructReal+NArrayReal, NStructInt+NArrayInt >
 The type of the "SuperParticle" which stored all components in AoS form. More...
 
using RealType = typename Particle< NStructReal, NStructInt >::RealType
 The type of the Real data. More...
 
using ParticleContainerType = ParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 
using ParticleTileType = ParticleTile< NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 
using ParticleInitData = ParticleInitType< NStructReal, NStructInt, NArrayReal, NArrayInt >
 
using ParticleLevel = std::map< std::pair< int, int >, ParticleTileType >
 
using AoS = typename ParticleTileType::AoS
 
using SoA = typename ParticleTileType::SoA
 
using RealVector = typename SoA::RealVector
 
using IntVector = typename SoA::IntVector
 
using ParticleVector = typename AoS::ParticleVector
 
using CharVector = Gpu::DeviceVector< char >
 
using SendBuffer = Gpu::PolymorphicVector< char >
 
using ParIterType = ParIter< NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 
using ParConstIterType = ParConstIter< NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 
template<template< class > class NewAllocator = amrex::DefaultAllocator>
using ContainerLike = amrex::ParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt, NewAllocator >
 

Public Member Functions

 ParticleContainer ()
 Default constructor - construct an empty particle container that has no concept of a level hierarchy. Must be properly initialized later. More...
 
 ParticleContainer (ParGDBBase *gdb)
 Construct a particle container using a ParGDB object. The container will track changes in the grid structure of the ParGDB automatically. More...
 
 ParticleContainer (const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
 Construct a particle container using a given Geometry, DistributionMapping, and BoxArray. Single level version. More...
 
 ParticleContainer (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
 Construct a particle container using a given Geometry, DistributionMapping, BoxArray and Vector of refinement ratios. Multi-level version. More...
 
 ParticleContainer (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
 Same as the above, but accepts different refinement ratios in each direction. More...
 
 ~ParticleContainer () override=default
 
 ParticleContainer (const ParticleContainer &)=delete
 
ParticleContaineroperator= (const ParticleContainer &)=delete
 
 ParticleContainer (ParticleContainer &&) noexcept=default
 
ParticleContaineroperator= (ParticleContainer &&) noexcept=default
 
void Define (ParGDBBase *gdb)
 Define a default-constructed ParticleContainer using a ParGDB object. The container will track changes in the grid structure of the ParGDB automatically. More...
 
void Define (const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
 Define a default-constructed ParticleContainer using a ParGDB object. Single-level version. More...
 
void Define (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
 Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version. More...
 
void Define (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
 Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version. More...
 
int numLocalTilesAtLevel (int lev) const
 The total number of tiles on this rank on this level. More...
 
void reserveData () override
 
void resizeData () override
 
void InitFromAsciiFile (const std::string &file, int extradata, const IntVect *Nrep=nullptr)
 
void InitFromBinaryFile (const std::string &file, int extradata)
 
void InitFromBinaryMetaFile (const std::string &file, int extradata)
 
void InitRandom (Long icount, ULong iseed, const ParticleInitData &pdata, bool serialize=false, RealBox bx=RealBox())
 This initializes the particle container with icount randomly distributed particles. If serialize is true, then the particles will all be generated on the IO Process, and the particle positions will be broadcast to all other process. If serialize is false, then the particle positions will be randomly generated in parallel, which each process using the random seed iseed + MyProc. The particles can be constrained to lie within the RealBox bx, if so desired. The default is the full domain. More...
 
void InitRandomPerBox (Long icount, ULong iseed, const ParticleInitData &pdata)
 This initializes the container with icount randomly distributed particles per box, using the random seed iseed. All the particles have the same data and attributes, which are passed using the pdata struct. More...
 
void InitOnePerCell (Real x_off, Real y_off, Real z_off, const ParticleInitData &pdata)
 This initializes the particle container with one particle per cell, where the other particle data and attributes are all contant. The coarsest level is used to generate the particle positions. The particle variable values are passed in through the pdata struct. The parameters x_off, y_off, and z_off represent offsets between 0 and 1 that show where inside the cells to place the particles. 0.5 means cell centered. More...
 
void InitNRandomPerCell (int n_per_cell, const ParticleInitData &pdata)
 This initializes the particle container with n_per_cell randomly distributed particles per cell, where the other particle data and and attributes are all constant. The cells on the coarsest level are used to generate the particle positions. The particle variable values are passed in through the pdata struct. More...
 
void Increment (MultiFab &mf, int level)
 
Long IncrementWithTotal (MultiFab &mf, int level, bool local=false)
 
void Redistribute (int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
 Redistribute puts all the particles back in the right places (for some value of right) More...
 
template<class index_type >
void ReorderParticles (int lev, const MFIter &mfi, const index_type *permutations)
 Reorder particles on the tile given by lev and mfi using a the permutations array. More...
 
void SortParticlesForDeposition (IntVect idx_type)
 Sort particles on each tile such that particles adjacent in memory are likely to map to adjacent cells. This ordering can be beneficial for performance on GPU when deposition quantities onto a grid. More...
 
void SortParticlesByCell ()
 Sort the particles on each tile by cell, using Fortran ordering. More...
 
void SortParticlesByBin (IntVect bin_size)
 Sort the particles on each tile by groups of cells, given an IntVect bin_size. More...
 
bool OK (int lev_min=0, int lev_max=-1, int nGrow=0) const
 OK checks that all particles are in the right places (for some value of right) More...
 
std::array< Long, 3 > ByteSpread () const
 
std::array< Long, 3 > PrintCapacity () const
 
void ShrinkToFit ()
 
Long NumberOfParticlesAtLevel (int level, bool only_valid=true, bool only_local=false) const
 Returns # of particles at specified the level. More...
 
Vector< Long > NumberOfParticlesInGrid (int level, bool only_valid=true, bool only_local=false) const
 
Long TotalNumberOfParticles (bool only_valid=true, bool only_local=false) const
 Returns # of particles at all levels. More...
 
void RemoveParticlesAtLevel (int level)
 The Following methods are for managing Virtual and Ghost Particles. More...
 
void RemoveParticlesNotAtFinestLevel ()
 
void CreateVirtualParticles (int level, AoS &virts) const
 Creates virtual particles for a given level that represent in some capacity all particles at finer levels. More...
 
void CreateGhostParticles (int level, int ngrow, AoS &ghosts) const
 Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in level-1. More...
 
void AddParticlesAtLevel (AoS &particles, int level, int nGrow=0)
 Add particles from a pbox to the grid at this level. More...
 
void CreateVirtualParticles (int level, ParticleTileType &virts) const
 Creates virtual particles for a given level that represent in some capacity all particles at finer levels. More...
 
void CreateGhostParticles (int level, int ngrow, ParticleTileType &ghosts) const
 Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in level-1. More...
 
void AddParticlesAtLevel (ParticleTileType &particles, int level, int nGrow=0)
 Add particles from a pbox to the grid at this level. More...
 
void clearParticles ()
 Clear all the particles in this container. This does not free memory. More...
 
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0>
void copyParticles (const PCType &other, bool local=false)
 Copy particles from other to this ParticleContainer. Will clear all the particles from this container first. local controls whether or not to call Redistribute() after adding the particles. More...
 
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0>
void addParticles (const PCType &other, bool local=false)
 Add particles from other to this ParticleContainer. local controls whether or not to call Redistribute after adding the particles. More...
 
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0, std::enable_if_t<! std::is_integral< F >::value, int > bar = 0>
void copyParticles (const PCType &other, F &&f, bool local=false)
 Copy particles from other to this ParticleContainer. Will clear all the particles from this container first. local controls whether or not to call Redistribute() after adding the particles. More...
 
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0, std::enable_if_t<! std::is_integral< F >::value, int > bar = 0>
void addParticles (const PCType &other, F &&f, bool local=false)
 Add particles from other to this ParticleContainer. local controls whether or not to call Redistribute after adding the particles. More...
 
void WriteParticleRealData (void *data, size_t size, std::ostream &os) const
 Write a contiguous chunk of real particle data to an ostream. More...
 
void ReadParticleRealData (void *data, size_t size, std::istream &is)
 Read a contiguous chunk of real particle data from an istream. More...
 
void Checkpoint (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names=Vector< std::string >(), const Vector< std::string > &int_comp_names=Vector< std::string >()) const
 Writes a particle checkpoint to file, suitable for restarting. More...
 
void Checkpoint (const std::string &dir, const std::string &name, bool is_checkpoint, const Vector< std::string > &real_comp_names=Vector< std::string >(), const Vector< std::string > &int_comp_names=Vector< std::string >()) const
 Writes a particle checkpoint to file, suitable for restarting. This version allows the particle component names to be passed in. This overload exists for backwards comptability. The is_checkpoint parameter is ignored. More...
 
template<class F >
void WriteBinaryParticleData (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f, bool is_checkpoint=false) const
 Writes particle data to disk in the AMReX native format. More...
 
void CheckpointPre ()
 
void CheckpointPost ()
 
void Restart (const std::string &dir, const std::string &file)
 Restart from checkpoint. More...
 
void Restart (const std::string &dir, const std::string &file, bool is_checkpoint)
 Older version, for backwards compatibility. More...
 
void WritePlotFile (const std::string &dir, const std::string &name) const
 This version of WritePlotFile writes all components and assigns component names. More...
 
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string > & >::value >::type * = nullptr>
void WritePlotFile (const std::string &dir, const std::string &name, F &&f) const
 This version of WritePlotFile writes all components and assigns component names. More...
 
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names) const
 This version of WritePlotFile writes all components and allows the user to specify the names of the components. More...
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f) const
 This version of WritePlotFile writes all components and allows the user to specify the names of the components. More...
 
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names) const
 This version of WritePlotFile writes all components and allows the user to specify the names of the components. More...
 
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string >>::value >::type * = nullptr>
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, F &&f) const
 This version of WritePlotFile writes all components and allows the user to specify the names of the components. More...
 
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp) const
 This version of WritePlotFile assigns component names, but allows the user to pass in a vector of ints that toggle on / off the writing of specific components. More...
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, F &&f) const
 This version of WritePlotFile assigns component names, but allows the user to pass in a vector of ints that toggle on / off the writing of specific components. More...
 
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names) const
 This is the most general version of WritePlotFile, which takes component names and flags for whether to write each variable as components. Note that the user should pass in vectors containing names of all the components, whether they are written or not. More...
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f) const
 This is the most general version of WritePlotFile, which takes component names and flags for whether to write each variable as components. Note that the user should pass in vectors containing names of all the components, whether they are written or not. More...
 
void WritePlotFilePre ()
 
void WritePlotFilePost ()
 
void WriteAsciiFile (const std::string &file)
 
const Vector< ParticleLevel > & GetParticles () const
 Return the underyling Vector (over AMR levels) of ParticleLevels. Const version. More...
 
Vector< ParticleLevel > & GetParticles ()
 Return the underyling Vector (over AMR levels) of ParticleLevels. Non-const version. More...
 
const ParticleLevelGetParticles (int lev) const
 Return the ParticleLevel for level "lev". Const version. More...
 
ParticleLevelGetParticles (int lev)
 Return the ParticleLevel for level "lev". Non-const version. More...
 
const ParticleTileTypeParticlesAt (int lev, int grid, int tile) const
 Return the ParticleTile for level "lev", grid "grid" and tile "tile." Const version. More...
 
ParticleTileTypeParticlesAt (int lev, int grid, int tile)
 Return the ParticleTile for level "lev", grid "grid" and tile "tile." Non-const version. More...
 
template<class Iterator >
const ParticleTileTypeParticlesAt (int lev, const Iterator &iter) const
 Return the ParticleTile for level "lev" and Iterator "iter". Const version. More...
 
template<class Iterator >
ParticleTileTypeParticlesAt (int lev, const Iterator &iter)
 Return the ParticleTile for level "lev" and Iterator "iter". Non-const version. More...
 
ParticleTileTypeDefineAndReturnParticleTile (int lev, int grid, int tile)
 Define and return the ParticleTile for level "lev", grid "grid" and tile "tile.". More...
 
template<class Iterator >
ParticleTileTypeDefineAndReturnParticleTile (int lev, const Iterator &iter)
 Define and return the ParticleTile for level "lev", and Iterator "iter". More...
 
void AssignDensity (int rho_index, Vector< std::unique_ptr< MultiFab > > &mf_to_be_filled, int lev_min, int ncomp, int finest_level, int ngrow=2) const
 Functions depending the layout of the data. Use with caution. More...
 
void Interpolate (Vector< std::unique_ptr< MultiFab > > &mesh_data, int lev_min, int lev_max)
 
void InterpolateSingleLevel (MultiFab &mesh_data, int lev)
 
void AssignCellDensitySingleLevel (int rho_index, MultiFab &mf, int level, int ncomp=1, int particle_lvl_offset=0) const
 
template<typename P >
IntVect Index (const P &p, int lev) const
 
ParticleLocData Reset (ParticleType &prt, bool update, bool verbose=true, ParticleLocData pld=ParticleLocData()) const
 Updates a particle's location (Where), tries to periodic shift any particles that have left the domain. May need work (see inline comments) More...
 
bool PeriodicShift (ParticleType &prt) const
 Returns true if the particle was shifted. More...
 
void SetLevelDirectoriesCreated (bool tf)
 
bool GetLevelDirectoriesCreated () const
 
void SetUsePrePost (bool tf) const
 
bool GetUsePrePost () const
 
int GetMaxNextIDPrePost () const
 
Long GetNParticlesPrePost () const
 
void SetUseUnlink (bool tf) const
 
bool GetUseUnlink () const
 
void RedistributeCPU (int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
 
void RedistributeGPU (int lev_min=0, int lev_max=-1, int nGrow=0, int local=0, bool remove_negative=true)
 
Long superParticleSize () const
 
template<typename T , typename std::enable_if< std::is_same< T, bool >::value, int >::type = 0>
void AddRealComp (T communicate=true)
 
template<typename T , typename std::enable_if< std::is_same< T, bool >::value, int >::type = 0>
void AddIntComp (T communicate=true)
 
int NumRuntimeRealComps () const
 
int NumRuntimeIntComps () const
 
int NumRealComps () const
 
int NumIntComps () const
 
void ResizeRuntimeRealComp (int new_size, bool communicate)
 
void ResizeRuntimeIntComp (int new_size, bool communicate)
 
template<template< class > class NewAllocator = amrex::DefaultAllocator>
ContainerLike< NewAllocator > make_alike () const
 
void WriteParticles (int level, std::ofstream &ofs, int fnum, Vector< int > &which, Vector< int > &count, Vector< Long > &where, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::map< std::pair< int, int >, IntVector >> &particle_io_flags, bool is_checkpoint) const
 
template<typename P >
IntVect Index (const P &p, int lev) const
 
template<typename P >
bool Where (const P &p, ParticleLocData &pld, int lev_min, int lev_max, int nGrow, int local_grid) const
 
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo>
void copyParticles (const PCType &other, bool local)
 
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo>
void addParticles (const PCType &other, bool local)
 
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo, std::enable_if_t<! std::is_integral< F >::value, int > bar>
void copyParticles (const PCType &other, F &&f, bool local)
 
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo, std::enable_if_t<! std::is_integral< F >::value, int > bar>
void addParticles (const PCType &other, F &&f, bool local)
 
template<class index_type >
void ReorderParticles (int lev, const MFIter &mfi, const index_type *permutations)
 
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string > & >::value >::type * >
void WritePlotFile (const std::string &dir, const std::string &name, F &&f) const
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f) const
 
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string >>::value >::type * >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< std::string > &real_comp_names, F &&f) const
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, F &&f) const
 
template<class F >
void WritePlotFile (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f) const
 
template<class F >
void WriteBinaryParticleData (const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F &&f, bool is_checkpoint) const
 
template<class RTYPE >
void ReadParticles (int cnt, int grd, int lev, std::ifstream &ifs, int finest_level_in_file, bool convert_ids)
 
- Public Member Functions inherited from amrex::ParticleContainerBase
 ParticleContainerBase ()=default
 
 ParticleContainerBase (ParGDBBase *gdb)
 
 ParticleContainerBase (const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
 
 ParticleContainerBase (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
 
 ParticleContainerBase (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
 
virtual ~ParticleContainerBase ()=default
 
 ParticleContainerBase (const ParticleContainerBase &)=delete
 
ParticleContainerBaseoperator= (const ParticleContainerBase &)=delete
 
 ParticleContainerBase (ParticleContainerBase &&)=default
 
ParticleContainerBaseoperator= (ParticleContainerBase &&)=default
 
void Define (ParGDBBase *gdb)
 
void Define (const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
 
void Define (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
 
void Define (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
 
bool isDefined () const
 
void RedefineDummyMF (int lev)
 
MFIter MakeMFIter (int lev, const MFItInfo &info) const
 
MFIter MakeMFIter (int lev) const
 
MFIter MakeMFIter (int lev, bool tile) const
 
void SetParGDB (const Geometry &geom, const DistributionMapping &dmap, const BoxArray &ba)
 Set the particle Geometry, DistributionMapping, and BoxArray. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. This is the single-level version. More...
 
void SetParGDB (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< IntVect > &rr)
 Set the particle Geometry, DistributionMapping, ref ratios, and BoxArray. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. This is the multi-level version. More...
 
void SetParGDB (const Vector< Geometry > &geom, const Vector< DistributionMapping > &dmap, const Vector< BoxArray > &ba, const Vector< int > &rr)
 Set the particle Geometry, DistributionMapping, ref ratios, and BoxArray. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. This is the multi-level version. More...
 
void SetParticleBoxArray (int lev, BoxArray new_ba)
 Set the particle BoxArray. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. More...
 
void SetParticleDistributionMap (int lev, DistributionMapping new_dmap)
 Set the particle DistributionMapping. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. More...
 
void SetParticleGeometry (int lev, Geometry new_geom)
 Set the particle Geometry. If the container was previously set to to track the AMR hierarchy of an AmrCore or AmrLevel object, that correspondence will be broken here. More...
 
const BoxArrayParticleBoxArray (int lev) const
 Get the BoxArray for a given level. More...
 
const DistributionMappingParticleDistributionMap (int lev) const
 Get the DistributionMapping for a given level. More...
 
const GeometryGeom (int lev) const
 Get the Geometry for a given level. More...
 
const GeometryParticleGeom (int lev) const
 Get the particle Geometry for a given level. More...
 
int finestLevel () const
 the finest level actually defined for the ParticleContainer More...
 
int maxLevel () const
 the finest allowed level in the ParticleContainer, whether it is defined or not. More...
 
int numLevels () const
 the number of defined levels in the ParticleContainer More...
 
const ParGDBBaseGetParGDB () const
 Get the ParGDB object used to define this container (const version) More...
 
ParGDBBaseGetParGDB ()
 Get the ParGDB object used to define this container. More...
 
int Verbose () const
 
void SetVerbose (int verbose)
 
const ParticleBufferMapBufferMap () const
 
Vector< intNeighborProcs (int ngrow) const
 
template<class MF >
bool OnSameGrids (int level, const MF &mf) const
 

Public Attributes

RealDescriptor ParticleRealDescriptor = FPC::Native64RealDescriptor()
 
Vector< inth_redistribute_real_comp
 
Vector< inth_redistribute_int_comp
 
bool levelDirectoriesCreated
 Variables for i/o optimization saved for pre and post checkpoint. More...
 
bool usePrePost
 
bool doUnlink
 
int maxnextidPrePost
 
int nOutFilesPrePost
 
Long nparticlesPrePost
 
Vector< Long > nParticlesAtLevelPrePost
 
Vector< Vector< int > > whichPrePost
 
Vector< Vector< int > > countPrePost
 
Vector< Vector< Long > > wherePrePost
 
std::string HdrFileNamePrePost
 
Vector< std::string > filePrefixPrePost
 
- Public Attributes inherited from amrex::ParticleContainerBase
AmrParticleLocator< DenseBins< Box > > m_particle_locator
 

Static Public Attributes

static constexpr int NStructReal = T_NStructReal
 Number of extra Real components in the particle struct. More...
 
static constexpr int NStructInt = T_NStructInt
 Number of extra integer components in the particle struct. More...
 
static constexpr int NArrayReal = T_NArrayReal
 Number of extra Real components stored in struct-of-array form. More...
 
static constexpr int NArrayInt = T_NArrayInt
 Number of extra integer components stored in struct-of-array form. More...
 
- Static Public Attributes inherited from amrex::ParticleContainerBase
static AMREX_EXPORT bool do_tiling = false
 
static AMREX_EXPORT IntVect tile_size { AMREX_D_DECL(1024000,8,8) }
 
static AMREX_EXPORT bool memEfficientSort = true
 

Protected Member Functions

template<typename P >
bool Where (const P &prt, ParticleLocData &pld, int lev_min=0, int lev_max=-1, int nGrow=0, int local_grid=-1) const
 Checks a particle's location on levels lev_min and higher. Returns false if the particle does not exist on that level. Only if lev_min == lev_max, nGrow can be > 0 (i.e., including nGrow ghost cells). More...
 
bool EnforcePeriodicWhere (ParticleType &prt, ParticleLocData &pld, int lev_min=0, int lev_max=-1, int local_grid=-1) const
 Checks whether the particle has crossed a periodic boundary in such a way that it is on levels lev_min and higher. More...
 
template<class RTYPE >
void ReadParticles (int cnt, int grd, int lev, std::ifstream &ifs, int finest_level_in_file, bool convert_ids)
 
void SetParticleSize ()
 
- Protected Member Functions inherited from amrex::ParticleContainerBase
void BuildRedistributeMask (int lev, int nghost=1) const
 
void defineBufferMap () const
 

Protected Attributes

DenseBins< ParticleTypem_bins
 
- Protected Attributes inherited from amrex::ParticleContainerBase
int m_verbose {0}
 
std::unique_ptr< ParGDBm_gdb_object = std::make_unique<ParGDB>()
 
ParGDBBasem_gdb {nullptr}
 
Vector< std::unique_ptr< MultiFab > > m_dummy_mf
 
std::unique_ptr< iMultiFabredistribute_mask_ptr
 
int redistribute_mask_nghost = std::numeric_limits<int>::min()
 
amrex::Vector< intneighbor_procs
 
ParticleBufferMap m_buffer_map
 

Private Member Functions

virtual void particlePostLocate (ParticleType &, const ParticleLocData &, const int)
 
virtual void correctCellVectors (int, int, int, const ParticleType &)
 
void RedistributeMPI (std::map< int, Vector< char > > &not_ours, int lev_min=0, int lev_max=0, int nGrow=0, int local=0)
 
void locateParticle (ParticleType &p, ParticleLocData &pld, int lev_min, int lev_max, int nGrow, int local_grid=-1) const
 
void Initialize ()
 

Private Attributes

bool m_runtime_comps_defined {false}
 
int m_num_runtime_real {0}
 
int m_num_runtime_int {0}
 
size_t particle_size
 
size_t superparticle_size
 
int num_real_comm_comps
 
int num_int_comm_comps
 
Vector< ParticleLevelm_particles
 

Friends

class ParIterBase< true, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 
class ParIterBase< false, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
 

Additional Inherited Members

- Static Public Member Functions inherited from amrex::ParticleContainerBase
static const std::string & CheckpointVersion ()
 
static const std::string & PlotfileVersion ()
 
static const std::string & DataPrefix ()
 
static int MaxReaders ()
 
static Long MaxParticlesPerRead ()
 
static const std::string & AggregationType ()
 
static int AggregationBuffer ()
 

Detailed Description

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
class amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >

A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured AMR hierarchy.

The data layout on a single tile is determined by the value of the following template parameters:

Template Parameters
T_NStructRealThe number of extra Real components in the particle struct
T_NStructIntThe number of extra integer components in the particle struct
T_NArrayRealThe number of extra Real components stored in struct-of-array form
T_NArrayIntThe number of extra integer components stored in struct-of-array form

Member Typedef Documentation

◆ AllocatorType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename T >
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AllocatorType = Allocator<T>

The memory allocator in use.

◆ AoS

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AoS = typename ParticleTileType::AoS

◆ CharVector

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::CharVector = Gpu::DeviceVector<char>

◆ ContainerLike

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<template< class > class NewAllocator = amrex::DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ContainerLike = amrex::ParticleContainer<NStructReal, NStructInt, NArrayReal, NArrayInt, NewAllocator>

type trait to translate one particle container to another, with changed allocator

◆ IntVector

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::IntVector = typename SoA::IntVector

◆ ParConstIterType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParConstIterType = ParConstIter<NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator>

◆ ParIterType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParIterType = ParIter<NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator>

◆ ParticleContainerType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainerType = ParticleContainer<NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator>

◆ ParticleInitData

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleInitData = ParticleInitType<NStructReal, NStructInt, NArrayReal, NArrayInt>

◆ ParticleLevel

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleLevel = std::map<std::pair<int, int>, ParticleTileType>

A single level worth of particles is indexed (grid id, tile id) for both SoA and AoS data.

◆ ParticleTileType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleTileType = ParticleTile<NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator>

◆ ParticleType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleType = Particle<NStructReal, NStructInt>

The type of Particles we hold.

◆ ParticleVector

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleVector = typename AoS::ParticleVector

◆ RealType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::RealType = typename Particle<NStructReal, NStructInt>::RealType

The type of the Real data.

◆ RealVector

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::RealVector = typename SoA::RealVector

◆ SendBuffer

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SendBuffer = Gpu::PolymorphicVector<char>

◆ SoA

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SoA = typename ParticleTileType::SoA

◆ SuperParticleType

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
using amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SuperParticleType = Particle<NStructReal+NArrayReal, NStructInt+NArrayInt>

The type of the "SuperParticle" which stored all components in AoS form.

Constructor & Destructor Documentation

◆ ParticleContainer() [1/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( )
inline

Default constructor - construct an empty particle container that has no concept of a level hierarchy. Must be properly initialized later.

◆ ParticleContainer() [2/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( ParGDBBase gdb)
inline

Construct a particle container using a ParGDB object. The container will track changes in the grid structure of the ParGDB automatically.

Parameters
gdbA pointer to a ParGDBBase, which contains pointers to the Geometry, DistributionMapping, and BoxArray objects that define the AMR hierarchy. Usually, this is generated by an AmrCore or AmrLevel object.

◆ ParticleContainer() [3/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( const Geometry geom,
const DistributionMapping dmap,
const BoxArray ba 
)
inline

Construct a particle container using a given Geometry, DistributionMapping, and BoxArray. Single level version.

Parameters
theGeometry object, which describes the problem domain
ADistributionMapping, which describes how the boxes are distributed onto MPI tasks
ABoxArray, which gives the set of grid boxes

◆ ParticleContainer() [4/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( const Vector< Geometry > &  geom,
const Vector< DistributionMapping > &  dmap,
const Vector< BoxArray > &  ba,
const Vector< int > &  rr 
)
inline

Construct a particle container using a given Geometry, DistributionMapping, BoxArray and Vector of refinement ratios. Multi-level version.

Parameters
geomA Vector of Geometry objects, one for each level
dmapA Vector of DistributionMappings, one for each level
baA Vector of BoxArrays, one for each level
rrA Vector of integer refinement ratios, of size num_levels - 1. rr[n] gives the refinement ratio between levels n and n+1

◆ ParticleContainer() [5/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( const Vector< Geometry > &  geom,
const Vector< DistributionMapping > &  dmap,
const Vector< BoxArray > &  ba,
const Vector< IntVect > &  rr 
)
inline

Same as the above, but accepts different refinement ratios in each direction.

Parameters
geomA Vector of Geometry objects, one for each level
dmapA Vector of DistributionMappings, one for each level
baA Vector of BoxArrays, one for each level
rrA Vector of IntVect refinement ratios, of size num_levels - 1. rr[n] gives the refinement ratio between levels n and n+1

◆ ~ParticleContainer()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::~ParticleContainer ( )
overridedefault

◆ ParticleContainer() [6/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( const ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator > &  )
delete

◆ ParticleContainer() [7/7]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleContainer ( ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator > &&  )
defaultnoexcept

Member Function Documentation

◆ AddIntComp()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename T , typename std::enable_if< std::is_same< T, bool >::value, int >::type = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AddIntComp ( communicate = true)
inline

◆ addParticles() [1/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::addParticles ( const PCType &  other,
bool  local 
)

◆ addParticles() [2/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::addParticles ( const PCType &  other,
bool  local = false 
)

Add particles from other to this ParticleContainer. local controls whether or not to call Redistribute after adding the particles.

Parameters
otherthe other pc to copy from
localwhether to call redistribute after

◆ addParticles() [3/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo, std::enable_if_t<! std::is_integral< F >::value, int > bar>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::addParticles ( const PCType &  other,
F &&  f,
bool  local 
)

◆ addParticles() [4/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0, std::enable_if_t<! std::is_integral< F >::value, int > bar = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::addParticles ( const PCType &  other,
F &&  f,
bool  local = false 
)

Add particles from other to this ParticleContainer. local controls whether or not to call Redistribute after adding the particles.

This version conditionally copies based on a predicate function applied to each particle.

Template Parameters
callablethat takes a SuperParticle and returns a bool
Parameters
otherthe other pc to copy from
ffunction to apply to each particle as a filter
localwhether to call redistribute after

◆ AddParticlesAtLevel() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AddParticlesAtLevel ( AoS particles,
int  level,
int  nGrow = 0 
)

Add particles from a pbox to the grid at this level.

Parameters
particles
level
nGrow

◆ AddParticlesAtLevel() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AddParticlesAtLevel ( ParticleTileType particles,
int  level,
int  nGrow = 0 
)

Add particles from a pbox to the grid at this level.

Parameters
particles
level
nGrow

◆ AddRealComp()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename T , typename std::enable_if< std::is_same< T, bool >::value, int >::type = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::AddRealComp ( communicate = true)
inline

◆ AssignCellDensitySingleLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::AssignCellDensitySingleLevel ( int  rho_index,
MultiFab mf,
int  level,
int  ncomp = 1,
int  particle_lvl_offset = 0 
) const

◆ AssignDensity()

template<int NStructReal, int NStructInt, int NArrayReal, int NArrayInt, template< class > class Allocator>
void amrex::ParticleContainer< NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >::AssignDensity ( int  rho_index,
Vector< std::unique_ptr< MultiFab > > &  mf_to_be_filled,
int  lev_min,
int  ncomp,
int  finest_level,
int  ngrow = 2 
) const

Functions depending the layout of the data. Use with caution.

Parameters
rho_index
mf_to_be_filled
lev_min
ncomp
finest_level
ngrow

◆ ByteSpread()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
std::array< Long, 3 > ParticleContainer::ByteSpread ( ) const

◆ Checkpoint() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Checkpoint ( const std::string &  dir,
const std::string &  name,
bool  is_checkpoint,
const Vector< std::string > &  real_comp_names = Vector< std::string >(),
const Vector< std::string > &  int_comp_names = Vector< std::string >() 
) const

Writes a particle checkpoint to file, suitable for restarting. This version allows the particle component names to be passed in. This overload exists for backwards comptability. The is_checkpoint parameter is ignored.

◆ Checkpoint() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Checkpoint ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names = Vector<std::string>(),
const Vector< std::string > &  int_comp_names = Vector<std::string>() 
) const
inline

Writes a particle checkpoint to file, suitable for restarting.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
nameThe name of the sub-directory for this particle type (i.e. "Tracer")
real_comp_namesvector of real component names, optional
int_comp_namesvector of int component names, optional

◆ CheckpointPost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::CheckpointPost ( )

◆ CheckpointPre()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::CheckpointPre ( )

◆ clearParticles()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::clearParticles ( )

Clear all the particles in this container. This does not free memory.

◆ copyParticles() [1/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::copyParticles ( const PCType &  other,
bool  local 
)

◆ copyParticles() [2/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::copyParticles ( const PCType &  other,
bool  local = false 
)

Copy particles from other to this ParticleContainer. Will clear all the particles from this container first. local controls whether or not to call Redistribute() after adding the particles.

Parameters
otherthe other pc to copy from
localwhether to call redistribute after

◆ copyParticles() [3/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo, std::enable_if_t<! std::is_integral< F >::value, int > bar>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::copyParticles ( const PCType &  other,
F &&  f,
bool  local 
)

◆ copyParticles() [4/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , class PCType , std::enable_if_t< IsParticleContainer< PCType >::value, int > foo = 0, std::enable_if_t<! std::is_integral< F >::value, int > bar = 0>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::copyParticles ( const PCType &  other,
F &&  f,
bool  local = false 
)

Copy particles from other to this ParticleContainer. Will clear all the particles from this container first. local controls whether or not to call Redistribute() after adding the particles.

This version conditionally copies based on a predicate function applied to each particle.

Template Parameters
callablethat takes a SuperParticle and returns a bool
Parameters
otherthe other pc to copy from
ffunction to apply to each particle as a filter
localwhether to call redistribute after

◆ correctCellVectors()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
virtual void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::correctCellVectors ( int  ,
int  ,
int  ,
const ParticleType  
)
inlineprivatevirtual

◆ CreateGhostParticles() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::CreateGhostParticles ( int  level,
int  ngrow,
AoS ghosts 
) const

Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in level-1.

Parameters
level
ngrow
ghosts

◆ CreateGhostParticles() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::CreateGhostParticles ( int  level,
int  ngrow,
ParticleTileType ghosts 
) const

Create ghost particles for a given level that are copies of particles near coarse->fine boundaries in level-1.

Parameters
level
ngrow
ghosts

◆ CreateVirtualParticles() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::CreateVirtualParticles ( int  level,
AoS virts 
) const

Creates virtual particles for a given level that represent in some capacity all particles at finer levels.

Parameters
level
virts

◆ CreateVirtualParticles() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::CreateVirtualParticles ( int  level,
ParticleTileType virts 
) const

Creates virtual particles for a given level that represent in some capacity all particles at finer levels.

Parameters
level
virts

◆ Define() [1/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Define ( const Geometry geom,
const DistributionMapping dmap,
const BoxArray ba 
)
inline

Define a default-constructed ParticleContainer using a ParGDB object. Single-level version.

Parameters
theGeometry object, which describes the problem domain
ADistributionMapping, which describes how the boxes are distributed onto MPI tasks
ABoxArray, which gives the set of grid boxes

◆ Define() [2/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Define ( const Vector< Geometry > &  geom,
const Vector< DistributionMapping > &  dmap,
const Vector< BoxArray > &  ba,
const Vector< int > &  rr 
)
inline

Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version.

Parameters
geomA Vector of Geometry objects, one for each level
dmapA Vector of DistributionMappings, one for each level
baA Vector of BoxArrays, one for each level
rrA Vector of integer refinement ratios, of size num_levels - 1. rr[n] gives the refinement ratio between levels n and n+1

◆ Define() [3/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Define ( const Vector< Geometry > &  geom,
const Vector< DistributionMapping > &  dmap,
const Vector< BoxArray > &  ba,
const Vector< IntVect > &  rr 
)
inline

Define a default-constructed ParticleContainer using a ParGDB object. Multi-level version.

Parameters
geomA Vector of Geometry objects, one for each level
dmapA Vector of DistributionMappings, one for each level
baA Vector of BoxArrays, one for each level
rrA Vector of integer refinement ratios, of size num_levels - 1. rr[n] gives the refinement ratio between levels n and n+1

◆ Define() [4/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Define ( ParGDBBase gdb)
inline

Define a default-constructed ParticleContainer using a ParGDB object. The container will track changes in the grid structure of the ParGDB automatically.

Parameters
gdbA pointer to a ParGDBBase, which contains pointers to the Geometry, DistributionMapping, and BoxArray objects that define the AMR hierarchy. Usually, this is generated by an AmrCore or AmrLevel object.

◆ DefineAndReturnParticleTile() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class Iterator >
ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::DefineAndReturnParticleTile ( int  lev,
const Iterator &  iter 
)
inline

Define and return the ParticleTile for level "lev", and Iterator "iter".

   Here, iter is either an MFIter or ParIter object pointing to the
   tile you want.

   If a ParticleTile at "grid" and "tile" has not been created yet,
   this function call will add it. This call will also allocate space
   for any runtime-added particle components.

   The ParticleLevel must already exist, meaning that the "resizeData()"
   method of this ParticleContainer has been called.

   Note that, when using a ParticleContainer that has been constructed
   with an AmrCore*, "resizeData()" must be called *after* the grids
   have been created, meaning after the call to AmrCore::InitFromScratch
   or AmrCore::InitFromCheckpoint has been made.

   \param lev  the level at which to get the particles
   \param iter MFIter or ParIter pointing to the tile to return

◆ DefineAndReturnParticleTile() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::DefineAndReturnParticleTile ( int  lev,
int  grid,
int  tile 
)
inline

Define and return the ParticleTile for level "lev", grid "grid" and tile "tile.".

   Here, grid and tile are integers that give the index and LocalTileIndex
   of the tile you want.

   If a ParticleTile at "grid" and "tile" has not been created yet,
   this function call will add it. This call will also allocate space
   for any runtime-added particle components.

   The ParticleLevel must already exist, meaning that the "resizeData()"
   method of this ParticleContainer has been called.

   Note that, when using a ParticleContainer that has been constructed
   with an AmrCore*, "resizeData()" must be called *after* the grids
   have been created, meaning after the call to AmrCore::InitFromScratch
   or AmrCore::InitFromCheckpoint has been made.

   \param lev  the level at which to get the particles
   \param grid the index of the grid at which to get the particles
   \param tile the LocalTileIndex of the tile at which to get the particles

◆ EnforcePeriodicWhere()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool ParticleContainer::EnforcePeriodicWhere ( ParticleType prt,
ParticleLocData pld,
int  lev_min = 0,
int  lev_max = -1,
int  local_grid = -1 
) const
protected

Checks whether the particle has crossed a periodic boundary in such a way that it is on levels lev_min and higher.

Parameters
prt
pld
lev_min
lev_max
local_grid

◆ GetLevelDirectoriesCreated()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetLevelDirectoriesCreated ( ) const
inline

◆ GetMaxNextIDPrePost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetMaxNextIDPrePost ( ) const
inline

◆ GetNParticlesPrePost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetNParticlesPrePost ( ) const
inline

◆ GetParticles() [1/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<ParticleLevel>& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetParticles ( )
inline

Return the underyling Vector (over AMR levels) of ParticleLevels. Non-const version.

◆ GetParticles() [2/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
const Vector<ParticleLevel>& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetParticles ( ) const
inline

Return the underyling Vector (over AMR levels) of ParticleLevels. Const version.

◆ GetParticles() [3/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleLevel& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetParticles ( int  lev)
inline

Return the ParticleLevel for level "lev". Non-const version.

   The ParticleLevel must already exist, meaning that the "resizeData()"
   method of this ParticleContainer has been called.

   Note that, when using a ParticleContainer that has been constructed
   with an AmrCore*, "resizeData()" must be called *after* the grids
   have been created, meaning after the call to AmrCore::InitFromScratch
   or AmrCore::InitFromCheckpoint has been made.

   \param lev the level at which to get the particles

◆ GetParticles() [4/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
const ParticleLevel& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetParticles ( int  lev) const
inline

Return the ParticleLevel for level "lev". Const version.

   The ParticleLevel must already exist, meaning that the "resizeData()"
   method of this ParticleContainer has been called.

   Note that, when using a ParticleContainer that has been constructed
   with an AmrCore*, "resizeData()" must be called *after* the grids
   have been created, meaning after the call to AmrCore::InitFromScratch
   or AmrCore::InitFromCheckpoint has been made.

   \param lev the level at which to get the particles

◆ GetUsePrePost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetUsePrePost ( ) const
inline

◆ GetUseUnlink()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::GetUseUnlink ( ) const
inline

◆ Increment()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::Increment ( MultiFab mf,
int  lev 
)

Adds the number of particles in each cell to the values currently located in the input MultiFab.

◆ IncrementWithTotal()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long ParticleContainer::IncrementWithTotal ( MultiFab mf,
int  level,
bool  local = false 
)

◆ Index() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename P >
IntVect amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Index ( const P &  p,
int  lev 
) const

◆ Index() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename P >
IntVect amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Index ( const P &  p,
int  lev 
) const

◆ InitFromAsciiFile()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitFromAsciiFile ( const std::string &  file,
int  extradata,
const IntVect Nrep = nullptr 
)

◆ InitFromBinaryFile()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitFromBinaryFile ( const std::string &  file,
int  extradata 
)

◆ InitFromBinaryMetaFile()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitFromBinaryMetaFile ( const std::string &  file,
int  extradata 
)

◆ Initialize()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::Initialize ( )
private

◆ InitNRandomPerCell()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitNRandomPerCell ( int  n_per_cell,
const ParticleInitData pdata 
)

This initializes the particle container with n_per_cell randomly distributed particles per cell, where the other particle data and and attributes are all constant. The cells on the coarsest level are used to generate the particle positions. The particle variable values are passed in through the pdata struct.

Parameters
n_per_cell
pdata

◆ InitOnePerCell()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitOnePerCell ( Real  x_off,
Real  y_off,
Real  z_off,
const ParticleInitData pdata 
)

This initializes the particle container with one particle per cell, where the other particle data and attributes are all contant. The coarsest level is used to generate the particle positions. The particle variable values are passed in through the pdata struct. The parameters x_off, y_off, and z_off represent offsets between 0 and 1 that show where inside the cells to place the particles. 0.5 means cell centered.

Parameters
x_off
y_off
z_off
pdata

◆ InitRandom()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitRandom ( Long  icount,
ULong  iseed,
const ParticleInitData pdata,
bool  serialize = false,
RealBox  bx = RealBox() 
)

This initializes the particle container with icount randomly distributed particles. If serialize is true, then the particles will all be generated on the IO Process, and the particle positions will be broadcast to all other process. If serialize is false, then the particle positions will be randomly generated in parallel, which each process using the random seed iseed + MyProc. The particles can be constrained to lie within the RealBox bx, if so desired. The default is the full domain.

Parameters
icount
iseed
mass
serialize
bx

◆ InitRandomPerBox()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InitRandomPerBox ( Long  icount,
ULong  iseed,
const ParticleInitData pdata 
)

This initializes the container with icount randomly distributed particles per box, using the random seed iseed. All the particles have the same data and attributes, which are passed using the pdata struct.

This routine is used when we want to replicate a box for a scaling study – within each box the distribution is random but the particle data is replicated across all boxes in the container. The boxes are assumed to be those on the coarsest level.

Parameters
icount
iseed
pdata

◆ Interpolate()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::Interpolate ( Vector< std::unique_ptr< MultiFab > > &  mesh_data,
int  lev_min,
int  lev_max 
)

◆ InterpolateSingleLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::InterpolateSingleLevel ( MultiFab mesh_data,
int  lev 
)

◆ locateParticle()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::locateParticle ( ParticleType p,
ParticleLocData pld,
int  lev_min,
int  lev_max,
int  nGrow,
int  local_grid = -1 
) const
private

◆ make_alike()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<template< class > class NewAllocator = amrex::DefaultAllocator>
ContainerLike<NewAllocator> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::make_alike ( ) const
inline

Create an empty particle container

This creates a new AMReX particle container type with same compile-time and run-time attributes. But, it can change its allocator. This is helpful when creating temporary particle buffers for filter operations and device-to-host copies.

Template Parameters
AllocatorAMReX allocator, e.g., amrex::PinnedArenaAllocator
Returns
an empty particle container

◆ NumberOfParticlesAtLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long ParticleContainer::NumberOfParticlesAtLevel ( int  level,
bool  only_valid = true,
bool  only_local = false 
) const

Returns # of particles at specified the level.

If "only_valid" is true it only counts valid particles.

Parameters
level
only_valid
only_local

◆ NumberOfParticlesInGrid()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector< Long > ParticleContainer::NumberOfParticlesInGrid ( int  level,
bool  only_valid = true,
bool  only_local = false 
) const

◆ NumIntComps()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NumIntComps ( ) const
inline

◆ numLocalTilesAtLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::numLocalTilesAtLevel ( int  lev) const
inline

The total number of tiles on this rank on this level.

◆ NumRealComps()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NumRealComps ( ) const
inline

◆ NumRuntimeIntComps()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NumRuntimeIntComps ( ) const
inline

◆ NumRuntimeRealComps()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NumRuntimeRealComps ( ) const
inline

◆ OK()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool ParticleContainer::OK ( int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0 
) const

OK checks that all particles are in the right places (for some value of right)

These flags are used to do proper checking for subcycling particles the default values are fine for non-subcycling methods

Parameters
lev_min
lev_max
nGrow

◆ operator=() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleContainer& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::operator= ( const ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator > &  )
delete

◆ operator=() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleContainer& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::operator= ( ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator > &&  )
defaultnoexcept

◆ particlePostLocate()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
virtual void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::particlePostLocate ( ParticleType ,
const ParticleLocData ,
const int   
)
inlineprivatevirtual

◆ ParticlesAt() [1/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class Iterator >
ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticlesAt ( int  lev,
const Iterator &  iter 
)
inline

Return the ParticleTile for level "lev" and Iterator "iter". Non-const version.

Here, iter is either an MFIter or ParIter object pointing to the tile you want.

This is a runtime error if a ParticleTile at "iter" has not been created yet.

The ParticleLevel must already exist, meaning that the "resizeData()" method of this ParticleContainer has been called.

Note that, when using a ParticleContainer that has been constructed with an AmrCore*, "resizeData()" must be called after the grids have been created, meaning after the call to AmrCore::InitFromScratch or AmrCore::InitFromCheckpoint has been made.

Parameters
levthe level at which to get the particles
iterMFIter or ParIter pointing to the tile to return

◆ ParticlesAt() [2/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class Iterator >
const ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticlesAt ( int  lev,
const Iterator &  iter 
) const
inline

Return the ParticleTile for level "lev" and Iterator "iter". Const version.

Here, iter is either an MFIter or ParIter object pointing to the tile you want.

This is a runtime error if a ParticleTile at "iter" has not been created yet.

The ParticleLevel must already exist, meaning that the "resizeData()" method of this ParticleContainer has been called.

Note that, when using a ParticleContainer that has been constructed with an AmrCore*, "resizeData()" must be called after the grids have been created, meaning after the call to AmrCore::InitFromScratch or AmrCore::InitFromCheckpoint has been made.

Parameters
levthe level at which to get the particles
iterMFIter or ParIter pointing to the tile to return

◆ ParticlesAt() [3/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticlesAt ( int  lev,
int  grid,
int  tile 
)
inline

Return the ParticleTile for level "lev", grid "grid" and tile "tile." Non-const version.

Here, grid and tile are integers that give the index and LocalTileIndex of the tile you want.

This is a runtime error if a ParticleTile at "grid" and "tile" has not been created yet.

The ParticleLevel must already exist, meaning that the "resizeData()" method of this ParticleContainer has been called.

Note that, when using a ParticleContainer that has been constructed with an AmrCore*, "resizeData()" must be called after the grids have been created, meaning after the call to AmrCore::InitFromScratch or AmrCore::InitFromCheckpoint has been made.

Parameters
levthe level at which to get the particles
gridthe index of the grid at which to get the particles
tilethe LocalTileIndex of the tile at which to get the particles

◆ ParticlesAt() [4/4]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
const ParticleTileType& amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticlesAt ( int  lev,
int  grid,
int  tile 
) const
inline

Return the ParticleTile for level "lev", grid "grid" and tile "tile." Const version.

Here, grid and tile are integers that give the index and LocalTileIndex of the tile you want.

This is a runtime error if a ParticleTile at "grid" and "tile" has not been created yet.

The ParticleLevel must already exist, meaning that the "resizeData()" method of this ParticleContainer has been called.

Note that, when using a ParticleContainer that has been constructed with an AmrCore*, "resizeData()" must be called after the grids have been created, meaning after the call to AmrCore::InitFromScratch or AmrCore::InitFromCheckpoint has been made.

Parameters
levthe level at which to get the particles
gridthe index of the grid at which to get the particles
tilethe LocalTileIndex of the tile at which to get the particles

◆ PeriodicShift()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool ParticleContainer::PeriodicShift ( ParticleType prt) const

Returns true if the particle was shifted.

Parameters
prt

◆ PrintCapacity()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
std::array< Long, 3 > ParticleContainer::PrintCapacity ( ) const

◆ ReadParticleRealData()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::ReadParticleRealData ( void *  data,
size_t  size,
std::istream &  is 
)

Read a contiguous chunk of real particle data from an istream.

Parameters
dataA pointer to the start of the buffer into which to read
sizeThe number of elements to read
osThe istream from which to read the data
rdA RealDescriptor describing the type of the floating point data

◆ ReadParticles() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class RTYPE >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ReadParticles ( int  cnt,
int  grd,
int  lev,
std::ifstream &  ifs,
int  finest_level_in_file,
bool  convert_ids 
)
protected

◆ ReadParticles() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class RTYPE >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ReadParticles ( int  cnt,
int  grd,
int  lev,
std::ifstream &  ifs,
int  finest_level_in_file,
bool  convert_ids 
)

◆ Redistribute()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::Redistribute ( int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0,
int  local = 0,
bool  remove_negative = true 
)

Redistribute puts all the particles back in the right places (for some value of right)

Assigns particles to the levels, grids, and tiles that contain their current positions. If periodic boundary conditions are used, those will be enforced here.

If Redistribute is called with default arguments, all particles will be placed on the finest level that covers their current positions.

The lev_min, lev_max, and nGrow flags are used to do proper checking for subcycling particles. The default values are fine for non-subcycling methods

The local flag controls whether this is local or global Redistribute. In a local Redistribute, particles can only have moved a certain distance since the last time Redistribute() was called. Thus, communication only needs to happen between neighboring ranks. In a global Redistribute, the particles can potentially go from any rank to any rank. This usually happens after initialization or when doing dynamic load balancing.

Parameters
lev_minThe minimum level consider. Particles on levels less than this will not be touched, and particles on finer levels will not be assigned to levels less than this, either. Default: 0.
lev_maxThe maximum level consider. Particles on levels greater than this will not be touched, and particles on coarser levels will not be assigned to levels greater than this, either. If negative, will use the finest level in the hierarchy. Default: -1.
nGrowIf particles are within nGrow cells of their current box, they will not moved. This is useful for subcycling methods, when fine level particles need to be redistributed but are not necessarily at the same time as those on the coarse level. Default: 0
localIf 0, this will be a non-local redistribute, meaning that particle can potentially go to any other box in the simulation. If > 0, this is the maximum number of cells a particle can have moved since the last Redistribute() call. Knowing this number allows an optimized MPI communication pattern to be used.

◆ RedistributeCPU()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::RedistributeCPU ( int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0,
int  local = 0,
bool  remove_negative = true 
)

◆ RedistributeGPU()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::RedistributeGPU ( int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0,
int  local = 0,
bool  remove_negative = true 
)

◆ RedistributeMPI()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::RedistributeMPI ( std::map< int, Vector< char > > &  not_ours,
int  lev_min = 0,
int  lev_max = 0,
int  nGrow = 0,
int  local = 0 
)
private

◆ RemoveParticlesAtLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::RemoveParticlesAtLevel ( int  level)

The Following methods are for managing Virtual and Ghost Particles.

Removes all particles at a given level

Parameters
level

◆ RemoveParticlesNotAtFinestLevel()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::RemoveParticlesNotAtFinestLevel ( )

◆ ReorderParticles() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class index_type >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ReorderParticles ( int  lev,
const MFIter mfi,
const index_type *  permutations 
)

Reorder particles on the tile given by lev and mfi using a the permutations array.

permutations is a pointer to an array on the GPU of size numParticles() with permutations[new index] = old index.

Parameters
lev
mfi
permutations

◆ ReorderParticles() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class index_type >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ReorderParticles ( int  lev,
const MFIter mfi,
const index_type *  permutations 
)

◆ reserveData()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::reserveData ( )
overridevirtual

Reimplemented from amrex::ParticleContainerBase.

◆ Reset()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
ParticleLocData ParticleContainer::Reset ( ParticleType prt,
bool  update,
bool  verbose = true,
ParticleLocData  pld = ParticleLocData() 
) const

Updates a particle's location (Where), tries to periodic shift any particles that have left the domain. May need work (see inline comments)

Parameters
prt
update
verbose
pld

◆ resizeData()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::resizeData ( )
overridevirtual

Reimplemented from amrex::ParticleContainerBase.

◆ ResizeRuntimeIntComp()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::ResizeRuntimeIntComp ( int  new_size,
bool  communicate 
)

◆ ResizeRuntimeRealComp()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::ResizeRuntimeRealComp ( int  new_size,
bool  communicate 
)

◆ Restart() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::Restart ( const std::string &  dir,
const std::string &  file 
)

Restart from checkpoint.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")

◆ Restart() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Restart ( const std::string &  dir,
const std::string &  file,
bool  is_checkpoint 
)

Older version, for backwards compatibility.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
is_checkpointWhether the particle id and cpu are included in the file.

◆ SetLevelDirectoriesCreated()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SetLevelDirectoriesCreated ( bool  tf)
inline

◆ SetParticleSize()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::SetParticleSize ( )
protected

◆ SetUsePrePost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SetUsePrePost ( bool  tf) const
inline

◆ SetUseUnlink()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::SetUseUnlink ( bool  tf) const
inline

◆ ShrinkToFit()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::ShrinkToFit ( )

◆ SortParticlesByBin()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::SortParticlesByBin ( IntVect  bin_size)

Sort the particles on each tile by groups of cells, given an IntVect bin_size.

If bin_size is the zero vector, this operation is a no-op.

◆ SortParticlesByCell()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::SortParticlesByCell ( )

Sort the particles on each tile by cell, using Fortran ordering.

◆ SortParticlesForDeposition()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::SortParticlesForDeposition ( IntVect  idx_type)

Sort particles on each tile such that particles adjacent in memory are likely to map to adjacent cells. This ordering can be beneficial for performance on GPU when deposition quantities onto a grid.

idx_type = {0, 0, 0}: Sort particles to a cell centered grid idx_type = {1, 1, 1}: Sort particles to a node centered grid idx_type = {2, 2, 2}: Compromise between a cell and node centered grid. This last option uses more memory than the fist two. Mixed versions are also possible.

Parameters
idx_type

◆ superParticleSize()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::superParticleSize ( ) const
inline

◆ TotalNumberOfParticles()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long ParticleContainer::TotalNumberOfParticles ( bool  only_valid = true,
bool  only_local = false 
) const

Returns # of particles at all levels.

If "only_valid" is true it only counts valid particles.

Parameters
only_valid
only_local

◆ Where() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename P >
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Where ( const P &  p,
ParticleLocData pld,
int  lev_min,
int  lev_max,
int  nGrow,
int  local_grid 
) const

◆ Where() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<typename P >
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::Where ( const P &  prt,
ParticleLocData pld,
int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0,
int  local_grid = -1 
) const
protected

Checks a particle's location on levels lev_min and higher. Returns false if the particle does not exist on that level. Only if lev_min == lev_max, nGrow can be > 0 (i.e., including nGrow ghost cells).

Parameters
prt
pld
lev_min
lev_max
nGrow
local_grid

◆ WriteAsciiFile()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WriteAsciiFile ( const std::string &  file)

◆ WriteBinaryParticleData() [1/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WriteBinaryParticleData ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f,
bool  is_checkpoint 
) const

◆ WriteBinaryParticleData() [2/2]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WriteBinaryParticleData ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f,
bool  is_checkpoint = false 
) const

Writes particle data to disk in the AMReX native format.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
nameThe name of the sub-directory for this particle type (i.e. "Tracer")
write_real_compfor each real component, whether or not we include that component in the file
write_int_compfor each integer component, whether or not we include that component in the file
real_comp_namesfor each real component, a name to label the data with
int_comp_namesfor each integer component, a name to label the data with
fcallable that returns whether a given particle should be written or not
is_checkpointwhether the data is written to a checkpoint or plotfile

◆ WriteParticleRealData()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WriteParticleRealData ( void *  data,
size_t  size,
std::ostream &  os 
) const

Write a contiguous chunk of real particle data to an ostream.

Parameters
dataA pointer to the start of the buffer to write
sizeThe number of elements to write
osThe ostream into which to write the data

◆ WriteParticles()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WriteParticles ( int  level,
std::ofstream &  ofs,
int  fnum,
Vector< int > &  which,
Vector< int > &  count,
Vector< Long > &  where,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::map< std::pair< int, int >, IntVector >> &  particle_io_flags,
bool  is_checkpoint 
) const

◆ WritePlotFile() [1/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFile ( const std::string &  dir,
const std::string &  name 
) const

This version of WritePlotFile writes all components and assigns component names.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
nameThe name of the sub-directory for this particle type (i.e. "Tracer")

◆ WritePlotFile() [2/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp 
) const

This version of WritePlotFile assigns component names, but allows the user to pass in a vector of ints that toggle on / off the writing of specific components.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
write_real_compfor each real component, whether to include that comp in the file
write_int_compfor each integer component, whether to include that comp in the file

◆ WritePlotFile() [3/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names 
) const

This is the most general version of WritePlotFile, which takes component names and flags for whether to write each variable as components. Note that the user should pass in vectors containing names of all the components, whether they are written or not.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
write_real_compfor each real component, whether to include that comp in the file
write_int_compfor each integer component, whether to include that comp in the file
real_comp_namesfor each real component, a name to label the data with
int_comp_namesfor each integer component, a name to label the data with

◆ WritePlotFile() [4/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f 
) const

This is the most general version of WritePlotFile, which takes component names and flags for whether to write each variable as components. Note that the user should pass in vectors containing names of all the components, whether they are written or not.

This version also lets you pass in a functor to toggle whether each particle gets output.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
write_real_compfor each real component, whether to include that comp in the file
write_int_compfor each integer component, whether to include that comp in the file
real_comp_namesfor each real component, a name to label the data with
int_comp_namesfor each integer component, a name to label the data with
fcallable that returns whether or not to write each particle

◆ WritePlotFile() [5/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f 
) const

◆ WritePlotFile() [6/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
F &&  f 
) const

This version of WritePlotFile assigns component names, but allows the user to pass in a vector of ints that toggle on / off the writing of specific components.

This version also lets you pass in a functor to toggle whether each particle gets output.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
write_real_compfor each real component, whether to include that comp in the file
write_int_compfor each integer component, whether to include that comp in the file
fcallable that returns whether or not to write each particle

◆ WritePlotFile() [7/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< int > &  write_real_comp,
const Vector< int > &  write_int_comp,
F &&  f 
) const

◆ WritePlotFile() [8/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names 
) const

This version of WritePlotFile writes all components and allows the user to specify the names of the components.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
real_comp_namesfor each real component, a name to label the data with

◆ WritePlotFile() [9/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names 
) const

This version of WritePlotFile writes all components and allows the user to specify the names of the components.

Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
real_comp_namesfor each real component, a name to label the data with
int_comp_namesfor each integer component, a name to label the data with

◆ WritePlotFile() [10/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f 
) const

This version of WritePlotFile writes all components and allows the user to specify the names of the components.

This version also lets you pass in a functor to toggle whether each particle gets output.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
real_comp_namesfor each real component, a name to label the data with
int_comp_namesfor each integer component, a name to label the data with
fcallable that returns whether or not to write each particle

◆ WritePlotFile() [11/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names,
const Vector< std::string > &  int_comp_names,
F &&  f 
) const

◆ WritePlotFile() [12/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string >>::value >::type * = nullptr>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names,
F &&  f 
) const

This version of WritePlotFile writes all components and allows the user to specify the names of the components.

This version also lets you pass in a functor to toggle whether each particle gets output.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
fileThe name of the sub-directory for this particle type (i.e. "Tracer")
real_comp_namesfor each real component, a name to label the data with
fcallable that returns whether or not to write each particle

◆ WritePlotFile() [13/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string >>::value >::type * >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
const Vector< std::string > &  real_comp_names,
F &&  f 
) const

◆ WritePlotFile() [14/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string > & >::value >::type * = nullptr>
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
F &&  f 
) const

This version of WritePlotFile writes all components and assigns component names.

This version also lets you pass in a functor to toggle whether each particle gets output.

Template Parameters
Ffunction type
Parameters
dirThe base directory into which to write (i.e. "plt00000")
nameThe name of the sub-directory for this particle type (i.e. "Tracer")
fcallable that returns whether or not to write each particle

◆ WritePlotFile() [15/15]

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
template<class F , typename std::enable_if<!std::is_same< F, Vector< std::string > & >::value >::type * >
void amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::WritePlotFile ( const std::string &  dir,
const std::string &  name,
F &&  f 
) const

◆ WritePlotFilePost()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFilePost ( )

◆ WritePlotFilePre()

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
void ParticleContainer::WritePlotFilePre ( )

Friends And Related Function Documentation

◆ ParIterBase< false, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
friend class ParIterBase< false, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
friend

◆ ParIterBase< true, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
friend class ParIterBase< true, NStructReal, NStructInt, NArrayReal, NArrayInt, Allocator >
friend

Member Data Documentation

◆ countPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<Vector<int> > amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::countPrePost
mutable

◆ doUnlink

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::doUnlink
mutable

◆ filePrefixPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<std::string> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::filePrefixPrePost
mutable

◆ h_redistribute_int_comp

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<int> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::h_redistribute_int_comp

◆ h_redistribute_real_comp

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<int> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::h_redistribute_real_comp

◆ HdrFileNamePrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
std::string amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::HdrFileNamePrePost
mutable

◆ levelDirectoriesCreated

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::levelDirectoriesCreated
mutable

Variables for i/o optimization saved for pre and post checkpoint.

◆ m_bins

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
DenseBins<ParticleType> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::m_bins
protected

◆ m_num_runtime_int

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::m_num_runtime_int {0}
private

◆ m_num_runtime_real

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::m_num_runtime_real {0}
private

◆ m_particles

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<ParticleLevel> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::m_particles
private

◆ m_runtime_comps_defined

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::m_runtime_comps_defined {false}
private

◆ maxnextidPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::maxnextidPrePost

◆ NArrayInt

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
constexpr int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NArrayInt = T_NArrayInt
staticconstexpr

Number of extra integer components stored in struct-of-array form.

◆ NArrayReal

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
constexpr int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NArrayReal = T_NArrayReal
staticconstexpr

Number of extra Real components stored in struct-of-array form.

◆ nOutFilesPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::nOutFilesPrePost
mutable

◆ nParticlesAtLevelPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<Long> amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::nParticlesAtLevelPrePost

◆ nparticlesPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Long amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::nparticlesPrePost

◆ NStructInt

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
constexpr int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NStructInt = T_NStructInt
staticconstexpr

Number of extra integer components in the particle struct.

◆ NStructReal

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
constexpr int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::NStructReal = T_NStructReal
staticconstexpr

Number of extra Real components in the particle struct.

◆ num_int_comm_comps

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::num_int_comm_comps
private

◆ num_real_comm_comps

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
int amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::num_real_comm_comps
private

◆ particle_size

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
size_t amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::particle_size
private

◆ ParticleRealDescriptor

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
RealDescriptor amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::ParticleRealDescriptor = FPC::Native64RealDescriptor()

◆ superparticle_size

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
size_t amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::superparticle_size
private

◆ usePrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
bool amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::usePrePost
mutable

◆ wherePrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<Vector<Long> > amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::wherePrePost
mutable

◆ whichPrePost

template<int T_NStructReal, int T_NStructInt = 0, int T_NArrayReal = 0, int T_NArrayInt = 0, template< class > class Allocator = DefaultAllocator>
Vector<Vector<int> > amrex::ParticleContainer< T_NStructReal, T_NStructInt, T_NArrayReal, T_NArrayInt, Allocator >::whichPrePost
mutable

The documentation for this class was generated from the following files: