MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ParticleSet Class Reference

ParticleSet initializes and manages data associated with particles. More...

#include <particleset.hpp>

Collaboration diagram for mfem::ParticleSet:
[legend]

Public Types

using IDType = unsigned long long
 

Public Member Functions

 ParticleSet (int num_particles, int dim, Ordering::Type coords_ordering=Ordering::byVDIM)
 Construct a serial ParticleSet.
 
 ParticleSet (int num_particles, int dim, const Array< int > &field_vdims, int num_tags, Ordering::Type all_ordering=Ordering::byVDIM)
 Construct a serial ParticleSet with specified fields and tags at construction.
 
 ParticleSet (int num_particles, int dim, const Array< int > &field_vdims, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_, Ordering::Type all_ordering=Ordering::byVDIM)
 Construct a serial ParticleSet with specified fields and tags at construction, with names.
 
 ParticleSet (int num_particles, int dim, Ordering::Type coords_ordering, const Array< int > &field_vdims, const Array< Ordering::Type > &field_orderings, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_)
 Comprehensive serial constructor of ParticleSet.
 
 ParticleSet (MPI_Comm comm_, int rank_num_particles, int dim, Ordering::Type coords_ordering=Ordering::byVDIM)
 Construct a parallel ParticleSet.
 
 ParticleSet (MPI_Comm comm_, int rank_num_particles, int dim, const Array< int > &field_vdims, int num_tags, Ordering::Type all_ordering=Ordering::byVDIM)
 Construct a parallel ParticleSet with specified fields and tags at construction.
 
 ParticleSet (MPI_Comm comm_, int rank_num_particles, int dim, const Array< int > &field_vdims, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_, Ordering::Type all_ordering=Ordering::byVDIM)
 Construct a parallel ParticleSet with specified fields and tags at construction, with names (for PrintCSV()).
 
 ParticleSet (MPI_Comm comm_, int rank_num_particles, int dim, Ordering::Type coords_ordering, const Array< int > &field_vdims, const Array< Ordering::Type > &field_orderings, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_)
 Comprehensive parallel constructor of ParticleSet.
 
MPI_Comm GetComm () const
 Get the MPI communicator for this ParticleSet.
 
IDType GetGlobalNParticles () const
 Get the global number of active particles across all ranks.
 
int GetDim () const
 Get the spatial dimension.
 
const Array< IDType > & GetIDs () const
 Get the global IDs of the active particles owned by this ParticleSet.
 
int AddField (int vdim, Ordering::Type field_ordering=Ordering::byVDIM, const char *field_name=nullptr)
 Add a field to the ParticleSet.
 
int AddNamedField (int vdim, const char *field_name, Ordering::Type field_ordering=Ordering::byVDIM)
 Add a field to the ParticleSet.
 
int AddTag (const char *tag_name=nullptr)
 Add a tag to the ParticleSet.
 
void Reserve (int res)
 Reserve memory for res particles.
 
int GetNParticles () const
 Get the number of active particles currently held by this ParticleSet.
 
int GetNFields () const
 Get the number of fields registered to particles.
 
const Array< int > GetFieldVDims () const
 Get an Array<int> of the field vector-dimensions registered to particles.
 
int FieldVDim (int f) const
 Get Field vector-dimension.
 
int GetNTags () const
 Get the number of tags registered to particles.
 
void AddParticle (const Particle &p)
 Add a particle using Particle .
 
void AddParticles (int num_particles, Array< int > *new_indices=nullptr)
 Add num_particles particles, and optionally get the local indices of new particles in new_indices .
 
void RemoveParticles (const Array< int > &list)
 Remove particle data specified by list of particle indices.
 
ParticleVectorCoords ()
 Get a reference to the coordinates ParticleVector.
 
const ParticleVectorCoords () const
 Get a const reference to the coordinates ParticleVector.
 
ParticleVectorField (int f)
 Get a reference to field f 's ParticleVector.
 
const ParticleVectorField (int f) const
 Get a const reference to field f 's ParticleVector.
 
