Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | Static Private Attributes | List of all members
Particles::Particle< dim, spacedim > Class Template Reference

#include <deal.II/particles/particle.h>

Public Member Functions

 Particle ()
 
 Particle (const Point< spacedim > &location, const Point< dim > &reference_location, const types::particle_index id)
 
 Particle (const Particle< dim, spacedim > &particle)
 
 Particle (const void *&begin_data, PropertyPool< dim, spacedim > *const property_pool=nullptr)
 
 Particle (Particle< dim, spacedim > &&particle) noexcept
 
Particle< dim, spacedim > & operator= (const Particle< dim, spacedim > &particle)
 
Particle< dim, spacedim > & operator= (Particle< dim, spacedim > &&particle) noexcept
 
 ~Particle ()
 
void * write_particle_data_to_memory (void *data) const
 
const void * read_particle_data_from_memory (const void *data)
 
void set_location (const Point< spacedim > &new_location)
 
const Point< spacedim > & get_location () const
 
Point< spacedim > & get_location ()
 
void set_reference_location (const Point< dim > &new_reference_location)
 
const Point< dim > & get_reference_location () const
 
types::particle_index get_id () const
 
void set_id (const types::particle_index &new_id)
 
void set_property_pool (PropertyPool< dim, spacedim > &property_pool)
 
bool has_properties () const
 
void set_properties (const ArrayView< const double > &new_properties)
 
ArrayView< double > get_properties ()
 
ArrayView< const double > get_properties () const
 
std::size_t serialized_size_in_bytes () const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
void free_properties ()
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 

Private Attributes

PropertyPool< dim, spacedim > * property_pool
 
PropertyPool< dim, spacedim >::Handle property_pool_handle
 

Static Private Attributes

static PropertyPool< dim, spacedim > global_property_pool
 

Detailed Description

template<int dim, int spacedim = dim>
class Particles::Particle< dim, spacedim >

A class that represents a particle in a domain that is meshed by a triangulation of some kind. The data this class stores is the position of the particle in the overall space, the position of the particle in the reference coordinate system of the cell it is currently in, an ID number that is unique among all particles, and a variable number of "properties".

The "properties" attached to each object of this class are stored by a PropertyPool object. These properties are stored as an array of double variables that can be accessed via an ArrayView object. For example, if one wanted to equip each particle with a "temperature" and "chemical composition" property that is advected along with the particle (and may change from time step to time step based on some differential equation, for example), then one would allocate two properties per particle in the PropertyPool object.

In practice, however, one often wants to associate properties with particles that are not just independent numbers as in the situation above. An example would be if one wanted to track the stress or strain that a particle is subjected to – a tensor-valued quantity. In these cases, one would interpret these scalar properties as the components of the stress or strain. In other words, one would first tell the PropertyPool to allocate as many properties per particle as there are components in the tensor one wants to track, and then write small conversion functions that take the ArrayView of scalar properties returned by the get_properties() function and convert it to a tensor of the appropriate type. This can then be evaluated and evolved in each time step. A second conversion function would convert back from a tensor to an ArrayView object to store the updated data back in the particle via the set_properties() function.

There are of course cases where the properties one cares about are not real (or, in computers, floating point) numbers but rather categorical: For example, one may want to mark some particles as "red", "blue", or "green". The property might then either be represented as an integer, or as an element of an enum. In these cases, one would need to come up with a way to represent these sorts of categorical fields in terms of floating point numbers. For example, one could map "red" to the floating point number 1.0, "blue" to 2.0, and "green" to 3.0. The conversion functions to translate between these two representations should then not be very difficult to write either.

Definition at line 95 of file particle.h.

Constructor & Destructor Documentation

◆ Particle() [1/5]

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::Particle ( )

Empty constructor for Particle, creates a particle at the origin.

Definition at line 29 of file particle.cc.

◆ Particle() [2/5]

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::Particle ( const Point< spacedim > &  location,
const Point< dim > &  reference_location,
const types::particle_index  id 
)

Constructor for Particle. This function creates a particle with the specified ID at the specified location. Note that there is no check for duplicate particle IDs so the user must make sure the IDs are unique over all processes. Data is stored in a global PropertyPool object (corresponding to the global "heap") but can later be transferred to another property pool by calling set_property_pool().

Parameters
[in]locationInitial location of particle.
[in]reference_locationInitial location of the particle in the coordinate system of the reference cell.
[in]idGlobally unique ID number of particle.

Definition at line 37 of file particle.cc.

◆ Particle() [3/5]

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::Particle ( const Particle< dim, spacedim > &  particle)

Copy-constructor for Particle. This function creates a particle with exactly the state of the input argument. The copied data is stored in a global PropertyPool object (corresponding to the global "heap") but can later be transferred to another property pool by calling set_property_pool().

Definition at line 51 of file particle.cc.

◆ Particle() [4/5]

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::Particle ( const void *&  begin_data,
PropertyPool< dim, spacedim > *const  property_pool = nullptr 
)

