Public Types

enum  KernelTypes
{
  GAUSSIAN_KERNEL
,
  EPANECHNIKOV_KERNEL
,
  LAPLACIAN_KERNEL
,
  SPHERICAL_KERNEL
,
  TRIANGULAR_KERNEL

}
 
enum  TreeTypes
{
  KD_TREE
,
  BALL_TREE
,
  COVER_TREE
,
  OCTREE
,
  R_TREE

}
 

Public Member Functions

 KDEModel (const double bandwidth=1.0, const double relError=KDEDefaultParams::relError, const double absError=KDEDefaultParams::absError, const KernelTypes kernelType=KernelTypes::GAUSSIAN_KERNEL, const TreeTypes treeType=TreeTypes::KD_TREE, const bool monteCarlo=KDEDefaultParams::mode, const double mcProb=KDEDefaultParams::mcProb, const size_t initialSampleSize=KDEDefaultParams::initialSampleSize, const double mcEntryCoef=KDEDefaultParams::mcEntryCoef, const double mcBreakCoef=KDEDefaultParams::mcBreakCoef)
 Initialize KDEModel. More...

 
 KDEModel (const KDEModel &other)
 Copy constructor of the given model. More...

 
 KDEModel (KDEModel &&other)
 Move constructor of the given model. Takes ownership of the model. More...

 
 ~KDEModel ()
 Destroy the KDEModel object. More...

 
double AbsoluteError () const
 Get the absolute error tolerance. More...

 
void AbsoluteError (const double newAbsError)
 Modify the absolute error tolerance. More...

 
double Bandwidth () const
 Get the bandwidth of the kernel. More...

 
void Bandwidth (const double newBandwidth)
 Modify the bandwidth of the kernel. More...

 
void BuildModel (arma::mat &&referenceSet)
 Build the KDE model with the given parameters and then trains it with the given reference data. More...

 
void Evaluate (arma::mat &&querySet, arma::vec &estimations)
 Perform kernel density estimation on the given query set. More...

 
void Evaluate (arma::vec &estimations)
 Perform kernel density estimation on the reference set. More...

 
KernelTypes KernelType () const
 Get the kernel type of the model. More...

 
KernelTypesKernelType ()
 Modify the kernel type of the model. More...

 
double MCBreakCoefficient () const
 Get Monte Carlo break coefficient. More...

 
void MCBreakCoefficient (const double newBreakCoef)
 Modify Monte Carlo break coefficient. More...

 
double MCEntryCoefficient () const
 Get Monte Carlo entry coefficient. More...

 
void MCEntryCoefficient (const double newEntryCoef)
 Modify Monte Carlo entry coefficient. More...

 
size_t MCInitialSampleSize () const
 Get the initial sample size for Monte Carlo estimations. More...

 
void MCInitialSampleSize (const size_t newSampleSize)
 Modify the initial sample size for Monte Carlo estimations. More...

 
double MCProbability () const
 Get Monte Carlo probability of error being bounded by relative error. More...

 
void MCProbability (const double newMCProb)
 Modify Monte Carlo probability of error being bounded by relative error. More...

 
KDEMode Mode () const
 Get the mode of the model. More...

 
KDEModeMode ()
 Modify the mode of the model. More...

 
bool MonteCarlo () const
 Get whether the model is using Monte Carlo estimations or not. More...

 
void MonteCarlo (const bool newMonteCarlo)
 Modify whether the model is using Monte Carlo estimations or not. More...

 
KDEModeloperator= (KDEModel other)
 Copy the given model. More...

 
double RelativeError () const
 Get the relative error tolerance. More...

 
void RelativeError (const double newRelError)
 Modify the relative error tolerance. More...

 
template
<
typename
Archive
>
void serialize (Archive &ar, const unsigned int version)
 Serialize the KDE model. More...

 
TreeTypes TreeType () const
 Get the tree type of the model. More...

 
