Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07: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
particle_handler.cc
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
17
19
21
22#include <limits>
23#include <memory>
24#include <utility>
25
27
28namespace Particles
29{
30 namespace
31 {
32 template <int dim, int spacedim>
33 std::vector<char>
34 pack_particles(std::vector<ParticleIterator<dim, spacedim>> &particles)
35 {
36 std::vector<char> buffer;
37
38 if (particles.empty())
39 return buffer;
40
41 buffer.resize(particles.size() *
42 particles.front()->serialized_size_in_bytes());
43 void *current_data = buffer.data();
44
45 for (const auto &particle : particles)
46 {
47 current_data = particle->write_particle_data_to_memory(current_data);
48 }
49
50 return buffer;
51 }
52 } // namespace
53
54
55
56 template <int dim, int spacedim>
59 , mapping()
60 , property_pool(std::make_unique<PropertyPool<dim, spacedim>>(0))
61 , global_number_of_particles(0)
62 , number_of_locally_owned_particles(0)
63 , global_max_particles_per_cell(0)
64 , next_free_particle_index(0)
65 , size_callback()
66 , store_callback()
67 , load_callback()
68 , handle(numbers::invalid_unsigned_int)
69 , tria_listeners()
70 {
72 }
73
74
75
76 template <int dim, int spacedim>
79 const Mapping<dim, spacedim> &mapping,
80 const unsigned int n_properties)
81 : triangulation(&triangulation, typeid(*this).name())
82 , mapping(&mapping, typeid(*this).name())
83 , property_pool(std::make_unique<PropertyPool<dim, spacedim>>(n_properties))
84 , cells_to_particle_cache(triangulation.n_active_cells(), particles.end())
85 , global_number_of_particles(0)
86 , number_of_locally_owned_particles(0)
87 , global_max_particles_per_cell(0)
88 , next_free_particle_index(0)
89 , size_callback()
90 , store_callback()
91 , load_callback()
92 , handle(numbers::invalid_unsigned_int)
93 , triangulation_cache(
94 std::make_unique<GridTools::Cache<dim, spacedim>>(triangulation,
95 mapping))
96 , tria_listeners()
97 {
100 }
101
102
103
104 template <int dim, int spacedim>
106 {
107 clear_particles();
108
109 for (const auto &connection : tria_listeners)
110 connection.disconnect();
111 }
112
113
114
115 template <int dim, int spacedim>
116 void
118 const Triangulation<dim, spacedim> &new_triangulation,
119 const Mapping<dim, spacedim> &new_mapping,
120 const unsigned int n_properties)
121 {
122 clear();
123
124 triangulation = &new_triangulation;
125 mapping = &new_mapping;
126
127 reset_particle_container(particles);
128
129 // Create the memory pool that will store all particle properties
130 property_pool = std::make_unique<PropertyPool<dim, spacedim>>(n_properties);
131
132 // Create the grid cache to cache the information about the triangulation
133 // that is used to locate the particles into subdomains and cells
134 triangulation_cache =
135 std::make_unique<GridTools::Cache<dim, spacedim>>(new_triangulation,
136 new_mapping);
137
138 cells_to_particle_cache.resize(triangulation->n_active_cells(),
139 particles.end());
140
141 connect_to_triangulation_signals();
142 }
143
144
145
146 template <int dim, int spacedim>
147 void
149 const ParticleHandler<dim, spacedim> &particle_handler)
150 {
151 const unsigned int n_properties =
152 particle_handler.property_pool->n_properties_per_slot();
153 initialize(*particle_handler.triangulation,
154 *particle_handler.mapping,
155 n_properties);
156
157 property_pool = std::make_unique<PropertyPool<dim, spacedim>>(
158 *(particle_handler.property_pool));
159
160 // copy static members
161 global_number_of_particles = particle_handler.global_number_of_particles;
162 number_of_locally_owned_particles =
163 particle_handler.number_of_locally_owned_particles;
164
165 global_max_particles_per_cell =
166 particle_handler.global_max_particles_per_cell;
167 next_free_particle_index = particle_handler.next_free_particle_index;
168
169 // Manually copy over the particles because we do not want to touch the
170 // anchor iterators set by initialize()
171 particles.insert(particle_container_owned_end(),
172 particle_handler.particle_container_owned_begin(),
173 particle_handler.particle_container_owned_end());
174 particles.insert(particle_container_ghost_end(),
175 particle_handler.particle_container_ghost_begin(),
176 particle_handler.particle_container_ghost_end());
177
178 for (auto it = particles.begin(); it != particles.end(); ++it)
179 if (!it->particles.empty())
180 cells_to_particle_cache[it->cell->active_cell_index()] = it;
181
182 ghost_particles_cache.ghost_particles_by_domain =
183 particle_handler.ghost_particles_cache.ghost_particles_by_domain;
184 handle = particle_handler.handle;
185 }
186
187
188
189 template <int dim, int spacedim>
190 void
192 {
193 clear_particles();
194 global_number_of_particles = 0;
195 number_of_locally_owned_particles = 0;
196 next_free_particle_index = 0;
197 global_max_particles_per_cell = 0;
198 }
199
200
201
202 template <int dim, int spacedim>
203 void
205 {
206 for (auto &particles_in_cell : particles)
207 for (auto &particle : particles_in_cell.particles)
209 property_pool->deregister_particle(particle);
210
211 cells_to_particle_cache.clear();
212 reset_particle_container(particles);
213 if (triangulation != nullptr)
214 cells_to_particle_cache.resize(triangulation->n_active_cells(),
215 particles.end());
216
217 // the particle properties have already been deleted by their destructor,
218 // but the memory is still allocated. Return the memory as well.
219 property_pool->clear();
220 }
221
222
223
224 template <int dim, int spacedim>
225 void
226 ParticleHandler<dim, spacedim>::reserve(const std::size_t n_particles)
227 {
228 property_pool->reserve(n_particles);
229 }
230
231
232
233 template <int dim, int spacedim>
234 void
236 particle_container &given_particles)
237 {
238 // Make sure to set a valid past-the-end iterator also in case we have no
239 // triangulation
241 past_the_end_iterator =
242 triangulation != nullptr ?
243 triangulation->end() :
244 typename Triangulation<dim, spacedim>::cell_iterator(nullptr, -1, -1);
245
246 given_particles.clear();
247 for (unsigned int i = 0; i < 3; ++i)
248 given_particles.emplace_back(
249 std::vector<typename PropertyPool<dim, spacedim>::Handle>(),
250 past_the_end_iterator);
251
252 // Set the end of owned particles to the middle of the three elements
253 const_cast<typename particle_container::iterator &>(owned_particles_end) =
254 ++given_particles.begin();
255 }
256
257
258
259 template <int dim, int spacedim>
260 void
262 {
263 // first sort the owned particles by the active cell index
264 bool sort_is_necessary = false;
265 {
266 auto previous = particle_container_owned_begin();
267 for (auto next = previous; next != particle_container_owned_end(); ++next)
268 {
269 if (previous->cell.state() == IteratorState::valid &&
270 next->cell.state() == IteratorState::valid &&
271 previous->cell > next->cell)
272 {
273 sort_is_necessary = true;
274 break;
275 }
276 previous = next;
277 }
278 }
279 if (sort_is_necessary)
280 {
281 // we could have tried to call std::list::sort with a custom
282 // comparator to get things sorted, but things are complicated by the
283 // three anchor entries that we do not want to move, and we would
284 // hence pay an O(N log(N)) algorithm with a large constant in front
285 // of it (on the order of 20+ instructions with many
286 // difficult-to-predict branches). Therefore, we simply copy the list
287 // into a new one (keeping alive the possibly large vectors with
288 // particles on cells) into a new container.
289 particle_container sorted_particles;
290
291 // note that this call updates owned_particles_end, so that
292 // particle_container_owned_end() below already points to the
293 // new container
294 reset_particle_container(sorted_particles);
295
296 // iterate over cells and insert the entries in the new order
297 for (const auto &cell : triangulation->active_cell_iterators())
298 if (!cell->is_artificial())
299 if (cells_to_particle_cache[cell->active_cell_index()] !=
300 particles.end())
301 {
302 // before we move the sorted_particles into particles
303 // particle_container_ghost_end() still points to the
304 // old particles container. Therefore this condition looks
305 // quirky.
306 typename particle_container::iterator insert_position =
307 cell->is_locally_owned() ? particle_container_owned_end() :
308 --sorted_particles.end();
309 typename particle_container::iterator new_entry =
310 sorted_particles.insert(
311 insert_position, typename particle_container::value_type());
312 new_entry->cell = cell;
313 new_entry->particles =
314 std::move(cells_to_particle_cache[cell->active_cell_index()]
315 ->particles);
316 }
317 particles = std::move(sorted_particles);
318
319 // refresh cells_to_particle_cache
320 cells_to_particle_cache.clear();
321 cells_to_particle_cache.resize(triangulation->n_active_cells(),
322 particles.end());
323 for (auto it = particles.begin(); it != particles.end(); ++it)
324 if (!it->particles.empty())
325 cells_to_particle_cache[it->cell->active_cell_index()] = it;
326 }
327
328 // Ensure that we did not accidentally modify the anchor entries with
329 // special purpose.
330 Assert(particles.front().cell.state() == IteratorState::past_the_end &&
331 particles.front().particles.empty() &&
332 particles.back().cell.state() == IteratorState::past_the_end &&
333 particles.back().particles.empty() &&
334 owned_particles_end->cell.state() == IteratorState::past_the_end &&
335 owned_particles_end->particles.empty(),
337
338#ifdef DEBUG
339 // check that no cache element hits the three anchor states in the list of
340 // particles
341 for (const auto &it : cells_to_particle_cache)
342 Assert(it != particles.begin() && it != owned_particles_end &&
343 it != --(particles.end()),
345
346 // check that we only have locally owned particles in the first region of
347 // cells; note that we skip the very first anchor element
348 for (auto it = particle_container_owned_begin();
349 it != particle_container_owned_end();
350 ++it)
351 Assert(it->cell->is_locally_owned(), ExcInternalError());
352
353 // check that the cache is consistent with the iterators
354 std::vector<typename particle_container::iterator> verify_cache(
355 triangulation->n_active_cells(), particles.end());
356 for (auto it = particles.begin(); it != particles.end(); ++it)
357 if (!it->particles.empty())
358 verify_cache[it->cell->active_cell_index()] = it;
359
360 for (unsigned int i = 0; i < verify_cache.size(); ++i)
361 Assert(verify_cache[i] == cells_to_particle_cache[i], ExcInternalError());
362#endif
363
364 // now compute local result with the function above and then compute the
365 // collective results
366 number_of_locally_owned_particles = 0;
367
368 types::particle_index result[2] = {};
369 for (const auto &particles_in_cell : particles)
370 {
371 const types::particle_index n_particles_in_cell =
372 particles_in_cell.particles.size();
373
374 // local_max_particles_per_cell
375 result[0] = std::max(result[0], n_particles_in_cell);
376
377 // number of locally owned particles
378 if (n_particles_in_cell > 0 &&
379 particles_in_cell.cell->is_locally_owned())
380 number_of_locally_owned_particles += n_particles_in_cell;
381
382 // local_max_particle_index
383 for (const auto &particle : particles_in_cell.particles)
384 result[1] = std::max(result[1], property_pool->get_id(particle));
385 }
386
387 global_number_of_particles =
388 ::Utilities::MPI::sum(number_of_locally_owned_particles,
389 triangulation->get_communicator());
390
391 if (global_number_of_particles == 0)
392 {
393 next_free_particle_index = 0;
394 global_max_particles_per_cell = 0;
395 }
396 else
397 {
398 Utilities::MPI::max(result, triangulation->get_communicator(), result);
399
400 next_free_particle_index = result[1] + 1;
401 global_max_particles_per_cell = result[0];
402 }
403 }
404
405
406
407 template <int dim, int spacedim>
411 const
412 {
413 if (cells_to_particle_cache.empty())
414 return 0;
415
416 if (cell->is_artificial() == false)
417 {
418 return cells_to_particle_cache[cell->active_cell_index()] !=
419 particles.end() ?
420 cells_to_particle_cache[cell->active_cell_index()]
421 ->particles.size() :
422 0;
423 }
424 else
425 AssertThrow(false,
426 ExcMessage("You can't ask for the particles on an artificial "
427 "cell since we don't know what exists on these "
428 "kinds of cells."));
429
431 }
432
433
434
435 template <int dim, int spacedim>
439 const
440 {
441 return (const_cast<ParticleHandler<dim, spacedim> *>(this))
442 ->particles_in_cell(cell);
443 }
444
445
446
447 template <int dim, int spacedim>
451 {
452 const unsigned int active_cell_index = cell->active_cell_index();
453
454 if (cell->is_artificial() == false)
455 {
456 if (cells_to_particle_cache[active_cell_index] == particles.end())
457 {
458 return boost::make_iterator_range(
459 particle_iterator(particles.begin(), *property_pool, 0),
460 particle_iterator(particles.begin(), *property_pool, 0));
461 }
462 else
463 {
464 const typename particle_container::iterator
465 particles_in_current_cell =
466 cells_to_particle_cache[active_cell_index];
467 typename particle_container::iterator particles_in_next_cell =
468 particles_in_current_cell;
469 ++particles_in_next_cell;
470 return boost::make_iterator_range(
471 particle_iterator(particles_in_current_cell, *property_pool, 0),
472 particle_iterator(particles_in_next_cell, *property_pool, 0));
473 }
474 }
475 else
476 AssertThrow(false,
477 ExcMessage("You can't ask for the particles on an artificial "
478 "cell since we don't know what exists on these "
479 "kinds of cells."));
480
481 return {};
482 }
483
484
485
486 template <int dim, int spacedim>
487 void
490 {
491 auto &particles_on_cell = particle->particles_in_cell->particles;
492
493 // if the particle has an invalid handle (e.g. because it has
494 // been duplicated before calling this function) do not try
495 // to deallocate its memory again
496 auto handle = particle->get_handle();
498 property_pool->deregister_particle(handle);
499
500 // need to reduce the cached number before deleting, because the iterator
501 // may be invalid after removing the particle even if only
502 // accessing the cell
503 const auto cell = particle->get_surrounding_cell();
504 const bool owned_cell = cell->is_locally_owned();
505 if (owned_cell)
506 --number_of_locally_owned_particles;
507
508 if (particles_on_cell.size() > 1)
509 {
510 particles_on_cell[particle->particle_index_within_cell] =
511 std::move(particles_on_cell.back());
512 particles_on_cell.resize(particles_on_cell.size() - 1);
513 }
514 else
515 {
516 particles.erase(particle->particles_in_cell);
517 cells_to_particle_cache[cell->active_cell_index()] = particles.end();
518 }
519 }
520
521
522
523 template <int dim, int spacedim>
524 void
527 &particles_to_remove)
528 {
529 // We need to remove particles backwards on each cell to keep the particle
530 // iterators alive as we keep removing particles on the same cell. To
531 // ensure that this is safe, we either check if we already have sorted
532 // iterators or if we need to manually sort
533 const auto check_greater = [](const particle_iterator &a,
534 const particle_iterator &b) {
535 return a->particles_in_cell->cell > b->particles_in_cell->cell ||
536 (a->particles_in_cell->cell == b->particles_in_cell->cell &&
537 a->particle_index_within_cell > b->particle_index_within_cell);
538 };
539
540 bool particles_are_sorted = true;
541 auto previous = particles_to_remove.begin();
542 for (auto next = previous; next != particles_to_remove.end(); ++next)
543 {
544 if (check_greater(*previous, *next))
545 {
546 particles_are_sorted = false;
547 break;
548 }
549 previous = next;
550 }
551 if (particles_are_sorted)
552 {
553 // pass along backwards in array
554 for (auto it = particles_to_remove.rbegin();
555 it != particles_to_remove.rend();
556 ++it)
557 remove_particle(*it);
558 }
559 else
560 {
561 std::vector<ParticleHandler<dim, spacedim>::particle_iterator>
562 sorted_particles(particles_to_remove);
563 std::sort(sorted_particles.begin(),
564 sorted_particles.end(),
565 check_greater);
566
567 for (const auto &particle : sorted_particles)
568 remove_particle(particle);
569 }
570
571 update_cached_numbers();
572 }
573
574
575
576 template <int dim, int spacedim>
579 const Particle<dim, spacedim> &particle,
581 {
582 return insert_particle(particle.get_location(),
583 particle.get_reference_location(),
584 particle.get_id(),
585 cell,
586 particle.get_properties());
587 }
588
589
590
591 template <int dim, int spacedim>
594 const typename PropertyPool<dim, spacedim>::Handle handle,
596 {
597 const unsigned int active_cell_index = cell->active_cell_index();
598 typename particle_container::iterator &cache =
599 cells_to_particle_cache[active_cell_index];
600 if (cache == particles.end())
601 {
602 const typename particle_container::iterator insert_position =
603 cell->is_locally_owned() ? particle_container_owned_end() :
604 particle_container_ghost_end();
605 cache = particles.emplace(
606 insert_position,
607 std::vector<typename PropertyPool<dim, spacedim>::Handle>{handle},
608 cell);
609 }
610 else
611 {
612 cache->particles.push_back(handle);
613 Assert(cache->cell == cell, ExcInternalError());
614 }
615 return particle_iterator(cache,
616 *property_pool,
617 cache->particles.size() - 1);
618 }
619
620
621
622 template <int dim, int spacedim>
625 const void *&data,
627 {
629 Assert(cells_to_particle_cache.size() == triangulation->n_active_cells(),
631 Assert(cell->is_locally_owned(),
632 ExcMessage("You tried to insert particles into a cell that is not "
633 "locally owned. This is not supported."));
634
635 particle_iterator particle_it =
636 insert_particle(property_pool->register_particle(), cell);
637
638 data = particle_it->read_particle_data_from_memory(data);
639
640 ++number_of_locally_owned_particles;
641
642 return particle_it;
643 }
644
645
646
647 template <int dim, int spacedim>
650 const Point<spacedim> &position,
651 const Point<dim> &reference_position,
652 const types::particle_index particle_index,
654 const ArrayView<const double> &properties)
655 {
657 Assert(cells_to_particle_cache.size() == triangulation->n_active_cells(),
660 Assert(cell->is_locally_owned(),
661 ExcMessage("You tried to insert particles into a cell that is not "
662 "locally owned. This is not supported."));
663
664 particle_iterator particle_it =
665 insert_particle(property_pool->register_particle(), cell);
666
667 particle_it->set_location(position);
668 particle_it->set_reference_location(reference_position);
669 particle_it->set_id(particle_index);
670
671 if (properties.size() != 0)
672 particle_it->set_properties(properties);
673
674 ++number_of_locally_owned_particles;
675
676 return particle_it;
677 }
678
679
680
681 template <int dim, int spacedim>
682 void
684 const std::multimap<
686 Particle<dim, spacedim>> &new_particles)
687 {
688 reserve(n_locally_owned_particles() + new_particles.size());
689 for (const auto &cell_and_particle : new_particles)
690 insert_particle(cell_and_particle.second, cell_and_particle.first);
691
692 update_cached_numbers();
693 }
694
695
696
697 template <int dim, int spacedim>
698 void
700 const std::vector<Point<spacedim>> &positions)
701 {
703
704 update_cached_numbers();
705 reserve(n_locally_owned_particles() + positions.size());
706
707 // Determine the starting particle index of this process, which
708 // is the highest currently existing particle index plus the sum
709 // of the number of newly generated particles of all
710 // processes with a lower rank if in a parallel computation.
711 const types::particle_index local_next_particle_index =
712 get_next_free_particle_index();
713 types::particle_index local_start_index = 0;
714
715#ifdef DEAL_II_WITH_MPI
716 if (const auto parallel_triangulation =
717 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
718 &*triangulation))
719 {
720 types::particle_index particles_to_add_locally = positions.size();
721 const int ierr = MPI_Scan(&particles_to_add_locally,
722 &local_start_index,
723 1,
725 MPI_SUM,
726 parallel_triangulation->get_communicator());
727 AssertThrowMPI(ierr);
728 local_start_index -= particles_to_add_locally;
729 }
730#endif
731
732 local_start_index += local_next_particle_index;
733
734 auto point_locations =
736 positions);
737
738 auto &cells = std::get<0>(point_locations);
739 auto &local_positions = std::get<1>(point_locations);
740 auto &index_map = std::get<2>(point_locations);
741 auto &missing_points = std::get<3>(point_locations);
742 // If a point was not found, throwing an error, as the old
743 // implementation of compute_point_locations would have done
744 AssertThrow(missing_points.empty(),
746
747 (void)missing_points;
748
749 for (unsigned int i = 0; i < cells.size(); ++i)
750 for (unsigned int p = 0; p < local_positions[i].size(); ++p)
751 insert_particle(positions[index_map[i][p]],
752 local_positions[i][p],
753 local_start_index + index_map[i][p],
754 cells[i]);
755
756 update_cached_numbers();
757 }
758
759
760
761 template <int dim, int spacedim>
762 std::map<unsigned int, IndexSet>
764 const std::vector<Point<spacedim>> &positions,
765 const std::vector<std::vector<BoundingBox<spacedim>>>
766 &global_bounding_boxes,
767 const std::vector<std::vector<double>> &properties,
768 const std::vector<types::particle_index> &ids)
769 {
770 if (!properties.empty())
771 {
772 AssertDimension(properties.size(), positions.size());
773#ifdef DEBUG
774 for (const auto &p : properties)
775 AssertDimension(p.size(), n_properties_per_particle());
776#endif
777 }
778
779 if (!ids.empty())
780 AssertDimension(ids.size(), positions.size());
781
782 const auto comm = triangulation->get_communicator();
783
784 const auto n_mpi_processes = Utilities::MPI::n_mpi_processes(comm);
785
786 // Compute the global number of properties
787 const auto n_global_properties =
788 Utilities::MPI::sum(properties.size(), comm);
789
790 // Gather the number of points per processor
791 const auto n_particles_per_proc =
792 Utilities::MPI::all_gather(comm, positions.size());
793
794 // Calculate all starting points locally
795 std::vector<unsigned int> particle_start_indices(n_mpi_processes);
796
797 unsigned int particle_start_index = get_next_free_particle_index();
798 for (unsigned int process = 0; process < particle_start_indices.size();
799 ++process)
800 {
801 particle_start_indices[process] = particle_start_index;
802 particle_start_index += n_particles_per_proc[process];
803 }
804
805 // Get all local information
806 const auto cells_positions_and_index_maps =
808 positions,
809 global_bounding_boxes);
810
811 // Unpack the information into several vectors:
812 // All cells that contain at least one particle
813 const auto &local_cells_containing_particles =
814 std::get<0>(cells_positions_and_index_maps);
815
816 // The reference position of every particle in the local part of the
817 // triangulation.
818 const auto &local_reference_positions =
819 std::get<1>(cells_positions_and_index_maps);
820 // The original index in the positions vector for each particle in the
821 // local part of the triangulation
822 const auto &original_indices_of_local_particles =
823 std::get<2>(cells_positions_and_index_maps);
824 // The real spatial position of every particle in the local part of the
825 // triangulation.
826 const auto &local_positions = std::get<3>(cells_positions_and_index_maps);
827 // The MPI process that inserted each particle
828 const auto &calling_process_indices =
829 std::get<4>(cells_positions_and_index_maps);
830
831 // Create the map of cpu to indices, indicating who sent us what particle
832 std::map<unsigned int, std::vector<unsigned int>>
833 original_process_to_local_particle_indices_tmp;
834 for (unsigned int i_cell = 0;
835 i_cell < local_cells_containing_particles.size();
836 ++i_cell)
837 {
838 for (unsigned int i_particle = 0;
839 i_particle < local_positions[i_cell].size();
840 ++i_particle)
841 {
842 const unsigned int local_id_on_calling_process =
843 original_indices_of_local_particles[i_cell][i_particle];
844 const unsigned int calling_process =
845 calling_process_indices[i_cell][i_particle];
846
847 original_process_to_local_particle_indices_tmp[calling_process]
848 .push_back(local_id_on_calling_process);
849 }
850 }
851 std::map<unsigned int, IndexSet> original_process_to_local_particle_indices;
852 for (auto &process_and_particle_indices :
853 original_process_to_local_particle_indices_tmp)
854 {
855 const unsigned int calling_process = process_and_particle_indices.first;
856 original_process_to_local_particle_indices.insert(
857 {calling_process, IndexSet(n_particles_per_proc[calling_process])});
858 std::sort(process_and_particle_indices.second.begin(),
859 process_and_particle_indices.second.end());
860 original_process_to_local_particle_indices[calling_process].add_indices(
861 process_and_particle_indices.second.begin(),
862 process_and_particle_indices.second.end());
863 original_process_to_local_particle_indices[calling_process].compress();
864 }
865
866 // A map from mpi process to properties, ordered as in the IndexSet.
867 // Notice that this ordering may be different from the ordering in the
868 // vectors above, since no local ordering is guaranteed by the
869 // distribute_compute_point_locations() call.
870 // This is only filled if n_global_properties is > 0
871 std::map<unsigned int, std::vector<std::vector<double>>>
872 locally_owned_properties_from_other_processes;
873
874 // A map from mpi process to ids, ordered as in the IndexSet.
875 // Notice that this ordering may be different from the ordering in the
876 // vectors above, since no local ordering is guaranteed by the
877 // distribute_compute_point_locations() call.
878 // This is only filled if ids.size() is > 0
879 std::map<unsigned int, std::vector<types::particle_index>>
880 locally_owned_ids_from_other_processes;
881
882 if (n_global_properties > 0 || !ids.empty())
883 {
884 // Gather whom I sent my own particles to, to decide whom to send
885 // the particle properties or the ids
886 auto send_to_cpu = Utilities::MPI::some_to_some(
887 comm, original_process_to_local_particle_indices);
888
889 // Prepare the vector of properties to send
890 if (n_global_properties > 0)
891 {
892 std::map<unsigned int, std::vector<std::vector<double>>>
893 non_locally_owned_properties;
894
895 for (const auto &it : send_to_cpu)
896 {
897 std::vector<std::vector<double>> properties_to_send(
898 it.second.n_elements(),
899 std::vector<double>(n_properties_per_particle()));
900 unsigned int index = 0;
901 for (const auto el : it.second)
902 properties_to_send[index++] = properties[el];
903 non_locally_owned_properties.insert(
904 {it.first, properties_to_send});
905 }
906
907 // Send the non locally owned properties to each mpi process
908 // that needs them
909 locally_owned_properties_from_other_processes =
910 Utilities::MPI::some_to_some(comm, non_locally_owned_properties);
911
913 locally_owned_properties_from_other_processes.size(),
914 original_process_to_local_particle_indices.size());
915 }
916
917 if (!ids.empty())
918 {
919 std::map<unsigned int, std::vector<types::particle_index>>
920 non_locally_owned_ids;
921 for (const auto &it : send_to_cpu)
922 {
923 std::vector<types::particle_index> ids_to_send(
924 it.second.n_elements());
925 unsigned int index = 0;
926 for (const auto el : it.second)
927 ids_to_send[index++] = ids[el];
928 non_locally_owned_ids.insert({it.first, ids_to_send});
929 }
930
931 // Send the non locally owned ids to each mpi process
932 // that needs them
933 locally_owned_ids_from_other_processes =
934 Utilities::MPI::some_to_some(comm, non_locally_owned_ids);
935
936 AssertDimension(locally_owned_ids_from_other_processes.size(),
937 original_process_to_local_particle_indices.size());
938 }
939 }
940
941 // Now fill up the actual particles
942 for (unsigned int i_cell = 0;
943 i_cell < local_cells_containing_particles.size();
944 ++i_cell)
945 {
946 for (unsigned int i_particle = 0;
947 i_particle < local_positions[i_cell].size();
948 ++i_particle)
949 {
950 const unsigned int local_id_on_calling_process =
951 original_indices_of_local_particles[i_cell][i_particle];
952
953 const unsigned int calling_process =
954 calling_process_indices[i_cell][i_particle];
955
956 const unsigned int index_within_set =
957 original_process_to_local_particle_indices[calling_process]
958 .index_within_set(local_id_on_calling_process);
959
960 const unsigned int particle_id =
961 ids.empty() ?
962 local_id_on_calling_process +
963 particle_start_indices[calling_process] :
964 locally_owned_ids_from_other_processes[calling_process]
965 [index_within_set];
966
967 auto particle_it =
968 insert_particle(local_positions[i_cell][i_particle],
969 local_reference_positions[i_cell][i_particle],
970 particle_id,
971 local_cells_containing_particles[i_cell]);
972
973 if (n_global_properties > 0)
974 {
975 particle_it->set_properties(
976 locally_owned_properties_from_other_processes
977 [calling_process][index_within_set]);
978 }
979 }
980 }
981
982 update_cached_numbers();
983
984 return original_process_to_local_particle_indices;
985 }
986
987
988
989 template <int dim, int spacedim>
990 std::map<unsigned int, IndexSet>
992 const std::vector<Particle<dim, spacedim>> &particles,
993 const std::vector<std::vector<BoundingBox<spacedim>>>
994 &global_bounding_boxes)
995 {
996 // Store the positions in a vector of points, the ids in a vector of ids,
997 // and the properties, if any, in a vector of vector of properties.
998 std::vector<Point<spacedim>> positions;
999 std::vector<std::vector<double>> properties;
1000 std::vector<types::particle_index> ids;
1001 positions.resize(particles.size());
1002 ids.resize(particles.size());
1003 if (n_properties_per_particle() > 0)
1004 properties.resize(particles.size(),
1005 std::vector<double>(n_properties_per_particle()));
1006
1007 unsigned int i = 0;
1008 for (const auto &p : particles)
1009 {
1010 positions[i] = p.get_location();
1011 ids[i] = p.get_id();
1012 if (p.has_properties())
1013 properties[i] = {p.get_properties().begin(),
1014 p.get_properties().end()};
1015 ++i;
1016 }
1017
1018 return insert_global_particles(positions,
1019 global_bounding_boxes,
1020 properties,
1021 ids);
1022 }
1023
1024
1025
1026 template <int dim, int spacedim>
1029 {
1030 return global_number_of_particles;
1031 }
1032
1033
1034
1035 template <int dim, int spacedim>
1038 {
1039 return global_max_particles_per_cell;
1040 }
1041
1042
1043
1044 template <int dim, int spacedim>
1047 {
1048 return number_of_locally_owned_particles;
1049 }
1050
1051
1052
1053 template <int dim, int spacedim>
1054 unsigned int
1056 {
1057 return property_pool->n_properties_per_slot();
1058 }
1059
1060
1061
1062 template <int dim, int spacedim>
1065 {
1066 return next_free_particle_index;
1067 }
1068
1069
1070
1071 template <int dim, int spacedim>
1072 IndexSet
1074 {
1075 IndexSet set(get_next_free_particle_index());
1076 std::vector<types::particle_index> indices;
1077 indices.reserve(n_locally_owned_particles());
1078 for (const auto &p : *this)
1079 indices.push_back(p.get_id());
1080 set.add_indices(indices.begin(), indices.end());
1081 set.compress();
1082 return set;
1083 }
1084
1085
1086
1087 template <int dim, int spacedim>
1090 {
1091 return property_pool->n_slots();
1092 }
1093
1094
1095
1096 template <int dim, int spacedim>
1097 void
1099 std::vector<Point<spacedim>> &positions,
1100 const bool add_to_output_vector)
1101 {
1102 // There should be one point per particle to gather
1103 AssertDimension(positions.size(), n_locally_owned_particles());
1104
1105 unsigned int i = 0;
1106 for (auto it = begin(); it != end(); ++it, ++i)
1107 {
1108 if (add_to_output_vector)
1109 positions[i] = positions[i] + it->get_location();
1110 else
1111 positions[i] = it->get_location();
1112 }
1113 }
1114
1115
1116
1117 template <int dim, int spacedim>
1118 void
1120 const std::vector<Point<spacedim>> &new_positions,
1121 const bool displace_particles)
1122 {
1123 // There should be one point per particle to fix the new position
1124 AssertDimension(new_positions.size(), n_locally_owned_particles());
1125
1126 unsigned int i = 0;
1127 for (auto it = begin(); it != end(); ++it, ++i)
1128 {
1129 Point<spacedim> &location = it->get_location();
1130 if (displace_particles)
1131 location += new_positions[i];
1132 else
1133 location = new_positions[i];
1134 }
1135 sort_particles_into_subdomains_and_cells();
1136 }
1137
1138
1139
1140 template <int dim, int spacedim>
1141 void
1143 const Function<spacedim> &function,
1144 const bool displace_particles)
1145 {
1146 // The function should have sufficient components to displace the
1147 // particles
1148 AssertDimension(function.n_components, spacedim);
1149
1150 Vector<double> new_position(spacedim);
1151 for (auto &particle : *this)
1152 {
1153 Point<spacedim> &particle_location = particle.get_location();
1154 function.vector_value(particle_location, new_position);
1155 if (displace_particles)
1156 for (unsigned int d = 0; d < spacedim; ++d)
1157 particle_location[d] += new_position[d];
1158 else
1159 for (unsigned int d = 0; d < spacedim; ++d)
1160 particle_location[d] = new_position[d];
1161 }
1162 sort_particles_into_subdomains_and_cells();
1163 }
1164
1165
1166
1167 template <int dim, int spacedim>
1170 {
1171 return *property_pool;
1172 }
1173
1174
1175
1176 namespace
1177 {
1187 template <int dim>
1188 bool
1189 compare_particle_association(
1190 const unsigned int a,
1191 const unsigned int b,
1192 const Tensor<1, dim> &particle_direction,
1193 const std::vector<Tensor<1, dim>> &center_directions)
1194 {
1195 const double scalar_product_a = center_directions[a] * particle_direction;
1196 const double scalar_product_b = center_directions[b] * particle_direction;
1197
1198 // The function is supposed to return if a is before b. We are looking
1199 // for the alignment of particle direction and center direction,
1200 // therefore return if the scalar product of a is larger.
1201 return (scalar_product_a > scalar_product_b);
1202 }
1203 } // namespace
1204
1205
1206
1207 template <int dim, int spacedim>
1208 void
1210 {
1211 Assert(triangulation != nullptr, ExcInternalError());
1212 Assert(cells_to_particle_cache.size() == triangulation->n_active_cells(),
1214
1215 // TODO: The current algorithm only works for particles that are in
1216 // the local domain or in ghost cells, because it only knows the
1217 // subdomain_id of ghost cells, but not of artificial cells. This
1218 // can be extended using the distributed version of compute point
1219 // locations.
1220 // TODO: Extend this function to allow keeping particles on other
1221 // processes around (with an invalid cell).
1222
1223 std::vector<particle_iterator> particles_out_of_cell;
1224
1225 // Reserve some space for particles that need sorting to avoid frequent
1226 // re-allocation. Guess 25% of particles need sorting. Balance memory
1227 // overhead and performance.
1228 particles_out_of_cell.reserve(n_locally_owned_particles() / 4);
1229
1230 // Now update the reference locations of the moved particles
1231 std::vector<Point<spacedim>> real_locations;
1232 std::vector<Point<dim>> reference_locations;
1233 real_locations.reserve(global_max_particles_per_cell);
1234 reference_locations.reserve(global_max_particles_per_cell);
1235
1236 for (const auto &cell : triangulation->active_cell_iterators())
1237 {
1238 // Particles can be inserted into arbitrary cells, e.g. if their cell is
1239 // not known. However, for artificial cells we can not evaluate
1240 // the reference position of particles. Do not sort particles that are
1241 // not locally owned, because they will be sorted by the process that
1242 // owns them.
1243 if (cell->is_locally_owned() == false)
1244 {
1245 continue;
1246 }
1247
1248 const unsigned int n_pic = n_particles_in_cell(cell);
1249 auto pic = particles_in_cell(cell);
1250
1251 real_locations.clear();
1252 for (const auto &particle : pic)
1253 real_locations.push_back(particle.get_location());
1254
1255 reference_locations.resize(n_pic);
1256 mapping->transform_points_real_to_unit_cell(cell,
1257 real_locations,
1258 reference_locations);
1259
1260 auto particle = pic.begin();
1261 for (const auto &p_unit : reference_locations)
1262 {
1263 if (numbers::is_finite(p_unit[0]) &&
1265 tolerance_inside_cell))
1266 particle->set_reference_location(p_unit);
1267 else
1268 particles_out_of_cell.push_back(particle);
1269
1270 ++particle;
1271 }
1272 }
1273
1274 // There are three reasons why a particle is not in its old cell:
1275 // It moved to another cell, to another subdomain or it left the mesh.
1276 // Particles that moved to another cell are updated and moved inside the
1277 // particles vector, particles that moved to another domain are
1278 // collected in the moved_particles_domain vector. Particles that left
1279 // the mesh completely are ignored and removed.
1280 std::map<types::subdomain_id, std::vector<particle_iterator>>
1281 moved_particles;
1282 std::map<
1284 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1285 moved_cells;
1286
1287 // We do not know exactly how many particles are lost, exchanged between
1288 // domains, or remain on this process. Therefore we pre-allocate
1289 // approximate sizes for these vectors. If more space is needed an
1290 // automatic and relatively fast (compared to other parts of this
1291 // algorithm) re-allocation will happen.
1292 using vector_size = typename std::vector<particle_iterator>::size_type;
1293
1294 std::set<types::subdomain_id> ghost_owners;
1295 if (const auto parallel_triangulation =
1296 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1297 &*triangulation))
1298 ghost_owners = parallel_triangulation->ghost_owners();
1299
1300 // Reserve some space for particles that need communication to avoid
1301 // frequent re-allocation. Guess 25% of particles out of their old cell need
1302 // communication. Balance memory overhead and performance.
1303 for (const auto &ghost_owner : ghost_owners)
1304 moved_particles[ghost_owner].reserve(particles_out_of_cell.size() / 4);
1305 for (const auto &ghost_owner : ghost_owners)
1306 moved_cells[ghost_owner].reserve(particles_out_of_cell.size() / 4);
1307
1308 {
1309 // Create a map from vertices to adjacent cells using grid cache
1310 const std::vector<
1311 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1312 &vertex_to_cells = triangulation_cache->get_vertex_to_cell_map();
1313
1314 // Create a corresponding map of vectors from vertex to cell center
1315 // using grid cache
1316 const std::vector<std::vector<Tensor<1, spacedim>>>
1317 &vertex_to_cell_centers =
1318 triangulation_cache->get_vertex_to_cell_centers_directions();
1319
1320 std::vector<unsigned int> neighbor_permutation;
1321
1322 // Reuse these vectors below, but only with a single element.
1323 // Avoid resizing for every particle.
1324 Point<dim> invalid_reference_point;
1325 Point<spacedim> invalid_point;
1326 invalid_reference_point[0] = std::numeric_limits<double>::infinity();
1327 invalid_point[0] = std::numeric_limits<double>::infinity();
1328 reference_locations.resize(1, invalid_reference_point);
1329 real_locations.resize(1, invalid_point);
1330
1331 // Find the cells that the particles moved to.
1332 for (auto &out_particle : particles_out_of_cell)
1333 {
1334 // make a copy of the current cell, since we will modify the
1335 // variable current_cell below, but we need the original in
1336 // the case the particle is not found
1337 auto current_cell = out_particle->get_surrounding_cell();
1338
1339 real_locations[0] = out_particle->get_location();
1340
1341 // Record if the new cell was found
1342 bool found_cell = false;
1343
1344 // Check if the particle is in one of the old cell's neighbors
1345 // that are adjacent to the closest vertex
1346 const unsigned int closest_vertex =
1347 GridTools::find_closest_vertex_of_cell<dim, spacedim>(
1348 current_cell, out_particle->get_location(), *mapping);
1349 Tensor<1, spacedim> vertex_to_particle =
1350 out_particle->get_location() - current_cell->vertex(closest_vertex);
1351 vertex_to_particle /= vertex_to_particle.norm();
1352
1353 const unsigned int closest_vertex_index =
1354 current_cell->vertex_index(closest_vertex);
1355 const unsigned int n_neighbor_cells =
1356 vertex_to_cells[closest_vertex_index].size();
1357
1358 neighbor_permutation.resize(n_neighbor_cells);
1359 for (unsigned int i = 0; i < n_neighbor_cells; ++i)
1360 neighbor_permutation[i] = i;
1361
1362 const auto &cell_centers =
1363 vertex_to_cell_centers[closest_vertex_index];
1364 std::sort(neighbor_permutation.begin(),
1365 neighbor_permutation.end(),
1366 [&vertex_to_particle, &cell_centers](const unsigned int a,
1367 const unsigned int b) {
1368 return compare_particle_association(a,
1369 b,
1370 vertex_to_particle,
1371 cell_centers);
1372 });
1373
1374 // Search all of the cells adjacent to the closest vertex of the
1375 // previous cell. Most likely we will find the particle in them.
1376 for (unsigned int i = 0; i < n_neighbor_cells; ++i)
1377 {
1378 typename std::set<typename Triangulation<dim, spacedim>::
1379 active_cell_iterator>::const_iterator cell =
1380 vertex_to_cells[closest_vertex_index].begin();
1381
1382 std::advance(cell, neighbor_permutation[i]);
1383 mapping->transform_points_real_to_unit_cell(*cell,
1384 real_locations,
1385 reference_locations);
1386
1387 if (GeometryInfo<dim>::is_inside_unit_cell(reference_locations[0],
1388 tolerance_inside_cell))
1389 {
1390 current_cell = *cell;
1391 found_cell = true;
1392 break;
1393 }
1394 }
1395
1396 if (!found_cell)
1397 {
1398 // For some clang-based compilers and boost versions the call to
1399 // RTree::query doesn't compile. We use a slower implementation as
1400 // workaround.
1401 // This is fixed in boost in
1402 // https://github.com/boostorg/numeric_conversion/commit/50a1eae942effb0a9b90724323ef8f2a67e7984a
1403#if defined(DEAL_II_WITH_BOOST_BUNDLED) || \
1404 !(defined(__clang_major__) && __clang_major__ >= 16) || \
1405 BOOST_VERSION >= 108100
1406 // The particle is not in a neighbor of the old cell.
1407 // Look for the new cell in the whole local domain.
1408 // This case is rare.
1409 std::vector<std::pair<Point<spacedim>, unsigned int>>
1410 closest_vertex_in_domain;
1411 triangulation_cache->get_used_vertices_rtree().query(
1412 boost::geometry::index::nearest(out_particle->get_location(),
1413 1),
1414 std::back_inserter(closest_vertex_in_domain));
1415
1416 // We should have one and only one result
1417 AssertDimension(closest_vertex_in_domain.size(), 1);
1418 const unsigned int closest_vertex_index_in_domain =
1419 closest_vertex_in_domain[0].second;
1420#else
1421 const unsigned int closest_vertex_index_in_domain =
1424 out_particle->get_location());
1425#endif
1426
1427 // Search all of the cells adjacent to the closest vertex of the
1428 // domain. Most likely we will find the particle in them.
1429 for (const auto &cell :
1430 vertex_to_cells[closest_vertex_index_in_domain])
1431 {
1432 mapping->transform_points_real_to_unit_cell(
1433 cell, real_locations, reference_locations);
1434
1436 reference_locations[0], tolerance_inside_cell))
1437 {
1438 current_cell = cell;
1439 found_cell = true;
1440 break;
1441 }
1442 }
1443 }
1444
1445 if (!found_cell)
1446 {
1447 // We can find no cell for this particle. It has left the
1448 // domain due to an integration error or an open boundary.
1449 // Signal the loss and move on.
1450 signals.particle_lost(out_particle,
1451 out_particle->get_surrounding_cell());
1452 continue;
1453 }
1454
1455 // If we are here, we found a cell and reference position for this
1456 // particle
1457 out_particle->set_reference_location(reference_locations[0]);
1458
1459 // Reinsert the particle into our domain if we own its cell.
1460 // Mark it for MPI transfer otherwise
1461 if (current_cell->is_locally_owned())
1462 {
1464 out_particle->particles_in_cell
1465 ->particles[out_particle->particle_index_within_cell];
1466
1467 // Avoid deallocating the memory of this particle
1468 const auto old_value = old;
1470
1471 // Allocate particle with the old handle
1472 insert_particle(old_value, current_cell);
1473 }
1474 else
1475 {
1476 moved_particles[current_cell->subdomain_id()].push_back(
1477 out_particle);
1478 moved_cells[current_cell->subdomain_id()].push_back(current_cell);
1479 }
1480 }
1481 }
1482
1483 // Exchange particles between processors if we have more than one process
1484#ifdef DEAL_II_WITH_MPI
1485 if (const auto parallel_triangulation =
1486 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1487 &*triangulation))
1488 {
1490 parallel_triangulation->get_communicator()) > 1)
1491 send_recv_particles(moved_particles, moved_cells);
1492 }
1493#endif
1494
1495 // remove_particles also calls update_cached_numbers()
1496 remove_particles(particles_out_of_cell);
1497
1498 // now make sure particle data is sorted in order of iteration
1499 std::vector<typename PropertyPool<dim, spacedim>::Handle> unsorted_handles;
1500 unsorted_handles.reserve(property_pool->n_registered_slots());
1501
1502 typename PropertyPool<dim, spacedim>::Handle sorted_handle = 0;
1503 for (auto &particles_in_cell : particles)
1504 for (auto &particle : particles_in_cell.particles)
1505 {
1506 unsorted_handles.push_back(particle);
1507 particle = sorted_handle++;
1508 }
1509
1510 property_pool->sort_memory_slots(unsorted_handles);
1511
1512 } // namespace Particles
1513
1514
1515
1516 template <int dim, int spacedim>
1517 void
1519 const bool enable_cache)
1520 {
1521 // Nothing to do in serial computations
1522 const auto parallel_triangulation =
1523 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1524 &*triangulation);
1525 if (parallel_triangulation != nullptr)
1526 {
1528 parallel_triangulation->get_communicator()) == 1)
1529 return;
1530 }
1531 else
1532 return;
1533
1534#ifndef DEAL_II_WITH_MPI
1535 (void)enable_cache;
1536#else
1537 // Clear ghost particles and their properties
1538 for (const auto &cell : triangulation->active_cell_iterators())
1539 if (cell->is_ghost() &&
1540 cells_to_particle_cache[cell->active_cell_index()] != particles.end())
1541 {
1542 Assert(cells_to_particle_cache[cell->active_cell_index()]->cell ==
1543 cell,
1545 // Clear particle properties
1546 for (auto &ghost_particle :
1547 cells_to_particle_cache[cell->active_cell_index()]->particles)
1548 property_pool->deregister_particle(ghost_particle);
1549
1550 // Clear particles themselves
1551 particles.erase(cells_to_particle_cache[cell->active_cell_index()]);
1552 cells_to_particle_cache[cell->active_cell_index()] = particles.end();
1553 }
1554
1555 // Clear ghost particles cache and invalidate it
1556 ghost_particles_cache.ghost_particles_by_domain.clear();
1557 ghost_particles_cache.valid = false;
1558
1559 // In the case of a parallel simulation with periodic boundary conditions
1560 // the vertices associated with periodic boundaries are not directly
1561 // connected to the ghost cells but they are connected to the ghost cells
1562 // through their coinciding vertices. We gather this information using the
1563 // vertices_with_ghost_neighbors map
1564 const std::map<unsigned int, std::set<types::subdomain_id>>
1566 triangulation_cache->get_vertices_with_ghost_neighbors();
1567
1568 const std::set<types::subdomain_id> ghost_owners =
1569 parallel_triangulation->ghost_owners();
1570 for (const auto ghost_owner : ghost_owners)
1571 ghost_particles_cache.ghost_particles_by_domain[ghost_owner].reserve(
1572 n_locally_owned_particles() / 4);
1573
1574 const std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain =
1575 triangulation_cache->get_vertex_to_neighbor_subdomain();
1576
1577 for (const auto &cell : triangulation->active_cell_iterators())
1578 {
1579 if (cell->is_locally_owned())
1580 {
1581 std::set<unsigned int> cell_to_neighbor_subdomain;
1582 for (const unsigned int v : cell->vertex_indices())
1583 {
1584 const auto vertex_ghost_neighbors =
1585 vertices_with_ghost_neighbors.find(cell->vertex_index(v));
1586 if (vertex_ghost_neighbors !=
1588 {
1589 cell_to_neighbor_subdomain.insert(
1590 vertex_ghost_neighbors->second.begin(),
1591 vertex_ghost_neighbors->second.end());
1592 }
1593 }
1594
1595 if (cell_to_neighbor_subdomain.size() > 0)
1596 {
1597 const particle_iterator_range particle_range =
1598 particles_in_cell(cell);
1599
1600 for (const auto domain : cell_to_neighbor_subdomain)
1601 {
1602 for (typename particle_iterator_range::iterator particle =
1603 particle_range.begin();
1604 particle != particle_range.end();
1605 ++particle)
1606 ghost_particles_cache.ghost_particles_by_domain[domain]
1607 .push_back(particle);
1608 }
1609 }
1610 }
1611 }
1612
1613 send_recv_particles(
1614 ghost_particles_cache.ghost_particles_by_domain,
1615 std::map<
1617 std::vector<
1619 enable_cache);
1620#endif
1621 }
1622
1623
1624
1625 template <int dim, int spacedim>
1626 void
1628 {
1629 // Nothing to do in serial computations
1630 const auto parallel_triangulation =
1631 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1632 &*triangulation);
1633 if (parallel_triangulation == nullptr ||
1635 parallel_triangulation->get_communicator()) == 1)
1636 {
1637 return;
1638 }
1639
1640
1641#ifdef DEAL_II_WITH_MPI
1642 // First clear the current ghost_particle information
1643 // ghost_particles.clear();
1644 Assert(ghost_particles_cache.valid,
1645 ExcMessage(
1646 "Ghost particles cannot be updated if they first have not been "
1647 "exchanged at least once with the cache enabled"));
1648
1649
1650 send_recv_particles_properties_and_location(
1651 ghost_particles_cache.ghost_particles_by_domain);
1652#endif
1653 }
1654
1655
1656
1657#ifdef DEAL_II_WITH_MPI
1658 template <int dim, int spacedim>
1659 void
1661 const std::map<types::subdomain_id, std::vector<particle_iterator>>
1662 &particles_to_send,
1663 const std::map<
1666 &send_cells,
1667 const bool build_cache)
1668 {
1669 Assert(triangulation != nullptr, ExcInternalError());
1670 Assert(cells_to_particle_cache.size() == triangulation->n_active_cells(),
1672
1673 ghost_particles_cache.valid = build_cache;
1674
1675 const auto parallel_triangulation =
1676 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1677 &*triangulation);
1678 Assert(parallel_triangulation,
1679 ExcMessage("This function is only implemented for "
1680 "parallel::TriangulationBase objects."));
1681
1682 // Determine the communication pattern
1683 const std::set<types::subdomain_id> ghost_owners =
1684 parallel_triangulation->ghost_owners();
1685 const std::vector<types::subdomain_id> neighbors(ghost_owners.begin(),
1686 ghost_owners.end());
1687 const unsigned int n_neighbors = neighbors.size();
1688
1689 if (send_cells.size() != 0)
1690 Assert(particles_to_send.size() == send_cells.size(), ExcInternalError());
1691
1692 // If we do not know the subdomain this particle needs to be send to,
1693 // throw an error
1694 Assert(particles_to_send.find(numbers::artificial_subdomain_id) ==
1695 particles_to_send.end(),
1697
1698 // TODO: Implement the shipping of particles to processes that are not
1699 // ghost owners of the local domain
1700 for (auto send_particles = particles_to_send.begin();
1701 send_particles != particles_to_send.end();
1702 ++send_particles)
1703 Assert(ghost_owners.find(send_particles->first) != ghost_owners.end(),
1705
1706 std::size_t n_send_particles = 0;
1707 for (auto send_particles = particles_to_send.begin();
1708 send_particles != particles_to_send.end();
1709 ++send_particles)
1710 n_send_particles += send_particles->second.size();
1711
1712 const unsigned int cellid_size = sizeof(CellId::binary_type);
1713
1714 // Containers for the amount and offsets of data we will send
1715 // to other processors and the data itself.
1716 std::vector<unsigned int> n_send_data(n_neighbors, 0);
1717 std::vector<unsigned int> send_offsets(n_neighbors, 0);
1718 std::vector<char> send_data;
1719
1720 Particle<dim, spacedim> test_particle;
1721 test_particle.set_property_pool(*property_pool);
1722
1723 const unsigned int individual_particle_data_size =
1724 test_particle.serialized_size_in_bytes() +
1725 (size_callback ? size_callback() : 0);
1726
1727 const unsigned int individual_total_particle_data_size =
1728 individual_particle_data_size + cellid_size;
1729
1730 // Only serialize things if there are particles to be send.
1731 // We can not return early even if no particles
1732 // are send, because we might receive particles from other processes
1733 if (n_send_particles > 0)
1734 {
1735 // Allocate space for sending particle data
1736 send_data.resize(n_send_particles *
1737 individual_total_particle_data_size);
1738
1739 void *data = static_cast<void *>(&send_data.front());
1740
1741 // Serialize the data sorted by receiving process
1742 for (unsigned int i = 0; i < n_neighbors; ++i)
1743 {
1744 send_offsets[i] = reinterpret_cast<std::size_t>(data) -
1745 reinterpret_cast<std::size_t>(&send_data.front());
1746
1747 const unsigned int n_particles_to_send =
1748 particles_to_send.at(neighbors[i]).size();
1749
1750 Assert(static_cast<std::size_t>(n_particles_to_send) *
1751 individual_total_particle_data_size ==
1752 static_cast<std::size_t>(
1753 n_particles_to_send *
1754 individual_total_particle_data_size),
1755 ExcMessage("Overflow when trying to send particle "
1756 "data"));
1757
1758 for (unsigned int j = 0; j < n_particles_to_send; ++j)
1759 {
1760 // If no target cells are given, use the iterator
1761 // information
1763 cell;
1764 if (send_cells.empty())
1765 cell = particles_to_send.at(neighbors[i])[j]
1766 ->get_surrounding_cell();
1767 else
1768 cell = send_cells.at(neighbors[i])[j];
1769
1770 const CellId::binary_type cellid =
1771 cell->id().template to_binary<dim>();
1772 memcpy(data, &cellid, cellid_size);
1773 data = static_cast<char *>(data) + cellid_size;
1774
1775 data = particles_to_send.at(neighbors[i])[j]
1776 ->write_particle_data_to_memory(data);
1777 if (store_callback)
1778 data =
1779 store_callback(particles_to_send.at(neighbors[i])[j], data);
1780 }
1781 n_send_data[i] = n_particles_to_send;
1782 }
1783 }
1784
1785 // Containers for the data we will receive from other processors
1786 std::vector<unsigned int> n_recv_data(n_neighbors);
1787 std::vector<unsigned int> recv_offsets(n_neighbors);
1788
1789 {
1790 const int mpi_tag = Utilities::MPI::internal::Tags::
1792
1793 std::vector<MPI_Request> n_requests(2 * n_neighbors);
1794 for (unsigned int i = 0; i < n_neighbors; ++i)
1795 {
1796 const int ierr = MPI_Irecv(&(n_recv_data[i]),
1797 1,
1798 MPI_UNSIGNED,
1799 neighbors[i],
1800 mpi_tag,
1801 parallel_triangulation->get_communicator(),
1802 &(n_requests[2 * i]));
1803 AssertThrowMPI(ierr);
1804 }
1805 for (unsigned int i = 0; i < n_neighbors; ++i)
1806 {
1807 const int ierr = MPI_Isend(&(n_send_data[i]),
1808 1,
1809 MPI_UNSIGNED,
1810 neighbors[i],
1811 mpi_tag,
1812 parallel_triangulation->get_communicator(),
1813 &(n_requests[2 * i + 1]));
1814 AssertThrowMPI(ierr);
1815 }
1816 const int ierr =
1817 MPI_Waitall(2 * n_neighbors, n_requests.data(), MPI_STATUSES_IGNORE);
1818 AssertThrowMPI(ierr);
1819 }
1820
1821 // Determine how many particles and data we will receive
1822 unsigned int total_recv_data = 0;
1823 for (unsigned int neighbor_id = 0; neighbor_id < n_neighbors; ++neighbor_id)
1824 {
1825 recv_offsets[neighbor_id] = total_recv_data;
1826 total_recv_data +=
1827 n_recv_data[neighbor_id] * individual_total_particle_data_size;
1828 }
1829
1830 // Set up the space for the received particle data
1831 std::vector<char> recv_data(total_recv_data);
1832
1833 // Exchange the particle data between domains
1834 {
1835 std::vector<MPI_Request> requests(2 * n_neighbors);
1836 unsigned int send_ops = 0;
1837 unsigned int recv_ops = 0;
1838
1839 const int mpi_tag = Utilities::MPI::internal::Tags::
1841
1842 for (unsigned int i = 0; i < n_neighbors; ++i)
1843 if (n_recv_data[i] > 0)
1844 {
1845 const int ierr =
1846 MPI_Irecv(&(recv_data[recv_offsets[i]]),
1847 n_recv_data[i] * individual_total_particle_data_size,
1848 MPI_CHAR,
1849 neighbors[i],
1850 mpi_tag,
1851 parallel_triangulation->get_communicator(),
1852 &(requests[send_ops]));
1853 AssertThrowMPI(ierr);
1854 ++send_ops;
1855 }
1856
1857 for (unsigned int i = 0; i < n_neighbors; ++i)
1858 if (n_send_data[i] > 0)
1859 {
1860 const int ierr =
1861 MPI_Isend(&(send_data[send_offsets[i]]),
1862 n_send_data[i] * individual_total_particle_data_size,
1863 MPI_CHAR,
1864 neighbors[i],
1865 mpi_tag,
1866 parallel_triangulation->get_communicator(),
1867 &(requests[send_ops + recv_ops]));
1868 AssertThrowMPI(ierr);
1869 ++recv_ops;
1870 }
1871 const int ierr =
1872 MPI_Waitall(send_ops + recv_ops, requests.data(), MPI_STATUSES_IGNORE);
1873 AssertThrowMPI(ierr);
1874 }
1875
1876 // Put the received particles into the domain if they are in the
1877 // triangulation
1878 const void *recv_data_it = static_cast<const void *>(recv_data.data());
1879
1880 // Store the particle iterators in the cache
1881 auto &ghost_particles_iterators =
1882 ghost_particles_cache.ghost_particles_iterators;
1883
1884 if (build_cache)
1885 {
1886 ghost_particles_iterators.clear();
1887
1888 auto &send_pointers_particles = ghost_particles_cache.send_pointers;
1889 send_pointers_particles.assign(n_neighbors + 1, 0);
1890
1891 for (unsigned int i = 0; i < n_neighbors; ++i)
1892 send_pointers_particles[i + 1] =
1893 send_pointers_particles[i] +
1894 n_send_data[i] * individual_particle_data_size;
1895
1896 auto &recv_pointers_particles = ghost_particles_cache.recv_pointers;
1897 recv_pointers_particles.assign(n_neighbors + 1, 0);
1898
1899 for (unsigned int i = 0; i < n_neighbors; ++i)
1900 recv_pointers_particles[i + 1] =
1901 recv_pointers_particles[i] +
1902 n_recv_data[i] * individual_particle_data_size;
1903
1904 ghost_particles_cache.neighbors = neighbors;
1905
1906 ghost_particles_cache.send_data.resize(
1907 ghost_particles_cache.send_pointers.back());
1908 ghost_particles_cache.recv_data.resize(
1909 ghost_particles_cache.recv_pointers.back());
1910 }
1911
1912 while (reinterpret_cast<std::size_t>(recv_data_it) -
1913 reinterpret_cast<std::size_t>(recv_data.data()) <
1914 total_recv_data)
1915 {
1916 CellId::binary_type binary_cellid;
1917 memcpy(&binary_cellid, recv_data_it, cellid_size);
1918 const CellId id(binary_cellid);
1919 recv_data_it = static_cast<const char *>(recv_data_it) + cellid_size;
1920
1922 triangulation->create_cell_iterator(id);
1923
1924 insert_particle(property_pool->register_particle(), cell);
1925 const typename particle_container::iterator &cache =
1926 cells_to_particle_cache[cell->active_cell_index()];
1927 Assert(cache->cell == cell, ExcInternalError());
1928
1929 particle_iterator particle_it(cache,
1930 *property_pool,
1931 cache->particles.size() - 1);
1932
1933 recv_data_it =
1934 particle_it->read_particle_data_from_memory(recv_data_it);
1935
1936 if (load_callback)
1937 recv_data_it = load_callback(particle_it, recv_data_it);
1938
1939 if (build_cache) // TODO: is this safe?
1940 ghost_particles_iterators.push_back(particle_it);
1941 }
1942
1943 AssertThrow(recv_data_it == recv_data.data() + recv_data.size(),
1944 ExcMessage(
1945 "The amount of data that was read into new particles "
1946 "does not match the amount of data sent around."));
1947 }
1948#endif
1949
1950
1951
1952#ifdef DEAL_II_WITH_MPI
1953 template <int dim, int spacedim>
1954 void
1956 const std::map<types::subdomain_id, std::vector<particle_iterator>>
1957 &particles_to_send)
1958 {
1959 const auto parallel_triangulation =
1960 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
1961 &*triangulation);
1962 Assert(
1963 parallel_triangulation,
1964 ExcMessage(
1965 "This function is only implemented for parallel::TriangulationBase "
1966 "objects."));
1967
1968 const auto &neighbors = ghost_particles_cache.neighbors;
1969 const auto &send_pointers = ghost_particles_cache.send_pointers;
1970 const auto &recv_pointers = ghost_particles_cache.recv_pointers;
1971
1972 std::vector<char> &send_data = ghost_particles_cache.send_data;
1973
1974 // Fill data to send
1975 if (send_pointers.back() > 0)
1976 {
1977 void *data = static_cast<void *>(&send_data.front());
1978
1979 // Serialize the data sorted by receiving process
1980 for (const auto i : neighbors)
1981 for (const auto &p : particles_to_send.at(i))
1982 {
1983 data = p->write_particle_data_to_memory(data);
1984 if (store_callback)
1985 data = store_callback(p, data);
1986 }
1987 }
1988
1989 std::vector<char> &recv_data = ghost_particles_cache.recv_data;
1990
1991 // Exchange the particle data between domains
1992 {
1993 std::vector<MPI_Request> requests(2 * neighbors.size());
1994 unsigned int send_ops = 0;
1995 unsigned int recv_ops = 0;
1996
1997 const int mpi_tag = Utilities::MPI::internal::Tags::
1999
2000 for (unsigned int i = 0; i < neighbors.size(); ++i)
2001 if ((recv_pointers[i + 1] - recv_pointers[i]) > 0)
2002 {
2003 const int ierr =
2004 MPI_Irecv(recv_data.data() + recv_pointers[i],
2005 recv_pointers[i + 1] - recv_pointers[i],
2006 MPI_CHAR,
2007 neighbors[i],
2008 mpi_tag,
2009 parallel_triangulation->get_communicator(),
2010 &(requests[send_ops]));
2011 AssertThrowMPI(ierr);
2012 ++send_ops;
2013 }
2014
2015 for (unsigned int i = 0; i < neighbors.size(); ++i)
2016 if ((send_pointers[i + 1] - send_pointers[i]) > 0)
2017 {
2018 const int ierr =
2019 MPI_Isend(send_data.data() + send_pointers[i],
2020 send_pointers[i + 1] - send_pointers[i],
2021 MPI_CHAR,
2022 neighbors[i],
2023 mpi_tag,
2024 parallel_triangulation->get_communicator(),
2025 &(requests[send_ops + recv_ops]));
2026 AssertThrowMPI(ierr);
2027 ++recv_ops;
2028 }
2029 const int ierr =
2030 MPI_Waitall(send_ops + recv_ops, requests.data(), MPI_STATUSES_IGNORE);
2031 AssertThrowMPI(ierr);
2032 }
2033
2034 // Put the received particles into the domain if they are in the
2035 // triangulation
2036 const void *recv_data_it = static_cast<const void *>(recv_data.data());
2037
2038 // Gather ghost particle iterators from the cache
2039 auto &ghost_particles_iterators =
2040 ghost_particles_cache.ghost_particles_iterators;
2041
2042 for (auto &recv_particle : ghost_particles_iterators)
2043 {
2044 // Update particle data using previously allocated memory space
2045 // for efficiency reasons
2046 recv_data_it =
2047 recv_particle->read_particle_data_from_memory(recv_data_it);
2048
2049 Assert(recv_particle->particles_in_cell->cell->is_ghost(),
2051
2052 if (load_callback)
2053 recv_data_it = load_callback(
2054 particle_iterator(recv_particle->particles_in_cell,
2055 *property_pool,
2056 recv_particle->particle_index_within_cell),
2057 recv_data_it);
2058 }
2059
2060 AssertThrow(recv_data_it == recv_data.data() + recv_data.size(),
2061 ExcMessage(
2062 "The amount of data that was read into new particles "
2063 "does not match the amount of data sent around."));
2064 }
2065#endif
2066
2067 template <int dim, int spacedim>
2068 void
2070 const std::function<std::size_t()> &size_callb,
2071 const std::function<void *(const particle_iterator &, void *)> &store_callb,
2072 const std::function<const void *(const particle_iterator &, const void *)>
2073 &load_callb)
2074 {
2075 size_callback = size_callb;
2076 store_callback = store_callb;
2077 load_callback = load_callb;
2078 }
2079
2080
2081 template <int dim, int spacedim>
2082 void
2084 {
2085 // First disconnect existing connections
2086 for (const auto &connection : tria_listeners)
2087 connection.disconnect();
2088
2089 tria_listeners.clear();
2090
2091 tria_listeners.push_back(triangulation->signals.create.connect([&]() {
2092 this->initialize(*(this->triangulation),
2093 *(this->mapping),
2094 this->property_pool->n_properties_per_slot());
2095 }));
2096
2097 this->tria_listeners.push_back(
2098 this->triangulation->signals.clear.connect([&]() { this->clear(); }));
2099
2100 // for distributed triangulations, connect to distributed signals
2102 *>(&(*triangulation)) != nullptr)
2103 {
2104 tria_listeners.push_back(
2105 triangulation->signals.post_distributed_refinement.connect(
2106 [&]() { this->post_mesh_change_action(); }));
2107 tria_listeners.push_back(
2108 triangulation->signals.post_distributed_repartition.connect(
2109 [&]() { this->post_mesh_change_action(); }));
2110 tria_listeners.push_back(
2111 triangulation->signals.post_distributed_load.connect(
2112 [&]() { this->post_mesh_change_action(); }));
2113 }
2114 else
2115 {
2116 tria_listeners.push_back(triangulation->signals.post_refinement.connect(
2117 [&]() { this->post_mesh_change_action(); }));
2118 }
2119 }
2120
2121
2122
2123 template <int dim, int spacedim>
2124 void
2126 {
2127 Assert(triangulation != nullptr, ExcInternalError());
2128
2129 const bool distributed_triangulation =
2130 dynamic_cast<
2132 &(*triangulation)) != nullptr;
2133 (void)distributed_triangulation;
2134
2135 Assert(
2136 distributed_triangulation || number_of_locally_owned_particles == 0,
2137 ExcMessage(
2138 "Mesh refinement in a non-distributed triangulation is not supported "
2139 "by the ParticleHandler class. Either insert particles after mesh "
2140 "creation, or use a distributed triangulation."));
2141
2142 // Resize the container if it is possible without
2143 // transferring particles
2144 if (number_of_locally_owned_particles == 0)
2145 cells_to_particle_cache.resize(triangulation->n_active_cells(),
2146 particles.end());
2147 }
2148
2149
2150
2151 template <int dim, int spacedim>
2152 void
2157
2158
2159
2160 template <int dim, int spacedim>
2161 void
2163 {
2164 register_data_attach();
2165 }
2166
2167
2168
2169 template <int dim, int spacedim>
2170 void
2172 {
2173 register_data_attach();
2174 }
2175
2176
2177
2178 template <int dim, int spacedim>
2179 void
2181 {
2183 *distributed_triangulation =
2185 dynamic_cast<
2187 &(*triangulation)));
2188 (void)distributed_triangulation;
2189
2190 Assert(
2191 distributed_triangulation != nullptr,
2192 ExcMessage(
2193 "Mesh refinement in a non-distributed triangulation is not supported "
2194 "by the ParticleHandler class. Either insert particles after mesh "
2195 "creation and do not refine afterwards, or use a distributed triangulation."));
2196
2197 const auto callback_function =
2199 &cell_iterator,
2200 const CellStatus cell_status) {
2201 return this->pack_callback(cell_iterator, cell_status);
2202 };
2203
2204 handle = distributed_triangulation->register_data_attach(
2205 callback_function, /*returns_variable_size_data=*/true);
2206 }
2207
2208
2209
2210 template <int dim, int spacedim>
2211 void
2213 {
2214 const bool serialization = false;
2215 notify_ready_to_unpack(serialization);
2216 }
2217
2218
2219
2220 template <int dim, int spacedim>
2221 void
2223 {
2224 const bool serialization = true;
2225 notify_ready_to_unpack(serialization);
2226 }
2227
2228
2229 template <int dim, int spacedim>
2230 void
2232 const bool serialization)
2233 {
2234 notify_ready_to_unpack(serialization);
2235 }
2236
2237
2238
2239 template <int dim, int spacedim>
2240 void
2242 const bool serialization)
2243 {
2245 *distributed_triangulation =
2247 dynamic_cast<
2249 &(*triangulation)));
2250 (void)distributed_triangulation;
2251
2252 Assert(
2253 distributed_triangulation != nullptr,
2254 ExcMessage(
2255 "Mesh refinement in a non-distributed triangulation is not supported "
2256 "by the ParticleHandler class. Either insert particles after mesh "
2257 "creation and do not refine afterwards, or use a distributed triangulation."));
2258
2259 // First prepare container for insertion
2260 clear();
2261
2262 // If we are resuming from a checkpoint, we first have to register the
2263 // store function again, to set the triangulation to the same state as
2264 // before the serialization. Only afterwards we know how to deserialize the
2265 // data correctly.
2266 if (serialization)
2267 register_data_attach();
2268
2269 // Check if something was stored and load it
2270 if (handle != numbers::invalid_unsigned_int)
2271 {
2272 const auto callback_function =
2274 &cell_iterator,
2275 const CellStatus cell_status,
2276 const boost::iterator_range<std::vector<char>::const_iterator>
2277 &range_iterator) {
2278 this->unpack_callback(cell_iterator, cell_status, range_iterator);
2279 };
2280
2281 distributed_triangulation->notify_ready_to_unpack(handle,
2282 callback_function);
2283
2284 // Reset handle and update global numbers.
2286 update_cached_numbers();
2287 }
2288 }
2289
2290
2291
2292 template <int dim, int spacedim>
2293 std::vector<char>
2296 const CellStatus status) const
2297 {
2298 std::vector<particle_iterator> stored_particles_on_cell;
2299
2300 switch (status)
2301 {
2304 // If the cell persist or is refined store all particles of the
2305 // current cell.
2306 {
2307 const unsigned int n_particles = n_particles_in_cell(cell);
2308 stored_particles_on_cell.reserve(n_particles);
2309
2310 for (unsigned int i = 0; i < n_particles; ++i)
2311 stored_particles_on_cell.push_back(particle_iterator(
2312 cells_to_particle_cache[cell->active_cell_index()],
2313 *property_pool,
2314 i));
2315 }
2316 break;
2317
2319 // If this cell is the parent of children that will be coarsened,
2320 // collect the particles of all children.
2321 {
2322 for (const auto &child : cell->child_iterators())
2323 {
2324 const unsigned int n_particles = n_particles_in_cell(child);
2325
2326 stored_particles_on_cell.reserve(
2327 stored_particles_on_cell.size() + n_particles);
2328
2329 const typename particle_container::iterator &cache =
2330 cells_to_particle_cache[child->active_cell_index()];
2331 for (unsigned int i = 0; i < n_particles; ++i)
2332 stored_particles_on_cell.push_back(
2333 particle_iterator(cache, *property_pool, i));
2334 }
2335 }
2336 break;
2337
2338 default:
2340 break;
2341 }
2342
2343 return pack_particles(stored_particles_on_cell);
2344 }
2345
2346
2347
2348 template <int dim, int spacedim>
2349 void
2352 const CellStatus status,
2353 const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
2354 {
2355 if (data_range.begin() == data_range.end())
2356 return;
2357
2358 const auto cell_to_store_particles =
2359 (status != CellStatus::cell_will_be_refined) ? cell : cell->child(0);
2360
2361 // deserialize particles and insert into local storage
2362 if (data_range.begin() != data_range.end())
2363 {
2364 const void *data = static_cast<const void *>(&(*data_range.begin()));
2365 const void *end = static_cast<const void *>(
2366 &(*data_range.begin()) + (data_range.end() - data_range.begin()));
2367
2368 while (data < end)
2369 {
2370 const void *old_data = data;
2371 const auto x = insert_particle(data, cell_to_store_particles);
2372
2373 // Ensure that the particle read exactly as much data as
2374 // it promised it needs to store its data
2375 const void *new_data = data;
2376 (void)old_data;
2377 (void)new_data;
2378 (void)x;
2379 AssertDimension((const char *)new_data - (const char *)old_data,
2380 x->serialized_size_in_bytes());
2381 }
2382
2383 Assert(data == end,
2384 ExcMessage(
2385 "The particle data could not be deserialized successfully. "
2386 "Check that when deserializing the particles you expect "
2387 "the same number of properties that were serialized."));
2388 }
2389
2390 auto loaded_particles_on_cell = particles_in_cell(cell_to_store_particles);
2391
2392 // now update particle storage location and properties if necessary
2393 switch (status)
2394 {
2396 {
2397 // all particles are correctly inserted
2398 }
2399 break;
2400
2402 {
2403 // all particles are in correct cell, but their reference location
2404 // has changed
2405 for (auto &particle : loaded_particles_on_cell)
2406 {
2407 const Point<dim> p_unit =
2408 mapping->transform_real_to_unit_cell(cell_to_store_particles,
2409 particle.get_location());
2410 particle.set_reference_location(p_unit);
2411 }
2412 }
2413 break;
2414
2416 {
2417 // we need to find the correct child to store the particles and
2418 // their reference location has changed
2419 typename particle_container::iterator &cache =
2420 cells_to_particle_cache[cell_to_store_particles
2421 ->active_cell_index()];
2422
2423 // make sure that the call above has inserted an entry
2424 Assert(cache != particles.end(), ExcInternalError());
2425
2426 // Cannot use range-based loop, because number of particles in cell
2427 // is going to change
2428 auto particle = loaded_particles_on_cell.begin();
2429 for (unsigned int i = 0; i < cache->particles.size();)
2430 {
2431 bool found_new_cell = false;
2432
2433 for (const auto &child : cell->child_iterators())
2434 {
2435 Assert(child->is_locally_owned(), ExcInternalError());
2436
2437 try
2438 {
2439 const Point<dim> p_unit =
2440 mapping->transform_real_to_unit_cell(
2441 child, particle->get_location());
2443 p_unit, tolerance_inside_cell))
2444 {
2445 found_new_cell = true;
2446 particle->set_reference_location(p_unit);
2447
2448 // if the particle is not in the cell we stored it
2449 // in above, its handle is in the wrong place
2450 if (child != cell_to_store_particles)
2451 {
2452 // move handle into correct cell
2453 insert_particle(cache->particles[i], child);
2454 // remove handle by replacing it with last one
2455 cache->particles[i] = cache->particles.back();
2456 cache->particles.pop_back();
2457 // no loop increment, we need to process
2458 // the new i-th particle.
2459 }
2460 else
2461 {
2462 // move on to next particle
2463 ++i;
2464 ++particle;
2465 }
2466 break;
2467 }
2468 }
2469 catch (typename Mapping<dim>::ExcTransformationFailed &)
2470 {}
2471 }
2472
2473 if (found_new_cell == false)
2474 {
2475 // If we get here, we did not find the particle in any
2476 // child. This case may happen for particles that are at the
2477 // boundary for strongly curved cells. We apply a tolerance
2478 // in the call to GeometryInfo<dim>::is_inside_unit_cell to
2479 // account for this, but if that is not enough, we still
2480 // need to prevent an endless loop here. Delete the particle
2481 // and move on.
2482 signals.particle_lost(particle,
2483 particle->get_surrounding_cell());
2484 if (cache->particles[i] !=
2486 property_pool->deregister_particle(cache->particles[i]);
2487 cache->particles[i] = cache->particles.back();
2488 cache->particles.pop_back();
2489 }
2490 }
2491 // clean up in case child 0 has no particle left
2492 if (cache->particles.empty())
2493 {
2494 particles.erase(cache);
2495 cache = particles.end();
2496 }
2497 }
2498 break;
2499
2500 default:
2502 break;
2503 }
2504 }
2505} // namespace Particles
2506
2507#include "particle_handler.inst"
2508
CellStatus
Definition cell_status.h:31
@ cell_will_be_refined
@ children_will_be_coarsened
std::size_t size() const
Definition array_view.h:684
std::array< unsigned int, 4 > binary_type
Definition cell_id.h:80
const unsigned int n_components
Definition function.h:163
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
Abstract base class for mapping classes.
Definition mapping.h:318
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
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)
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
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)
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)
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::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
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
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)
void copy_from(const ParticleHandler< dim, spacedim > &particle_handler)
SmartPointer< const Mapping< dim, spacedim >, ParticleHandler< dim, spacedim > > mapping
types::particle_index get_next_free_particle_index() const
IndexSet locally_owned_particle_ids() const
void set_property_pool(PropertyPool< dim, spacedim > &property_pool)
Definition particle.h:610
const Point< dim > & get_reference_location() const
Definition particle.h:583
const Point< spacedim > & get_location() const
Definition particle.h:556
std::size_t serialized_size_in_bytes() const
Definition particle.cc:286
types::particle_index get_id() const
Definition particle.h:592
ArrayView< double > get_properties()
Definition particle.cc:330
Definition point.h:111
numbers::NumberTraits< Number >::real_type norm() const
IteratorState::IteratorStates state() const
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const ::CellStatus)> &pack_callback, const bool returns_variable_size_data)
void notify_ready_to_unpack(const unsigned int handle, const std::function< void(const cell_iterator &, const ::CellStatus, const boost::iterator_range< std::vector< char >::const_iterator > &)> &unpack_callback)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcPointNotAvailableHere()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1550
#define DEAL_II_ASSERT_UNREACHABLE()
unsigned int find_closest_vertex(const std::map< unsigned int, Point< spacedim > > &vertices, const Point< spacedim > &p)
return_type compute_point_locations_try_all(const Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &points, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator())
return_type distributed_compute_point_locations(const GridTools::Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &local_points, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bboxes, const double tolerance=1e-10, const std::vector< bool > &marked_vertices={}, const bool enforce_unique_mapping=true)
@ past_the_end
Iterator reached end of container.
@ valid
Iterator points to a valid object.
@ particle_handler_send_recv_particles_send
ParticleHandler<dim, spacedim>::send_recv_particles.
Definition mpi_tags.h:113
@ particle_handler_send_recv_particles_setup
ParticleHandler<dim, spacedim>::send_recv_particles.
Definition mpi_tags.h:109
T sum(const T &t, const MPI_Comm mpi_communicator)
std::map< unsigned int, T > some_to_some(const MPI_Comm comm, const std::map< unsigned int, T > &objects_to_send)
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
Definition mpi.cc:92
T max(const T &t, const MPI_Comm mpi_communicator)
std::vector< T > all_gather(const MPI_Comm comm, const T &object_to_send)
const types::subdomain_id artificial_subdomain_id
Definition types.h:362
static const unsigned int invalid_unsigned_int
Definition types.h:220
bool is_finite(const double x)
Definition numbers.h:538
STL namespace.
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned int subdomain_id
Definition types.h:43
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors
*braid_SplitCommworld & comm
#define DEAL_II_PARTICLE_INDEX_MPI_TYPE
static bool is_inside_unit_cell(const Point< dim > &p)