Constructor for Particle. This function creates a particle from a data vector. Data is stored in a global PropertyPool object (corresponding to the global "heap") but can later be transferred to another property pool by calling set_property_pool(). This constructor is usually called after serializing a particle by calling the write_data() function.

Parameters
[in,out]begin_dataA pointer to a memory location from which to read the information that completely describes a particle. This class then de-serializes its data from this memory location and advances the pointer beyond the data that has been read to initialize the particle information.
[in,out]property_poolAn optional pointer to a property pool that is used to manage the property data used by this particle. If this argument is not provided, then a global property pool is used; on the other hand, if a non-null pointer is provided, this constructor assumes begin_data contains serialized data of the same length and type that is allocated by property_pool. If the data pointer provided here corresponds to data for a particle that has properties, then this function will only succeed if a property pool is provided as second argument that is able to store the correct number of properties per particle.

Definition at line 75 of file particle.cc.

◆ Particle() [5/5]

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::Particle ( Particle< dim, spacedim > &&  particle)
noexcept

Move constructor for Particle, creates a particle from an existing one by stealing its state.

Definition at line 112 of file particle.cc.

◆ ~Particle()

template<int dim, int spacedim>
Particles::Particle< dim, spacedim >::~Particle ( )

Destructor. Releases the property handle if it is valid, and therefore frees that memory space for other particles. (Note: the memory is managed by the property pool, and the pool is responsible for what happens to the memory.

Definition at line 190 of file particle.cc.

Member Function Documentation

◆ operator=() [1/2]

template<int dim, int spacedim>
Particle< dim, spacedim > & Particles::Particle< dim, spacedim >::operator= ( const Particle< dim, spacedim > &  particle)

Copy assignment operator.

Definition at line 129 of file particle.cc.

◆ operator=() [2/2]

template<int dim, int spacedim>
Particle< dim, spacedim > & Particles::Particle< dim, spacedim >::operator= ( Particle< dim, spacedim > &&  particle)
noexcept

Move assignment operator.

Definition at line 161 of file particle.cc.

◆ write_particle_data_to_memory()

template<int dim, int spacedim>
void * Particles::Particle< dim, spacedim >::write_particle_data_to_memory ( void *  data) const

Write particle data into a data array. The array is expected to be large enough to take the data, and the void pointer should point to the first entry of the array to which the data should be written. This function is meant for serializing all particle properties and later de-serializing the properties by calling the appropriate constructor Particle(void *&data, PropertyPool *property_pool = nullptr);

Parameters
[in]dataThe memory location to write particle data into.
Returns
A pointer to the next byte after the array to which data has been written.

Definition at line 218 of file particle.cc.

◆ read_particle_data_from_memory()

template<int dim, int spacedim>
const void * Particles::Particle< dim, spacedim >::read_particle_data_from_memory ( const void *  data)

Update all of the data associated with a particle: id, location, reference location and, if any, properties by using a data array. The array is expected to be large enough to take the data, and the void pointer should point to the first entry of the array to which the data should be written. This function is meant for de-serializing the particle data without requiring that a new Particle class be built. This is used in the ParticleHandler to update the ghost particles without de-allocating and re-allocating memory.

Parameters
[in]dataA pointer to a memory location from which to read the information that completely describes a particle. This class then de-serializes its data from this memory location.
Returns
A pointer to the next byte after the array from which data has been read.

Definition at line 251 of file particle.cc.

◆ set_location()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::set_location ( const Point< spacedim > &  new_location)
inline

Set the location of this particle. Note that this does not check whether this is a valid location in the simulation domain.

Parameters
[in]new_locationThe new location for this particle.
Note
In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as const objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.

Definition at line 547 of file particle.h.

◆ get_location() [1/2]

template<int dim, int spacedim>
const Point< spacedim > & Particles::Particle< dim, spacedim >::get_location ( ) const
inline

Get the location of this particle.

Returns
The location of this particle.

Definition at line 556 of file particle.h.

◆ get_location() [2/2]

template<int dim, int spacedim>
Point< spacedim > & Particles::Particle< dim, spacedim >::get_location ( )
inline

Get read- and write-access to the location of this particle. Note that changing the location does not check whether this is a valid location in the simulation domain.

Note
In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as const objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Returns
The location of this particle.

Definition at line 565 of file particle.h.

◆ set_reference_location()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::set_reference_location ( const Point< dim > &  new_reference_location)
inline

Set the reference location of this particle.

Parameters
[in]new_reference_locationThe new reference location for this particle.
Note
In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as const objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.

Definition at line 574 of file particle.h.

◆ get_reference_location()

template<int dim, int spacedim>
const Point< dim > & Particles::Particle< dim, spacedim >::get_reference_location ( ) const
inline

Return the reference location of this particle in its current cell.

Definition at line 583 of file particle.h.

◆ get_id()

template<int dim, int spacedim>
types::particle_index Particles::Particle< dim, spacedim >::get_id ( ) const
inline

Return the ID number of this particle. The ID of a particle is intended to be a property that is globally unique even in parallel computations and is transferred along with other properties of a particle if it moves from a cell owned by the current processor to a cell owned by a different processor, or if ownership of the cell it is on is transferred to a different processor.

Definition at line 592 of file particle.h.

◆ set_id()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::set_id ( const types::particle_index new_id)
inline

Set the ID number of this particle. The ID of a particle is intended to be a property that is globally unique even in parallel computations and is transferred along with other properties of a particle if it moves from a cell owned by the current processor to a cell owned by a different processor, or if ownership of the cell it is on is transferred to a different processor. As a consequence, when setting the ID of a particle, care needs to be taken to ensure that particles have globally unique IDs. (The ParticleHandler does not itself check whether particle IDs so set are globally unique in a parallel setting since this would be a very expensive operation.)

Parameters
[in]new_idThe new ID number for this particle.
Note
In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as const objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.

Definition at line 601 of file particle.h.

◆ set_property_pool()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::set_property_pool ( PropertyPool< dim, spacedim > &  property_pool)
inline

Tell the particle where to store its properties (even if it does not own properties). Usually this is only done once per particle, but since the particle does not know about the properties, we want to do it not at construction time. Another use for this function is after particle transfer to a new process.

If a particle already stores properties in a property pool, then their values are saved, the memory is released in the previous property pool, and a copy of the particle's properties will be allocated in the new property pool.

Definition at line 610 of file particle.h.

◆ has_properties()

template<int dim, int spacedim>
bool Particles::Particle< dim, spacedim >::has_properties ( ) const
inline

Return whether this particle has a valid property pool and a valid handle to properties.

Definition at line 672 of file particle.h.

◆ set_properties()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::set_properties ( const ArrayView< const double > &  new_properties)

Set the properties of this particle.

Parameters
[in]new_propertiesAn ArrayView containing the new properties for this particle.
Note
In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as const objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.

Definition at line 304 of file particle.cc.

◆ get_properties() [1/2]

template<int dim, int spacedim>
ArrayView< double > Particles::Particle< dim, spacedim >::get_properties ( )

Get write-access to properties of this particle. If the particle has no properties yet, but has access to a PropertyPool object it will allocate properties to allow writing into them. If it has no properties and has no access to a PropertyPool this function will throw an exception.

Returns
An ArrayView of the properties of this particle.

Definition at line 329 of file particle.cc.

◆ get_properties() [2/2]

template<int dim, int spacedim>
ArrayView< const double > Particles::Particle< dim, spacedim >::get_properties ( ) const
inline

Get read-access to properties of this particle. If the particle has no properties this function throws an exception.

Returns
An ArrayView of the properties of this particle.

Definition at line 660 of file particle.h.

◆ serialized_size_in_bytes()

template<int dim, int spacedim>
std::size_t Particles::Particle< dim, spacedim >::serialized_size_in_bytes ( ) const

Return the size in bytes this particle occupies if all of its data is serialized (i.e. the number of bytes that is written by the write_data function of this class).

Definition at line 286 of file particle.cc.

◆ save()

template<int dim, int spacedim>
template<class Archive >
void Particles::Particle< dim, spacedim >::save ( Archive &  ar,
const unsigned int  version 
) const
inline

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 525 of file particle.h.

◆ load()

template<int dim, int spacedim>
template<class Archive >
void Particles::Particle< dim, spacedim >::load ( Archive &  ar,
const unsigned int  version 
)
inline

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library. Note that in order to store the properties correctly, the property pool of this particle has to be known at the time of reading, i.e. set_property_pool() has to have been called, before this function is called.

Definition at line 491 of file particle.h.

◆ free_properties()

template<int dim, int spacedim>
void Particles::Particle< dim, spacedim >::free_properties ( )

Free the memory of the property pool

Definition at line 203 of file particle.cc.

◆ serialize()

template<int dim, int spacedim = dim>
template<class Archive >
void Particles::Particle< dim, spacedim >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

Member Data Documentation

◆ global_property_pool

template<int dim, int spacedim>
PropertyPool< dim, spacedim > Particles::Particle< dim, spacedim >::global_property_pool
staticprivate
Initial value:
= {
0}

A global property pool used when a particle is not associated with a property pool that belongs to, for example, a ParticleHandler.

Definition at line 470 of file particle.h.

◆ property_pool

template<int dim, int spacedim = dim>
PropertyPool<dim, spacedim>* Particles::Particle< dim, spacedim >::property_pool
private

A pointer to the property pool. Necessary to translate from the handle to the actual memory locations.

Definition at line 476 of file particle.h.

◆ property_pool_handle

template<int dim, int spacedim = dim>
PropertyPool<dim,spacedim>::Handle Particles::Particle< dim, spacedim >::property_pool_handle
private

A handle to all particle properties

Definition at line 481 of file particle.h.


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