ASPECT
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
aspect::Postprocess::ParticleDistributionStatistics< dim > Class Template Reference
Inheritance diagram for aspect::Postprocess::ParticleDistributionStatistics< dim >:
Inheritance graph
[legend]

Public Member Functions

std::pair< std::string, std::string > execute (TableHandler &statistics) override
 
std::list< std::string > required_other_postprocessors () const override
 
void parse_parameters (ParameterHandler &prm) override
 
- Public Member Functions inherited from aspect::Postprocess::Interface< dim >
virtual void save (std::map< std::string, std::string > &status_strings) const
 
virtual void load (const std::map< std::string, std::string > &status_strings)
 
- Public Member Functions inherited from aspect::Plugins::InterfaceBase
virtual ~InterfaceBase ()=default
 
virtual void initialize ()
 
virtual void update ()
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
const Introspection< dim > & introspection () const
 
const Simulator< dim > & get_simulator () const
 
const Parameters< dim > & get_parameters () const
 
SimulatorSignals< dim > & get_signals () const
 
MPI_Comm get_mpi_communicator () const
 
TimerOutput & get_computing_timer () const
 
const ConditionalOStream & get_pcout () const
 
double get_time () const
 
double get_timestep () const
 
double get_old_timestep () const
 
unsigned int get_timestep_number () const
 
const TimeStepping::Manager< dim > & get_timestepping_manager () const
 
unsigned int get_nonlinear_iteration () const
 
const parallel::distributed::Triangulation< dim > & get_triangulation () const
 
double get_volume () const
 
const Mapping< dim > & get_mapping () const
 
std::string get_output_directory () const
 
bool include_adiabatic_heating () const
 
bool include_latent_heat () const
 
bool include_melt_transport () const
 
int get_stokes_velocity_degree () const
 
double get_adiabatic_surface_temperature () const
 
double get_surface_pressure () const
 
bool convert_output_to_years () const
 
unsigned int get_pre_refinement_step () const
 
unsigned int n_compositional_fields () const
 
double get_end_time () const
 
void get_refinement_criteria (Vector< float > &estimated_error_per_cell) const
 
void get_artificial_viscosity (Vector< float > &viscosity_per_cell, const bool skip_interior_cells=false) const
 
void get_artificial_viscosity_composition (Vector< float > &viscosity_per_cell, const unsigned int compositional_variable) const
 
const LinearAlgebra::BlockVectorget_current_linearization_point () const
 
const LinearAlgebra::BlockVectorget_solution () const
 
const LinearAlgebra::BlockVectorget_old_solution () const
 
const LinearAlgebra::BlockVectorget_old_old_solution () const
 
const LinearAlgebra::BlockVectorget_reaction_vector () const
 
const LinearAlgebra::BlockVectorget_mesh_velocity () const
 
const DoFHandler< dim > & get_dof_handler () const
 
const FiniteElement< dim > & get_fe () const
 
const LinearAlgebra::BlockSparseMatrixget_system_matrix () const
 
const LinearAlgebra::BlockSparseMatrixget_system_preconditioner_matrix () const
 
const MaterialModel::Interface< dim > & get_material_model () const
 
const GravityModel::Interface< dim > & get_gravity_model () const
 
const InitialTopographyModel::Interface< dim > & get_initial_topography_model () const
 
const std::shared_ptr< const InitialTopographyModel::Interface< dim > > get_initial_topography_model_pointer () const
 
const GeometryModel::Interface< dim > & get_geometry_model () const
 
const AdiabaticConditions::Interface< dim > & get_adiabatic_conditions () const
 
bool has_boundary_temperature () const
 
const BoundaryTemperature::Manager< dim > & get_boundary_temperature_manager () const
 
const BoundaryConvectiveHeating::Manager< dim > & get_boundary_convective_heating_manager () const
 
const BoundaryHeatFlux::Interface< dim > & get_boundary_heat_flux () const
 
bool has_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const BoundaryTraction::Manager< dim > & get_boundary_traction_manager () const
 
std::shared_ptr< const InitialTemperature::Manager< dim > > get_initial_temperature_manager_pointer () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
std::shared_ptr< const InitialComposition::Manager< dim > > get_initial_composition_manager_pointer () const
 
const InitialComposition::Manager< dim > & get_initial_composition_manager () const
 
const std::set< types::boundary_id > & get_fixed_temperature_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_heat_flux_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_convective_heating_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_composition_boundary_indicators () const
 
const std::set< types::boundary_id > & get_mesh_deformation_boundary_indicators () const
 
const BoundaryVelocity::Manager< dim > & get_boundary_velocity_manager () const
 
const HeatingModel::Manager< dim > & get_heating_model_manager () const
 
const MeshRefinement::Manager< dim > & get_mesh_refinement_manager () const
 
const MeltHandler< dim > & get_melt_handler () const
 
const VolumeOfFluidHandler< dim > & get_volume_of_fluid_handler () const
 