Array< int > & Tag (int t)
 Get a reference to tag t 's Array<int>.
 
const Array< int > & Tag (int t) const
 Get a const reference to tag t 's Array<int>.
 
Particle GetParticle (int i) const
 Get new Particle object with copy of data associated with particle i .
 
Particle GetParticleRef (int i)
 Get Particle object whose members reference the actual data associated with particle i in this ParticleSet.
 
bool IsParticleRefValid () const
 Determine if GetParticleRef is valid.
 
void SetParticle (int i, const Particle &p)
 Set data for particle at index i with data from provided particle p.
 
void PrintCSV (const char *fname, int precision=16)
 Print all particle data to a comma-delimited CSV file.
 
void PrintCSV (const char *fname, const Array< int > &field_idxs, const Array< int > &tag_idxs, int precision=16)
 Print only particle field and tags given by field_idxs and tag_idxs respectively to a CSV file.
 
void Redistribute (const Array< unsigned int > &rank_list)
 Redistribute particle data to rank_list.
 
 ~ParticleSet ()
 Destructor.
 
 ParticleSet (const ParticleSet &)=delete
 
ParticleSetoperator= (const ParticleSet &)=delete
 

Protected Member Functions

void AddParticles (const Array< IDType > &new_ids, Array< int > *new_indices=nullptr)
 Add particles with global identifiers new_ids and optionally get the local indices of new particles in new_indices .
 
void UpdateID (int local_idx, IDType new_global_id)
 Update global ID of a particle.
 
Particle CreateParticle () const
 Create a Particle object with the same spatial dimension, number of fields and field vdims, and number of tags as this ParticleSet.
 
void WriteToFile (const char *fname, const std::stringstream &ss_header, const std::stringstream &ss_data)
 Write string in ss_header , followed by ss_data , to a single file; compatible in parallel.
 
bool IsValidParticle (const Particle &p) const
 Check if a particle could belong in this ParticleSet by comparing field and tag dimension.
 
 ParticleSet (int id_stride_, IDType id_counter_, int num_particles, int dim, Ordering::Type coords_ordering, const Array< int > &field_vdims, const Array< Ordering::Type > &field_orderings, const Array< const char * > &field_names_, int num_tags, const Array< const char * > &tag_names_)
 Hidden main constructor of ParticleSet.
 

Protected Attributes

const int id_stride
 Stride for IDs (used internally when new particles are added).
 
IDType id_counter
 Current globally unique ID to be assigned to the next particle added.
 
Array< IDTypeids
 Global unique IDs of particles owned by this rank.
 
ParticleVector coords
 Spatial coordinates of particles owned by this rank.
 
std::vector< std::unique_ptr< ParticleVector > > fields
 All particle fields for particles owned by this rank.
 
std::vector< std::unique_ptr< Array< int > > > tags
 All particle tags for particles owned by this rank.
 
std::vector< std::string > field_names
 Field names, to be written when PrintCSV() is called.
 
std::vector< std::string > tag_names
 Tag names, to be written when PrintCSV() is called.
 
MPI_Comm comm
 
struct gslib::crystal * cr = nullptr
 
struct gslib::comm * gsl_comm = nullptr
 

Detailed Description

ParticleSet initializes and manages data associated with particles.

Particles are inherently initialized to have a position and an ID, and optionally can have any number of Vector (of arbitrary vdim) and scalar integer data in the form of fields and tags respectively. All particle data are internally stored in a Struct-of-Arrays fashion, as elaborated on below.

