ASPECT
particle_distribution_statistics.h
Go to the documentation of this file.
1 
2 /*
3  Copyright (C) 2025 - by the authors of the ASPECT code.
4 
5  This file is part of ASPECT.
6 
7  ASPECT is free software; you can redistribute it and/or modify
8  it under the terms of the GNU General Public License as published by
9  the Free Software Foundation; either version 2, or (at your option)
10  any later version.
11 
12  ASPECT is distributed in the hope that it will be useful,
13  but WITHOUT ANY WARRANTY; without even the implied warranty of
14  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  GNU General Public License for more details.
16 
17  You should have received a copy of the GNU General Public License
18  along with ASPECT; see the file LICENSE. If not see
19  <http://www.gnu.org/licenses/>.
20 */
21 
22 #ifndef _aspect_postprocess_particle_distribution_statistics_h
23 #define _aspect_postprocess_particle_distribution_statistics_h
24 
26 #include <deal.II/base/function_lib.h>
28 
29 namespace aspect
30 {
31  namespace Postprocess
32  {
33 
42  template <int dim>
44  {
45  public:
49  void initialize() override;
50 
54  std::pair<std::string,std::string>
55  execute (TableHandler &statistics) override;
56 
61  std::list<std::string>
62  required_other_postprocessors() const override;
63 
67  static
68  void
69  declare_parameters (ParameterHandler &prm);
70 
74  void
75  parse_parameters (ParameterHandler &prm) override;
76 
77  private:
84 
92  unsigned int granularity;
93 
102  double bandwidth;
103 
109 
117  std::vector<typename Particles::ParticleHandler<dim>::particle_iterator_range>
118  get_neighboring_particle_ranges(const typename Triangulation<dim>::active_cell_iterator &cell,
119  const typename Particles::ParticleHandler<dim> &particle_handler);
120 
127  void
128  fill_PDF_from_cell(const typename Triangulation<dim>::active_cell_iterator &cell,
130 
143  void insert_kernel_sum_into_pdf(const typename Triangulation<dim>::active_cell_iterator &cell,
144  const Point<dim> reference_point,
145  const std::array<unsigned int,dim> table_index,
146  const unsigned int particles_in_cell,
148 
154  double
155  apply_selected_kernel_function(const double distance) const;
156 
163  double
164  kernelfunction_uniform(const double distance) const;
165 
172  double
173  kernelfunction_triangular(const double distance) const;
174 
181  double
182  kernelfunction_gaussian(const double distance) const;
183 
188  std::unique_ptr<GridTools::Cache<dim>> grid_cache;
189  };
190  }
191 }
192 
193 
194 #endif
void insert_kernel_sum_into_pdf(const typename Triangulation< dim >::active_cell_iterator &cell, const Point< dim > reference_point, const std::array< unsigned int, dim > table_index, const unsigned int particles_in_cell, Particle::ParticlePDF< dim > &pdf)
std::vector< typename Particles::ParticleHandler< dim >::particle_iterator_range > get_neighboring_particle_ranges(const typename Triangulation< dim >::active_cell_iterator &cell, const typename Particles::ParticleHandler< dim > &particle_handler)
double apply_selected_kernel_function(const double distance) const
void fill_PDF_from_cell(const typename Triangulation< dim >::active_cell_iterator &cell, Particle::ParticlePDF< dim > &pdf)
double kernelfunction_triangular(const double distance) const
Particle::ParticlePDF< dim >::KernelFunction kernel_function
double kernelfunction_gaussian(const double distance) const
double kernelfunction_uniform(const double distance) const
std::pair< std::string, std::string > execute(TableHandler &statistics) override
static void declare_parameters(ParameterHandler &prm)
void parse_parameters(ParameterHandler &prm) override
std::list< std::string > required_other_postprocessors() const override