33 template <
int dim,
int spacedim>
35 pack_particles(std::vector<ParticleIterator<dim, spacedim>> &particles)
37 std::vector<char> buffer;
39 if (particles.size() == 0)
42 buffer.resize(particles.size() *
43 particles.front()->serialized_size_in_bytes());
44 void *current_data = buffer.data();
46 for (
const auto &particle : particles)
48 current_data = particle->write_particle_data_to_memory(current_data);
57 template <
int dim,
int spacedim>
62 , global_number_of_particles(0)
63 , number_of_locally_owned_particles(0)
64 , global_max_particles_per_cell(0)
65 , next_free_particle_index(0)
69 , handle(
numbers::invalid_unsigned_int)
77 template <
int dim,
int spacedim>
81 const unsigned int n_properties)
83 , mapping(&mapping, typeid(*this).name())
84 , property_pool(
std::make_unique<
PropertyPool<dim, spacedim>>(n_properties))
85 , cells_to_particle_cache(
triangulation.n_active_cells(), particles.end())
86 , global_number_of_particles(0)
87 , number_of_locally_owned_particles(0)
88 , global_max_particles_per_cell(0)
89 , next_free_particle_index(0)
93 , handle(
numbers::invalid_unsigned_int)
94 , triangulation_cache(
105 template <
int dim,
int spacedim>
110 for (
const auto &connection : tria_listeners)
111 connection.disconnect();
116 template <
int dim,
int spacedim>
121 const unsigned int n_properties)
126 mapping = &new_mapping;
128 reset_particle_container(particles);
131 property_pool = std::make_unique<PropertyPool<dim, spacedim>>(n_properties);
135 triangulation_cache =
136 std::make_unique<GridTools::Cache<dim, spacedim>>(new_triangulation,
139 cells_to_particle_cache.resize(
triangulation->n_active_cells(),
142 connect_to_triangulation_signals();
147 template <
int dim,
int spacedim>
152 const unsigned int n_properties =
158 property_pool = std::make_unique<PropertyPool<dim, spacedim>>(
163 number_of_locally_owned_particles =
166 global_max_particles_per_cell =
172 particles.insert(particle_container_owned_end(),
175 particles.insert(particle_container_ghost_end(),
179 for (
auto it = particles.begin(); it != particles.end(); ++it)
180 if (!it->particles.empty())
181 cells_to_particle_cache[it->cell->active_cell_index()] = it;
183 ghost_particles_cache.ghost_particles_by_domain =
185 handle = particle_handler.
handle;
190 template <
int dim,
int spacedim>
195 global_number_of_particles = 0;
196 number_of_locally_owned_particles = 0;
197 next_free_particle_index = 0;
198 global_max_particles_per_cell = 0;
203 template <
int dim,
int spacedim>
207 for (
auto &particles_in_cell : particles)
208 for (
auto &particle : particles_in_cell.particles)
210 property_pool->deregister_particle(particle);
212 cells_to_particle_cache.clear();
213 reset_particle_container(particles);
215 cells_to_particle_cache.resize(
triangulation->n_active_cells(),
220 property_pool->clear();
225 template <
int dim,
int spacedim>
229 property_pool->reserve(n_particles);
234 template <
int dim,
int spacedim>
242 past_the_end_iterator =
247 given_particles.clear();
248 for (
unsigned int i = 0; i < 3; ++i)
249 given_particles.emplace_back(
251 past_the_end_iterator);
254 const_cast<typename particle_container::iterator &
>(owned_particles_end) =
255 ++given_particles.begin();
260 template <
int dim,
int spacedim>
265 bool sort_is_necessary =
false;
267 auto previous = particle_container_owned_begin();
268 for (
auto next = previous; next != particle_container_owned_end(); ++next)
272 previous->cell > next->cell)
274 sort_is_necessary =
true;
280 if (sort_is_necessary)
295 reset_particle_container(sorted_particles);
298 for (
const auto &cell :
triangulation->active_cell_iterators())
299 if (!cell->is_artificial())
300 if (cells_to_particle_cache[cell->active_cell_index()] !=
307 typename particle_container::iterator insert_position =
308 cell->is_locally_owned() ? particle_container_owned_end() :
309 --sorted_particles.end();
310 typename particle_container::iterator new_entry =
311 sorted_particles.insert(
312 insert_position,
typename particle_container::value_type());
313 new_entry->cell = cell;
314 new_entry->particles =
315 std::move(cells_to_particle_cache[cell->active_cell_index()]
318 particles = std::move(sorted_particles);
321 cells_to_particle_cache.clear();
322 cells_to_particle_cache.resize(
triangulation->n_active_cells(),
324 for (
auto it = particles.begin(); it != particles.end(); ++it)
325 if (!it->particles.empty())
326 cells_to_particle_cache[it->cell->active_cell_index()] = it;
332 particles.front().particles.empty() &&
334 particles.back().particles.empty() &&
336 owned_particles_end->particles.empty(),
342 for (
const auto &it : cells_to_particle_cache)
343 Assert(it != particles.begin() && it != owned_particles_end &&
344 it != --(particles.end()),
349 for (
auto it = particle_container_owned_begin();
350 it != particle_container_owned_end();
355 std::vector<typename particle_container::iterator> verify_cache(
357 for (
auto it = particles.begin(); it != particles.end(); ++it)
358 if (!it->particles.empty())
359 verify_cache[it->cell->active_cell_index()] = it;
361 for (
unsigned int i = 0; i < verify_cache.size(); ++i)
367 number_of_locally_owned_particles = 0;
370 for (
const auto &particles_in_cell : particles)
373 particles_in_cell.particles.size();
376 result[0] =
std::max(result[0], n_particles_in_cell);
379 if (n_particles_in_cell > 0 &&
380 particles_in_cell.cell->is_locally_owned())
381 number_of_locally_owned_particles += n_particles_in_cell;
384 for (
const auto &particle : particles_in_cell.particles)
385 result[1] =
std::max(result[1], property_pool->get_id(particle));
388 global_number_of_particles =
392 if (global_number_of_particles == 0)
394 next_free_particle_index = 0;
395 global_max_particles_per_cell = 0;
401 next_free_particle_index = result[1] + 1;
402 global_max_particles_per_cell = result[0];
408 template <
int dim,
int spacedim>
414 if (cells_to_particle_cache.size() == 0)
417 if (cell->is_artificial() ==
false)
419 return cells_to_particle_cache[cell->active_cell_index()] !=
421 cells_to_particle_cache[cell->active_cell_index()]
427 ExcMessage(
"You can't ask for the particles on an artificial "
428 "cell since we don't know what exists on these "
436 template <
int dim,
int spacedim>
443 ->particles_in_cell(cell);
448 template <
int dim,
int spacedim>
453 const unsigned int active_cell_index = cell->active_cell_index();
455 if (cell->is_artificial() ==
false)
457 if (cells_to_particle_cache[active_cell_index] == particles.end())
459 return boost::make_iterator_range(
465 const typename particle_container::iterator
466 particles_in_current_cell =
467 cells_to_particle_cache[active_cell_index];
468 typename particle_container::iterator particles_in_next_cell =
469 particles_in_current_cell;
470 ++particles_in_next_cell;
471 return boost::make_iterator_range(
478 ExcMessage(
"You can't ask for the particles on an artificial "
479 "cell since we don't know what exists on these "
487 template <
int dim,
int spacedim>
492 auto &particles_on_cell = particle->particles_in_cell->particles;
497 auto handle = particle->get_handle();
499 property_pool->deregister_particle(handle);
504 const auto cell = particle->get_surrounding_cell();
505 const bool owned_cell = cell->is_locally_owned();
507 --number_of_locally_owned_particles;
509 if (particles_on_cell.size() > 1)
511 particles_on_cell[particle->particle_index_within_cell] =
512 std::move(particles_on_cell.back());
513 particles_on_cell.resize(particles_on_cell.size() - 1);
517 particles.erase(particle->particles_in_cell);
518 cells_to_particle_cache[cell->active_cell_index()] = particles.end();
524 template <
int dim,
int spacedim>
528 &particles_to_remove)
536 return a->particles_in_cell->cell > b->particles_in_cell->cell ||
537 (a->particles_in_cell->cell == b->particles_in_cell->cell &&
538 a->particle_index_within_cell > b->particle_index_within_cell);
541 bool particles_are_sorted =
true;
542 auto previous = particles_to_remove.begin();
543 for (
auto next = previous; next != particles_to_remove.end(); ++next)
545 if (check_greater(*previous, *next))
547 particles_are_sorted =
false;
552 if (particles_are_sorted)
555 for (
auto it = particles_to_remove.rbegin();
556 it != particles_to_remove.rend();
558 remove_particle(*it);
562 std::vector<ParticleHandler<dim, spacedim>::particle_iterator>
563 sorted_particles(particles_to_remove);
564 std::sort(sorted_particles.begin(),
565 sorted_particles.end(),
568 for (
const auto &particle : sorted_particles)
569 remove_particle(particle);
572 update_cached_numbers();
577 template <
int dim,
int spacedim>
592 template <
int dim,
int spacedim>
598 const unsigned int active_cell_index = cell->active_cell_index();
599 typename particle_container::iterator &cache =
600 cells_to_particle_cache[active_cell_index];
601 if (cache == particles.end())
603 const typename particle_container::iterator insert_position =
604 cell->is_locally_owned() ? particle_container_owned_end() :
605 particle_container_ghost_end();
606 cache = particles.emplace(
613 cache->particles.push_back(handle);
618 cache->particles.size() - 1);
623 template <
int dim,
int spacedim>
632 Assert(cell->is_locally_owned(),
633 ExcMessage(
"You tried to insert particles into a cell that is not "
634 "locally owned. This is not supported."));
637 insert_particle(property_pool->register_particle(), cell);
639 data = particle_it->read_particle_data_from_memory(data);
641 ++number_of_locally_owned_particles;
648 template <
int dim,
int spacedim>
661 Assert(cell->is_locally_owned(),
662 ExcMessage(
"You tried to insert particles into a cell that is not "
663 "locally owned. This is not supported."));
666 insert_particle(property_pool->register_particle(), cell);
668 particle_it->set_location(position);
669 particle_it->set_reference_location(reference_position);
670 particle_it->set_id(particle_index);
672 if (properties.
size() != 0)
673 particle_it->set_properties(properties);
675 ++number_of_locally_owned_particles;
682 template <
int dim,
int spacedim>
689 reserve(n_locally_owned_particles() + new_particles.size());
690 for (
const auto &cell_and_particle : new_particles)
691 insert_particle(cell_and_particle.second, cell_and_particle.first);
693 update_cached_numbers();
698 template <
int dim,
int spacedim>
705 update_cached_numbers();
706 reserve(n_locally_owned_particles() + positions.size());
713 get_next_free_particle_index();
716#ifdef DEAL_II_WITH_MPI
717 if (
const auto parallel_triangulation =
722 const int ierr = MPI_Scan(&particles_to_add_locally,
727 parallel_triangulation->get_communicator());
729 local_start_index -= particles_to_add_locally;
733 local_start_index += local_next_particle_index;
735 auto point_locations =
739 auto &cells = std::get<0>(point_locations);
740 auto &local_positions = std::get<1>(point_locations);
741 auto &index_map = std::get<2>(point_locations);
742 auto &missing_points = std::get<3>(point_locations);
748 (void)missing_points;
750 for (
unsigned int i = 0; i < cells.size(); ++i)
751 for (
unsigned int p = 0; p < local_positions[i].size(); ++p)
752 insert_particle(positions[index_map[i][p]],
753 local_positions[i][p],
754 local_start_index + index_map[i][p],
757 update_cached_numbers();
762 template <
int dim,
int spacedim>
763 std::map<unsigned int, IndexSet>
767 & global_bounding_boxes,
768 const std::vector<std::vector<double>> & properties,
769 const std::vector<types::particle_index> &ids)
771 if (!properties.empty())
775 for (
const auto &p : properties)
788 const auto n_global_properties =
792 const auto n_particles_per_proc =
796 std::vector<unsigned int> particle_start_indices(n_mpi_processes);
798 unsigned int particle_start_index = get_next_free_particle_index();
799 for (
unsigned int process = 0; process < particle_start_indices.size();
802 particle_start_indices[process] = particle_start_index;
803 particle_start_index += n_particles_per_proc[process];
807 const auto cells_positions_and_index_maps =
810 global_bounding_boxes);
814 const auto &local_cells_containing_particles =
815 std::get<0>(cells_positions_and_index_maps);
819 const auto &local_reference_positions =
820 std::get<1>(cells_positions_and_index_maps);
823 const auto &original_indices_of_local_particles =
824 std::get<2>(cells_positions_and_index_maps);
827 const auto &local_positions = std::get<3>(cells_positions_and_index_maps);
829 const auto &calling_process_indices =
830 std::get<4>(cells_positions_and_index_maps);
833 std::map<unsigned int, std::vector<unsigned int>>
834 original_process_to_local_particle_indices_tmp;
835 for (
unsigned int i_cell = 0;
836 i_cell < local_cells_containing_particles.size();
839 for (
unsigned int i_particle = 0;
840 i_particle < local_positions[i_cell].size();
843 const unsigned int local_id_on_calling_process =
844 original_indices_of_local_particles[i_cell][i_particle];
845 const unsigned int calling_process =
846 calling_process_indices[i_cell][i_particle];
848 original_process_to_local_particle_indices_tmp[calling_process]
849 .push_back(local_id_on_calling_process);
852 std::map<unsigned int, IndexSet> original_process_to_local_particle_indices;
853 for (
auto &process_and_particle_indices :
854 original_process_to_local_particle_indices_tmp)
856 const unsigned int calling_process = process_and_particle_indices.first;
857 original_process_to_local_particle_indices.insert(
858 {calling_process,
IndexSet(n_particles_per_proc[calling_process])});
859 std::sort(process_and_particle_indices.second.begin(),
860 process_and_particle_indices.second.end());
861 original_process_to_local_particle_indices[calling_process].add_indices(
862 process_and_particle_indices.second.begin(),
863 process_and_particle_indices.second.end());
864 original_process_to_local_particle_indices[calling_process].compress();
872 std::map<unsigned int, std::vector<std::vector<double>>>
873 locally_owned_properties_from_other_processes;
880 std::map<unsigned int, std::vector<types::particle_index>>
881 locally_owned_ids_from_other_processes;
883 if (n_global_properties > 0 || !ids.empty())
888 comm, original_process_to_local_particle_indices);
891 if (n_global_properties > 0)
893 std::map<unsigned int, std::vector<std::vector<double>>>
894 non_locally_owned_properties;
896 for (
const auto &it : send_to_cpu)
898 std::vector<std::vector<double>> properties_to_send(
899 it.second.n_elements(),
900 std::vector<double>(n_properties_per_particle()));
901 unsigned int index = 0;
902 for (
const auto el : it.second)
903 properties_to_send[index++] = properties[el];
904 non_locally_owned_properties.insert(
905 {it.first, properties_to_send});
910 locally_owned_properties_from_other_processes =
914 locally_owned_properties_from_other_processes.size(),
915 original_process_to_local_particle_indices.size());
920 std::map<unsigned int, std::vector<types::particle_index>>
921 non_locally_owned_ids;
922 for (
const auto &it : send_to_cpu)
924 std::vector<types::particle_index> ids_to_send(
925 it.second.n_elements());
926 unsigned int index = 0;
927 for (
const auto el : it.second)
928 ids_to_send[index++] = ids[el];
929 non_locally_owned_ids.insert({it.first, ids_to_send});
934 locally_owned_ids_from_other_processes =
938 original_process_to_local_particle_indices.size());
943 for (
unsigned int i_cell = 0;
944 i_cell < local_cells_containing_particles.size();
947 for (
unsigned int i_particle = 0;
948 i_particle < local_positions[i_cell].size();
951 const unsigned int local_id_on_calling_process =
952 original_indices_of_local_particles[i_cell][i_particle];
954 const unsigned int calling_process =
955 calling_process_indices[i_cell][i_particle];
957 const unsigned int index_within_set =
958 original_process_to_local_particle_indices[calling_process]
959 .index_within_set(local_id_on_calling_process);
961 const unsigned int particle_id =
963 local_id_on_calling_process +
964 particle_start_indices[calling_process] :
965 locally_owned_ids_from_other_processes[calling_process]
969 insert_particle(local_positions[i_cell][i_particle],
970 local_reference_positions[i_cell][i_particle],
972 local_cells_containing_particles[i_cell]);
974 if (n_global_properties > 0)
976 particle_it->set_properties(
977 locally_owned_properties_from_other_processes
978 [calling_process][index_within_set]);
983 update_cached_numbers();
985 return original_process_to_local_particle_indices;
990 template <
int dim,
int spacedim>
991 std::map<unsigned int, IndexSet>
995 &global_bounding_boxes)
999 std::vector<Point<spacedim>> positions;
1000 std::vector<std::vector<double>> properties;
1001 std::vector<types::particle_index> ids;
1002 positions.resize(particles.size());
1003 ids.resize(particles.size());
1004 if (n_properties_per_particle() > 0)
1005 properties.resize(particles.size(),
1006 std::vector<double>(n_properties_per_particle()));
1009 for (
const auto &p : particles)
1011 positions[i] = p.get_location();
1012 ids[i] = p.get_id();
1013 if (p.has_properties())
1014 properties[i] = {p.get_properties().begin(),
1015 p.get_properties().end()};
1019 return insert_global_particles(positions,
1020 global_bounding_boxes,
1027 template <
int dim,
int spacedim>
1031 return global_number_of_particles;
1036 template <
int dim,
int spacedim>
1040 return global_max_particles_per_cell;
1045 template <
int dim,
int spacedim>
1049 return number_of_locally_owned_particles;
1054 template <
int dim,
int spacedim>
1058 return property_pool->n_properties_per_slot();
1063 template <
int dim,
int spacedim>
1067 return next_free_particle_index;
1072 template <
int dim,
int spacedim>
1076 IndexSet set(get_next_free_particle_index());
1077 std::vector<types::particle_index> indices;
1078 indices.reserve(n_locally_owned_particles());
1079 for (
const auto &p : *
this)
1080 indices.push_back(p.get_id());
1081 set.add_indices(indices.begin(), indices.end());
1088 template <
int dim,
int spacedim>
1092 return property_pool->n_slots();
1097 template <
int dim,
int spacedim>
1101 const bool add_to_output_vector)
1107 for (
auto it = begin(); it != end(); ++it, ++i)
1109 if (add_to_output_vector)
1110 positions[i] = positions[i] + it->get_location();
1112 positions[i] = it->get_location();
1118 template <
int dim,
int spacedim>
1122 const bool displace_particles)
1128 for (
auto it = begin(); it != end(); ++it, ++i)
1130 if (displace_particles)
1131 it->set_location(it->get_location() + new_positions[i]);
1133 it->set_location(new_positions[i]);
1135 sort_particles_into_subdomains_and_cells();
1140 template <
int dim,
int spacedim>
1144 const bool displace_particles)
1151 for (
auto &particle : *
this)
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];
1159 for (
unsigned int d = 0; d < spacedim; ++d)
1160 particle_location[d] = new_position[d];
1161 particle.set_location(particle_location);
1163 sort_particles_into_subdomains_and_cells();
1168 template <
int dim,
int spacedim>
1172 return *property_pool;
1190 compare_particle_association(
1191 const unsigned int a,
1192 const unsigned int b,
1196 const double scalar_product_a = center_directions[a] * particle_direction;
1197 const double scalar_product_b = center_directions[b] * particle_direction;
1202 return (scalar_product_a > scalar_product_b);
1208 template <
int dim,
int spacedim>
1224 std::vector<particle_iterator> particles_out_of_cell;
1229 particles_out_of_cell.reserve(n_locally_owned_particles() / 4);
1232 std::vector<Point<spacedim>> real_locations;
1233 std::vector<Point<dim>> reference_locations;
1234 real_locations.reserve(global_max_particles_per_cell);
1235 reference_locations.reserve(global_max_particles_per_cell);
1237 for (
const auto &cell :
triangulation->active_cell_iterators())
1244 if (cell->is_locally_owned() ==
false)
1249 const unsigned int n_pic = n_particles_in_cell(cell);
1250 auto pic = particles_in_cell(cell);
1252 real_locations.clear();
1253 for (
const auto &particle : pic)
1254 real_locations.push_back(particle.get_location());
1256 reference_locations.resize(n_pic);
1257 mapping->transform_points_real_to_unit_cell(cell,
1259 reference_locations);
1261 auto particle = pic.begin();
1262 for (
const auto &p_unit : reference_locations)
1266 particle->set_reference_location(p_unit);
1268 particles_out_of_cell.push_back(particle);
1280 std::map<types::subdomain_id, std::vector<particle_iterator>>
1284 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1292 using vector_size =
typename std::vector<particle_iterator>::size_type;
1294 std::set<types::subdomain_id> ghost_owners;
1295 if (
const auto parallel_triangulation =
1298 ghost_owners = parallel_triangulation->ghost_owners();
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);
1311 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1312 &vertex_to_cells = triangulation_cache->get_vertex_to_cell_map();
1316 const std::vector<std::vector<Tensor<1, spacedim>>>
1317 &vertex_to_cell_centers =
1318 triangulation_cache->get_vertex_to_cell_centers_directions();
1320 std::vector<unsigned int> neighbor_permutation;
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);
1332 for (
auto &out_particle : particles_out_of_cell)
1337 auto current_cell = out_particle->get_surrounding_cell();
1339 real_locations[0] = out_particle->get_location();
1342 bool found_cell =
false;
1346 const unsigned int closest_vertex =
1347 GridTools::find_closest_vertex_of_cell<dim, spacedim>(
1348 current_cell, out_particle->get_location(), *mapping);
1350 out_particle->get_location() - current_cell->vertex(closest_vertex);
1351 vertex_to_particle /= vertex_to_particle.
norm();
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();
1358 neighbor_permutation.resize(n_neighbor_cells);
1359 for (
unsigned int i = 0; i < n_neighbor_cells; ++i)
1360 neighbor_permutation[i] = i;
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,
1376 for (
unsigned int i = 0; i < n_neighbor_cells; ++i)
1378 typename std::set<typename Triangulation<dim, spacedim>::
1379 active_cell_iterator>::const_iterator cell =
1380 vertex_to_cells[closest_vertex_index].begin();
1382 std::advance(cell, neighbor_permutation[i]);
1383 mapping->transform_points_real_to_unit_cell(*cell,
1385 reference_locations);
1388 reference_locations[0]))
1390 current_cell = *cell;
1403#if defined(DEAL_II_WITH_BOOST_BUNDLED) || \
1404 !(defined(__clang_major__) && __clang_major__ >= 16) || \
1405 BOOST_VERSION >= 108100
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(),
1414 std::back_inserter(closest_vertex_in_domain));
1418 const unsigned int closest_vertex_index_in_domain =
1419 closest_vertex_in_domain[0].second;
1421 const unsigned int closest_vertex_index_in_domain =
1424 out_particle->get_location());
1429 for (
const auto &cell :
1430 vertex_to_cells[closest_vertex_index_in_domain])
1432 mapping->transform_points_real_to_unit_cell(
1433 cell, real_locations, reference_locations);
1436 reference_locations[0]))
1438 current_cell = cell;
1450 signals.particle_lost(out_particle,
1451 out_particle->get_surrounding_cell());
1457 out_particle->set_reference_location(reference_locations[0]);
1461 if (current_cell->is_locally_owned())
1464 out_particle->particles_in_cell
1465 ->particles[out_particle->particle_index_within_cell];
1468 const auto old_value = old;
1472 insert_particle(old_value, current_cell);
1476 moved_particles[current_cell->subdomain_id()].push_back(
1478 moved_cells[current_cell->subdomain_id()].push_back(current_cell);
1484#ifdef DEAL_II_WITH_MPI
1485 if (
const auto parallel_triangulation =
1490 parallel_triangulation->get_communicator()) > 1)
1491 send_recv_particles(moved_particles, moved_cells);
1496 remove_particles(particles_out_of_cell);
1499 std::vector<typename PropertyPool<dim, spacedim>::Handle> unsorted_handles;
1500 unsorted_handles.reserve(property_pool->n_registered_slots());
1503 for (
auto &particles_in_cell : particles)
1504 for (
auto &particle : particles_in_cell.particles)
1506 unsorted_handles.push_back(particle);
1507 particle = sorted_handle++;
1510 property_pool->sort_memory_slots(unsorted_handles);
1516 template <
int dim,
int spacedim>
1519 const bool enable_cache)
1522 const auto parallel_triangulation =
1525 if (parallel_triangulation !=
nullptr)
1528 parallel_triangulation->get_communicator()) == 1)
1534#ifndef DEAL_II_WITH_MPI
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())
1542 Assert(cells_to_particle_cache[cell->active_cell_index()]->cell ==
1546 for (
auto &ghost_particle :
1547 cells_to_particle_cache[cell->active_cell_index()]->particles)
1548 property_pool->deregister_particle(ghost_particle);
1551 particles.erase(cells_to_particle_cache[cell->active_cell_index()]);
1552 cells_to_particle_cache[cell->active_cell_index()] = particles.end();
1556 ghost_particles_cache.ghost_particles_by_domain.clear();
1557 ghost_particles_cache.valid =
false;
1564 const std::map<unsigned int, std::set<types::subdomain_id>>
1566 triangulation_cache->get_vertices_with_ghost_neighbors();
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);
1574 const std::vector<std::set<unsigned int>> vertex_to_neighbor_subdomain =
1575 triangulation_cache->get_vertex_to_neighbor_subdomain();
1577 for (
const auto &cell :
triangulation->active_cell_iterators())
1579 if (cell->is_locally_owned())
1581 std::set<unsigned int> cell_to_neighbor_subdomain;
1582 for (
const unsigned int v : cell->vertex_indices())
1584 const auto vertex_ghost_neighbors =
1586 if (vertex_ghost_neighbors !=
1589 cell_to_neighbor_subdomain.insert(
1590 vertex_ghost_neighbors->second.begin(),
1591 vertex_ghost_neighbors->second.end());
1595 if (cell_to_neighbor_subdomain.size() > 0)
1598 particles_in_cell(cell);
1600 for (
const auto domain : cell_to_neighbor_subdomain)
1602 for (
typename particle_iterator_range::iterator particle =
1603 particle_range.begin();
1604 particle != particle_range.end();
1606 ghost_particles_cache.ghost_particles_by_domain[domain]
1607 .push_back(particle);
1613 send_recv_particles(
1614 ghost_particles_cache.ghost_particles_by_domain,
1625 template <
int dim,
int spacedim>
1630 const auto parallel_triangulation =
1633 if (parallel_triangulation ==
nullptr ||
1635 parallel_triangulation->get_communicator()) == 1)
1641#ifdef DEAL_II_WITH_MPI
1644 Assert(ghost_particles_cache.valid,
1646 "Ghost particles cannot be updated if they first have not been "
1647 "exchanged at least once with the cache enabled"));
1650 send_recv_particles_properties_and_location(
1651 ghost_particles_cache.ghost_particles_by_domain);
1657#ifdef DEAL_II_WITH_MPI
1658 template <
int dim,
int spacedim>
1667 const bool build_cache)
1673 ghost_particles_cache.valid = build_cache;
1675 const auto parallel_triangulation =
1678 Assert(parallel_triangulation,
1679 ExcMessage(
"This function is only implemented for "
1680 "parallel::TriangulationBase objects."));
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();
1689 if (send_cells.size() != 0)
1695 particles_to_send.end(),
1700 for (
auto send_particles = particles_to_send.begin();
1701 send_particles != particles_to_send.end();
1703 Assert(ghost_owners.find(send_particles->first) != ghost_owners.end(),
1706 std::size_t n_send_particles = 0;
1707 for (
auto send_particles = particles_to_send.begin();
1708 send_particles != particles_to_send.end();
1710 n_send_particles += send_particles->second.size();
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;
1723 const unsigned int individual_particle_data_size =
1725 (size_callback ? size_callback() : 0);
1727 const unsigned int individual_total_particle_data_size =
1728 individual_particle_data_size + cellid_size;
1733 if (n_send_particles > 0)
1736 send_data.resize(n_send_particles *
1737 individual_total_particle_data_size);
1739 void *data =
static_cast<void *
>(&send_data.front());
1742 for (
unsigned int i = 0; i < n_neighbors; ++i)
1744 send_offsets[i] =
reinterpret_cast<std::size_t
>(data) -
1745 reinterpret_cast<std::size_t
>(&send_data.front());
1747 const unsigned int n_particles_to_send =
1748 particles_to_send.at(neighbors[i]).size();
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 "
1758 for (
unsigned int j = 0; j < n_particles_to_send; ++j)
1764 if (send_cells.size() == 0)
1765 cell = particles_to_send.at(neighbors[i])[j]
1766 ->get_surrounding_cell();
1768 cell = send_cells.at(neighbors[i])[j];
1771 cell->id().template to_binary<dim>();
1772 memcpy(data, &cellid, cellid_size);
1773 data =
static_cast<char *
>(data) + cellid_size;
1775 data = particles_to_send.at(neighbors[i])[j]
1776 ->write_particle_data_to_memory(data);
1779 store_callback(particles_to_send.at(neighbors[i])[j], data);
1781 n_send_data[i] = n_particles_to_send;
1786 std::vector<unsigned int> n_recv_data(n_neighbors);
1787 std::vector<unsigned int> recv_offsets(n_neighbors);
1793 std::vector<MPI_Request> n_requests(2 * n_neighbors);
1794 for (
unsigned int i = 0; i < n_neighbors; ++i)
1796 const int ierr = MPI_Irecv(&(n_recv_data[i]),
1801 parallel_triangulation->get_communicator(),
1802 &(n_requests[2 * i]));
1805 for (
unsigned int i = 0; i < n_neighbors; ++i)
1807 const int ierr = MPI_Isend(&(n_send_data[i]),
1812 parallel_triangulation->get_communicator(),
1813 &(n_requests[2 * i + 1]));
1817 MPI_Waitall(2 * n_neighbors, n_requests.data(), MPI_STATUSES_IGNORE);
1822 unsigned int total_recv_data = 0;
1823 for (
unsigned int neighbor_id = 0; neighbor_id < n_neighbors; ++neighbor_id)
1825 recv_offsets[neighbor_id] = total_recv_data;
1827 n_recv_data[neighbor_id] * individual_total_particle_data_size;
1831 std::vector<char> recv_data(total_recv_data);
1835 std::vector<MPI_Request> requests(2 * n_neighbors);
1836 unsigned int send_ops = 0;
1837 unsigned int recv_ops = 0;
1842 for (
unsigned int i = 0; i < n_neighbors; ++i)
1843 if (n_recv_data[i] > 0)
1846 MPI_Irecv(&(recv_data[recv_offsets[i]]),
1847 n_recv_data[i] * individual_total_particle_data_size,
1851 parallel_triangulation->get_communicator(),
1852 &(requests[send_ops]));
1857 for (
unsigned int i = 0; i < n_neighbors; ++i)
1858 if (n_send_data[i] > 0)
1861 MPI_Isend(&(send_data[send_offsets[i]]),
1862 n_send_data[i] * individual_total_particle_data_size,
1866 parallel_triangulation->get_communicator(),
1867 &(requests[send_ops + recv_ops]));
1872 MPI_Waitall(send_ops + recv_ops, requests.data(), MPI_STATUSES_IGNORE);
1878 const void *recv_data_it =
static_cast<const void *
>(recv_data.data());
1881 auto &ghost_particles_iterators =
1882 ghost_particles_cache.ghost_particles_iterators;
1886 ghost_particles_iterators.clear();
1888 auto &send_pointers_particles = ghost_particles_cache.send_pointers;
1889 send_pointers_particles.assign(n_neighbors + 1, 0);
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;
1896 auto &recv_pointers_particles = ghost_particles_cache.recv_pointers;
1897 recv_pointers_particles.assign(n_neighbors + 1, 0);
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;
1904 ghost_particles_cache.neighbors = neighbors;
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());
1912 while (
reinterpret_cast<std::size_t
>(recv_data_it) -
1913 reinterpret_cast<std::size_t
>(recv_data.data()) <
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;
1924 insert_particle(property_pool->register_particle(), cell);
1925 const typename particle_container::iterator &cache =
1926 cells_to_particle_cache[cell->active_cell_index()];
1931 cache->particles.size() - 1);
1934 particle_it->read_particle_data_from_memory(recv_data_it);
1937 recv_data_it = load_callback(particle_it, recv_data_it);
1940 ghost_particles_iterators.push_back(particle_it);
1943 AssertThrow(recv_data_it == recv_data.data() + recv_data.size(),
1945 "The amount of data that was read into new particles "
1946 "does not match the amount of data sent around."));
1952#ifdef DEAL_II_WITH_MPI
1953 template <
int dim,
int spacedim>
1959 const auto parallel_triangulation =
1963 parallel_triangulation,
1965 "This function is only implemented for parallel::TriangulationBase "
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;
1972 std::vector<char> &send_data = ghost_particles_cache.send_data;
1975 if (send_pointers.back() > 0)
1977 void *data =
static_cast<void *
>(&send_data.front());
1980 for (
const auto i : neighbors)
1981 for (
const auto &p : particles_to_send.at(i))
1983 data = p->write_particle_data_to_memory(data);
1985 data = store_callback(p, data);
1989 std::vector<char> &recv_data = ghost_particles_cache.recv_data;
1993 std::vector<MPI_Request> requests(2 * neighbors.size());
1994 unsigned int send_ops = 0;
1995 unsigned int recv_ops = 0;
2000 for (
unsigned int i = 0; i < neighbors.size(); ++i)
2001 if ((recv_pointers[i + 1] - recv_pointers[i]) > 0)
2004 MPI_Irecv(recv_data.data() + recv_pointers[i],
2005 recv_pointers[i + 1] - recv_pointers[i],
2009 parallel_triangulation->get_communicator(),
2010 &(requests[send_ops]));
2015 for (
unsigned int i = 0; i < neighbors.size(); ++i)
2016 if ((send_pointers[i + 1] - send_pointers[i]) > 0)
2019 MPI_Isend(send_data.data() + send_pointers[i],
2020 send_pointers[i + 1] - send_pointers[i],
2024 parallel_triangulation->get_communicator(),
2025 &(requests[send_ops + recv_ops]));
2030 MPI_Waitall(send_ops + recv_ops, requests.data(), MPI_STATUSES_IGNORE);
2036 const void *recv_data_it =
static_cast<const void *
>(recv_data.data());
2039 auto &ghost_particles_iterators =
2040 ghost_particles_cache.ghost_particles_iterators;
2042 for (
auto &recv_particle : ghost_particles_iterators)
2047 recv_particle->read_particle_data_from_memory(recv_data_it);
2049 Assert(recv_particle->particles_in_cell->cell->is_ghost(),
2053 recv_data_it = load_callback(
2056 recv_particle->particle_index_within_cell),
2060 AssertThrow(recv_data_it == recv_data.data() + recv_data.size(),
2062 "The amount of data that was read into new particles "
2063 "does not match the amount of data sent around."));
2067 template <
int dim,
int spacedim>
2070 const std::function<std::size_t()> & size_callb,
2075 size_callback = size_callb;
2076 store_callback = store_callb;
2077 load_callback = load_callb;
2081 template <
int dim,
int spacedim>
2086 for (
const auto &connection : tria_listeners)
2087 connection.disconnect();
2089 tria_listeners.clear();
2091 tria_listeners.push_back(
triangulation->signals.create.connect([&]() {
2092 this->initialize(*(this->triangulation),
2094 this->property_pool->n_properties_per_slot());
2097 this->tria_listeners.push_back(
2098 this->
triangulation->signals.clear.connect([&]() { this->clear(); }));
2104 tria_listeners.push_back(
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(
2112 [&]() { this->post_mesh_change_action(); }));
2116 tria_listeners.push_back(
triangulation->signals.post_refinement.connect(
2117 [&]() { this->post_mesh_change_action(); }));
2123 template <
int dim,
int spacedim>
2129 const bool distributed_triangulation =
2132 &(*triangulation)) !=
nullptr;
2133 (void)distributed_triangulation;
2136 distributed_triangulation || number_of_locally_owned_particles == 0,
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."));
2144 if (number_of_locally_owned_particles == 0)
2145 cells_to_particle_cache.resize(
triangulation->n_active_cells(),
2151 template <
int dim,
int spacedim>
2155 register_data_attach();
2160 template <
int dim,
int spacedim>
2164 register_data_attach();
2169 template <
int dim,
int spacedim>
2173 register_data_attach();
2178 template <
int dim,
int spacedim>
2183 *distributed_triangulation =
2186 *
>(&(*triangulation)));
2187 (void)distributed_triangulation;
2190 distributed_triangulation !=
nullptr,
2192 "Mesh refinement in a non-distributed triangulation is not supported "
2193 "by the ParticleHandler class. Either insert particles after mesh "
2194 "creation and do not refine afterwards, or use a distributed triangulation."));
2196#ifdef DEAL_II_WITH_P4EST
2197 const auto callback_function =
2202 return this->pack_callback(cell_iterator, cell_status);
2206 callback_function,
true);
2212 template <
int dim,
int spacedim>
2216 const bool serialization =
false;
2217 notify_ready_to_unpack(serialization);
2222 template <
int dim,
int spacedim>
2226 const bool serialization =
true;
2227 notify_ready_to_unpack(serialization);
2231 template <
int dim,
int spacedim>
2234 const bool serialization)
2236 notify_ready_to_unpack(serialization);
2241 template <
int dim,
int spacedim>
2244 const bool serialization)
2247 *distributed_triangulation =
2250 *
>(&(*triangulation)));
2251 (void)distributed_triangulation;
2254 distributed_triangulation !=
nullptr,
2256 "Mesh refinement in a non-distributed triangulation is not supported "
2257 "by the ParticleHandler class. Either insert particles after mesh "
2258 "creation and do not refine afterwards, or use a distributed triangulation."));
2263#ifdef DEAL_II_WITH_P4EST
2269 register_data_attach();
2274 const auto callback_function =
2279 const boost::iterator_range<std::vector<char>::const_iterator>
2281 this->unpack_callback(cell_iterator, cell_status, range_iterator);
2289 update_cached_numbers();
2292 (void)serialization;
2298 template <
int dim,
int spacedim>
2304 std::vector<particle_iterator> stored_particles_on_cell;
2313 const unsigned int n_particles = n_particles_in_cell(cell);
2314 stored_particles_on_cell.reserve(n_particles);
2316 for (
unsigned int i = 0; i < n_particles; ++i)
2318 cells_to_particle_cache[cell->active_cell_index()],
2328 for (
const auto &child : cell->child_iterators())
2330 const unsigned int n_particles = n_particles_in_cell(child);
2332 stored_particles_on_cell.reserve(
2333 stored_particles_on_cell.size() + n_particles);
2335 const typename particle_container::iterator &cache =
2336 cells_to_particle_cache[child->active_cell_index()];
2337 for (
unsigned int i = 0; i < n_particles; ++i)
2338 stored_particles_on_cell.push_back(
2349 return pack_particles(stored_particles_on_cell);
2354 template <
int dim,
int spacedim>
2359 const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
2361 if (data_range.begin() == data_range.end())
2364 const auto cell_to_store_particles =
2370 if (data_range.begin() != data_range.end())
2372 const void *data =
static_cast<const void *
>(&(*data_range.begin()));
2373 const void *end =
static_cast<const void *
>(
2374 &(*data_range.begin()) + (data_range.end() - data_range.begin()));
2377 insert_particle(data, cell_to_store_particles);
2381 "The particle data could not be deserialized successfully. "
2382 "Check that when deserializing the particles you expect "
2383 "the same number of properties that were serialized."));
2386 auto loaded_particles_on_cell = particles_in_cell(cell_to_store_particles);
2401 for (
auto &particle : loaded_particles_on_cell)
2404 mapping->transform_real_to_unit_cell(cell_to_store_particles,
2405 particle.get_location());
2406 particle.set_reference_location(p_unit);
2415 typename particle_container::iterator &cache =
2416 cells_to_particle_cache[cell_to_store_particles
2417 ->active_cell_index()];
2424 auto particle = loaded_particles_on_cell.begin();
2425 for (
unsigned int i = 0; i < cache->particles.size();)
2427 bool found_new_cell =
false;
2429 for (
unsigned int child_index = 0;
2430 child_index < GeometryInfo<dim>::max_children_per_cell;
2434 child = cell->child(child_index);
2440 mapping->transform_real_to_unit_cell(
2441 child, particle->get_location());
2445 found_new_cell =
true;
2446 particle->set_reference_location(p_unit);
2451 if (child_index != 0)
2453 insert_particle(cache->particles[i], child);
2455 cache->particles[i] = cache->particles.back();
2456 cache->particles.resize(
2457 cache->particles.size() - 1);
2471 if (found_new_cell ==
false)
2480 signals.particle_lost(particle,
2481 particle->get_surrounding_cell());
2482 cache->particles[i] = cache->particles.back();
2483 cache->particles.resize(cache->particles.size() - 1);
2487 if (cache->particles.empty())
2489 particles.erase(cache);
2490 cache = particles.end();
2502#include "particle_handler.inst"
std::array< unsigned int, 4 > binary_type
const unsigned int n_components
virtual void vector_value(const Point< dim > &p, Vector< RangeNumberType > &values) const
Abstract base class for mapping classes.
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
void unpack_after_coarsening_and_refinement()
void register_data_attach()
internal::GhostParticlePartitioner< dim, spacedim > ghost_particles_cache
particle_container::iterator particle_container_ghost_begin() const
void prepare_for_serialization()
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
void update_cached_numbers()
void prepare_for_coarsening_and_refinement()
particle_container::iterator particle_container_ghost_end() const
boost::iterator_range< particle_iterator > particle_iterator_range
void post_mesh_change_action()
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 connect_to_triangulation_signals()
void register_load_callback_function(const bool serialization)
std::enable_if_t< std::is_convertible< VectorType *, Function< spacedim > * >::value==false > set_particle_positions(const VectorType &input_vector, const bool displace_particles=true)
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 unpack_callback(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)
void register_store_callback_function()
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
void update_ghost_particles()
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 sort_particles_into_subdomains_and_cells()
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
virtual ~ParticleHandler()
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
particle_container particles
std::vector< char > pack_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status) const
IndexSet locally_owned_particle_ids() const
void set_property_pool(PropertyPool< dim, spacedim > &property_pool)
const Point< dim > & get_reference_location() const
const Point< spacedim > & get_location() const
std::size_t serialized_size_in_bytes() const
const ArrayView< double > get_properties()
types::particle_index get_id() const
numbers::NumberTraits< Number >::real_type norm() const
IteratorState::IteratorStates state() const
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)
unsigned int register_data_attach(const std::function< std::vector< char >(const cell_iterator &, const CellStatus)> &pack_callback, const bool returns_variable_size_data)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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
@ past_the_end
Iterator reached end of container.
@ valid
Iterator points to a valid object.
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)
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
static const unsigned int invalid_unsigned_int
bool is_finite(const double x)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned int subdomain_id
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::map< unsigned int, std::set<::types::subdomain_id > > * vertices_with_ghost_neighbors
#define DEAL_II_PARTICLE_INDEX_MPI_TYPE
static bool is_inside_unit_cell(const Point< dim > &p)