Coordinates:
All particle coordinates are stored in a ParticleVector with vector dimension equal to the spatial dimension, ordered either byNODES or byVDIM. The ParticleVector coords contains the coordinates of all particles.
IDs:
Each particle is assigned a unique global ID of type IDType. In parallel, IDs are initialized starting with rank and striding by size. The IDs of all particles owned by this rank are stored in ids.
Fields:
Fields represent scalar or vector real_t data to be associated with each particles, such as mass, momentum, or moment. For a given field, all particle data is stored in a single ParticleVector with a given vector dimension (1 for scalar data) and Ordering::Type (byNODES or byVDIM). The unique_ptrs to all the ParticleVectors are stored in the std::vector fields.
Tags:
Tags represent integers associated with each particle. For a given tag, all particle data are stored in a single Array<int>. The unique_ptrs to all the Array<int> is stored in the std::vector tags.
Names:
Each field and tag can optionally be given a name (string) to be used when printing particle data in CSV format using PrintCSV(). The names of all fields and tags are stored in the std::vectors field_names and tag_names, respectively.
Note
We assume that all particles in a ParticleSet have the same number of fields and tags.

Following the example in the Particle class, we will use the particles below to illustrate the data layout for coords, ids, fields, tags, field_names, and tag_names. In each case, the name of the field and tag is enclosed in '...' for clarity. Additionally, we assume for this example that the particle coordinates and the 'vel' field are ordered byVDIM in their respective ParticleVector.

Particle_0: id = id0, coords = (x0, y0),
fields = {'mass'=m0, 'vel' = (vx0, vy0)},
tags = {'type'=t0, 'color'=c0}
Particle_1: id = id1, coords = (x1, y1),
fields = {'mass'=m1, 'vel' = (vx1, vy1)},
tags = {'type'=t1, 'color'=c1}
Particle_2: id = id2, coords = (x2, y2),
fields = {'mass'=m2, 'vel' = (vx2, vy2)},
tags = {'type'=t2, 'color'=c2}
ParticleVector coords
Spatial coordinates of particles owned by this rank.
std::vector< std::unique_ptr< ParticleVector > > fields
All particle fields for particles owned by this rank.
std::vector< std::unique_ptr< Array< int > > > tags
All particle tags for particles owned by this rank.

Definition at line 248 of file particleset.hpp.

Member Typedef Documentation

◆ IDType

using mfem::ParticleSet::IDType = unsigned long long

Definition at line 251 of file particleset.hpp.

Constructor & Destructor Documentation

◆ ParticleSet() [1/10]

mfem::ParticleSet::ParticleSet ( int id_stride_,
IDType id_counter_,
int num_particles,
int dim,
Ordering::Type coords_ordering,
const Array< int > & field_vdims,
const Array< Ordering::Type > & field_orderings,
const Array< const char * > & field_names_,
int num_tags,
const Array< const char * > & tag_names_ )
protected

Hidden main constructor of ParticleSet.

Parameters
[in]id_stride_ID stride.
[in]id_counter_Starting ID counter.
[in]num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]coords_orderingOrdering of coordinates
[in]field_vdimsArray of field vector dimensions
[in]field_orderingsArray of field ordering types.
[in]field_names_Array of field names.
[in]num_tagsNumber of tags to register.
[in]tag_names_Array of tag names.

Definition at line 525 of file particleset.cpp.

◆ ParticleSet() [2/10]

mfem::ParticleSet::ParticleSet ( int num_particles,
int dim,
Ordering::Type coords_ordering = Ordering::byVDIM )

Construct a serial ParticleSet.

Parameters
[in]num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]coords_orderingOrdering of coordinates.

Definition at line 582 of file particleset.cpp.

◆ ParticleSet() [3/10]

mfem::ParticleSet::ParticleSet ( int num_particles,
int dim,
const Array< int > & field_vdims,
int num_tags,
Ordering::Type all_ordering = Ordering::byVDIM )

Construct a serial ParticleSet with specified fields and tags at construction.

Parameters
[in]num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]field_vdimsArray of field vector dimensions.
[in]num_tagsNumber of tags to register.
[in]all_ordering(Optional) Ordering of coordinates and field ParticleVector.

Definition at line 591 of file particleset.cpp.

◆ ParticleSet() [4/10]

mfem::ParticleSet::ParticleSet ( int num_particles,
int dim,
const Array< int > & field_vdims,
const Array< const char * > & field_names_,
int num_tags,
const Array< const char * > & tag_names_,
Ordering::Type all_ordering = Ordering::byVDIM )

