Reference documentation for deal.II version GIT relicensing-1321-g19b0133728 2024-07-26 07:50:01+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
particle_handler.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_particles_particle_handler_h
16#define dealii_particles_particle_handler_h
17
18#include <deal.II/base/config.h>
19
25
27
28#include <deal.II/fe/mapping.h>
29
31
36
37#include <boost/range/iterator_range.hpp>
38#include <boost/serialization/map.hpp>
39
41
42namespace Particles
43{
61 template <int dim, int spacedim = dim>
63 {
64 public:
69
73 using particle_iterator_range = boost::iterator_range<particle_iterator>;
74
80
85
96 const unsigned int n_properties = 0);
97
101 virtual ~ParticleHandler();
102
108 void
111 const unsigned int n_properties = 0);
112
133 void
134 copy_from(const ParticleHandler<dim, spacedim> &particle_handler);
135
139 void
140 clear();
141
147 void
149
164 void
165 reserve(const std::size_t n_particles);
166
175 void
177
182 begin() const;
183
188 begin();
189
194 end() const;
195
200 end();
201
206 begin_ghost() const;
207
212 begin_ghost();
213
218 end_ghost() const;
219
224 end_ghost();
225
232 const;
233
245
257 const;
258
264 void
265 remove_particle(const particle_iterator &particle);
266
272 void
273 remove_particles(const std::vector<particle_iterator> &particles);
274
283 const Particle<dim, spacedim> &particle,
285
303 const Point<spacedim> &position,
304 const Point<dim> &reference_position,
305 const types::particle_index particle_index,
307 const ArrayView<const double> &properties = {});
308
315 void
317 const std::multimap<
320
330 void
331 insert_particles(const std::vector<Point<spacedim>> &positions);
332
393 std::map<unsigned int, IndexSet>
395 const std::vector<Point<spacedim>> &positions,
396 const std::vector<std::vector<BoundingBox<spacedim>>>
397 &global_bounding_boxes,
398 const std::vector<std::vector<double>> &properties = {},
399 const std::vector<types::particle_index> &ids = {});
400
431 std::map<unsigned int, IndexSet>
433 const std::vector<Particle<dim, spacedim>> &particles,
434 const std::vector<std::vector<BoundingBox<spacedim>>>
435 &global_bounding_boxes);
436
480 template <typename VectorType>
481 std::enable_if_t<
482 std::is_convertible_v<VectorType *, Function<spacedim> *> == false>
483 set_particle_positions(const VectorType &input_vector,
484 const bool displace_particles = true);
485
506 void
507 set_particle_positions(const std::vector<Point<spacedim>> &new_positions,
508 const bool displace_particles = true);
509
510
529 void
531 const bool displace_particles = true);
532
560 template <typename VectorType>
561 void
562 get_particle_positions(VectorType &output_vector,
563 const bool add_to_output_vector = false);
564
579 void
580 get_particle_positions(std::vector<Point<spacedim>> &positions,
581 const bool add_to_output_vector = false);
582
609 void
611 const std::function<std::size_t()> &size_callback,
612 const std::function<void *(const particle_iterator &, void *)>
614 const std::function<const void *(const particle_iterator &, const void *)>
616
626 n_global_particles() const;
627
636
643
651
670
674 unsigned int
676
692
697 PropertyPool<dim, spacedim> &
698 get_property_pool() const;
699
714 void
716
722 void
723 exchange_ghost_particles(const bool enable_ghost_cache = false);
724
732 void
734
744 void
746
765 void
767
779 void
781
790 void
791 deserialize();
792
804 void
806
818 void
819 register_load_callback_function(const bool serialization);
820
827 template <class Archive>
828 void
829 serialize(Archive &ar, const unsigned int version);
830
843 struct Signals
844 {
858 boost::signals2::signal<void(
859 const typename Particles::ParticleIterator<dim, spacedim> &particle,
861 &cell)>
863 };
864
870
871 private:
880 const void *&data,
882
892
902 void
904
911
917
926 std::unique_ptr<PropertyPool<dim, spacedim>> property_pool;
927
932
938 const typename particle_container::iterator owned_particles_end;
939
945 std::vector<typename particle_container::iterator> cells_to_particle_cache;
946
952
958
968
974
986 std::function<std::size_t()> size_callback;
987
998 std::function<void *(const particle_iterator &, void *)> store_callback;
999
1010 std::function<const void *(const particle_iterator &, const void *)>
1012
1020
1024 const double tolerance_inside_cell = 1e-12;
1025
1035 std::unique_ptr<GridTools::Cache<dim, spacedim>> triangulation_cache;
1036
1037#ifdef DEAL_II_WITH_MPI
1061 void
1063 const std::map<types::subdomain_id, std::vector<particle_iterator>>
1064 &particles_to_send,
1065 const std::map<
1067 std::vector<
1069 &new_cells_for_particles = std::map<
1071 std::vector<
1073 const bool enable_cache = false);
1074
1086 void
1088 const std::map<types::subdomain_id, std::vector<particle_iterator>>
1089 &particles_to_send);
1090
1091#endif
1092
1100
1105 void
1107
1112 std::vector<boost::signals2::connection> tria_listeners;
1113
1122 void
1124
1130 void
1132
1140 void
1141 notify_ready_to_unpack(const bool serialization);
1142
1148 std::vector<char>
1151 const CellStatus status) const;
1152
1157 void
1160 const CellStatus status,
1161 const boost::iterator_range<std::vector<char>::const_iterator>
1162 &data_range);
1163
1168 typename particle_container::iterator
1170
1175 typename particle_container::iterator
1177
1182 typename particle_container::iterator
1184
1189 typename particle_container::iterator
1191 };
1192
1193
1194
1195 /* ---------------------- inline and template functions ------------------
1196 */
1197
1198 template <int dim, int spacedim>
1201 {
1202 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin();
1203 }
1204
1205
1206
1207 template <int dim, int spacedim>
1210 {
1211 return particle_iterator(particle_container_owned_begin(),
1212 *property_pool,
1213 0);
1214 }
1215
1216
1217
1218 template <int dim, int spacedim>
1221 {
1222 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end();
1223 }
1224
1225
1226
1227 template <int dim, int spacedim>
1230 {
1231 return particle_iterator(particle_container_owned_end(), *property_pool, 0);
1232 }
1233
1234
1235
1236 template <int dim, int spacedim>
1239 {
1240 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin_ghost();
1241 }
1242
1243
1244
1245 template <int dim, int spacedim>
1248 {
1249 return particle_iterator(particle_container_ghost_begin(),
1250 *property_pool,
1251 0);
1252 }
1253
1254
1255
1256 template <int dim, int spacedim>
1259 {
1260 return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end_ghost();
1261 }
1262
1263
1264
1265 template <int dim, int spacedim>
1268 {
1269 return particle_iterator(particle_container_ghost_end(), *property_pool, 0);
1270 }
1271
1272
1273
1274 template <int dim, int spacedim>
1277 {
1278 // We should always have at least the three anchor entries in the list of
1279 // particles
1280 Assert(!particles.empty(), ExcInternalError());
1281 typename particle_container::iterator begin =
1282 const_cast<particle_container &>(particles).begin();
1283 return ++begin;
1284 }
1285
1286
1287
1288 template <int dim, int spacedim>
1291 {
1292 Assert(!particles.empty(), ExcInternalError());
1293 return owned_particles_end;
1294 }
1295
1296
1297
1298 template <int dim, int spacedim>
1301 {
1302 Assert(!particles.empty(), ExcInternalError());
1303 typename particle_container::iterator begin = owned_particles_end;
1304 return ++begin;
1305 }
1306
1307
1308
1309 template <int dim, int spacedim>
1312 {
1313 Assert(!particles.empty(), ExcInternalError());
1314 typename particle_container::iterator end =
1315 const_cast<particle_container &>(particles).end();
1316 return --end;
1317 }
1318
1319
1320
1321 template <int dim, int spacedim>
1322 template <class Archive>
1323 inline void
1324 ParticleHandler<dim, spacedim>::serialize(Archive &ar, const unsigned int)
1325 {
1326 // Note that we do not serialize the particle data itself (i.e., the
1327 // 'particles' member variable). Instead we
1328 // use the serialization functionality of the triangulation class, because
1329 // this guarantees that data is immediately shipped to new processes if
1330 // the domain is distributed differently after resuming from a checkpoint.
1331 //
1332 // See @ref step_83 "step-83" for how to serialize ParticleHandler objects.
1333 ar &global_number_of_particles &global_max_particles_per_cell
1334 &next_free_particle_index;
1335 }
1336
1337
1338
1339 template <int dim, int spacedim>
1340 template <typename VectorType>
1341 inline std::enable_if_t<
1342 std::is_convertible_v<VectorType *, Function<spacedim> *> == false>
1344 const VectorType &input_vector,
1345 const bool displace_particles)
1346 {
1347 AssertDimension(input_vector.size(),
1348 get_next_free_particle_index() * spacedim);
1349 for (auto &p : *this)
1350 {
1351 Point<spacedim> &position = p.get_location();
1352 const auto id = p.get_id();
1353
1354 if (displace_particles)
1355 for (unsigned int i = 0; i < spacedim; ++i)
1356 position[i] += input_vector[id * spacedim + i];
1357 else
1358 for (unsigned int i = 0; i < spacedim; ++i)
1359 position[i] = input_vector[id * spacedim + i];
1360 }
1361 sort_particles_into_subdomains_and_cells();
1362 }
1363
1364
1365
1366 template <int dim, int spacedim>
1367 template <typename VectorType>
1368 inline void
1370 VectorType &output_vector,
1371 const bool add_to_output_vector)
1372 {
1373 AssertDimension(output_vector.size(),
1374 get_next_free_particle_index() * spacedim);
1375 for (const auto &p : *this)
1376 {
1377 auto point = p.get_location();
1378 const auto id = p.get_id();
1379 if (add_to_output_vector)
1380 for (unsigned int i = 0; i < spacedim; ++i)
1381 output_vector[id * spacedim + i] += point[i];
1382 else
1383 for (unsigned int i = 0; i < spacedim; ++i)
1384 output_vector[id * spacedim + i] = point[i];
1385 }
1386 if (add_to_output_vector)
1387 output_vector.compress(VectorOperation::add);
1388 else
1389 output_vector.compress(VectorOperation::insert);
1390 }
1391
1392} // namespace Particles
1393
1395
1396#endif
CellStatus
Definition cell_status.h:31
Abstract base class for mapping classes.
Definition mapping.h:318
std::list< ParticlesInCell > particle_container
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
particle_iterator end_ghost() const
unsigned int global_max_particles_per_cell
internal::GhostParticlePartitioner< dim, spacedim > ghost_particles_cache
particle_container::iterator particle_container_ghost_begin() 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)
particle_container::iterator particle_container_owned_end() const
particle_container::iterator particle_container_ghost_end() const
boost::iterator_range< particle_iterator > particle_iterator_range
void send_recv_particles(const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send, 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)
ParticleIterator< dim, spacedim > particle_iterator
void send_recv_particles_properties_and_location(const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send)
types::particle_index number_of_locally_owned_particles
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
types::particle_index get_max_local_particle_index() const
void reserve(const std::size_t n_particles)
particle_iterator begin() const
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 notify_ready_to_unpack(const bool serialization)
void register_load_callback_function(const bool serialization)
std::enable_if_t< std::is_convertible_v< VectorType *, Function< spacedim > * >==false > set_particle_positions(const VectorType &input_vector, const bool displace_particles=true)
std::function< std::size_t()> size_callback
std::vector< char > pack_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const 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
const particle_container::iterator owned_particles_end
particle_iterator end() const
void remove_particle(const particle_iterator &particle)
types::particle_index n_locally_owned_particles() const
void reset_particle_container(particle_container &particles)
types::particle_index n_particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
typename ParticleAccessor< dim, spacedim >::particle_container particle_container
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
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)
void remove_particles(const std::vector< particle_iterator > &particles)
types::particle_index n_global_max_particles_per_cell() const
std::vector< boost::signals2::connection > tria_listeners
particle_container::iterator particle_container_owned_begin() const
void unpack_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
std::vector< typename particle_container::iterator > cells_to_particle_cache
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
types::particle_index get_next_free_particle_index() const
IndexSet locally_owned_particle_ids() const
Definition point.h:111
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
boost::signals2::signal< void(const typename Particles::ParticleIterator< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)> particle_lost