Reference documentation for deal.II version 9.3.3
\(\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\}}\)
particle_handler.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_particles_particle_handler_h
17#define dealii_particles_particle_handler_h
18
19#include <deal.II/base/config.h>
20
24#include <deal.II/base/mpi.h>
27
29
30#include <deal.II/fe/mapping.h>
31
33
38
39#include <boost/range/iterator_range.hpp>
40#include <boost/serialization/map.hpp>
41
43
44namespace Particles
45{
60 template <int dim, int spacedim = dim>
62 {
63 public:
68
72 using particle_iterator_range = boost::iterator_range<particle_iterator>;
73
78
89 const unsigned int n_properties = 0);
90
94 virtual ~ParticleHandler() override = default;
95
101 void
104 const unsigned int n_properties = 0);
105
126 void
127 copy_from(const ParticleHandler<dim, spacedim> &particle_handler);
128
132 void
133 clear();
134
140 void
142
151 void
153
158 begin() const;
159
164 begin();
165
170 end() const;
171
176 end();
177
182 begin_ghost() const;
183
188 begin_ghost();
189
194 end_ghost() const;
195
200 end_ghost();
201
228 const;
229
260
291 const;
292
296 void
297 remove_particle(const particle_iterator &particle);
298
307 const Particle<dim, spacedim> &particle,
309
316 void
318 const std::multimap<
321
331 void
332 insert_particles(const std::vector<Point<spacedim>> &positions);
333
394 std::map<unsigned int, IndexSet>
396 const std::vector<Point<spacedim>> &positions,
397 const std::vector<std::vector<BoundingBox<spacedim>>>
398 & global_bounding_boxes,
399 const std::vector<std::vector<double>> & properties = {},
400 const std::vector<types::particle_index> &ids = {});
401
432 std::map<unsigned int, IndexSet>
434 const std::vector<Particle<dim, spacedim>> &particles,
435 const std::vector<std::vector<BoundingBox<spacedim>>>
436 &global_bounding_boxes);
437
481 template <class VectorType>
482 typename std::enable_if<
483 std::is_convertible<VectorType *, Function<spacedim> *>::value ==
484 false>::type
485 set_particle_positions(const VectorType &input_vector,
486 const bool displace_particles = true);
487
508 void
509 set_particle_positions(const std::vector<Point<spacedim>> &new_positions,
510 const bool displace_particles = true);
511
512
531 void
533 const bool displace_particles = true);
534
562 template <class VectorType>
563 void
564 get_particle_positions(VectorType &output_vector,
565 const bool add_to_output_vector = false);
566
581 void
582 get_particle_positions(std::vector<Point<spacedim>> &positions,
583 const bool add_to_output_vector = false);
584
611 void
613 const std::function<std::size_t()> &size_callback,
614 const std::function<void *(const particle_iterator &, void *)>
616 const std::function<const void *(const particle_iterator &, const void *)>
618
628 n_global_particles() const;
629
638
645
653
673 locally_relevant_ids() const;
674
693
697 unsigned int
699
704 PropertyPool<dim, spacedim> &
705 get_property_pool() const;
706
721 void
723
729 void
730 exchange_ghost_particles(const bool enable_ghost_cache = false);
731
739 void
741
748 void
750
757 void
758 register_load_callback_function(const bool serialization);
759
764 template <class Archive>
765 void
766 serialize(Archive &ar, const unsigned int version);
767
780 struct Signals
781 {
795 boost::signals2::signal<void(
796 const typename Particles::ParticleIterator<dim, spacedim> &particle,
798 &cell)>
800 };
801
807
808 private:
815
821
830 std::unique_ptr<PropertyPool<dim, spacedim>> property_pool;
831
836 std::multimap<internal::LevelInd, Particle<dim, spacedim>> particles;
837
843 std::multimap<internal::LevelInd, Particle<dim, spacedim>> ghost_particles;
844
850
860
866
878 std::function<std::size_t()> size_callback;
879
890 std::function<void *(const particle_iterator &, void *)> store_callback;
891
902 std::function<const void *(const particle_iterator &, const void *)>
904
911 unsigned int handle;
912
922 std::unique_ptr<GridTools::Cache<dim, spacedim>> triangulation_cache;
923
924#ifdef DEAL_II_WITH_MPI
953 void
955 const std::map<types::subdomain_id, std::vector<particle_iterator>>
956 &particles_to_send,
958 &received_particles,
959 const std::map<
961 std::vector<
963 &new_cells_for_particles = std::map<
965 std::vector<
967 const bool enable_cache = false);
968
986 void
988 const std::map<types::subdomain_id, std::vector<particle_iterator>>
989 &particles_to_send,
991 &received_particles);
992
993
994#endif
995
1003
1009 std::vector<char>
1012 const typename Triangulation<dim, spacedim>::CellStatus status) const;
1013
1018 void
1021 const typename Triangulation<dim, spacedim>::CellStatus status,
1022 const boost::iterator_range<std::vector<char>::const_iterator>
1023 &data_range);
1024 };
1025
1026
1027
1028 /* ---------------------- inline and template functions ------------------
1029 */
1030
1031 template <int dim, int spacedim>
1034 {
1035 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin();
1036 }
1037
1038
1039
1040 template <int dim, int spacedim>
1043 {
1044 return particle_iterator(particles, particles.begin());
1045 }
1046
1047
1048
1049 template <int dim, int spacedim>
1052 {
1053 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end();
1054 }
1055
1056
1057
1058 template <int dim, int spacedim>
1061 {
1062 return particle_iterator(particles, particles.end());
1063 }
1064
1065
1066
1067 template <int dim, int spacedim>
1070 {
1071 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin_ghost();
1072 }
1073
1074
1075
1076 template <int dim, int spacedim>
1079 {
1080 return particle_iterator(ghost_particles, ghost_particles.begin());
1081 }
1082
1083
1084
1085 template <int dim, int spacedim>
1088 {
1089 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end_ghost();
1090 }
1091
1092
1093
1094 template <int dim, int spacedim>
1097 {
1098 return particle_iterator(ghost_particles, ghost_particles.end());
1099 }
1100
1101
1102
1103 template <int dim, int spacedim>
1104 template <class Archive>
1105 inline void
1106 ParticleHandler<dim, spacedim>::serialize(Archive &ar, const unsigned int)
1107 {
1108 // Note that we do not serialize the particle data itself. Instead we
1109 // use the serialization functionality of the triangulation class, because
1110 // this guarantees that data is immediately shipped to new processes if
1111 // the domain is distributed differently after resuming from a checkpoint.
1112 ar //&particles
1113 &global_number_of_particles &global_max_particles_per_cell
1114 & next_free_particle_index;
1115 }
1116
1117
1118
1119 template <int dim, int spacedim>
1120 template <class VectorType>
1121 inline typename std::enable_if<
1122 std::is_convertible<VectorType *, Function<spacedim> *>::value ==
1123 false>::type
1125 const VectorType &input_vector,
1126 const bool displace_particles)
1127 {
1128 AssertDimension(input_vector.size(),
1129 get_next_free_particle_index() * spacedim);
1130 for (auto &p : *this)
1131 {
1132 auto new_point(displace_particles ? p.get_location() :
1133 Point<spacedim>());
1134 const auto id = p.get_id();
1135 for (unsigned int i = 0; i < spacedim; ++i)
1136 new_point[i] += input_vector[id * spacedim + i];
1137 p.set_location(new_point);
1138 }
1139 sort_particles_into_subdomains_and_cells();
1140 }
1141
1142
1143
1144 template <int dim, int spacedim>
1145 template <class VectorType>
1146 inline void
1148 VectorType &output_vector,
1149 const bool add_to_output_vector)
1150 {
1151 AssertDimension(output_vector.size(),
1152 get_next_free_particle_index() * spacedim);
1153 for (const auto &p : *this)
1154 {
1155 auto point = p.get_location();
1156 const auto id = p.get_id();
1157 if (add_to_output_vector)
1158 for (unsigned int i = 0; i < spacedim; ++i)
1159 output_vector[id * spacedim + i] += point[i];
1160 else
1161 for (unsigned int i = 0; i < spacedim; ++i)
1162 output_vector[id * spacedim + i] = point[i];
1163 }
1164 if (add_to_output_vector)
1165 output_vector.compress(VectorOperation::add);
1166 else
1167 output_vector.compress(VectorOperation::insert);
1168 }
1169
1170
1171
1172 template <int dim, int spacedim>
1173 inline IndexSet
1175 {
1176 return this->locally_owned_particle_ids();
1177 }
1178
1179} // namespace Particles
1180
1182
1183#endif
std::function< void *(const particle_iterator &, void *)> store_callback
void register_additional_store_load_functions(const std::function< std::size_t()> &size_callback, const std::function< void *(const particle_iterator &, void *)> &store_callback, const std::function< const void *(const particle_iterator &, const void *)> &load_callback)
void exchange_ghost_particles(const bool enable_ghost_cache=false)
types::particle_index n_global_particles() const
std::multimap< internal::LevelInd, Particle< dim, spacedim > > particles
particle_iterator end_ghost() const
unsigned int global_max_particles_per_cell
internal::GhostParticlePartitioner< dim, spacedim > ghost_particles_cache
IndexSet locally_relevant_ids() const
unsigned int n_properties_per_particle() const
types::particle_index global_number_of_particles
void get_particle_positions(VectorType &output_vector, const bool add_to_output_vector=false)
boost::iterator_range< particle_iterator > particle_iterator_range
std::multimap< internal::LevelInd, Particle< dim, spacedim > > ghost_particles
ParticleIterator< dim, spacedim > particle_iterator
particle_iterator begin_ghost() const
PropertyPool< dim, spacedim > & get_property_pool() const
void initialize(const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0)
void serialize(Archive &ar, const unsigned int version)
std::unique_ptr< PropertyPool< dim, spacedim > > property_pool
particle_iterator begin() const
void send_recv_particles_properties_and_location(const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send, std::multimap< internal::LevelInd, Particle< dim, spacedim > > &received_particles)
std::map< unsigned int, IndexSet > insert_global_particles(const std::vector< Point< spacedim > > &positions, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, const std::vector< std::vector< double > > &properties={}, const std::vector< types::particle_index > &ids={})
void register_load_callback_function(const bool serialization)
virtual ~ParticleHandler() override=default
std::function< std::size_t()> size_callback
std::vector< char > store_particles(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status) const
void insert_particles(const std::multimap< typename Triangulation< dim, spacedim >::active_cell_iterator, Particle< dim, spacedim > > &particles)
types::particle_index next_free_particle_index
particle_iterator end() const
void remove_particle(const particle_iterator &particle)
types::particle_index n_locally_owned_particles() const
types::particle_index n_particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
void load_particles(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
SmartPointer< const Triangulation< dim, spacedim >, ParticleHandler< dim, spacedim > > triangulation
particle_iterator insert_particle(const Particle< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
types::particle_index n_global_max_particles_per_cell() const
void copy_from(const ParticleHandler< dim, spacedim > &particle_handler)
std::function< const void *(const particle_iterator &, const void *)> load_callback
SmartPointer< const Mapping< dim, spacedim >, ParticleHandler< dim, spacedim > > mapping
std::unique_ptr< GridTools::Cache< dim, spacedim > > triangulation_cache
std::enable_if< std::is_convertible< VectorType *, Function< spacedim > * >::value==false >::type set_particle_positions(const VectorType &input_vector, const bool displace_particles=true)
types::particle_index get_next_free_particle_index() const
void send_recv_particles(const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send, std::multimap< internal::LevelInd, Particle< dim, spacedim > > &received_particles, const std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > > &new_cells_for_particles=std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > >(), const bool enable_cache=false)
IndexSet locally_owned_particle_ids() const
#define DEAL_II_DEPRECATED
Definition: config.h:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
std::pair< int, int > LevelInd
Definition: particle.h:42
boost::signals2::signal< void(const typename Particles::ParticleIterator< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)> particle_lost