Construct a serial ParticleSet with specified fields and tags at construction, with names.

Parameters
[in]num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]field_vdimsArray of field vector dimensions.
[in]field_names_Array of field names.
[in]num_tagsNumber of tags to register.
[in]tag_names_Array of tag names.
[in]all_ordering(Optional) Ordering of coordinates and field ParticleVector.

Definition at line 601 of file particleset.cpp.

◆ ParticleSet() [5/10]

mfem::ParticleSet::ParticleSet ( int num_particles,
int dim,
Ordering::Type coords_ordering,
const Array< int > & field_vdims,
const Array< Ordering::Type > & field_orderings,
const Array< const char * > & field_names_,
int num_tags,
const Array< const char * > & tag_names_ )

Comprehensive serial constructor of ParticleSet.

Parameters
[in]num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]coords_orderingOrdering of coordinates.
[in]field_vdimsArray of field vector dimensions.
[in]field_orderingsArray of field ordering types.
[in]field_names_Array of field names.
[in]num_tagsNumber of tags to register.
[in]tag_names_Array of tag names.

Definition at line 614 of file particleset.cpp.

◆ ParticleSet() [6/10]

mfem::ParticleSet::ParticleSet ( MPI_Comm comm_,
int rank_num_particles,
int dim,
Ordering::Type coords_ordering = Ordering::byVDIM )

Construct a parallel ParticleSet.

Parameters
[in]comm_MPI communicator.
[in]rank_num_particlesNumber of particles to initialize.
[in]dimParticle spatial dimension.
[in]coords_ordering(Optional) Ordering of coordinates.

Definition at line 629 of file particleset.cpp.

◆ ParticleSet() [7/10]

mfem::ParticleSet::ParticleSet ( MPI_Comm comm_,
int rank_num_particles,
int dim,
const Array< int > & field_vdims,
int num_tags,
Ordering::Type all_ordering = Ordering::byVDIM )

Construct a parallel ParticleSet with specified fields and tags at construction.

Parameters
[in]comm_MPI communicator.
[in]rank_num_particles# of particles to initialize on this rank.
[in]dimParticle spatial dimension.
[in]field_vdimsArray of field vector dimensions.
[in]num_tagsNumber of tags to register.
[in]all_ordering(Optional) Ordering of coordinates and field ParticleVector.

Definition at line 638 of file particleset.cpp.

◆ ParticleSet() [8/10]

mfem::ParticleSet::ParticleSet ( MPI_Comm comm_,
int rank_num_particles,
int dim,
const Array< int > & field_vdims,
const Array< const char * > & field_names_,
int num_tags,
const Array< const char * > & tag_names_,
Ordering::Type all_ordering = Ordering::byVDIM )

Construct a parallel ParticleSet with specified fields and tags at construction, with names (for PrintCSV()).

Parameters
[in]comm_MPI communicator.
[in]rank_num_particles# of particles to initialize on this rank.
[in]dimParticle spatial dimension.
[in]field_vdimsArray of field vector dimension.
[in]field_names_Array of field names.
[in]num_tagsNumber of tags to register.
[in]tag_names_Array of tag names.
[in]all_ordering(Optional) Ordering of coordinates and field ParticleVector.

Definition at line 649 of file particleset.cpp.

◆ ParticleSet() [9/10]

mfem::ParticleSet::ParticleSet ( MPI_Comm comm_,
int rank_num_particles,
int dim,
Ordering::Type coords_ordering,
const Array< int > & field_vdims,
const Array< Ordering::Type > & field_orderings,
const Array< const char * > & field_names_,
int num_tags,
const Array< const char * > & tag_names_ )

Comprehensive parallel constructor of ParticleSet.

Parameters
[in]comm_MPI communicator.
[in]rank_num_particles# of particles to initialize on this rank.
[in]dimParticle spatial dimension.
[in]coords_orderingOrdering of coordinates.
[in]field_vdimsArray of field vector dimensions.
[in]field_orderingsArray of field ordering types.
[in]field_names_Array of field names.
[in]num_tagsNumber of tags to register.
[in]tag_names_Array of tag names.

