deal.II version GIT relicensing-2055-gd703fda1ad 2024-10-30 20:00:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
generators.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
19
21
24
26
28
29#include <limits>
30
32
33namespace Particles
34{
35 namespace Generators
36 {
37 namespace
38 {
42 template <int dim>
44 ProbabilityFunctionNegative,
46 << "Your probability density function in the particle generator "
47 "returned a negative value for the following position: "
48 << arg1 << ". Please check your function expression.");
49
50
51 // This function integrates the given probability density
52 // function over all cells of the triangulation. For each
53 // cell it stores the cumulative sum (of all previous
54 // cells including the current cell) in a vector that is then
55 // returned. Therefore the returned vector has as many entries
56 // as active cells, the first entry being the integral over the
57 // first cell, and the last entry the integral over the whole
58 // locally owned domain. Cells that are not locally owned
59 // simply store the same value as the cell before (equivalent
60 // to assuming a probability density function value of 0).
61 template <int dim, int spacedim = dim>
62 std::vector<double>
63 compute_local_cumulative_cell_weights(
65 const Mapping<dim, spacedim> &mapping,
66 const Function<spacedim> &probability_density_function)
67 {
68 std::vector<double> cumulative_cell_weights(
69 triangulation.n_active_cells());
70 double cumulative_weight = 0.0;
71
72 // Evaluate function at all cell midpoints
73 const QMidpoint<dim> quadrature_formula;
74
75 // In the simplest case we do not even need a FEValues object, because
76 // using cell->center() and cell->measure() would be equivalent. This
77 // fails however for higher-order mappings.
78 FE_Nothing<dim, spacedim> alibi_finite_element;
79 FEValues<dim, spacedim> fe_values(mapping,
80 alibi_finite_element,
81 quadrature_formula,
84
85 // compute the integral weight
86 for (const auto &cell : triangulation.active_cell_iterators())
87 {
88 if (cell->is_locally_owned())
89 {
90 fe_values.reinit(cell);
91 const double quadrature_point_weight =
92 probability_density_function.value(
93 fe_values.get_quadrature_points()[0]);
94
95 AssertThrow(quadrature_point_weight >= 0.0,
96 ProbabilityFunctionNegative<dim>(
97 quadrature_formula.point(0)));
98
99 // get_cell_weight makes sure to return positive values
100 cumulative_weight += quadrature_point_weight * fe_values.JxW(0);
101 }
102 cumulative_cell_weights[cell->active_cell_index()] =
103 cumulative_weight;
104 }
105
106 return cumulative_cell_weights;
107 }
108
109
110
111 // This function generates a random position in the given cell and
112 // returns the position and its coordinates in the reference cell. It
113 // first tries to generate a random and uniformly distributed point in
114 // real space, but if that fails (e.g. because the cell has a bad aspect
115 // ratio) it reverts to generating a random point in the unit cell.
116 template <int dim, int spacedim>
117 std::pair<Point<spacedim>, Point<dim>>
118 random_location_in_cell(
120 const Mapping<dim, spacedim> &mapping,
121 std::mt19937 &random_number_generator)
122 {
123 Assert(cell->reference_cell().is_hyper_cube() == true,
125
126 // Uniform distribution on the interval [0,1]. This
127 // will be used to generate random particle locations.
128 std::uniform_real_distribution<double> uniform_distribution_01(0, 1);
129
130 const BoundingBox<spacedim> cell_bounding_box(cell->bounding_box());
131 const std::pair<Point<spacedim>, Point<spacedim>> &cell_bounds(
132 cell_bounding_box.get_boundary_points());
133
134 // Generate random points in these bounds until one is within the cell
135 // or we exceed the maximum number of attempts.
136 const unsigned int n_attempts = 100;
137 for (unsigned int i = 0; i < n_attempts; ++i)
138 {
139 Point<spacedim> position;
140 for (unsigned int d = 0; d < spacedim; ++d)
141 {
142 position[d] = uniform_distribution_01(random_number_generator) *
143 (cell_bounds.second[d] - cell_bounds.first[d]) +
144 cell_bounds.first[d];
145 }
146
147 try
148 {
149 const Point<dim> reference_position =
150 mapping.transform_real_to_unit_cell(cell, position);
151
152 if (cell->reference_cell().contains_point(reference_position))
153 return {position, reference_position};
154 }
156 {
157 // The point is not in this cell. Do nothing, just try again.
158 }
159 }
160
161 // If the above algorithm has not worked (e.g. because of badly
162 // deformed cells), retry generating particles
163 // randomly within the reference cell and then mapping it to to
164 // real space. This is not generating a
165 // uniform distribution in real space, but will always succeed.
166 Point<dim> reference_position;
167 for (unsigned int d = 0; d < dim; ++d)
168 reference_position[d] =
169 uniform_distribution_01(random_number_generator);
170
171 const Point<spacedim> position =
172 mapping.transform_unit_to_real_cell(cell, reference_position);
173
174 return {position, reference_position};
175 }
176 } // namespace
177
178
179
180 template <int dim, int spacedim>
181 void
184 const std::vector<Point<dim>> &particle_reference_locations,
185 ParticleHandler<dim, spacedim> &particle_handler,
186 const Mapping<dim, spacedim> &mapping)
187 {
188 types::particle_index particle_index = 0;
189 types::particle_index n_particles_to_generate =
190 triangulation.n_active_cells() * particle_reference_locations.size();
191
192#ifdef DEAL_II_WITH_MPI
193 if (const auto tria =
194 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
196 {
197 n_particles_to_generate = tria->n_locally_owned_active_cells() *
198 particle_reference_locations.size();
199
200 // The local particle start index is the number of all particles
201 // generated on lower MPI ranks.
202 const int ierr = MPI_Exscan(&n_particles_to_generate,
203 &particle_index,
204 1,
206 MPI_SUM,
207 tria->get_mpi_communicator());
208 AssertThrowMPI(ierr);
209 }
210#endif
211
212 particle_handler.reserve(particle_handler.n_locally_owned_particles() +
213 n_particles_to_generate);
214
215 for (const auto &cell : triangulation.active_cell_iterators() |
217 {
218 for (const auto &reference_location : particle_reference_locations)
219 {
220 const Point<spacedim> position_real =
221 mapping.transform_unit_to_real_cell(cell, reference_location);
222
223 particle_handler.insert_particle(position_real,
224 reference_location,
225 particle_index,
226 cell);
227 ++particle_index;
228 }
229 }
230
231 particle_handler.update_cached_numbers();
232 }
233
234
235
236 template <int dim, int spacedim>
240 const types::particle_index id,
241 std::mt19937 &random_number_generator,
242 const Mapping<dim, spacedim> &mapping)
243 {
244 const auto position_and_reference_position =
245 random_location_in_cell(cell, mapping, random_number_generator);
246 return Particle<dim, spacedim>(position_and_reference_position.first,
247 position_and_reference_position.second,
248 id);
249 }
250
251
252
253 template <int dim, int spacedim>
257 const types::particle_index id,
258 std::mt19937 &random_number_generator,
259 ParticleHandler<dim, spacedim> &particle_handler,
260 const Mapping<dim, spacedim> &mapping)
261 {
262 const auto position_and_reference_position =
263 random_location_in_cell(cell, mapping, random_number_generator);
264 return particle_handler.insert_particle(
265 position_and_reference_position.first,
266 position_and_reference_position.second,
267 id,
268 cell);
269 }
270
271
272
273 template <int dim, int spacedim>
274 void
277 const Function<spacedim> &probability_density_function,
278 const bool random_cell_selection,
279 const types::particle_index n_particles_to_create,
280 ParticleHandler<dim, spacedim> &particle_handler,
281 const Mapping<dim, spacedim> &mapping,
282 const unsigned int random_number_seed)
283 {
284 unsigned int combined_seed = random_number_seed;
285 if (const auto tria =
286 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
288 {
289 const unsigned int my_rank =
290 Utilities::MPI::this_mpi_process(tria->get_mpi_communicator());
291 combined_seed += my_rank;
292 }
293 std::mt19937 random_number_generator(combined_seed);
294
295 types::particle_index start_particle_id(0);
296 types::particle_index n_local_particles(0);
297
298 std::vector<types::particle_index> particles_per_cell(
299 triangulation.n_active_cells(), 0);
300
301 // First determine how many particles to generate in which cell
302 {
303 // Get the local accumulated probabilities for every cell
304 const std::vector<double> cumulative_cell_weights =
305 compute_local_cumulative_cell_weights(triangulation,
306 mapping,
307 probability_density_function);
308
309 // Sum the local integrals over all nodes
310 double local_weight_integral = (cumulative_cell_weights.size() > 0) ?
311 cumulative_cell_weights.back() :
312 0.0;
313
314
315 double local_start_weight = numbers::signaling_nan<double>();
316 double global_weight_integral = numbers::signaling_nan<double>();
317
318 if (const auto tria =
319 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
321 {
322 std::tie(local_start_weight, global_weight_integral) =
324 local_weight_integral, tria->get_mpi_communicator());
325 }
326 else
327 {
328 local_start_weight = 0;
329 global_weight_integral = local_weight_integral;
330 }
331
332 AssertThrow(global_weight_integral > std::numeric_limits<double>::min(),
334 "The integral of the user prescribed probability "
335 "density function over the domain equals zero, "
336 "deal.II has no way to determine the cell of "
337 "generated particles. Please ensure that the "
338 "provided function is positive in at least a "
339 "part of the domain; also check the syntax of "
340 "the function."));
341
342 // Calculate start id
343 start_particle_id =
344 std::llround(static_cast<double>(n_particles_to_create) *
345 local_start_weight / global_weight_integral);
346
347 // Calculate number of local particles
348 const types::particle_index end_particle_id =
349 std::llround(static_cast<double>(n_particles_to_create) *
350 ((local_start_weight + local_weight_integral) /
351 global_weight_integral));
352 n_local_particles = end_particle_id - start_particle_id;
353
354 if (random_cell_selection)
355 {
356 // Uniform distribution on the interval [0,local_weight_integral).
357 // This will be used to randomly select cells for all local
358 // particles.
359 std::uniform_real_distribution<double> uniform_distribution(
360 0.0, local_weight_integral);
361
362 // Loop over all particles to create locally and pick their cells
363 for (types::particle_index current_particle_index = 0;
364 current_particle_index < n_local_particles;
365 ++current_particle_index)
366 {
367 // Draw the random number that determines the cell of the
368 // particle
369 const double random_weight =
370 uniform_distribution(random_number_generator);
371
372 const std::vector<double>::const_iterator selected_cell =
373 std::lower_bound(cumulative_cell_weights.begin(),
374 cumulative_cell_weights.end(),
375 random_weight);
376 const unsigned int cell_index =
377 std::distance(cumulative_cell_weights.begin(), selected_cell);
378
379 ++particles_per_cell[cell_index];
380 }
381 }
382 else
383 {
384 // Compute number of particles per cell according to the ratio
385 // between their weight and the local weight integral
386 types::particle_index particles_created = 0;
387
388 for (const auto &cell : triangulation.active_cell_iterators() |
390 {
391 const types::particle_index cumulative_particles_to_create =
392 std::llround(
393 static_cast<double>(n_local_particles) *
394 cumulative_cell_weights[cell->active_cell_index()] /
395 local_weight_integral);
396
397 // Compute particles for this cell as difference between
398 // number of particles that should be created including this
399 // cell minus the number of particles already created.
400 particles_per_cell[cell->active_cell_index()] =
401 cumulative_particles_to_create - particles_created;
402 particles_created +=
403 particles_per_cell[cell->active_cell_index()];
404 }
405 }
406 }
407
408 // Now generate as many particles per cell as determined above
409 {
410 particle_handler.reserve(particle_handler.n_locally_owned_particles() +
411 n_local_particles);
412 unsigned int current_particle_index = start_particle_id;
413
414 for (const auto &cell : triangulation.active_cell_iterators() |
416 {
417 for (unsigned int i = 0;
418 i < particles_per_cell[cell->active_cell_index()];
419 ++i)
420 {
422 current_particle_index,
423 random_number_generator,
424 particle_handler,
425 mapping);
426
427 ++current_particle_index;
428 }
429 }
430
431 particle_handler.update_cached_numbers();
432 }
433 }
434
435
436
437 template <int dim, int spacedim>
438 void
440 const std::vector<std::vector<BoundingBox<spacedim>>>
441 &global_bounding_boxes,
442 ParticleHandler<dim, spacedim> &particle_handler,
443 const Mapping<dim, spacedim> &mapping,
444 const ComponentMask &components,
445 const std::vector<std::vector<double>> &properties)
446 {
447 const auto &fe = dof_handler.get_fe();
448
449 // Take care of components
450 const ComponentMask mask =
451 (components.size() == 0 ? ComponentMask(fe.n_components(), true) :
452 components);
453
454 const std::map<types::global_dof_index, Point<spacedim>>
455 support_points_map =
456 DoFTools::map_dofs_to_support_points(mapping, dof_handler, mask);
457
458 // Generate the vector of points from the map
459 // Memory is reserved for efficiency reasons
460 std::vector<Point<spacedim>> support_points_vec;
461 support_points_vec.reserve(support_points_map.size());
462 for (const auto &element : support_points_map)
463 support_points_vec.push_back(element.second);
464
465 particle_handler.insert_global_particles(support_points_vec,
466 global_bounding_boxes,
467 properties);
468 }
469
470
471 template <int dim, int spacedim>
472 void
475 const Quadrature<dim> &quadrature,
476 // const std::vector<Point<dim>> & particle_reference_locations,
477 const std::vector<std::vector<BoundingBox<spacedim>>>
478 &global_bounding_boxes,
479 ParticleHandler<dim, spacedim> &particle_handler,
480 const Mapping<dim, spacedim> &mapping,
481 const std::vector<std::vector<double>> &properties)
482 {
483 const std::vector<Point<dim>> &particle_reference_locations =
484 quadrature.get_points();
485 std::vector<Point<spacedim>> points_to_generate;
486 points_to_generate.reserve(particle_reference_locations.size() *
487 triangulation.n_active_cells());
488
489 // Loop through cells and gather gauss points
490 for (const auto &cell : triangulation.active_cell_iterators() |
492 {
493 for (const auto &reference_location : particle_reference_locations)
494 {
495 const Point<spacedim> position_real =
496 mapping.transform_unit_to_real_cell(cell, reference_location);
497 points_to_generate.push_back(position_real);
498 }
499 }
500 particle_handler.insert_global_particles(points_to_generate,
501 global_bounding_boxes,
502 properties);
503 }
504 } // namespace Generators
505} // namespace Particles
506
507#include "generators.inst"
508
unsigned int size() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
virtual RangeNumberType value(const Point< dim > &p, const unsigned int component=0) const
Abstract base class for mapping classes.
Definition mapping.h:318
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
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={})
types::particle_index n_locally_owned_particles() const
particle_iterator insert_particle(const Particle< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
Definition point.h:111
const std::vector< Point< dim > > & get_points() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int cell_index
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertThrowMPI(error_code)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
@ update_JxW_values
Transformed quadrature weights.
@ update_quadrature_points
Transformed quadrature points.
void map_dofs_to_support_points(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, std::vector< Point< spacedim > > &support_points, const ComponentMask &mask={})
void probabilistic_locations(const Triangulation< dim, spacedim > &triangulation, const Function< spacedim > &probability_density_function, const bool random_cell_selection, const types::particle_index n_particles_to_create, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const unsigned int random_number_seed=5432)
ParticleIterator< dim, spacedim > random_particle_in_cell_insert(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const types::particle_index id, std::mt19937 &random_number_generator, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
void quadrature_points(const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
Particle< dim, spacedim > random_particle_in_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const types::particle_index id, std::mt19937 &random_number_generator, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
void dof_support_points(const DoFHandler< dim, spacedim > &dof_handler, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const ComponentMask &components={}, const std::vector< std::vector< double > > &properties={})
void regular_reference_locations(const Triangulation< dim, spacedim > &triangulation, const std::vector< Point< dim > > &particle_reference_locations, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::pair< T, T > partial_and_total_sum(const T &value, const MPI_Comm comm)
unsigned int this_mpi_process(const MPI_Comm mpi_communicator)
Definition mpi.cc:107
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
#define DEAL_II_PARTICLE_INDEX_MPI_TYPE