const NewtonHandler< dim > & get_newton_handler () const
 
const MeshDeformation::MeshDeformationHandler< dim > & get_mesh_deformation_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const AffineConstraints< double > & get_current_constraints () const
 
bool simulator_is_past_initialization () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandler & get_statistics_object () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 
unsigned int n_particle_managers () const
 
const Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index) const
 
Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index)
 
bool is_stokes_matrix_free ()
 
const StokesMatrixFreeHandler< dim > & get_stokes_matrix_free () const
 
RotationProperties< dim > compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const
 
void remove_nullspace (LinearAlgebra::BlockVector &solution, LinearAlgebra::BlockVector &distributed_stokes_solution) const
 
double normalize_pressure (LinearAlgebra::BlockVector &vector) const
 
void denormalize_pressure (const double pressure_adjustment, LinearAlgebra::BlockVector &vector) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::Plugins::InterfaceBase
static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::SimulatorAccess< dim >
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
 

Detailed Description

template<int dim>
class aspect::Postprocess::ParticleDistributionStatistics< dim >

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.

Member Function Documentation

§ execute()

template<int dim>
std::pair<std::string,std::string> aspect::Postprocess::ParticleDistributionStatistics< dim >::execute ( TableHandler &  statistics)
overridevirtual

Evaluate the solution for some particle statistics.

Implements aspect::Postprocess::Interface< dim >.

§ required_other_postprocessors()

template<int dim>
std::list<std::string> aspect::Postprocess::ParticleDistributionStatistics< dim >::required_other_postprocessors ( ) const
overridevirtual

Let the postprocessor manager know about the other postprocessors this one depends on. Specifically, the particles postprocessor.

Reimplemented from aspect::Postprocess::Interface< dim >.

§ declare_parameters()

template<int dim>
static void aspect::Postprocess::ParticleDistributionStatistics< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

template<int dim>
void aspect::Postprocess::ParticleDistributionStatistics< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Read the parameters this class declares from the parameter file.

Reimplemented from aspect::Plugins::InterfaceBase.

§ fill_PDF_from_cell()

template<int dim>
void aspect::Postprocess::ParticleDistributionStatistics< dim >::fill_PDF_from_cell ( const typename Triangulation< dim >::active_cell_iterator &  cell,
Particle::ParticlePDF< dim > &  pdf 
)
private

Fills the supplied PDF instance with values from the particles in the given cell.

Parameters
cellThe cell whose particles to populate the PDF with.
pdfThe instance of ParticlePDF to fill with data from the cell's particles.

§ insert_kernel_sum_into_pdf()

template<int dim>
void aspect::Postprocess::ParticleDistributionStatistics< dim >::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 
)
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.

Parameters
cellThe cell whose particles to populate the PDF with.
reference_pointThe point from which to get the value of the kernel function.
table_indexThe index in the PDF to insert the data into.
particles_in_cellThe number of particles in the cell.
pdfThe instance of ParticlePDF to fill with data from the cell's particles.

§ apply_selected_kernel_function()

template<int dim>
double aspect::Postprocess::ParticleDistributionStatistics< dim >::apply_selected_kernel_function ( const double  distance) const
private

Returns the value of the selected kernel function.

Parameters
distancethe distance to pass to the selected kernel function.

§ kernelfunction_uniform()

template<int dim>
double aspect::Postprocess::ParticleDistributionStatistics< dim >::kernelfunction_uniform ( const double  distance) const
private

The Uniform kernel function returns a value of 1.0 as long as the distance is less than the KDE's bandwidth.

Parameters
distancethe output of the kernel function depends on the distance between the reference point and the center of the kernel function.

§ kernelfunction_triangular()

template<int dim>
double aspect::Postprocess::ParticleDistributionStatistics< dim >::kernelfunction_triangular ( const double  distance) const
private

The Triangular kernel function returns a value of 1.0 minus the distance variable.

Parameters
distancethe output of the kernel function depends on the distance between the reference point and the center of the kernel function.

§ kernelfunction_gaussian()

template<int dim>
double aspect::Postprocess::ParticleDistributionStatistics< dim >::kernelfunction_gaussian ( const double  distance) const
private

The gaussian function returns the value of a gaussian distribution at the specified distance.

Parameters
distancethe output of the kernel function depends on the distance between the reference point and the center of the kernel function.

Member Data Documentation

§ KDE_per_particle

template<int dim>
bool aspect::Postprocess::ParticleDistributionStatistics< dim >::KDE_per_particle
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.

§ granularity

template<int dim>
unsigned int aspect::Postprocess::ParticleDistributionStatistics< dim >::granularity
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.

§ bandwidth

template<int dim>
double aspect::Postprocess::ParticleDistributionStatistics< dim >::bandwidth
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.

§ kernel_function

template<int dim>
Particle::ParticlePDF<dim>::KernelFunction aspect::Postprocess::ParticleDistributionStatistics< dim >::kernel_function
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.


The documentation for this class was generated from the following file: