ASPECT
Public Member Functions | Static Public Member Functions | List of all members
aspect::Particle::Interpolator::Interface< dim > Class Template Referenceabstract
Inheritance diagram for aspect::Particle::Interpolator::Interface< dim >:
Inheritance graph
[legend]

Public Member Functions

virtual ~Interface ()=default
 
virtual DEAL_II_DEPRECATED std::vector< std::vector< double > > properties_at_points (const ParticleHandler< dim > &particle_handler, const std::vector< Point< dim >> &positions, const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &cell=typename parallel::distributed::Triangulation< dim >::active_cell_iterator()) const
 
virtual std::vector< std::vector< double > > properties_at_points (const ParticleHandler< dim > &particle_handler, const std::vector< Point< dim >> &positions, const ComponentMask &selected_properties, const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &cell=typename parallel::distributed::Triangulation< dim >::active_cell_iterator()) const =0
 
virtual void parse_parameters (ParameterHandler &prm)
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 

Detailed Description

template<int dim>
class aspect::Particle::Interpolator::Interface< dim >

An abstract class defining virtual methods for performing interpolation of particle properties to arbitrary points.

Definition at line 50 of file interface.h.

Constructor & Destructor Documentation

§ ~Interface()

template<int dim>
virtual aspect::Particle::Interpolator::Interface< dim >::~Interface ( )
virtualdefault

Destructor. Made virtual so that derived classes can be created and destroyed through pointers to the base class.

Member Function Documentation

§ properties_at_points() [1/2]

template<int dim>
virtual DEAL_II_DEPRECATED std::vector<std::vector<double> > aspect::Particle::Interpolator::Interface< dim >::properties_at_points ( const ParticleHandler< dim > &  particle_handler,
const std::vector< Point< dim >> &  positions,
const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &  cell = typename parallel::distributed::Triangulation< dim >::active_cell_iterator() 
) const
virtual

Perform an interpolation of the properties of the particles in this cell onto a vector of positions in this cell. Implementations of this function must return a vector of a vector of doubles which contains a somehow computed value of all particle properties at all given positions.

Parameters
[in]particle_handlerReference to the particle handler that allows accessing the particles in the domain.
[in]positionsThe vector of positions where the properties should be evaluated.
[in]cellAn optional iterator to the cell containing the particles. Not all callers will know the cell of the particles, but providing the cell when known speeds up the interpolation significantly.
Returns
A vector with as many entries as positions. Every entry is a vector of interpolated particle properties at this position.

§ properties_at_points() [2/2]

template<int dim>
virtual std::vector<std::vector<double> > aspect::Particle::Interpolator::Interface< dim >::properties_at_points ( const ParticleHandler< dim > &  particle_handler,
const std::vector< Point< dim >> &  positions,
const ComponentMask &  selected_properties,
const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &  cell = typename parallel::distributed::Triangulation< dim >::active_cell_iterator() 
) const
pure virtual

Perform an interpolation of the properties of the particles in this cell onto a vector of positions in this cell. Implementations of this function must return a vector of a vector of doubles with as many entries as positions in positions. Each entry is a vector with as many entries as there are particle properties in this computation. All in selected_properties selected components will be filled with computed properties, all other components are not filled (or filled with invalid values).

Parameters
[in]particle_handlerReference to the particle handler that allows accessing the particles in the domain.
[in]positionsThe vector of positions where the properties should be evaluated.
[in]selected_propertiesA component mask that determines which particle properties are interpolated in this function.
[in]cellAn optional iterator to the cell containing the particles. Not all callers will know the cell of the particles, but providing the cell when known speeds up the interpolation significantly.
Returns
A vector with as many entries as positions. Every entry is a vector of interpolated particle properties at this position.

Implemented in aspect::Particle::Interpolator::QuadraticLeastSquares< dim >, aspect::Particle::Interpolator::BilinearLeastSquares< dim >, aspect::Particle::Interpolator::NearestNeighbor< dim >, aspect::Particle::Interpolator::CellAverage< dim >, and aspect::Particle::Interpolator::HarmonicAverage< dim >.

§ declare_parameters()

template<int dim>
static void aspect::Particle::Interpolator::Interface< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files. The default implementation of this function does not describe any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.

§ parse_parameters()

template<int dim>
virtual void aspect::Particle::Interpolator::Interface< dim >::parse_parameters ( ParameterHandler &  prm)
virtual

Read the parameters this class declares from the parameter file. The default implementation of this function does not read any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.

Reimplemented in aspect::Particle::Interpolator::QuadraticLeastSquares< dim >, aspect::Particle::Interpolator::BilinearLeastSquares< dim >, aspect::Particle::Interpolator::NearestNeighbor< dim >, aspect::Particle::Interpolator::CellAverage< dim >, and aspect::Particle::Interpolator::HarmonicAverage< dim >.


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