TreeTypesTreeType ()
 Modify the tree type of the model. More...

 

Detailed Description

Definition at line 362 of file kde_model.hpp.

Member Enumeration Documentation

◆ KernelTypes

Enumerator
GAUSSIAN_KERNEL 
EPANECHNIKOV_KERNEL 
LAPLACIAN_KERNEL 
SPHERICAL_KERNEL 
TRIANGULAR_KERNEL 

Definition at line 374 of file kde_model.hpp.

◆ TreeTypes

enum TreeTypes
Enumerator
KD_TREE 
BALL_TREE 
COVER_TREE 
OCTREE 
R_TREE 

Definition at line 365 of file kde_model.hpp.

Constructor & Destructor Documentation

◆ KDEModel() [1/3]

KDEModel ( const double  bandwidth = 1.0,
const double  relError = KDEDefaultParams::relError,
const double  absError = KDEDefaultParams::absError,
const KernelTypes  kernelType = KernelTypes::GAUSSIAN_KERNEL,
const TreeTypes  treeType = TreeTypes::KD_TREE,
const bool  monteCarlo = KDEDefaultParams::mode,
const double  mcProb = KDEDefaultParams::mcProb,
const size_t  initialSampleSize = KDEDefaultParams::initialSampleSize,
const double  mcEntryCoef = KDEDefaultParams::mcEntryCoef,
const double  mcBreakCoef = KDEDefaultParams::mcBreakCoef 
)

Initialize KDEModel.

Parameters
bandwidthBandwidth to use for the kernel.
relErrorMaximum relative error tolerance for each point in the model. For example, 0.05 means that each value must be within 5% of the true KDE value.
absErrorMaximum absolute error tolerance for each point in the model. For example, 0.1 means that for each point the value can have a maximum error of 0.1 units.
kernelTypeType of kernel to use.
treeTypeType of tree to use.
monteCarloWhether to use Monte Carlo estimations when possible.
mcProbof a Monte Carlo estimation to be bounded by relative error tolerance.
initialSampleSizeInitial sample size for Monte Carlo estimations.
mcEntryCoefCoefficient to control how much larger does the amount of node descendants has to be compared to the initial sample size in order for it to be a candidate for Monte Carlo estimations.
mcBreakCoefCoefficient to control what fraction of the node's descendants evaluated is the limit before Monte Carlo estimation recurses.

◆ KDEModel() [2/3]

KDEModel ( const KDEModel other)

Copy constructor of the given model.

◆ KDEModel() [3/3]

KDEModel ( KDEModel &&  other)

Move constructor of the given model. Takes ownership of the model.

◆ ~KDEModel()

~KDEModel ( )

Destroy the KDEModel object.

Member Function Documentation

◆ AbsoluteError() [1/2]

double AbsoluteError ( ) const
inline

Get the absolute error tolerance.

Definition at line 516 of file kde_model.hpp.

◆ AbsoluteError() [2/2]

void AbsoluteError ( const double  newAbsError)

Modify the absolute error tolerance.

◆ Bandwidth() [1/2]

double Bandwidth ( ) const
inline

Get the bandwidth of the kernel.

Definition at line 504 of file kde_model.hpp.

◆ Bandwidth() [2/2]

void Bandwidth ( const double  newBandwidth)

Modify the bandwidth of the kernel.

◆ BuildModel()

void BuildModel ( arma::mat &&  referenceSet)

Build the KDE model with the given parameters and then trains it with the given reference data.

Takes possession of the reference set to avoid a copy, so the reference set will not be usable after this.

Parameters
referenceSetSet of reference points.

◆ Evaluate() [1/2]

void Evaluate ( arma::mat &&  querySet,
arma::vec &  estimations 
)

Perform kernel density estimation on the given query set.

Takes possession of the query set to avoid a copy, so the query set will not be usable after this. If possible, it returns normalized estimations.

