16#ifndef dealii_particles_particle_accessor_h
17#define dealii_particles_particle_accessor_h
36 class ParticleIterator;
38 class ParticleHandler;
44 template <
int dim,
int spacedim = dim>
102 std::vector<typename PropertyPool<dim, spacedim>::Handle>
particles;
362 template <
class Archive>
364 save(Archive &ar,
const unsigned int version)
const;
375 template <
class Archive>
377 load(Archive &ar,
const unsigned int version);
385 template <
class Archive>
387 serialize(Archive &archive,
const unsigned int version);
391 BOOST_SERIALIZATION_SPLIT_MEMBER()
481 template <
int dim,
int spacedim>
482 template <
class Archive>
486 unsigned int n_properties = 0;
491 ar &location &reference_location &
id &n_properties;
493 set_location(location);
494 set_reference_location(reference_location);
497 if (n_properties > 0)
501 properties.
size() == n_properties,
503 "This particle was serialized with " +
504 std::to_string(n_properties) +
505 " properties, but the new property handler provides space for " +
506 std::to_string(properties.
size()) +
507 " properties. Deserializing a particle only works for matching property sizes."));
509 ar &boost::serialization::make_array(properties.
data(), n_properties);
515 template <
int dim,
int spacedim>
516 template <
class Archive>
520 unsigned int n_properties = 0;
521 if ((property_pool !=
nullptr) &&
523 n_properties = get_properties().size();
526 Point<dim> reference_location = get_reference_location();
529 ar &location &reference_location &
id &n_properties;
531 if (n_properties > 0)
532 ar &boost::serialization::make_array(get_properties().data(),
539 template <
int dim,
int spacedim>
542 , property_pool(nullptr)
543 , particle_index_within_cell(
numbers::invalid_unsigned_int)
548 template <
int dim,
int spacedim>
550 const typename particle_container::iterator particles_in_cell,
552 const unsigned int particle_index_within_cell)
553 : particles_in_cell(particles_in_cell)
554 , property_pool(const_cast<
PropertyPool<dim, spacedim> *>(&property_pool))
555 , particle_index_within_cell(particle_index_within_cell)
560 template <
int dim,
int spacedim>
570 const double *pdata =
reinterpret_cast<const double *
>(id_data);
573 for (
unsigned int i = 0; i < spacedim; ++i)
574 location(i) = *pdata++;
575 set_location(location);
578 for (
unsigned int i = 0; i < dim; ++i)
579 reference_location(i) = *pdata++;
580 set_reference_location(reference_location);
583 if (has_properties())
586 property_pool->get_properties(get_handle());
587 const unsigned int size = particle_properties.
size();
588 for (
unsigned int i = 0; i < size; ++i)
589 particle_properties[i] = *pdata++;
592 return static_cast<const void *
>(pdata);
597 template <
int dim,
int spacedim>
607 double *pdata =
reinterpret_cast<double *
>(id_data);
610 for (
unsigned int i = 0; i < spacedim; ++i, ++pdata)
611 *pdata = get_location()[i];
614 for (
unsigned int i = 0; i < dim; ++i, ++pdata)
615 *pdata = get_reference_location()[i];
618 if (has_properties())
621 property_pool->get_properties(get_handle());
622 for (
unsigned int i = 0; i < particle_properties.
size(); ++i, ++pdata)
623 *pdata = particle_properties[i];
626 return static_cast<void *
>(pdata);
631 template <
int dim,
int spacedim>
637 property_pool->set_location(get_handle(), new_loc);
642 template <
int dim,
int spacedim>
648 return property_pool->get_location(get_handle());
653 template <
int dim,
int spacedim>
660 property_pool->set_reference_location(get_handle(), new_loc);
665 template <
int dim,
int spacedim>
671 return property_pool->get_reference_location(get_handle());
676 template <
int dim,
int spacedim>
682 return property_pool->get_id(get_handle());
687 template <
int dim,
int spacedim>
698 template <
int dim,
int spacedim>
704 property_pool->set_id(get_handle(), new_id);
709 template <
int dim,
int spacedim>
722 return (property_pool->n_properties_per_slot() > 0);
727 template <
int dim,
int spacedim>
730 const std::vector<double> &new_properties)
740 template <
int dim,
int spacedim>
748 property_pool->get_properties(get_handle());
752 "You are trying to assign properties with an incompatible length. "
753 "The particle has space to store " +
754 std::to_string(property_values.
size()) +
755 " properties, but you are trying to assign " +
756 std::to_string(new_properties.
size()) +
757 " properties. This is not allowed."));
759 if (property_values.
size() > 0)
760 std::copy(new_properties.
begin(),
761 new_properties.
end(),
762 property_values.
begin());
767 template <
int dim,
int spacedim>
778 for (
unsigned int d = 0; d < dim; ++d)
779 array[d] = new_properties[d];
786 template <
int dim,
int spacedim>
792 return property_pool->get_properties(get_handle());
797 template <
int dim,
int spacedim>
803 (void)new_property_pool;
808 template <
int dim,
int spacedim>
816 return particles_in_cell->cell;
821 template <
int dim,
int spacedim>
830 return particles_in_cell->cell;
835 template <
int dim,
int spacedim>
841 return property_pool->get_properties(get_handle());
846 template <
int dim,
int spacedim>
852 std::size_t size =
sizeof(get_id()) +
sizeof(get_location()) +
853 sizeof(get_reference_location());
855 if (has_properties())
857 size +=
sizeof(double) * get_properties().size();
864 template <
int dim,
int spacedim>
870 ++particle_index_within_cell;
872 if (particle_index_within_cell >= particles_in_cell->particles.size())
874 particle_index_within_cell = 0;
881 template <
int dim,
int spacedim>
887 if (particle_index_within_cell > 0)
888 --particle_index_within_cell;
892 particle_index_within_cell = particles_in_cell->particles.empty() ?
894 particles_in_cell->particles.size() - 1;
900 template <
int dim,
int spacedim>
905 return !(*
this == other);
910 template <
int dim,
int spacedim>
922 template <
int dim,
int spacedim>
926 if (property_pool !=
nullptr &&
928 particle_index_within_cell < particles_in_cell->particles.size())
930 else if (property_pool !=
nullptr &&
932 particle_index_within_cell == 0)
942 template <
int dim,
int spacedim>
946 return particles_in_cell->particles[particle_index_within_cell];
951 template <
int dim,
int spacedim>
955 return particles_in_cell->particles[particle_index_within_cell];
976 template <
int dim,
int spacedim>
986 operator()(const ::Particles::ParticleAccessor<dim, spacedim>
989 return accessor.get_location();
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
value_type * data() const noexcept
const Triangulation< dim, spacedim >::cell_iterator & get_surrounding_cell(const Triangulation< dim, spacedim > &triangulation) const
particle_container::iterator particles_in_cell
const Point< spacedim > & get_location() const
types::particle_index get_id() const
void save(Archive &ar, const unsigned int version) const
unsigned int particle_index_within_cell
types::particle_index get_local_index() const
std::list< ParticlesInCell > particle_container
const void * read_particle_data_from_memory(const void *data)
void load(Archive &ar, const unsigned int version)
const Triangulation< dim, spacedim >::cell_iterator & get_surrounding_cell() const
IteratorState::IteratorStates state() const
bool operator!=(const ParticleAccessor< dim, spacedim > &other) const
void set_properties(const Tensor< 1, dim > &new_properties)
const ArrayView< double > get_properties()
void set_location(const Point< spacedim > &new_location)
const ArrayView< const double > get_properties() const
PropertyPool< dim, spacedim >::Handle & get_handle()
void serialize(Archive &archive, const unsigned int version)
void set_property_pool(PropertyPool< dim, spacedim > &property_pool)
void set_reference_location(const Point< dim > &new_reference_location)
const Point< dim > & get_reference_location() const
std::size_t serialized_size_in_bytes() const
const PropertyPool< dim, spacedim >::Handle & get_handle() const
void set_properties(const std::vector< double > &new_properties)
ParticleAccessor(const typename particle_container::iterator particles_in_cell, const PropertyPool< dim, spacedim > &property_pool, const unsigned int particle_index_within_cell)
bool has_properties() const
PropertyPool< dim, spacedim > * property_pool
void * write_particle_data_to_memory(void *data) const
void set_properties(const ArrayView< const double > &new_properties)
bool operator==(const ParticleAccessor< dim, spacedim > &other) const
void set_id(const types::particle_index &new_id)
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
@ past_the_end
Iterator reached end of container.
@ valid
Iterator points to a valid object.
@ invalid
Iterator is invalid, probably due to an error.
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Triangulation< dim, spacedim >::active_cell_iterator cell
ParticlesInCell()=default
ParticlesInCell(const std::vector< typename PropertyPool< dim, spacedim >::Handle > &particles, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
std::vector< typename PropertyPool< dim, spacedim >::Handle > particles
result_type operator()(const ::Particles::ParticleAccessor< dim, spacedim > &accessor) const