ASPECT
distribution.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2025 - by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef _aspect_particle_distribution_h
22 #define _aspect_particle_distribution_h
23 
26 #include <algorithm>
27 #include <limits>
28 #include <deal.II/base/table.h>
29 #include <deal.II/particles/property_pool.h>
30 #include <deal.II/particles/particle_handler.h>
31 #include <deal.II/base/function_lib.h>
32 #include <vector>
33 
34 namespace aspect
35 {
36  namespace Particle
37  {
47  template <int dim>
49  {
50  public:
51 
57  enum class KernelFunction
58  {
59  gaussian,
60  triangular,
61  uniform,
64  };
65 
76  ParticlePDF(const unsigned int granularity,
77  const double bandwidth,
79 
88  ParticlePDF(const double bandwidth,
89  const KernelFunction kernel_function);
90 
97  void
98  fill_from_particle_range(const typename Particles::ParticleHandler<dim>::particle_iterator_range particle_range,
99  const std::vector<typename Particles::ParticleHandler<dim>::particle_iterator_range> particle_ranges_to_sum_over,
100  const unsigned int n_particles_in_cell);
101 
112  void
113  insert_kernel_sum_from_particle_range(const Point<dim> reference_point,
114  std::array<unsigned int,dim> table_index,
115  const unsigned int n_particles_in_cell,
116  const typename Particles::ParticleHandler<dim>::particle_iterator_range particle_range);
117 
125  void
126  add_value_to_function_table(const unsigned int x_index,
127  const unsigned int y_index,
128  const unsigned int z_index,
129  const double input_value);
130 
137  void
138  add_value_to_function_table(const std::array<unsigned int,dim> &index_point,
139  const double input_value);
140 
148  void
149  add_value_to_function_table(const double input_value,
150  const types::particle_index reference_particle_id);
151 
159  double
160  evaluate_function_at_index(const unsigned int x_index,
161  const unsigned int y_index,
162  const unsigned int z_index) const;
163 
167  void
169 
173  double
174  get_max() const;
175 
179  double
180  get_min() const;
181 
185  double
186  get_standard_deviation() const;
187 
193  types::particle_index
194  get_max_particle() const;
195 
201  types::particle_index
202  get_min_particle() const;
203 
204  private:
211  Table<dim,double> function_output_table;
212 
221  double bandwidth;
222 
228 
234  unsigned int granularity;
235 
239  double max;
240 
244  double min;
245 
250 
254  double mean;
255 
262 
269 
274  types::particle_index max_particle_index;
275 
280  types::particle_index min_particle_index;
281 
287  double apply_selected_kernel_function(const double distance) const;
288 
295  double kernelfunction_uniform(const double distance) const;
296 
303  double kernelfunction_triangular(const double distance) const;
304 
311  double kernelfunction_gaussian(const double distance) const;
312  };
313  }
314 }
315 
316 #endif
double kernelfunction_triangular(const double distance) const
void fill_from_particle_range(const typename Particles::ParticleHandler< dim >::particle_iterator_range particle_range, const std::vector< typename Particles::ParticleHandler< dim >::particle_iterator_range > particle_ranges_to_sum_over, const unsigned int n_particles_in_cell)
ParticlePDF(const unsigned int granularity, const double bandwidth, const KernelFunction kernel_function)
types::particle_index get_min_particle() const
double kernelfunction_uniform(const double distance) const
double get_standard_deviation() const
small_vector< double > function_output_vector
Definition: distribution.h:261
types::particle_index get_max_particle() const
double kernelfunction_gaussian(const double distance) const
types::particle_index min_particle_index
Definition: distribution.h:280
void insert_kernel_sum_from_particle_range(const Point< dim > reference_point, std::array< unsigned int, dim > table_index, const unsigned int n_particles_in_cell, const typename Particles::ParticleHandler< dim >::particle_iterator_range particle_range)
double apply_selected_kernel_function(const double distance) const
void add_value_to_function_table(const unsigned int x_index, const unsigned int y_index, const unsigned int z_index, const double input_value)
boost::container::small_vector< T, N > small_vector
Definition: global.h:252
Table< dim, double > function_output_table
Definition: distribution.h:211
types::particle_index max_particle_index
Definition: distribution.h:274
double evaluate_function_at_index(const unsigned int x_index, const unsigned int y_index, const unsigned int z_index) const