Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
particle_handler.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 
20 
22 
23 #include <utility>
24 
26 
27 namespace Particles
28 {
29  namespace
30  {
31  template <int dim, int spacedim>
32  std::vector<char>
33  pack_particles(std::vector<Particle<dim, spacedim>> &particles)
34  {
35  std::vector<char> buffer;
36 
37  if (particles.size() == 0)
38  return buffer;
39 
40  buffer.resize(particles.size() *
41  particles.front().serialized_size_in_bytes());
42  void *current_data = buffer.data();
43 
44  for (const auto &particle : particles)
45  {
46  particle.write_data(current_data);
47  }
48 
49  return buffer;
50  }
51 
52  template <int dim, int spacedim>
53  std::vector<Particle<dim, spacedim>>
54  unpack_particles(
55  const boost::iterator_range<std::vector<char>::const_iterator>
56  & data_range,
57  PropertyPool &property_pool)
58  {
59  std::vector<Particle<dim, spacedim>> particles;
60 
61  if (data_range.empty())
62  return particles;
63 
64  Particle<dim, spacedim> particle;
65  particle.set_property_pool(property_pool);
66  const unsigned int particle_size = particle.serialized_size_in_bytes();
67 
68  particles.reserve(data_range.size() / particle_size);
69 
70  const void *data = static_cast<const void *>(&(*data_range.begin()));
71 
72  while (data < &(*data_range.end()))
73  {
74  particles.emplace_back(data, &property_pool);
75  }
76 
77  Assert(
78  data == &(*data_range.end()),
79  ExcMessage(
80  "The particle data could not be deserialized successfully. "
81  "Check that when deserializing the particles you expect the same "
82  "number of properties that were serialized."));
83 
84  return particles;
85  }
86  } // namespace
87 
88  template <int dim, int spacedim>
90  : triangulation()
91  , particles()
92  , ghost_particles()
93  , global_number_of_particles(0)
94  , global_max_particles_per_cell(0)
95  , next_free_particle_index(0)
96  , property_pool(new PropertyPool(0))
97  , size_callback()
98  , store_callback()
99  , load_callback()
100  , handle(numbers::invalid_unsigned_int)
101  {}
102 
103 
104 
105  template <int dim, int spacedim>
108  const Mapping<dim, spacedim> & mapping,
109  const unsigned int n_properties)
110  : triangulation(&triangulation, typeid(*this).name())
111  , mapping(&mapping, typeid(*this).name())
112  , particles()
113  , ghost_particles()
114  , global_number_of_particles(0)
115  , global_max_particles_per_cell(0)
116  , next_free_particle_index(0)
117  , property_pool(new PropertyPool(n_properties))
118  , size_callback()
119  , store_callback()
120  , load_callback()
121  , handle(numbers::invalid_unsigned_int)
122  {
124  std_cxx14::make_unique<GridTools::Cache<dim, spacedim>>(triangulation,
125  mapping);
126  }
127 
128 
129 
130  template <int dim, int spacedim>
131  void
133  const Triangulation<dim, spacedim> &new_triangulation,
134  const Mapping<dim, spacedim> & new_mapping,
135  const unsigned int n_properties)
136  {
137  triangulation = &new_triangulation;
138  mapping = &new_mapping;
139 
140  // Create the memory pool that will store all particle properties
141  property_pool = std_cxx14::make_unique<PropertyPool>(n_properties);
142 
143  // Create the grid cache to cache the information about the triangulation
144  // that is used to locate the particles into subdomains and cells
145  triangulation_cache =
146  std_cxx14::make_unique<GridTools::Cache<dim, spacedim>>(new_triangulation,
147  new_mapping);
148  }
149 
150 
151 
152  template <int dim, int spacedim>
153  void
155  {
156  clear_particles();
157  global_number_of_particles = 0;
158  next_free_particle_index = 0;
159  global_max_particles_per_cell = 0;
160  }
161 
162 
163 
164  template <int dim, int spacedim>
165  void
167  {
168  particles.clear();
169  }
170 
171 
172 
173  template <int dim, int spacedim>
174  void
176  {
177  types::particle_index locally_highest_index = 0;
178  unsigned int local_max_particles_per_cell = 0;
179  unsigned int current_particles_per_cell = 0;
181  triangulation->begin_active();
182 
183  for (particle_iterator particle = begin(); particle != end(); ++particle)
184  {
185  locally_highest_index =
186  std::max(locally_highest_index, particle->get_id());
187 
188  if (particle->get_surrounding_cell(*triangulation) != current_cell)
189  {
190  current_particles_per_cell = 0;
191  current_cell = particle->get_surrounding_cell(*triangulation);
192  }
193 
194  ++current_particles_per_cell;
195  local_max_particles_per_cell =
196  std::max(local_max_particles_per_cell, current_particles_per_cell);
197  }
198 
199  if (const auto parallel_triangulation =
200  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
201  &*triangulation))
202  {
203  global_number_of_particles = ::Utilities::MPI::sum(
204  particles.size(), parallel_triangulation->get_communicator());
205  next_free_particle_index =
206  global_number_of_particles == 0 ?
207  0 :
209  locally_highest_index,
210  parallel_triangulation->get_communicator()) +
211  1;
212  global_max_particles_per_cell = ::Utilities::MPI::max(
213  local_max_particles_per_cell,
214  parallel_triangulation->get_communicator());
215  }
216  else
217  {
218  global_number_of_particles = particles.size();
219  next_free_particle_index =
220  global_number_of_particles == 0 ? 0 : locally_highest_index + 1;
221  global_max_particles_per_cell = local_max_particles_per_cell;
222  }
223  }
224 
225 
226 
227  template <int dim, int spacedim>
230  {
231  return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin();
232  }
233 
234 
235 
236  template <int dim, int spacedim>
239  {
240  return particle_iterator(particles, particles.begin());
241  }
242 
243 
244 
245  template <int dim, int spacedim>
248  {
249  return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end();
250  }
251 
252 
253 
254  template <int dim, int spacedim>
257  {
258  return particle_iterator(particles, particles.end());
259  }
260 
261 
262 
263  template <int dim, int spacedim>
266  {
267  return (const_cast<ParticleHandler<dim, spacedim> *>(this))->begin_ghost();
268  }
269 
270 
271 
272  template <int dim, int spacedim>
275  {
276  return particle_iterator(ghost_particles, ghost_particles.begin());
277  }
278 
279 
280 
281  template <int dim, int spacedim>
284  {
285  return (const_cast<ParticleHandler<dim, spacedim> *>(this))->end_ghost();
286  }
287 
288 
289 
290  template <int dim, int spacedim>
293  {
294  return particle_iterator(ghost_particles, ghost_particles.end());
295  }
296 
297 
298 
299  template <int dim, int spacedim>
303  const
304  {
305  return (const_cast<ParticleHandler<dim, spacedim> *>(this))
306  ->particles_in_cell(cell);
307  }
308 
309 
310 
311  template <int dim, int spacedim>
315  {
316  const internal::LevelInd level_index =
317  std::make_pair(cell->level(), cell->index());
318 
319  if (cell->is_ghost())
320  {
321  const auto particles_in_cell = ghost_particles.equal_range(level_index);
323  particle_iterator(ghost_particles, particles_in_cell.first),
324  particle_iterator(ghost_particles, particles_in_cell.second));
325  }
326 
327  const auto particles_in_cell = particles.equal_range(level_index);
329  particle_iterator(particles, particles_in_cell.first),
330  particle_iterator(particles, particles_in_cell.second));
331  }
332 
333 
334 
335  template <int dim, int spacedim>
336  void
339  {
340  particles.erase(particle->particle);
341  }
342 
343 
344 
345  template <int dim, int spacedim>
348  const Particle<dim, spacedim> & particle,
350  {
351  typename std::multimap<internal::LevelInd,
352  Particle<dim, spacedim>>::iterator it =
353  particles.insert(
354  std::make_pair(internal::LevelInd(cell->level(), cell->index()),
355  particle));
356 
357  particle_iterator particle_it(particles, it);
358  particle_it->set_property_pool(*property_pool);
359 
360  if (particle.has_properties())
361  for (unsigned int n = 0; n < particle.get_properties().size(); ++n)
362  particle_it->get_properties()[n] = particle.get_properties()[n];
363 
364  return particle_it;
365  }
366 
367 
368 
369  template <int dim, int spacedim>
370  void
372  const std::multimap<
374  Particle<dim, spacedim>> &new_particles)
375  {
376  for (auto particle = new_particles.begin(); particle != new_particles.end();
377  ++particle)
378  particles.insert(
379  particles.end(),
380  std::make_pair(internal::LevelInd(particle->first->level(),
381  particle->first->index()),
382  particle->second));
383 
384  update_cached_numbers();
385  }
386 
387 
388 
389  template <int dim, int spacedim>
390  void
392  const std::vector<Point<spacedim>> &positions)
393  {
394  update_cached_numbers();
395 
396  // Determine the starting particle index of this process, which
397  // is the highest currently existing particle index plus the sum
398  // of the number of newly generated particles of all
399  // processes with a lower rank if in a parallel computation.
400  const types::particle_index local_next_particle_index =
401  get_next_free_particle_index();
402  types::particle_index local_start_index = 0;
403 
404 #ifdef DEAL_II_WITH_MPI
405  if (const auto parallel_triangulation =
406  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
407  &*triangulation))
408  {
409  types::particle_index particles_to_add_locally = positions.size();
410  const int ierr = MPI_Scan(&particles_to_add_locally,
411  &local_start_index,
412  1,
414  MPI_SUM,
415  parallel_triangulation->get_communicator());
416  AssertThrowMPI(ierr);
417  local_start_index -= particles_to_add_locally;
418  }
419 #endif
420 
421  local_start_index += local_next_particle_index;
422 
424  auto point_locations = GridTools::compute_point_locations(cache, positions);
425 
426  auto &cells = std::get<0>(point_locations);
427  auto &local_positions = std::get<1>(point_locations);
428  auto &index_map = std::get<2>(point_locations);
429 
430  if (cells.size() == 0)
431  return;
432 
433  auto hint =
434  particles.find(std::make_pair(cells[0]->level(), cells[0]->index()));
435  for (unsigned int i = 0; i < cells.size(); ++i)
436  {
437  internal::LevelInd current_cell(cells[i]->level(), cells[i]->index());
438  for (unsigned int p = 0; p < local_positions[i].size(); ++p)
439  {
440  hint = particles.insert(
441  hint,
442  std::make_pair(current_cell,
443  Particle<dim, spacedim>(positions[index_map[i][p]],
444  local_positions[i][p],
445  local_start_index +
446  index_map[i][p])));
447  }
448  }
449 
450  update_cached_numbers();
451  }
452 
453 
454 
455  template <int dim, int spacedim>
456  std::map<unsigned int, IndexSet>
458  const std::vector<Point<spacedim>> &positions,
459  const std::vector<std::vector<BoundingBox<spacedim>>>
460  & global_bounding_boxes,
461  const std::vector<std::vector<double>> &properties)
462  {
463  if (!properties.empty())
464  {
465  AssertDimension(properties.size(), positions.size());
466 #ifdef DEBUG
467  for (const auto &p : properties)
468  AssertDimension(p.size(), n_properties_per_particle());
469 #endif
470  }
471 
472  const auto tria =
474  &(*triangulation));
475  const auto comm =
476  (tria != nullptr ? tria->get_communicator() : MPI_COMM_WORLD);
477 
479 
481 
482  // Compute the global number of properties
483  const auto n_global_properties =
484  Utilities::MPI::sum(properties.size(), comm);
485 
486  // Gather the number of points per processor
487  const auto n_particles_per_proc =
488  Utilities::MPI::all_gather(comm, positions.size());
489 
490  // Calculate all starting points locally
491  std::vector<unsigned int> particle_start_indices(n_mpi_processes);
492 
493  unsigned int particle_start_index = 0;
494  for (unsigned int process = 0; process < particle_start_indices.size();
495  ++process)
496  {
497  particle_start_indices[process] = particle_start_index;
498  particle_start_index += n_particles_per_proc[process];
499  }
500 
501  // Get all local information
502  const auto cells_positions_and_index_maps =
504  positions,
505  global_bounding_boxes);
506 
507  // Unpack the information into several vectors:
508  // All cells that contain at least one particle
509  const auto &local_cells_containing_particles =
510  std::get<0>(cells_positions_and_index_maps);
511 
512  // The reference position of every particle in the local part of the
513  // triangulation.
514  const auto &local_reference_positions =
515  std::get<1>(cells_positions_and_index_maps);
516  // The original index in the positions vector for each particle in the
517  // local part of the triangulation
518  const auto &original_indices_of_local_particles =
519  std::get<2>(cells_positions_and_index_maps);
520  // The real spatial position of every particle in the local part of the
521  // triangulation.
522  const auto &local_positions = std::get<3>(cells_positions_and_index_maps);
523  // The MPI process that inserted each particle
524  const auto &calling_process_indices =
525  std::get<4>(cells_positions_and_index_maps);
526 
527  // Create the map of cpu to indices, indicating who sent us what
528  // particle.
529  std::map<unsigned int, IndexSet> original_process_to_local_particle_indices;
530  for (unsigned int i_cell = 0;
531  i_cell < local_cells_containing_particles.size();
532  ++i_cell)
533  {
534  for (unsigned int i_particle = 0;
535  i_particle < local_positions[i_cell].size();
536  ++i_particle)
537  {
538  const unsigned int local_id_on_calling_process =
539  original_indices_of_local_particles[i_cell][i_particle];
540  const unsigned int calling_process =
541  calling_process_indices[i_cell][i_particle];
542 
543  if (original_process_to_local_particle_indices.find(
544  calling_process) ==
545  original_process_to_local_particle_indices.end())
546  original_process_to_local_particle_indices.insert(
547  {calling_process,
548  IndexSet(n_particles_per_proc[calling_process])});
549 
550  original_process_to_local_particle_indices[calling_process]
551  .add_index(local_id_on_calling_process);
552  }
553  }
554  for (auto &process_and_particle_indices :
555  original_process_to_local_particle_indices)
556  process_and_particle_indices.second.compress();
557 
558 
559  // A map from mpi process to properties, ordered as in the IndexSet.
560  // Notice that this ordering maybe different from the ordering in the
561  // vectors above, since no local ordering is guaranteed by the
562  // distribute_compute_point_locations() call.
563  // This is only filled if n_global_properties is > 0
564  std::map<unsigned int, std::vector<std::vector<double>>>
565  locally_owned_properties_from_other_processes;
566 
567  if (n_global_properties > 0)
568  {
569  // Gather whom I sent my own particles to, to decide whom to send
570  // the particle properties
571  auto send_to_cpu = Utilities::MPI::some_to_some(
572  comm, original_process_to_local_particle_indices);
573 
574  std::map<unsigned int, std::vector<std::vector<double>>>
575  non_locally_owned_properties;
576 
577  // Prepare the vector of properties to send
578  for (const auto &it : send_to_cpu)
579  {
580  std::vector<std::vector<double>> properties_to_send(
581  it.second.n_elements(),
582  std::vector<double>(n_properties_per_particle()));
583  unsigned int index = 0;
584  for (const auto &el : it.second)
585  properties_to_send[index++] = properties[el];
586  non_locally_owned_properties.insert({it.first, properties_to_send});
587  }
588 
589  // Send the non locally owned properties to each mpi process
590  // that needs them
591  locally_owned_properties_from_other_processes =
592  Utilities::MPI::some_to_some(comm, non_locally_owned_properties);
593 
594  AssertDimension(locally_owned_properties_from_other_processes.size(),
595  original_process_to_local_particle_indices.size());
596  }
597 
598 
599  // Create the multimap of local particles
600  std::multimap<typename Triangulation<dim, spacedim>::active_cell_iterator,
602  particles;
603 
604  // Now fill up the actual particles
605  for (unsigned int i_cell = 0;
606  i_cell < local_cells_containing_particles.size();
607  ++i_cell)
608  {
609  for (unsigned int i_particle = 0;
610  i_particle < local_positions[i_cell].size();
611  ++i_particle)
612  {
613  const unsigned int local_id_on_calling_process =
614  original_indices_of_local_particles[i_cell][i_particle];
615  const unsigned int calling_process =
616  calling_process_indices[i_cell][i_particle];
617 
618  const unsigned int particle_id =
619  local_id_on_calling_process +
620  particle_start_indices[calling_process];
621 
622  Particle<dim, spacedim> particle(
623  local_positions[i_cell][i_particle],
624  local_reference_positions[i_cell][i_particle],
625  particle_id);
626 
627  if (n_global_properties > 0)
628  {
629  const unsigned int index_within_set =
630  original_process_to_local_particle_indices[calling_process]
631  .index_within_set(local_id_on_calling_process);
632 
633  const auto &this_particle_properties =
634  locally_owned_properties_from_other_processes
635  [calling_process][index_within_set];
636 
637  particle.set_property_pool(get_property_pool());
638  particle.set_properties(this_particle_properties);
639  }
640 
641  particles.emplace(local_cells_containing_particles[i_cell],
642  particle);
643  }
644  }
645 
646  this->insert_particles(particles);
647 
648  return original_process_to_local_particle_indices;
649  }
650 
651 
652 
653  template <int dim, int spacedim>
656  {
657  return global_number_of_particles;
658  }
659 
660 
661 
662  template <int dim, int spacedim>
665  {
666  return global_max_particles_per_cell;
667  }
668 
669 
670 
671  template <int dim, int spacedim>
674  {
675  return particles.size();
676  }
677 
678 
679 
680  template <int dim, int spacedim>
681  unsigned int
683  {
684  return property_pool->n_properties_per_slot();
685  }
686 
687 
688 
689  template <int dim, int spacedim>
690  unsigned int
693  const
694  {
695  const internal::LevelInd found_cell =
696  std::make_pair(cell->level(), cell->index());
697 
698  if (cell->is_locally_owned())
699  return particles.count(found_cell);
700  else if (cell->is_ghost())
701  return ghost_particles.count(found_cell);
702  else if (cell->is_artificial())
703  AssertThrow(false, ExcInternalError());
704 
705  return 0;
706  }
707 
708 
709 
710  template <int dim, int spacedim>
713  {
714  return next_free_particle_index;
715  }
716 
717 
718 
719  template <int dim, int spacedim>
720  IndexSet
722  {
723  IndexSet set(get_next_free_particle_index());
724  for (const auto &p : *this)
725  set.add_index(p.get_id());
726  set.compress();
727  return set;
728  }
729 
730 
731 
732  template <int dim, int spacedim>
733  void
735  std::vector<Point<spacedim>> &positions,
736  const bool add_to_output_vector)
737  {
738  // There should be one point per particle to gather
739  AssertDimension(positions.size(), n_locally_owned_particles());
740 
741  unsigned int i = 0;
742  for (auto it = begin(); it != end(); ++it, ++i)
743  {
744  if (add_to_output_vector)
745  positions[i] = positions[i] + it->get_location();
746  else
747  positions[i] = it->get_location();
748  }
749  }
750 
751 
752 
753  template <int dim, int spacedim>
754  void
756  const std::vector<Point<spacedim>> &new_positions,
757  const bool displace_particles)
758  {
759  // There should be one point per particle to fix the new position
760  AssertDimension(new_positions.size(), n_locally_owned_particles());
761 
762  unsigned int i = 0;
763  for (auto it = begin(); it != end(); ++it, ++i)
764  {
765  if (displace_particles)
766  it->set_location(it->get_location() + new_positions[i]);
767  else
768  it->set_location(new_positions[i]);
769  }
770  sort_particles_into_subdomains_and_cells();
771  }
772 
773 
774 
775  template <int dim, int spacedim>
776  void
778  const Function<spacedim> &function,
779  const bool displace_particles)
780  {
781  // The function should have sufficient components to displace the
782  // particles
783  AssertDimension(function.n_components, spacedim);
784 
785  Vector<double> new_position(spacedim);
786  for (auto &particle : *this)
787  {
788  Point<spacedim> particle_location = particle.get_location();
789  function.vector_value(particle_location, new_position);
790  if (displace_particles)
791  for (unsigned int d = 0; d < spacedim; ++d)
792  particle_location[d] += new_position[d];
793  else
794  for (unsigned int d = 0; d < spacedim; ++d)
795  particle_location[d] = new_position[d];
796  particle.set_location(particle_location);
797  }
798  sort_particles_into_subdomains_and_cells();
799  }
800 
801 
802 
803  template <int dim, int spacedim>
804  PropertyPool &
806  {
807  return *property_pool;
808  }
809 
810 
811 
812  namespace
813  {
823  template <int dim>
824  bool
825  compare_particle_association(
826  const unsigned int a,
827  const unsigned int b,
828  const Tensor<1, dim> & particle_direction,
829  const std::vector<Tensor<1, dim>> &center_directions)
830  {
831  const double scalar_product_a = center_directions[a] * particle_direction;
832  const double scalar_product_b = center_directions[b] * particle_direction;
833 
834  // The function is supposed to return if a is before b. We are looking
835  // for the alignment of particle direction and center direction,
836  // therefore return if the scalar product of a is larger.
837  return (scalar_product_a > scalar_product_b);
838  }
839  } // namespace
840 
841 
842 
843  template <int dim, int spacedim>
844  void
846  {
847  // TODO: The current algorithm only works for particles that are in
848  // the local domain or in ghost cells, because it only knows the
849  // subdomain_id of ghost cells, but not of artificial cells. This
850  // can be extended using the distributed version of compute point
851  // locations.
852  // TODO: Extend this function to allow keeping particles on other
853  // processes around (with an invalid cell).
854 
855  std::vector<particle_iterator> particles_out_of_cell;
856  particles_out_of_cell.reserve(n_locally_owned_particles());
857 
858  // Now update the reference locations of the moved particles
859  for (particle_iterator it = begin(); it != end(); ++it)
860  {
861  const typename Triangulation<dim, spacedim>::cell_iterator cell =
862  it->get_surrounding_cell(*triangulation);
863 
864  try
865  {
866  const Point<dim> p_unit =
867  mapping->transform_real_to_unit_cell(cell, it->get_location());
869  {
870  it->set_reference_location(p_unit);
871  }
872  else
873  {
874  // The particle has left the cell
875  particles_out_of_cell.push_back(it);
876  }
877  }
878  catch (typename Mapping<dim>::ExcTransformationFailed &)
879  {
880  // The particle has left the cell
881  particles_out_of_cell.push_back(it);
882  }
883  }
884 
885  // There are three reasons why a particle is not in its old cell:
886  // It moved to another cell, to another subdomain or it left the mesh.
887  // Particles that moved to another cell are updated and stored inside the
888  // sorted_particles vector, particles that moved to another domain are
889  // collected in the moved_particles_domain vector. Particles that left
890  // the mesh completely are ignored and removed.
891  std::vector<std::pair<internal::LevelInd, Particle<dim, spacedim>>>
892  sorted_particles;
893  std::map<types::subdomain_id, std::vector<particle_iterator>>
894  moved_particles;
895  std::map<
897  std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
898  moved_cells;
899 
900  // We do not know exactly how many particles are lost, exchanged between
901  // domains, or remain on this process. Therefore we pre-allocate
902  // approximate sizes for these vectors. If more space is needed an
903  // automatic and relatively fast (compared to other parts of this
904  // algorithm) re-allocation will happen.
905  using vector_size = typename std::vector<particle_iterator>::size_type;
906  sorted_particles.reserve(
907  static_cast<vector_size>(particles_out_of_cell.size() * 1.25));
908 
909  std::set<types::subdomain_id> ghost_owners;
910  if (const auto parallel_triangulation =
911  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
912  &*triangulation))
913  ghost_owners = parallel_triangulation->ghost_owners();
914 
915  for (const auto ghost_owner : ghost_owners)
916  moved_particles[ghost_owner].reserve(
917  static_cast<vector_size>(particles_out_of_cell.size() * 0.25));
918  for (const auto ghost_owner : ghost_owners)
919  moved_cells[ghost_owner].reserve(
920  static_cast<vector_size>(particles_out_of_cell.size() * 0.25));
921 
922  {
923  // Create a map from vertices to adjacent cells using grid cache
924  std::vector<
925  std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
926  vertex_to_cells = triangulation_cache->get_vertex_to_cell_map();
927 
928  // Create a corresponding map of vectors from vertex to cell center using
929  // grid cache
930  std::vector<std::vector<Tensor<1, spacedim>>> vertex_to_cell_centers =
931  triangulation_cache->get_vertex_to_cell_centers_directions();
932 
933  std::vector<unsigned int> neighbor_permutation;
934 
935  // Find the cells that the particles moved to.
936  typename std::vector<particle_iterator>::iterator
937  it = particles_out_of_cell.begin(),
938  end_particle = particles_out_of_cell.end();
939 
940  for (; it != end_particle; ++it)
941  {
942  // The cell the particle is in
943  Point<dim> current_reference_position;
944  bool found_cell = false;
945 
946  // Check if the particle is in one of the old cell's neighbors
947  // that are adjacent to the closest vertex
949  current_cell = (*it)->get_surrounding_cell(*triangulation);
950 
951  const unsigned int closest_vertex =
952  GridTools::find_closest_vertex_of_cell<dim, spacedim>(
953  current_cell, (*it)->get_location());
954  Tensor<1, spacedim> vertex_to_particle =
955  (*it)->get_location() - current_cell->vertex(closest_vertex);
956  vertex_to_particle /= vertex_to_particle.norm();
957 
958  const unsigned int closest_vertex_index =
959  current_cell->vertex_index(closest_vertex);
960  const unsigned int n_neighbor_cells =
961  vertex_to_cells[closest_vertex_index].size();
962 
963  neighbor_permutation.resize(n_neighbor_cells);
964  for (unsigned int i = 0; i < n_neighbor_cells; ++i)
965  neighbor_permutation[i] = i;
966 
967  const auto cell_centers =
968  vertex_to_cell_centers[closest_vertex_index];
969  std::sort(neighbor_permutation.begin(),
970  neighbor_permutation.end(),
971  [&vertex_to_particle, &cell_centers](const unsigned int a,
972  const unsigned int b) {
973  return compare_particle_association(a,
974  b,
975  vertex_to_particle,
976  cell_centers);
977  });
978 
979  // Search all of the cells adjacent to the closest vertex of the
980  // previous cell Most likely we will find the particle in them.
981  for (unsigned int i = 0; i < n_neighbor_cells; ++i)
982  {
983  try
984  {
985  typename std::set<typename Triangulation<dim, spacedim>::
986  active_cell_iterator>::const_iterator
987  cell = vertex_to_cells[closest_vertex_index].begin();
988 
989  std::advance(cell, neighbor_permutation[i]);
990  const Point<dim> p_unit =
991  mapping->transform_real_to_unit_cell(*cell,
992  (*it)->get_location());
994  {
995  current_cell = *cell;
996  current_reference_position = p_unit;
997  found_cell = true;
998  break;
999  }
1000  }
1001  catch (typename Mapping<dim>::ExcTransformationFailed &)
1002  {}
1003  }
1004 
1005  if (!found_cell)
1006  {
1007  // The particle is not in a neighbor of the old cell.
1008  // Look for the new cell in the whole local domain.
1009  // This case is rare.
1010  try
1011  {
1012  const std::pair<const typename Triangulation<dim, spacedim>::
1013  active_cell_iterator,
1014  Point<dim>>
1015  current_cell_and_position =
1016  GridTools::find_active_cell_around_point<>(
1017  *mapping, *triangulation, (*it)->get_location());
1018  current_cell = current_cell_and_position.first;
1019  current_reference_position = current_cell_and_position.second;
1020  }
1021  catch (GridTools::ExcPointNotFound<spacedim> &)
1022  {
1023  // We can find no cell for this particle. It has left the
1024  // domain due to an integration error or an open boundary.
1025  continue;
1026  }
1027  }
1028 
1029  // If we are here, we found a cell and reference position for this
1030  // particle
1031  (*it)->set_reference_location(current_reference_position);
1032 
1033  // Reinsert the particle into our domain if we own its cell.
1034  // Mark it for MPI transfer otherwise
1035  if (current_cell->is_locally_owned())
1036  {
1037  sorted_particles.push_back(
1038  std::make_pair(internal::LevelInd(current_cell->level(),
1039  current_cell->index()),
1040  (*it)->particle->second));
1041  }
1042  else
1043  {
1044  moved_particles[current_cell->subdomain_id()].push_back(*it);
1045  moved_cells[current_cell->subdomain_id()].push_back(current_cell);
1046  }
1047  }
1048  }
1049 
1050  // Sort the updated particles. This pre-sort speeds up inserting
1051  // them into particles to O(N) complexity.
1052  std::multimap<internal::LevelInd, Particle<dim, spacedim>>
1053  sorted_particles_map;
1054 
1055  // Exchange particles between processors if we have more than one process
1056 #ifdef DEAL_II_WITH_MPI
1057  if (const auto parallel_triangulation =
1058  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
1059  &*triangulation))
1060  {
1062  parallel_triangulation->get_communicator()) > 1)
1063  send_recv_particles(moved_particles,
1064  sorted_particles_map,
1065  moved_cells);
1066  }
1067 #endif
1068 
1069  sorted_particles_map.insert(sorted_particles.begin(),
1070  sorted_particles.end());
1071 
1072  for (unsigned int i = 0; i < particles_out_of_cell.size(); ++i)
1073  remove_particle(particles_out_of_cell[i]);
1074 
1075  particles.insert(sorted_particles_map.begin(), sorted_particles_map.end());
1076  update_cached_numbers();
1077  }
1078 
1079 
1080 
1081  template <int dim, int spacedim>
1082  void
1084  {
1085  // Nothing to do in serial computations
1086  const auto parallel_triangulation =
1087  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
1088  &*triangulation);
1089  if (parallel_triangulation != nullptr)
1090  {
1092  parallel_triangulation->get_communicator()) == 1)
1093  return;
1094  }
1095  else
1096  return;
1097 
1098 #ifdef DEAL_II_WITH_MPI
1099  // First clear the current ghost_particle information
1100  ghost_particles.clear();
1101 
1102  std::map<types::subdomain_id, std::vector<particle_iterator>>
1103  ghost_particles_by_domain;
1104 
1105  const std::set<types::subdomain_id> ghost_owners =
1106  parallel_triangulation->ghost_owners();
1107  for (const auto ghost_owner : ghost_owners)
1108  ghost_particles_by_domain[ghost_owner].reserve(
1109  static_cast<typename std::vector<particle_iterator>::size_type>(
1110  particles.size() * 0.25));
1111 
1112  std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain(
1113  triangulation->n_vertices());
1114 
1115  for (const auto &cell : triangulation->active_cell_iterators())
1116  {
1117  if (cell->is_ghost())
1118  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
1119  vertex_to_neighbor_subdomain[cell->vertex_index(v)].insert(
1120  cell->subdomain_id());
1121  }
1122 
1123  for (const auto &cell : triangulation->active_cell_iterators())
1124  {
1125  if (!cell->is_ghost())
1126  {
1127  std::set<unsigned int> cell_to_neighbor_subdomain;
1128  for (const unsigned int v : GeometryInfo<dim>::vertex_indices())
1129  {
1130  cell_to_neighbor_subdomain.insert(
1131  vertex_to_neighbor_subdomain[cell->vertex_index(v)].begin(),
1132  vertex_to_neighbor_subdomain[cell->vertex_index(v)].end());
1133  }
1134 
1135  if (cell_to_neighbor_subdomain.size() > 0)
1136  {
1137  const particle_iterator_range particle_range =
1138  particles_in_cell(cell);
1139 
1140  for (const auto domain : cell_to_neighbor_subdomain)
1141  {
1142  for (typename particle_iterator_range::iterator particle =
1143  particle_range.begin();
1144  particle != particle_range.end();
1145  ++particle)
1146  ghost_particles_by_domain[domain].push_back(particle);
1147  }
1148  }
1149  }
1150  }
1151 
1152  send_recv_particles(ghost_particles_by_domain, ghost_particles);
1153 #endif
1154  }
1155 
1156 
1157 
1158 #ifdef DEAL_II_WITH_MPI
1159  template <int dim, int spacedim>
1160  void
1162  const std::map<types::subdomain_id, std::vector<particle_iterator>>
1163  &particles_to_send,
1165  &received_particles,
1166  const std::map<
1169  &send_cells)
1170  {
1171  const auto parallel_triangulation =
1172  dynamic_cast<const parallel::Triangulation<dim, spacedim> *>(
1173  &*triangulation);
1174  Assert(
1175  parallel_triangulation,
1176  ExcMessage(
1177  "This function is only implemented for parallel::Triangulation objects."));
1178 
1179  // Determine the communication pattern
1180  const std::set<types::subdomain_id> ghost_owners =
1181  parallel_triangulation->ghost_owners();
1182  const std::vector<types::subdomain_id> neighbors(ghost_owners.begin(),
1183  ghost_owners.end());
1184  const unsigned int n_neighbors = neighbors.size();
1185 
1186  if (send_cells.size() != 0)
1187  Assert(particles_to_send.size() == send_cells.size(), ExcInternalError());
1188 
1189  // If we do not know the subdomain this particle needs to be send to,
1190  // throw an error
1191  Assert(particles_to_send.find(numbers::artificial_subdomain_id) ==
1192  particles_to_send.end(),
1193  ExcInternalError());
1194 
1195  // TODO: Implement the shipping of particles to processes that are not
1196  // ghost owners of the local domain
1197  for (auto send_particles = particles_to_send.begin();
1198  send_particles != particles_to_send.end();
1199  ++send_particles)
1200  Assert(ghost_owners.find(send_particles->first) != ghost_owners.end(),
1201  ExcNotImplemented());
1202 
1203  unsigned int n_send_particles = 0;
1204  for (auto send_particles = particles_to_send.begin();
1205  send_particles != particles_to_send.end();
1206  ++send_particles)
1207  n_send_particles += send_particles->second.size();
1208 
1209  const unsigned int cellid_size = sizeof(CellId::binary_type);
1210 
1211  // Containers for the amount and offsets of data we will send
1212  // to other processors and the data itself.
1213  std::vector<unsigned int> n_send_data(n_neighbors, 0);
1214  std::vector<unsigned int> send_offsets(n_neighbors, 0);
1215  std::vector<char> send_data;
1216 
1217  // Only serialize things if there are particles to be send.
1218  // We can not return early even if no particles
1219  // are send, because we might receive particles from other processes
1220  if (n_send_particles > 0)
1221  {
1222  // Allocate space for sending particle data
1223  const unsigned int particle_size =
1224  begin()->serialized_size_in_bytes() + cellid_size +
1225  (size_callback ? size_callback() : 0);
1226  send_data.resize(n_send_particles * particle_size);
1227  void *data = static_cast<void *>(&send_data.front());
1228 
1229  // Serialize the data sorted by receiving process
1230  for (unsigned int i = 0; i < n_neighbors; ++i)
1231  {
1232  send_offsets[i] = reinterpret_cast<std::size_t>(data) -
1233  reinterpret_cast<std::size_t>(&send_data.front());
1234 
1235  for (unsigned int j = 0;
1236  j < particles_to_send.at(neighbors[i]).size();
1237  ++j)
1238  {
1239  // If no target cells are given, use the iterator information
1241  cell;
1242  if (send_cells.size() == 0)
1243  cell =
1244  particles_to_send.at(neighbors[i])[j]->get_surrounding_cell(
1245  *triangulation);
1246  else
1247  cell = send_cells.at(neighbors[i])[j];
1248 
1249  const CellId::binary_type cellid =
1250  cell->id().template to_binary<dim>();
1251  memcpy(data, &cellid, cellid_size);
1252  data = static_cast<char *>(data) + cellid_size;
1253 
1254  particles_to_send.at(neighbors[i])[j]->write_data(data);
1255  if (store_callback)
1256  data =
1257  store_callback(particles_to_send.at(neighbors[i])[j], data);
1258  }
1259  n_send_data[i] = reinterpret_cast<std::size_t>(data) -
1260  send_offsets[i] -
1261  reinterpret_cast<std::size_t>(&send_data.front());
1262  }
1263  }
1264 
1265  // Containers for the data we will receive from other processors
1266  std::vector<unsigned int> n_recv_data(n_neighbors);
1267  std::vector<unsigned int> recv_offsets(n_neighbors);
1268 
1269  // Notify other processors how many particles we will send
1270  {
1271  const int mpi_tag = Utilities::MPI::internal::Tags::
1273 
1274  std::vector<MPI_Request> n_requests(2 * n_neighbors);
1275  for (unsigned int i = 0; i < n_neighbors; ++i)
1276  {
1277  const int ierr = MPI_Irecv(&(n_recv_data[i]),
1278  1,
1279  MPI_UNSIGNED,
1280  neighbors[i],
1281  mpi_tag,
1282  parallel_triangulation->get_communicator(),
1283  &(n_requests[2 * i]));
1284  AssertThrowMPI(ierr);
1285  }
1286  for (unsigned int i = 0; i < n_neighbors; ++i)
1287  {
1288  const int ierr = MPI_Isend(&(n_send_data[i]),
1289  1,
1290  MPI_UNSIGNED,
1291  neighbors[i],
1292  mpi_tag,
1293  parallel_triangulation->get_communicator(),
1294  &(n_requests[2 * i + 1]));
1295  AssertThrowMPI(ierr);
1296  }
1297  const int ierr =
1298  MPI_Waitall(2 * n_neighbors, n_requests.data(), MPI_STATUSES_IGNORE);
1299  AssertThrowMPI(ierr);
1300  }
1301 
1302  // Determine how many particles and data we will receive
1303  unsigned int total_recv_data = 0;
1304  for (unsigned int neighbor_id = 0; neighbor_id < n_neighbors; ++neighbor_id)
1305  {
1306  recv_offsets[neighbor_id] = total_recv_data;
1307  total_recv_data += n_recv_data[neighbor_id];
1308  }
1309 
1310  // Set up the space for the received particle data
1311  std::vector<char> recv_data(total_recv_data);
1312 
1313  // Exchange the particle data between domains
1314  {
1315  std::vector<MPI_Request> requests(2 * n_neighbors);
1316  unsigned int send_ops = 0;
1317  unsigned int recv_ops = 0;
1318 
1319  const int mpi_tag = Utilities::MPI::internal::Tags::
1321 
1322  for (unsigned int i = 0; i < n_neighbors; ++i)
1323  if (n_recv_data[i] > 0)
1324  {
1325  const int ierr =
1326  MPI_Irecv(&(recv_data[recv_offsets[i]]),
1327  n_recv_data[i],
1328  MPI_CHAR,
1329  neighbors[i],
1330  mpi_tag,
1331  parallel_triangulation->get_communicator(),
1332  &(requests[send_ops]));
1333  AssertThrowMPI(ierr);
1334  send_ops++;
1335  }
1336 
1337  for (unsigned int i = 0; i < n_neighbors; ++i)
1338  if (n_send_data[i] > 0)
1339  {
1340  const int ierr =
1341  MPI_Isend(&(send_data[send_offsets[i]]),
1342  n_send_data[i],
1343  MPI_CHAR,
1344  neighbors[i],
1345  mpi_tag,
1346  parallel_triangulation->get_communicator(),
1347  &(requests[send_ops + recv_ops]));
1348  AssertThrowMPI(ierr);
1349  recv_ops++;
1350  }
1351  const int ierr =
1352  MPI_Waitall(send_ops + recv_ops, requests.data(), MPI_STATUSES_IGNORE);
1353  AssertThrowMPI(ierr);
1354  }
1355 
1356  // Put the received particles into the domain if they are in the
1357  // triangulation
1358  const void *recv_data_it = static_cast<const void *>(recv_data.data());
1359 
1360  while (reinterpret_cast<std::size_t>(recv_data_it) -
1361  reinterpret_cast<std::size_t>(recv_data.data()) <
1362  total_recv_data)
1363  {
1364  CellId::binary_type binary_cellid;
1365  memcpy(&binary_cellid, recv_data_it, cellid_size);
1366  const CellId id(binary_cellid);
1367  recv_data_it = static_cast<const char *>(recv_data_it) + cellid_size;
1368 
1370  id.to_cell(*triangulation);
1371 
1372  typename std::multimap<internal::LevelInd,
1373  Particle<dim, spacedim>>::iterator
1374  recv_particle = received_particles.insert(std::make_pair(
1375  internal::LevelInd(cell->level(), cell->index()),
1376  Particle<dim, spacedim>(recv_data_it, property_pool.get())));
1377 
1378  if (load_callback)
1379  recv_data_it =
1380  load_callback(particle_iterator(received_particles, recv_particle),
1381  recv_data_it);
1382  }
1383 
1384  AssertThrow(recv_data_it == recv_data.data() + recv_data.size(),
1385  ExcMessage(
1386  "The amount of data that was read into new particles "
1387  "does not match the amount of data sent around."));
1388  }
1389 #endif
1390 
1391 
1392 
1393  template <int dim, int spacedim>
1394  void
1396  const std::function<std::size_t()> & size_callb,
1397  const std::function<void *(const particle_iterator &, void *)> &store_callb,
1398  const std::function<const void *(const particle_iterator &, const void *)>
1399  &load_callb)
1400  {
1401  size_callback = size_callb;
1402  store_callback = store_callb;
1403  load_callback = load_callb;
1404  }
1405 
1406 
1407 
1408  template <int dim, int spacedim>
1409  void
1411  {
1413  *non_const_triangulation =
1416  *>(&(*triangulation)));
1417  (void)non_const_triangulation;
1418 
1419  Assert(non_const_triangulation != nullptr, ::ExcNotImplemented());
1420 
1421 #ifdef DEAL_II_WITH_P4EST
1422  // Only save and load particles if there are any, we might get here for
1423  // example if somebody created a ParticleHandler but generated 0
1424  // particles.
1425  update_cached_numbers();
1426 
1427  if (global_max_particles_per_cell > 0)
1428  {
1429  const auto callback_function =
1430  [this](const typename Triangulation<dim, spacedim>::cell_iterator
1431  &cell_iterator,
1433  cell_status) {
1434  return this->store_particles(cell_iterator, cell_status);
1435  };
1436 
1437  handle = non_const_triangulation->register_data_attach(
1438  callback_function, /*returns_variable_size_data=*/true);
1439  }
1440 #endif
1441  }
1442 
1443 
1444 
1445  template <int dim, int spacedim>
1446  void
1448  const bool serialization)
1449  {
1450  // All particles have been stored, when we reach this point. Empty the
1451  // particle data.
1452  clear_particles();
1453 
1455  *non_const_triangulation =
1458  *>(&(*triangulation)));
1459  (void)non_const_triangulation;
1460 
1461  Assert(non_const_triangulation != nullptr, ::ExcNotImplemented());
1462 
1463 #ifdef DEAL_II_WITH_P4EST
1464  // If we are resuming from a checkpoint, we first have to register the
1465  // store function again, to set the triangulation in the same state as
1466  // before the serialization. Only by this it knows how to deserialize the
1467  // data correctly. Only do this if something was actually stored.
1468  if (serialization && (global_max_particles_per_cell > 0))
1469  {
1470  const auto callback_function =
1471  [this](const typename Triangulation<dim, spacedim>::cell_iterator
1472  &cell_iterator,
1474  cell_status) {
1475  return this->store_particles(cell_iterator, cell_status);
1476  };
1477 
1478  handle = non_const_triangulation->register_data_attach(
1479  callback_function, /*returns_variable_size_data=*/true);
1480  }
1481 
1482  // Check if something was stored and load it
1483  if (handle != numbers::invalid_unsigned_int)
1484  {
1485  const auto callback_function =
1486  [this](
1488  &cell_iterator,
1489  const typename Triangulation<dim, spacedim>::CellStatus cell_status,
1490  const boost::iterator_range<std::vector<char>::const_iterator>
1491  &range_iterator) {
1492  this->load_particles(cell_iterator, cell_status, range_iterator);
1493  };
1494 
1495  non_const_triangulation->notify_ready_to_unpack(handle,
1496  callback_function);
1497 
1498  // Reset handle and update global number of particles. The number
1499  // can change because of discarded or newly generated particles
1501  update_cached_numbers();
1502  }
1503 #else
1504  (void)serialization;
1505 #endif
1506  }
1507 
1508 
1509 
1510  template <int dim, int spacedim>
1511  std::vector<char>
1513  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
1514  const typename Triangulation<dim, spacedim>::CellStatus status) const
1515  {
1516  std::vector<Particle<dim, spacedim>> stored_particles_on_cell;
1517 
1518  switch (status)
1519  {
1522  // If the cell persist or is refined store all particles of the
1523  // current cell.
1524  {
1525  unsigned int n_particles = 0;
1526 
1527  const internal::LevelInd level_index = {cell->level(),
1528  cell->index()};
1529  const auto particles_in_cell =
1530  (cell->is_ghost() ? ghost_particles.equal_range(level_index) :
1531  particles.equal_range(level_index));
1532 
1533  n_particles = n_particles_in_cell(cell);
1534  stored_particles_on_cell.reserve(n_particles);
1535 
1536  std::for_each(
1537  particles_in_cell.first,
1538  particles_in_cell.second,
1539  [&stored_particles_on_cell](
1541  &particle) {
1542  stored_particles_on_cell.push_back(particle.second);
1543  });
1544 
1545  AssertDimension(n_particles, stored_particles_on_cell.size());
1546  }
1547  break;
1548 
1550  // If this cell is the parent of children that will be coarsened,
1551  // collect the particles of all children.
1552  {
1553  unsigned int n_particles = 0;
1554 
1555  for (unsigned int child_index = 0;
1556  child_index < GeometryInfo<dim>::max_children_per_cell;
1557  ++child_index)
1558  {
1560  child = cell->child(child_index);
1561  n_particles += n_particles_in_cell(child);
1562  }
1563 
1564  stored_particles_on_cell.reserve(n_particles);
1565 
1566  for (unsigned int child_index = 0;
1567  child_index < GeometryInfo<dim>::max_children_per_cell;
1568  ++child_index)
1569  {
1571  child = cell->child(child_index);
1572  const internal::LevelInd level_index = {child->level(),
1573  child->index()};
1574  const auto particles_in_cell =
1575  (child->is_ghost() ?
1576  ghost_particles.equal_range(level_index) :
1577  particles.equal_range(level_index));
1578 
1579  std::for_each(
1580  particles_in_cell.first,
1581  particles_in_cell.second,
1582  [&stored_particles_on_cell](
1584  &particle) {
1585  stored_particles_on_cell.push_back(particle.second);
1586  });
1587  }
1588 
1589  AssertDimension(n_particles, stored_particles_on_cell.size());
1590  }
1591  break;
1592 
1593  default:
1594  Assert(false, ExcInternalError());
1595  break;
1596  }
1597 
1598  return pack_particles(stored_particles_on_cell);
1599  }
1600 
1601  template <int dim, int spacedim>
1602  void
1604  const typename Triangulation<dim, spacedim>::cell_iterator & cell,
1605  const typename Triangulation<dim, spacedim>::CellStatus status,
1606  const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
1607  {
1608  // We leave this container non-const to be able to `std::move`
1609  // its contents directly into the particles multimap later.
1610  std::vector<Particle<dim, spacedim>> loaded_particles_on_cell =
1611  unpack_particles<dim, spacedim>(data_range, *property_pool);
1612 
1613  // Update the reference to the current property pool for all particles.
1614  // This was not stored, because they might be transported across process
1615  // domains.
1616  for (auto &particle : loaded_particles_on_cell)
1617  particle.set_property_pool(*property_pool);
1618 
1619  switch (status)
1620  {
1622  {
1623  auto position_hint = particles.end();
1624  for (const auto &particle : loaded_particles_on_cell)
1625  {
1626  // Use std::multimap::emplace_hint to speed up insertion of
1627  // particles. This is a C++11 function, but not all compilers
1628  // that report a -std=c++11 (like gcc 4.6) implement it, so
1629  // require C++14 instead.
1630 #ifdef DEAL_II_WITH_CXX14
1631  position_hint =
1632  particles.emplace_hint(position_hint,
1633  std::make_pair(cell->level(),
1634  cell->index()),
1635  std::move(particle));
1636 #else
1637  position_hint =
1638  particles.insert(position_hint,
1639  std::make_pair(std::make_pair(cell->level(),
1640  cell->index()),
1641  std::move(particle)));
1642 #endif
1643  // Move the hint position forward by one, i.e., for the next
1644  // particle. The 'hint' position will thus be right after the
1645  // one just inserted.
1646  ++position_hint;
1647  }
1648  }
1649  break;
1650 
1652  {
1653  typename std::multimap<internal::LevelInd,
1654  Particle<dim, spacedim>>::iterator
1655  position_hint = particles.end();
1656  for (auto &particle : loaded_particles_on_cell)
1657  {
1658  const Point<dim> p_unit =
1659  mapping->transform_real_to_unit_cell(cell,
1660  particle.get_location());
1661  particle.set_reference_location(p_unit);
1662  // Use std::multimap::emplace_hint to speed up insertion of
1663  // particles. This is a C++11 function, but not all compilers
1664  // that report a -std=c++11 (like gcc 4.6) implement it, so
1665  // require C++14 instead.
1666 #ifdef DEAL_II_WITH_CXX14
1667  position_hint =
1668  particles.emplace_hint(position_hint,
1669  std::make_pair(cell->level(),
1670  cell->index()),
1671  std::move(particle));
1672 #else
1673  position_hint =
1674  particles.insert(position_hint,
1675  std::make_pair(std::make_pair(cell->level(),
1676  cell->index()),
1677  std::move(particle)));
1678 #endif
1679  // Move the hint position forward by one, i.e., for the next
1680  // particle. The 'hint' position will thus be right after the
1681  // one just inserted.
1682  ++position_hint;
1683  }
1684  }
1685  break;
1686 
1688  {
1689  std::vector<
1690  typename std::multimap<internal::LevelInd,
1691  Particle<dim, spacedim>>::iterator>
1693  for (unsigned int child_index = 0;
1694  child_index < GeometryInfo<dim>::max_children_per_cell;
1695  ++child_index)
1696  {
1698  child = cell->child(child_index);
1699  position_hints[child_index] = particles.upper_bound(
1700  std::make_pair(child->level(), child->index()));
1701  }
1702 
1703  for (auto &particle : loaded_particles_on_cell)
1704  {
1705  for (unsigned int child_index = 0;
1706  child_index < GeometryInfo<dim>::max_children_per_cell;
1707  ++child_index)
1708  {
1710  child = cell->child(child_index);
1711 
1712  try
1713  {
1714  const Point<dim> p_unit =
1715  mapping->transform_real_to_unit_cell(
1716  child, particle.get_location());
1718  {
1719  particle.set_reference_location(p_unit);
1720  // Use std::multimap::emplace_hint to speed up
1721  // insertion of particles. This is a C++11
1722  // function, but not all compilers that report a
1723  // -std=c++11 (like gcc 4.6) implement it, so
1724  // require C++14 instead.
1725 #ifdef DEAL_II_WITH_CXX14
1726  position_hints[child_index] =
1727  particles.emplace_hint(
1728  position_hints[child_index],
1729  std::make_pair(child->level(), child->index()),
1730  std::move(particle));
1731 #else
1732  position_hints[child_index] = particles.insert(
1733  position_hints[child_index],
1734  std::make_pair(std::make_pair(child->level(),
1735  child->index()),
1736  std::move(particle)));
1737 #endif
1738  // Move the hint position forward by one, i.e.,
1739  // for the next particle. The 'hint' position will
1740  // thus be right after the one just inserted.
1741  ++position_hints[child_index];
1742  break;
1743  }
1744  }
1745  catch (typename Mapping<dim>::ExcTransformationFailed &)
1746  {}
1747  }
1748  }
1749  }
1750  break;
1751 
1752  default:
1753  Assert(false, ExcInternalError());
1754  break;
1755  }
1756  }
1757 } // namespace Particles
1758 
1759 #include "particle_handler.inst"
1760 
IndexSet
Definition: index_set.h:74
Particles::ParticleHandler::particles_in_cell
particle_iterator_range particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
Definition: particle_handler.cc:313
Particles::ParticleHandler::ParticleHandler
ParticleHandler()
Definition: particle_handler.cc:89
Particles::ParticleHandler::register_store_callback_function
void register_store_callback_function()
Definition: particle_handler.cc:1410
Particles::ParticleHandler::set_particle_positions
std::enable_if< std::is_convertible< VectorType *, Function< spacedim > * >::value==false >::type set_particle_positions(const VectorType &input_vector, const bool displace_particles=true)
Definition: particle_handler.h:807
Particles::ParticleHandler::n_properties_per_particle
unsigned int n_properties_per_particle() const
Definition: particle_handler.cc:682
DEAL_II_PARTICLE_INDEX_MPI_TYPE
#define DEAL_II_PARTICLE_INDEX_MPI_TYPE
Definition: particle.h:76
Particles::ParticleHandler::insert_particles
void insert_particles(const std::multimap< typename Triangulation< dim, spacedim >::active_cell_iterator, Particle< dim, spacedim >> &particles)
Definition: particle_handler.cc:371
Particles::ParticleHandler::begin_ghost
particle_iterator begin_ghost() const
Definition: particle_handler.cc:265
Utilities::MPI::internal::Tags::particle_handler_send_recv_particles_setup
@ particle_handler_send_recv_particles_setup
ParticleHandler<dim, spacedim>::send_recv_particles.
Definition: mpi_tags.h:110
particle_handler.h
Particles::ParticleHandler::get_next_free_particle_index
types::particle_index get_next_free_particle_index() const
Definition: particle_handler.cc:712
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
Triangulation
Definition: tria.h:1109
Particles::ParticleHandler::particle_iterator
ParticleIterator< dim, spacedim > particle_iterator
Definition: particle_handler.h:66
GeometryInfo
Definition: geometry_info.h:1224
Particles::PropertyPool
Definition: property_pool.h:43
Particles
Definition: data_out.h:27
Particles::ParticleHandler::particle_iterator_range
boost::iterator_range< particle_iterator > particle_iterator_range
Definition: particle_handler.h:71
Particles::ParticleHandler::send_recv_particles
void send_recv_particles(const std::map< types::subdomain_id, std::vector< particle_iterator >> &particles_to_send, std::multimap< internal::LevelInd, Particle< dim, spacedim >> &received_particles, const std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator >> &new_cells_for_particles=std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator >>())
Definition: particle_handler.cc:1161
parallel::distributed::Triangulation::notify_ready_to_unpack
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)
Definition: tria.cc:4378
Particles::Particle::set_properties
void set_properties(const ArrayView< const double > &new_properties)
Definition: particle.cc:295
Particles::ParticleIterator
Definition: particle_iterator.h:39
Particles::ParticleHandler::load_particles
void load_particles(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
Definition: particle_handler.cc:1603
DoFTools::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
Particles::ParticleHandler::triangulation_cache
std::unique_ptr< GridTools::Cache< dim, spacedim > > triangulation_cache
Definition: particle_handler.h:717
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Particles::ParticleHandler::insert_global_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={})
Definition: particle_handler.cc:457
Particles::ParticleHandler::register_additional_store_load_functions
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)
Definition: particle_handler.cc:1395
LinearAlgebra::CUDAWrappers::kernel::set
__global__ void set(Number *val, const Number s, const size_type N)
BoundingBox
Definition: bounding_box.h:128
level
unsigned int level
Definition: grid_out.cc:4355
Triangulation::CellStatus
CellStatus
Definition: tria.h:1969
Particles::Particle::set_property_pool
void set_property_pool(PropertyPool &property_pool)
Definition: particle.cc:286
Particles::ParticleHandler::get_property_pool
PropertyPool & get_property_pool() const
Definition: particle_handler.cc:805
Particles::Particle
Definition: particle.h:146
Particles::ParticleHandler::end_ghost
particle_iterator end_ghost() const
Definition: particle_handler.cc:283
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
Tensor::norm
numbers::NumberTraits< Number >::real_type norm() const
Particles::ParticleHandler::n_particles_in_cell
unsigned int n_particles_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
Definition: particle_handler.cc:691
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
Utilities::MPI::internal::Tags::particle_handler_send_recv_particles_send
@ particle_handler_send_recv_particles_send
ParticleHandler<dim, spacedim>::send_recv_particles.
Definition: mpi_tags.h:112
Particles::ParticleHandler::mapping
SmartPointer< const Mapping< dim, spacedim >, ParticleHandler< dim, spacedim > > mapping
Definition: particle_handler.h:619
CellId
Definition: cell_id.h:69
grid_tools.h
Tensor< 1, dim >
GridTools::compute_point_locations
return_type compute_point_locations(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())
Definition: grid_tools.cc:4193
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
parallel::distributed::Triangulation
Definition: tria.h:242
Particles::ParticleHandler::get_particle_positions
void get_particle_positions(VectorType &output_vector, const bool add_to_output_vector=false)
Definition: particle_handler.h:830
Particles::ParticleHandler::n_global_particles
types::particle_index n_global_particles() const
Definition: particle_handler.cc:655
Particles::ParticleHandler::sort_particles_into_subdomains_and_cells
void sort_particles_into_subdomains_and_cells()
Definition: particle_handler.cc:845
TrilinosWrappers::internal::end
VectorType::value_type * end(VectorType &V)
Definition: trilinos_sparse_matrix.cc:65
Utilities::MPI::sum
T sum(const T &t, const MPI_Comm &mpi_communicator)
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
AssertThrowMPI
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1707
GridTools::Cache
Definition: grid_tools.h:128
Particles::ParticleHandler::end
particle_iterator end() const
Definition: particle_handler.cc:247
numbers
Definition: numbers.h:207
Particles::ParticleHandler::n_locally_owned_particles
types::particle_index n_locally_owned_particles() const
Definition: particle_handler.cc:673
TriaActiveIterator
Definition: tria_iterator.h:759
types::subdomain_id
unsigned int subdomain_id
Definition: types.h:43
Particles::ParticleHandler::store_particles
std::vector< char > store_particles(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status) const
Definition: particle_handler.cc:1512
advance
void advance(std::tuple< I1, I2 > &t, const unsigned int n)
Definition: synchronous_iterator.h:146
Particles::ParticleHandler::remove_particle
void remove_particle(const particle_iterator &particle)
Definition: particle_handler.cc:337
unsigned int
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
LinearAlgebra::CUDAWrappers::kernel::size_type
types::global_dof_index size_type
Definition: cuda_kernels.h:45
Particles::ParticleHandler::locally_relevant_ids
IndexSet locally_relevant_ids() const
Definition: particle_handler.cc:721
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
Particles::ParticleHandler::insert_particle
particle_iterator insert_particle(const Particle< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
Definition: particle_handler.cc:347
Particles::internal::LevelInd
std::pair< int, int > LevelInd
Definition: particle.h:92
Utilities::MPI::n_mpi_processes
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:117
grid_tools_cache.h
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
Particles::ParticleHandler
Definition: data_out.h:30
Particles::ParticleHandler::triangulation
SmartPointer< const Triangulation< dim, spacedim >, ParticleHandler< dim, spacedim > > triangulation
Definition: particle_handler.h:613
Particles::ParticleHandler::update_cached_numbers
void update_cached_numbers()
Definition: particle_handler.cc:175
Particles::ParticleHandler::clear_particles
void clear_particles()
Definition: particle_handler.cc:166
Point< spacedim >
CellId::binary_type
std::array< unsigned int, 4 > binary_type
Definition: cell_id.h:79
parallel::TriangulationBase
Definition: tria_base.h:78
Particles::ParticleHandler::n_global_max_particles_per_cell
types::particle_index n_global_max_particles_per_cell() const
Definition: particle_handler.cc:664
Particles::Particle::has_properties
bool has_properties() const
Definition: particle.cc:360
Utilities::MPI::all_gather
std::vector< T > all_gather(const MPI_Comm &comm, const T &object_to_send)
numbers::artificial_subdomain_id
const types::subdomain_id artificial_subdomain_id
Definition: types.h:302
Function
Definition: function.h:151
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
memory.h
Particles::ParticleHandler::register_load_callback_function
void register_load_callback_function(const bool serialization)
Definition: particle_handler.cc:1447
parallel::distributed::Triangulation::register_data_attach
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const CellStatus)> &pack_callback, const bool returns_variable_size_data)
Definition: tria.cc:4348
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
GridTools::distributed_compute_point_locations
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)
Definition: grid_tools.cc:4940
Utilities::MPI::some_to_some
std::map< unsigned int, T > some_to_some(const MPI_Comm &comm, const std::map< unsigned int, T > &objects_to_send)
Vector< double >
make_iterator_range
IteratorRange< BaseIterator > make_iterator_range(const BaseIterator &begin, const typename identity< BaseIterator >::type &end)
Definition: iterator_range.h:296
Particles::ParticleHandler::begin
particle_iterator begin() const
Definition: particle_handler.cc:229
Particles::ParticleHandler::exchange_ghost_particles
void exchange_ghost_particles()
Definition: particle_handler.cc:1083
TriaIterator
Definition: tria_iterator.h:578
Particles::Particle::get_properties
const ArrayView< double > get_properties()
Definition: particle.cc:327
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
Triangulation::clear
virtual void clear()
Definition: tria.cc:10209
Particles::ParticleHandler::initialize
void initialize(const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0)
Definition: particle_handler.cc:132
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
Particles::ParticleHandler::clear
void clear()
Definition: particle_handler.cc:154