From: Wolfgang Bangerth Date: Mon, 29 Jun 2020 01:03:39 +0000 (-0600) Subject: Warn against the use of certain functions. X-Git-Tag: v9.3.0-rc1~1105^2~5 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=78a63e790a8dcf7570f6745aebfdfedd6c230f5e;p=dealii.git Warn against the use of certain functions. --- diff --git a/include/deal.II/particles/particle_handler.h b/include/deal.II/particles/particle_handler.h index a06a03b27c..4038da69bd 100644 --- a/include/deal.II/particles/particle_handler.h +++ b/include/deal.II/particles/particle_handler.h @@ -177,6 +177,26 @@ namespace Particles /** * Return the number of particles that live on the given cell. + * + * @note While this function is used in step-19, it is not an efficient + * function to use if the number of particles is large. That is because + * to find the particles that are located in one cell costs + * ${\cal O)(\log N)$ where $N$ is the number of overall particles. Since + * you will likely do this for every cell, and assuming that the number + * of particles and the number of cells are roughly proportional, + * you end up with an ${\cal O)(N \log N)$ algorithm. A better approach + * is to use the fact that internally, particles are arranged in the + * order of the active cells they are in. In other words, if you iterate + * over all particles, you will encounter them in the same order as + * you walk over the active cells. You can exploit this by keeping an + * iterator to the first particle of the first cell, and when you move + * to the next cell, you increment the particle iterator as well until + * you find a particle located on that next cell. Counting how many + * steps this took will then give you the number you are looking for, + * at a cost of ${\cal O)(\log N)$ when accumulated over all cells. + * This is the approach used in step-70, for example. The approach is + * also detailed in the "Possibilities for extensions section" + * of step-19. */ types::particle_index n_particles_in_cell( @@ -190,6 +210,25 @@ namespace Particles * * The number of elements in the returned range equals what the * n_particles_in_cell() function returns. + * + * @note While this function is used in step-19, it is not an efficient + * function to use if the number of particles is large. That is because + * to find the particles that are located in one cell costs + * ${\cal O)(\log N)$ where $N$ is the number of overall particles. Since + * you will likely do this for every cell, and assuming that the number + * of particles and the number of cells are roughly proportional, + * you end up with an ${\cal O)(N \log N)$ algorithm. A better approach + * is to use the fact that internally, particles are arranged in the + * order of the active cells they are in. In other words, if you iterate + * over all particles, you will encounter them in the same order as + * you walk over the active cells. You can exploit this by keeping an + * iterator to the first particle of the first cell, and when you move + * to the next cell, you increment the particle iterator as well until + * you find a particle located on that next cell. This is the approach + * used in step-70, for example, and has an overall cost of + * ${\cal O)(\log N)$ when accumulated over all cells. The approach is + * also detailed in the "Possibilities for extensions section" + * of step-19. */ particle_iterator_range particles_in_cell( @@ -202,6 +241,25 @@ namespace Particles * * The number of elements in the returned range equals what the * n_particles_in_cell() function returns. + * + * @note While this function is used in step-19, it is not an efficient + * function to use if the number of particles is large. That is because + * to find the particles that are located in one cell costs + * ${\cal O)(\log N)$ where $N$ is the number of overall particles. Since + * you will likely do this for every cell, and assuming that the number + * of particles and the number of cells are roughly proportional, + * you end up with an ${\cal O)(N \log N)$ algorithm. A better approach + * is to use the fact that internally, particles are arranged in the + * order of the active cells they are in. In other words, if you iterate + * over all particles, you will encounter them in the same order as + * you walk over the active cells. You can exploit this by keeping an + * iterator to the first particle of the first cell, and when you move + * to the next cell, you increment the particle iterator as well until + * you find a particle located on that next cell. This is the approach + * used in step-70, for example, and has an overall cost of + * ${\cal O)(\log N)$ when accumulated over all cells. The approach is + * also detailed in the "Possibilities for extensions section" + * of step-19. */ particle_iterator_range particles_in_cell(