Definition at line 662 of file particleset.cpp.

◆ ~ParticleSet()

mfem::ParticleSet::~ParticleSet ( )

Destructor.

Definition at line 933 of file particleset.cpp.

◆ ParticleSet() [10/10]

mfem::ParticleSet::ParticleSet ( const ParticleSet & )
delete

Member Function Documentation

◆ AddField()

int mfem::ParticleSet::AddField ( int vdim,
Ordering::Type field_ordering = Ordering::byVDIM,
const char * field_name = nullptr )

Add a field to the ParticleSet.

Parameters
[in]vdimVector dimension of the field.
[in]field_ordering(Optional) Ordering::Type of the field.
[in]field_name(Optional) Name of the field.
Returns
Index of the newly-added field.

Definition at line 698 of file particleset.cpp.

◆ AddNamedField()

int mfem::ParticleSet::AddNamedField ( int vdim,
const char * field_name,
Ordering::Type field_ordering = Ordering::byVDIM )
inline

Add a field to the ParticleSet.

Same as AddField() but with different parameter order for convenience

Definition at line 558 of file particleset.hpp.

◆ AddParticle()

void mfem::ParticleSet::AddParticle ( const Particle & p)

Add a particle using Particle .

Definition at line 726 of file particleset.cpp.

◆ AddParticles() [1/2]

void mfem::ParticleSet::AddParticles ( const Array< IDType > & new_ids,
Array< int > * new_indices = nullptr )
protected

Add particles with global identifiers new_ids and optionally get the local indices of new particles in new_indices .

Note the data of new particles is uninitialized and must be set.

Definition at line 211 of file particleset.cpp.

◆ AddParticles() [2/2]

void mfem::ParticleSet::AddParticles ( int num_particles,
Array< int > * new_indices = nullptr )

Add num_particles particles, and optionally get the local indices of new particles in new_indices .

The data of new particles is uninitialized and must be set.

Definition at line 741 of file particleset.cpp.

◆ AddTag()

int mfem::ParticleSet::AddTag ( const char * tag_name = nullptr)

Add a tag to the ParticleSet.

Parameters
[in]tag_name(Optional) Name of the tag.
Returns
Index of the newly-added tag.

Definition at line 713 of file particleset.cpp.

◆ Coords() [1/2]

ParticleVector & mfem::ParticleSet::Coords ( )
inline

Get a reference to the coordinates ParticleVector.

Definition at line 606 of file particleset.hpp.

◆ Coords() [2/2]

const ParticleVector & mfem::ParticleSet::Coords ( ) const
inline

Get a const reference to the coordinates ParticleVector.

Definition at line 609 of file particleset.hpp.

◆ CreateParticle()

Particle mfem::ParticleSet::CreateParticle ( ) const
protected

Create a Particle object with the same spatial dimension, number of fields and field vdims, and number of tags as this ParticleSet.

Definition at line 469 of file particleset.cpp.

◆ Field() [1/2]

ParticleVector & mfem::ParticleSet::Field ( int f)
inline

Get a reference to field f 's ParticleVector.

Definition at line 612 of file particleset.hpp.

◆ Field() [2/2]

const ParticleVector & mfem::ParticleSet::Field ( int f) const
inline

Get a const reference to field f 's ParticleVector.

Definition at line 615 of file particleset.hpp.

◆ FieldVDim()

int mfem::ParticleSet::FieldVDim ( int f) const
inline

Get Field vector-dimension.

Definition at line 586 of file particleset.hpp.

◆ GetComm()

MPI_Comm mfem::ParticleSet::GetComm ( ) const
inline

Get the MPI communicator for this ParticleSet.

Definition at line 531 of file particleset.hpp.

◆ GetDim()

int mfem::ParticleSet::GetDim ( ) const
inline

Get the spatial dimension.

Definition at line 537 of file particleset.hpp.

