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 
100  void
101  fill_from_particle_range(const typename Particles::ParticleHandler<dim>::particle_iterator_range &particle_range,
102  const std::vector<typename Particles::ParticleHandler<dim>::particle_iterator_range>
103  &particle_ranges_to_sum_over,
104  const unsigned int n_particles_in_cell,
105  const typename ::Mapping<dim> &mapping,
106  const typename Triangulation<dim>::active_cell_iterator &cell);
107 
119  void
120  insert_kernel_sum_from_particle_range(const Point<dim> &reference_point,
121  const std::array<unsigned int,dim> &table_index,
122  const typename Triangulation<dim>::active_cell_iterator &cell,
123  const std::vector<typename Particles::ParticleHandler<dim>::particle_iterator_range>
124  &particle_ranges_to_sum_over,
125  const typename ::Mapping<dim> &mapping);
126 
134  void
135  add_value_to_function_table(const unsigned int x_index,
136  const unsigned int y_index,
137  const unsigned int z_index,
138  const double input_value);
139 
146  void
147  add_value_to_function_table(const std::array<unsigned int,dim> &index_point,
148  const double input_value);
149 
157  void
158  add_value_to_function_table(const double input_value,
159  const types::particle_index reference_particle_id);
160 
168  double
169  evaluate_function_at_index(const unsigned int x_index,
170  const unsigned int y_index,
171  const unsigned int z_index) const;
172 
176  void
178 
182  double
183  get_max() const;
184 
188  double
189  get_min() const;
190 
194  Point<dim>
195  get_max_position() const;
196 
200  std::vector<Point<dim>>
201  get_min_positions() const;
202 
206  double
207  get_standard_deviation() const;
208 
214  types::particle_index
215  get_max_particle() const;
216 
222  types::particle_index
223  get_min_particle() const;
224 
225  private:
232  Table<dim,double> function_output_table;
233 
242  double bandwidth;
243 
249 
255  unsigned int granularity;
256 
260  double max;
261 
265  double min;
266 
270  Point<dim> max_position;
271 
277  std::vector<Point<dim>> min_positions;
278 
283 
287  double mean;
288 
295 
302 
307  types::particle_index max_particle_index;
308 
313  types::particle_index min_particle_index;
314 
320  double apply_selected_kernel_function(const double distance) const;
321 
328  double kernelfunction_uniform(const double distance) const;
329 
336  double kernelfunction_triangular(const double distance) const;
337 
344  double kernelfunction_gaussian(const double distance) const;
345  };
346  }
347 }
348 
349 #endif
double kernelfunction_triangular(const double distance) const
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:294
void insert_kernel_sum_from_particle_range(const Point< dim > &reference_point, const std::array< unsigned int, dim > &table_index, const typename Triangulation< dim >::active_cell_iterator &cell, const std::vector< typename Particles::ParticleHandler< dim >::particle_iterator_range > &particle_ranges_to_sum_over, const typename ::Mapping< dim > &mapping)
std::vector< Point< dim > > min_positions
Definition: distribution.h:277
types::particle_index get_max_particle() const
std::vector< Point< dim > > get_min_positions() const
double kernelfunction_gaussian(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, const typename ::Mapping< dim > &mapping, const typename Triangulation< dim >::active_cell_iterator &cell)
types::particle_index min_particle_index
Definition: distribution.h:313
Point< dim > get_max_position() const
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:232
types::particle_index max_particle_index
Definition: distribution.h:307
double evaluate_function_at_index(const unsigned int x_index, const unsigned int y_index, const unsigned int z_index) const