21 #ifndef _aspect_particle_distribution_h 22 #define _aspect_particle_distribution_h 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> 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);
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);
136 const unsigned int y_index,
137 const unsigned int z_index,
138 const double input_value);
148 const double input_value);
159 const types::particle_index reference_particle_id);
170 const unsigned int y_index,
171 const unsigned int z_index)
const;
200 std::vector<Point<dim>>
214 types::particle_index
222 types::particle_index
double kernelfunction_triangular(const double distance) const
double standard_deviation
void compute_statistical_values()
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
bool is_defined_per_particle
small_vector< double > function_output_vector
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
KernelFunction kernel_function
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
Point< dim > get_max_position() const
Point< dim > max_position
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
Table< dim, double > function_output_table
types::particle_index max_particle_index
double evaluate_function_at_index(const unsigned int x_index, const unsigned int y_index, const unsigned int z_index) const