![]() |
ASPECT
|
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | get_composition_values_at_q_point (const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point) |
Private Member Functions | |
void | fill_PDF_from_cell (const typename Triangulation< dim >::active_cell_iterator &cell, Particle::ParticlePDF< dim > &pdf) |
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) |
double | apply_selected_kernel_function (const double distance) const |
double | kernelfunction_uniform (const double distance) const |
double | kernelfunction_triangular (const double distance) const |
double | kernelfunction_gaussian (const double distance) const |
Private Attributes | |
bool | KDE_per_particle |
unsigned int | granularity |
double | bandwidth |
Particle::ParticlePDF< dim >::KernelFunction | kernel_function |
A postprocessor that generates a point-density function describing the distribution of particles within cells. The point-density function is generated using a Kernel Density Estimator. The postprocessor calculates standard deviation and max/min of the point-density functions per cell.
Definition at line 43 of file particle_distribution_statistics.h.
|
overridevirtual |
Evaluate the solution for some particle statistics.
Implements aspect::Postprocess::Interface< dim >.
|
overridevirtual |
Let the postprocessor manager know about the other postprocessors this one depends on. Specifically, the particles postprocessor.
Reimplemented from aspect::Postprocess::Interface< dim >.
|
static |
Declare the parameters this class takes through input files.
|
overridevirtual |
Read the parameters this class declares from the parameter file.
Reimplemented from aspect::Plugins::InterfaceBase.
|
private |
Fills the supplied PDF instance with values from the particles in the given cell.
cell | The cell whose particles to populate the PDF with. |
The instance of ParticlePDF to fill with data from the cell's particles. |
|
private |
This function is only called from fill_PDF_from_cell
. It iterates through every particle in the cell and sums the value of the kernel function between the reference point and the position of the cell.
cell | The cell whose particles to populate the PDF with. |
reference_point | The point from which to get the value of the kernel function. |
table_index | The index in the PDF to insert the data into. |
particles_in_cell | The number of particles in the cell. |
The instance of ParticlePDF to fill with data from the cell's particles. |
|
private |
Returns the value of the selected kernel function.
distance | the distance to pass to the selected kernel function. |
|
private |
The Uniform kernel function returns a value of 1.0 as long as the distance is less than the KDE's bandwidth.
distance | the output of the kernel function depends on the distance between the reference point and the center of the kernel function. |
|
private |
The Triangular kernel function returns a value of 1.0 minus the distance variable.
distance | the output of the kernel function depends on the distance between the reference point and the center of the kernel function. |
|
private |
The gaussian function returns the value of a gaussian distribution at the specified distance.
distance | the output of the kernel function depends on the distance between the reference point and the center of the kernel function. |
|
private |
If KDE_per_particle
is true, the point-density function is defined at the position of every particle in the cell. If it is false, the point-density-function is defined on a regular grid throughout the cell.
Definition at line 78 of file particle_distribution_statistics.h.
|
private |
The granularity
variable determines the number of points at which the point-density function is defined within each cell. For example, a value of 2 means that the point-density function is defined at \(2\times 2=4\) points in 2D. This variable only applies if KDE_per_particle
is false.
Definition at line 87 of file particle_distribution_statistics.h.
|
private |
The bandwidth
variable scales the point-density function. Choosing an appropriate bandwidth is important because a bandwidth value which is too low or too high can result in oversmoothing or undersmoothing of the point-density function. Oversmoothing or undersmoothing results in a function which represents the underlying data less accurately.
Definition at line 97 of file particle_distribution_statistics.h.
|
private |
kernel_function
is an internal variable to keep track of which kernel function was read from the .prm file.
Definition at line 103 of file particle_distribution_statistics.h.