12#ifndef MFEM_PARTICLESET
13#define MFEM_PARTICLESET
86 std::vector<Array<int>>
tags;
123 MFEM_ASSERT(
f >= 0 &&
static_cast<std::size_t
>(
f) <
fields.size(),
124 "Invalid field index");
125 MFEM_ASSERT(c >= 0 && c <
fields[
f].Size(),
126 "Invalid component index");
133 MFEM_ASSERT(
f >= 0 &&
static_cast<std::size_t
>(
f) <
fields.size(),
134 "invalid field index");
135 MFEM_ASSERT(c >= 0 && c <
fields[
f].Size(),
136 "invalid component index");
143 MFEM_ASSERT(
f >= 0 &&
static_cast<std::size_t
>(
f) <
fields.size(),
144 "invalid field index");
151 MFEM_ASSERT(
f >= 0 &&
static_cast<std::size_t
>(
f) <
fields.size(),
152 "invalid field index");
159 MFEM_ASSERT(t >= 0 &&
static_cast<std::size_t
>(t) <
tags.size(),
160 "invalid tag index");
165 const int&
Tag(
int t)
const
167 MFEM_ASSERT(t >= 0 &&
static_cast<std::size_t
>(t) <
tags.size(),
168 "invalid tag index");
257 static std::string GetDefaultFieldName(
int i);
260 static std::string GetDefaultTagName(
int i);
266 static int GetRank(MPI_Comm comm_);
267 static int GetSize(MPI_Comm comm_);
302 std::vector<std::unique_ptr<ParticleVector>>
fields;
309 std::vector<std::unique_ptr<Array<int>>>
tags;
338#if defined(MFEM_USE_MPI) && defined(MFEM_USE_GSLIB)
339 struct gslib::crystal *
cr =
nullptr;
343 template<std::
size_t NBytes>
344 static void TransferParticlesImpl(
ParticleSet &pset,
348 using TransferParticlesType = void (*)(
ParticleSet &pset,
353 MFEM_REGISTER_KERNELS(TransferParticles, TransferParticlesType, (
size_t));
354 friend TransferParticles;
372 {
ids[local_idx] = new_global_id; }
382 void WriteToFile(
const char *fname,
const std::stringstream &ss_header,
383 const std::stringstream &ss_data);
534 IDType GetGlobalNParticles() const;
551 const char* field_name=
nullptr);
561 return AddField(vdim, field_ordering, field_name);
570 int AddTag(
const char* tag_name=
nullptr);
657 void PrintCSV(
const char *fname,
int precision=16);
662 const Array<int> &tag_idxs,
int precision=16);
664#if defined(MFEM_USE_MPI) && defined(MFEM_USE_GSLIB)
int Size() const
Return the logical size of the array.
ParticleSet initializes and manages data associated with particles.
MPI_Comm GetComm() const
Get the MPI communicator for this ParticleSet.
const ParticleVector & Coords() const
Get a const reference to the coordinates ParticleVector.
Array< IDType > ids
Global unique IDs of particles owned by this rank.
ParticleVector coords
Spatial coordinates of particles owned by this rank.
~ParticleSet()
Destructor.
ParticleSet & operator=(const ParticleSet &)=delete
unsigned long long IDType
void Redistribute(const Array< unsigned int > &rank_list)
Redistribute particle data to rank_list.
int GetNTags() const
Get the number of tags registered to particles.
void AddParticle(const Particle &p)
Add a particle using Particle .
bool IsValidParticle(const Particle &p) const
Check if a particle could belong in this ParticleSet by comparing field and tag dimension.
Particle GetParticleRef(int i)
Get Particle object whose members reference the actual data associated with particle i in this Partic...
std::vector< std::string > tag_names
Tag names, to be written when PrintCSV() is called.
std::vector< std::unique_ptr< ParticleVector > > fields
All particle fields for particles owned by this rank.
Particle GetParticle(int i) const
Get new Particle object with copy of data associated with particle i .
const Array< int > GetFieldVDims() const
Get an Array<int> of the field vector-dimensions registered to particles.
ParticleVector & Coords()
Get a reference to the coordinates ParticleVector.
const Array< IDType > & GetIDs() const
Get the global IDs of the active particles owned by this ParticleSet.
int GetNFields() const
Get the number of fields registered to particles.
int AddTag(const char *tag_name=nullptr)
Add a tag to the ParticleSet.
const int id_stride
Stride for IDs (used internally when new particles are added).
void UpdateID(int local_idx, IDType new_global_id)
Update global ID of a particle.
bool IsParticleRefValid() const
Determine if GetParticleRef is valid.
ParticleVector & Field(int f)
Get a reference to field f 's ParticleVector.
int GetDim() const
Get the spatial dimension.
int GetNParticles() const
Get the number of active particles currently held by 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.
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.
Array< int > & Tag(int t)
Get a reference to tag t 's Array<int>.
int AddField(int vdim, Ordering::Type field_ordering=Ordering::byVDIM, const char *field_name=nullptr)
Add a field to the ParticleSet.
struct gslib::crystal * cr
void PrintCSV(const char *fname, int precision=16)
Print all particle data to a comma-delimited CSV file.
void Reserve(int res)
Reserve memory for res particles.
std::vector< std::unique_ptr< Array< int > > > tags
All particle tags for particles owned by this rank.
IDType id_counter
Current globally unique ID to be assigned to the next particle added.
const ParticleVector & Field(int f) const
Get a const reference to field f 's ParticleVector.
void RemoveParticles(const Array< int > &list)
Remove particle data specified by list of particle indices.
struct gslib::comm * gsl_comm
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 i...
std::vector< std::string > field_names
Field names, to be written when PrintCSV() is called.
int AddNamedField(int vdim, const char *field_name, Ordering::Type field_ordering=Ordering::byVDIM)
Add a field to the ParticleSet.
int FieldVDim(int f) const
Get Field vector-dimension.
const Array< int > & Tag(int t) const
Get a const reference to tag t 's Array<int>.
Particle CreateParticle() const
Create a Particle object with the same spatial dimension, number of fields and field vdims,...
void SetParticle(int i, const Particle &p)
Set data for particle at index i with data from provided particle p.
ParticleSet(const ParticleSet &)=delete
ParticleVector carries vector data (of a given vector dimension) for an arbitrary number of particles...
int GetVDim() const
Get the Vector dimension of the ParticleVector.
Container for data associated with a single particle.
std::vector< Array< int > > tags
A std::vector of Array<int> where each Array<int> holds data for a given tag.
int GetFieldVDim(int f) const
Get the vector dimension of field f .
int GetDim() const
Get the spatial dimension of this particle.
Particle(const Particle &)=default
Particle & operator=(Particle &&)=default
bool operator==(const Particle &rhs) const
Particle equality operator.
int & Tag(int t)
Get reference to tag t .
int GetNTags() const
Get the number of tags associated with this particle.
Particle(Particle &&)=default
Vector & Field(int f)
Get reference to field f Vector.
std::vector< Vector > fields
A std::vector of Vector where each Vector holds data for a given field (e.g., mass,...
const real_t & FieldValue(int f, int c=0) const
Get const reference to field f , component c value.
const Vector & Coords() const
Get const reference to particle coordinates Vector.
const Vector & Field(int f) const
Get const reference to field f Vector.
Vector coords
Spatial coordinates.
real_t & FieldValue(int f, int c=0)
Get reference to field f , component c value.
Particle(int dim, const Array< int > &field_vdims, int num_tags)
Construct a Particle instance.
Particle & operator=(const Particle &)=default
bool operator!=(const Particle &rhs) const
Particle inequality operator.
Vector & Coords()
Get reference to particle coordinates Vector.
void Print(std::ostream &os=mfem::out) const
Print all particle data to os.
const int & Tag(int t) const
Get const reference to tag t .
int GetNFields() const
Get the number of fields associated with this particle.
void SetTagRef(int t, int *tag_data)
Set tag t to reference external data.
void SetFieldRef(int f, real_t *field_data)
Set field f to reference external data.
int Size() const
Returns the size of the vector.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
std::function< real_t(const Vector &)> f(real_t mass_coeff)
real_t p(const Vector &x, real_t t)
Base class for Schrodinger solver kernels.