Precondition
The model has to be previously created with BuildModel.
Parameters
querySetSet of query points.
estimationsVector where the results will be stored in the same order as the query points.

◆ Evaluate() [2/2]

void Evaluate ( arma::vec &  estimations)

Perform kernel density estimation on the reference set.

If possible, it returns normalized estimations.

Precondition
The model has to be previously created with BuildModel.
Parameters
estimationsVector where the results will be stored in the same order as the query points.

◆ KernelType() [1/2]

KernelTypes KernelType ( ) const
inline

Get the kernel type of the model.

Definition at line 528 of file kde_model.hpp.

◆ KernelType() [2/2]

KernelTypes& KernelType ( )
inline

Modify the kernel type of the model.

Definition at line 531 of file kde_model.hpp.

◆ MCBreakCoefficient() [1/2]

double MCBreakCoefficient ( ) const
inline

Get Monte Carlo break coefficient.

Definition at line 558 of file kde_model.hpp.

References BOOST_TEMPLATE_CLASS_VERSION(), and mlpack::bindings::tests::CleanMemory().

◆ MCBreakCoefficient() [2/2]

void MCBreakCoefficient ( const double  newBreakCoef)

Modify Monte Carlo break coefficient.

◆ MCEntryCoefficient() [1/2]

double MCEntryCoefficient ( ) const
inline

Get Monte Carlo entry coefficient.

Definition at line 552 of file kde_model.hpp.

◆ MCEntryCoefficient() [2/2]

void MCEntryCoefficient ( const double  newEntryCoef)

Modify Monte Carlo entry coefficient.

◆ MCInitialSampleSize() [1/2]

size_t MCInitialSampleSize ( ) const
inline

Get the initial sample size for Monte Carlo estimations.

Definition at line 546 of file kde_model.hpp.

◆ MCInitialSampleSize() [2/2]

void MCInitialSampleSize ( const size_t  newSampleSize)

Modify the initial sample size for Monte Carlo estimations.

◆ MCProbability() [1/2]

double MCProbability ( ) const
inline

Get Monte Carlo probability of error being bounded by relative error.

Definition at line 540 of file kde_model.hpp.

◆ MCProbability() [2/2]

void MCProbability ( const double  newMCProb)

Modify Monte Carlo probability of error being bounded by relative error.

◆ Mode() [1/2]

KDEMode Mode ( ) const

Get the mode of the model.

◆ Mode() [2/2]

KDEMode& Mode ( )

Modify the mode of the model.

◆ MonteCarlo() [1/2]

bool MonteCarlo ( ) const
inline

Get whether the model is using Monte Carlo estimations or not.

Definition at line 534 of file kde_model.hpp.

◆ MonteCarlo() [2/2]

void MonteCarlo ( const bool  newMonteCarlo)

Modify whether the model is using Monte Carlo estimations or not.

◆ operator=()

KDEModel& operator= ( KDEModel  other)

Copy the given model.

Use std::move if the object to copy is no longer needed.

Parameters
otherKDEModel to copy.

◆ RelativeError() [1/2]

double RelativeError ( ) const
inline

Get the relative error tolerance.

Definition at line 510 of file kde_model.hpp.

◆ RelativeError() [2/2]

void RelativeError ( const double  newRelError)

Modify the relative error tolerance.

◆ serialize()

void serialize ( Archive &  ar,
const unsigned int  version 
)

Serialize the KDE model.

◆ TreeType() [1/2]

TreeTypes TreeType ( ) const
inline

Get the tree type of the model.

Definition at line 522 of file kde_model.hpp.

◆ TreeType() [2/2]

TreeTypes& TreeType ( )
inline

Modify the tree type of the model.

Definition at line 525 of file kde_model.hpp.


The documentation for this class was generated from the following file:
  • /home/ryan/src/mlpack.org/_src/mlpack-3.3.2/src/mlpack/methods/kde/kde_model.hpp