◆ GetFieldVDims()

const Array< int > mfem::ParticleSet::GetFieldVDims ( ) const

Get an Array<int> of the field vector-dimensions registered to particles.

Definition at line 201 of file particleset.cpp.

◆ GetGlobalNParticles()

ParticleSet::IDType mfem::ParticleSet::GetGlobalNParticles ( ) const

Get the global number of active particles across all ranks.

Definition at line 688 of file particleset.cpp.

◆ GetIDs()

const Array< IDType > & mfem::ParticleSet::GetIDs ( ) const
inline

Get the global IDs of the active particles owned by this ParticleSet.

Definition at line 540 of file particleset.hpp.

◆ GetNFields()

int mfem::ParticleSet::GetNFields ( ) const
inline

Get the number of fields registered to particles.

Definition at line 580 of file particleset.hpp.

◆ GetNParticles()

int mfem::ParticleSet::GetNParticles ( ) const
inline

Get the number of active particles currently held by this ParticleSet.

Definition at line 577 of file particleset.hpp.

◆ GetNTags()

int mfem::ParticleSet::GetNTags ( ) const
inline

Get the number of tags registered to particles.

Definition at line 589 of file particleset.hpp.

◆ GetParticle()

Particle mfem::ParticleSet::GetParticle ( int i) const

Get new Particle object with copy of data associated with particle i .

Definition at line 772 of file particleset.cpp.

◆ GetParticleRef()

Particle mfem::ParticleSet::GetParticleRef ( int i)

Get Particle object whose members reference the actual data associated with particle i in this ParticleSet.

See also
IsParticleRefValid for when this method can be used.
Warning
If particles are added, removed, or redistributed after invoking this, the returned Particle member references may be invalidated.

Definition at line 807 of file particleset.cpp.

◆ IsParticleRefValid()

bool mfem::ParticleSet::IsParticleRefValid ( ) const

Determine if GetParticleRef is valid.

If coordinates and all fields are ordered byVDIM, then returns true. Otherwise, false.

Definition at line 791 of file particleset.cpp.

◆ IsValidParticle()

bool mfem::ParticleSet::IsValidParticle ( const Particle & p) const
protected

Check if a particle could belong in this ParticleSet by comparing field and tag dimension.

Definition at line 556 of file particleset.cpp.

◆ operator=()

ParticleSet & mfem::ParticleSet::operator= ( const ParticleSet & )
delete

◆ PrintCSV() [1/2]

void mfem::ParticleSet::PrintCSV ( const char * fname,
const Array< int > & field_idxs,
const Array< int > & tag_idxs,
int precision = 16 )

Print only particle field and tags given by field_idxs and tag_idxs respectively to a CSV file.

Definition at line 863 of file particleset.cpp.

◆ PrintCSV() [2/2]

void mfem::ParticleSet::PrintCSV ( const char * fname,
int precision = 16 )

Print all particle data to a comma-delimited CSV file.

The first row contains the header. We include the particle ID, owning rank (in parallel), coordinates, followed by all fields and tags.

The output can be visualized in Paraview by loading the csv files, and applying the "Table To Points" filter.

Definition at line 846 of file particleset.cpp.

◆ Redistribute()

DO_NOT_DOCUMENT void mfem::ParticleSet::Redistribute ( const Array< unsigned int > & rank_list)

Redistribute particle data to rank_list.

Parameters
[in]rank_listArray of size GetNParticles() denoting ultimate destination of particle data. Index = this rank means no data is moved.

Definition at line 436 of file particleset.cpp.

◆ RemoveParticles()

void mfem::ParticleSet::RemoveParticles ( const Array< int > & list)

Remove particle data specified by list of particle indices.

Definition at line 753 of file particleset.cpp.

◆ Reserve()

void mfem::ParticleSet::Reserve ( int res)

Reserve memory for res particles.

Can help to avoid re-allocation for adding + removing particles.

Definition at line 182 of file particleset.cpp.

◆ SetParticle()

void mfem::ParticleSet::SetParticle ( int i,
const Particle & p )

Set data for particle at index i with data from provided particle p.

Definition at line 828 of file particleset.cpp.

◆ Tag() [1/2]

Array< int > & mfem::ParticleSet::Tag ( int t)
inline

Get a reference to tag t 's Array<int>.

Definition at line 618 of file particleset.hpp.

◆ Tag() [2/2]

const Array< int > & mfem::ParticleSet::Tag ( int t) const
inline

Get a const reference to tag t 's Array<int>.

Definition at line 621 of file particleset.hpp.

◆ UpdateID()

void mfem::ParticleSet::UpdateID ( int local_idx,
IDType new_global_id )
inlineprotected

Update global ID of a particle.

This method updates the global ID of the particle at given local index after Redistribute().

Note
This method must be used very carefully as it updates global ID of a particle.

Definition at line 371 of file particleset.hpp.

◆ WriteToFile()

void mfem::ParticleSet::WriteToFile ( const char * fname,
const std::stringstream & ss_header,
const std::stringstream & ss_data )
protected

Write string in ss_header , followed by ss_data , to a single file; compatible in parallel.

Definition at line 474 of file particleset.cpp.

Member Data Documentation

◆ comm

MPI_Comm mfem::ParticleSet::comm
protected

Definition at line 335 of file particleset.hpp.

◆ coords

ParticleVector mfem::ParticleSet::coords
protected

Spatial coordinates of particles owned by this rank.

For the sample_particleset_data, coords would be coords=(x0,y0,x1,y1,x2,y2) assuming coords ordering is byVDIM.

Definition at line 294 of file particleset.hpp.

◆ cr

struct gslib::crystal* mfem::ParticleSet::cr = nullptr
protected

Definition at line 339 of file particleset.hpp.

◆ field_names

std::vector<std::string> mfem::ParticleSet::field_names
protected

Field names, to be written when PrintCSV() is called.

For the sample_particleset_data, field_names would be field_names[0]='mass', field_names[1]='vel'.

Definition at line 316 of file particleset.hpp.

◆ fields

std::vector<std::unique_ptr<ParticleVector> > mfem::ParticleSet::fields
protected

All particle fields for particles owned by this rank.

For the sample_particleset_data, fields would be *fields[0]=(m0,m1,m2), *fields[1]=(vx0,vy0,vx1,vy1,vx2,vy2) assuming fields[1] ordering is byVDIM.

Definition at line 302 of file particleset.hpp.

◆ gsl_comm

struct gslib::comm* mfem::ParticleSet::gsl_comm = nullptr
protected

Definition at line 340 of file particleset.hpp.

◆ id_counter

IDType mfem::ParticleSet::id_counter
protected

Current globally unique ID to be assigned to the next particle added.

In parallel, this starts locally as the rank and increments with id_stride, ensuring a global unique identifier whenever a particle is added.

Definition at line 280 of file particleset.hpp.

◆ id_stride

const int mfem::ParticleSet::id_stride
protected

Stride for IDs (used internally when new particles are added).

In parallel, this defaults to the number of MPI ranks.

Definition at line 273 of file particleset.hpp.

◆ ids

Array<IDType> mfem::ParticleSet::ids
protected

Global unique IDs of particles owned by this rank.

For the sample_particleset_data, ids would be ids[0]=id0, ids[1]=id1, ids[2]=id2.

Definition at line 287 of file particleset.hpp.

◆ tag_names

std::vector<std::string> mfem::ParticleSet::tag_names
protected

Tag names, to be written when PrintCSV() is called.

For the sample_particleset_data, tag_names would be tag_names[0]='type', tag_names[1]='color'.

Definition at line 323 of file particleset.hpp.

◆ tags

std::vector<std::unique_ptr<Array<int> > > mfem::ParticleSet::tags
protected

All particle tags for particles owned by this rank.

For the sample_particleset_data, tags would be *tags[0]=(t0,t1,t2), *tags[1]=(c0,c1,c2).

Definition at line 309 of file particleset.hpp.


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