mlpack

mlpack core class documentation

Underlying the implementations of mlpack’s machine learning algorithms are mlpack core support classes, each of which are documented on this page.

πŸ”— Core math utilities

mlpack provides a number of additional mathematical utility classes and functions on top of Armadillo.


πŸ”— Aliases

Aliases are matrix, vector, or cube objects that share memory with another matrix, vector, or cube. They are often used internally inside of mlpack to avoid copies.

Important caveats about aliases:





πŸ”— Range

The Range class represents a simple mathematical range (i.e. [0, 3]), with the bounds represented as doubles.


Constructors


Accessing and modifying range properties


Working with ranges



Example:

mlpack::Range r1(5.0, 6.0); // [5, 6]
mlpack::Range r2(7.0, 8.0); // [7, 8]

mlpack::Range r3 = r1 | r2; // [5, 8]
mlpack::Range r4 = r1 & r2; // empty range

bool b1 = r1.Contains(r2); // false
bool b2 = r1.Contains(5.5); // true
bool b3 = r1.Contains(r3); // true
bool b4 = r3.Contains(r4); // false

// Create a range of `float`s and a range of `int`s.
mlpack::RangeType<float> r5(1.0f, 1.5f); // [1.0, 1.5]
mlpack::RangeType<int> r6(3, 4); // [3, 4]

Range is used by:


πŸ”— ColumnCovariance()

Example:

// Generate a random data matrix with 100 points in 5 dimensions.
arma::mat data(5, 100, arma::fill::randu);

// Compute the covariance matrix of the column-major matrix.
arma::mat cov = mlpack::ColumnCovariance(data);
cov.print("Covariance of random matrix:");

πŸ”— ColumnsToBlocks

The ColumnsToBlocks class provides a way to transform data points (e.g. columns in a matrix) into a block matrix format, primarily useful for visualization as an image.

As a simple example, given a matrix with four columns A, B, C, and D, ColumnsToBlocks can transform this matrix into the form

[[m m m m m]
 [m A m B m]
 [m m m m m]
 [m C m D m]
 [m m m m m]]

where m is a margin, and where each column may itself be reshaped into a block.


Constructors


Properties


πŸ”— Scaling values

ColumnsToBlocks also has the capability of linearly scaling values of the inputs to a given range.

Note: the margin element (ctb.BufValue()) is considered during the scaling process.


Transforming into block format


Examples

Reshape two 4-element vectors into one row of two blocks.

// This matrix has two columns.
arma::mat input;
input = { { -1.0000, 0.1429 },
          { -0.7143, 0.4286 },
          { -0.4286, 0.7143 },
          { -0.1429, 1.0000 } };
input.print("Input columns:");

arma::mat output;
mlpack::ColumnsToBlocks ctb(1, 2);
ctb.Transform(input, output);

// The columns of the input will be reshaped as a square which is
// surrounded by padding value -1 (this value could be changed with the
// BufValue() method):
// -1.0000  -1.0000  -1.0000  -1.0000  -1.0000  -1.0000  -1.0000
// -1.0000  -1.0000  -0.4286  -1.0000   0.1429   0.7143  -1.0000
// -1.0000  -0.7143  -0.1429  -1.0000   0.4286   1.0000  -1.0000
// -1.0000  -1.0000  -1.0000  -1.0000  -1.0000  -1.0000  -1.0000
output.print("Output using 2x2 block size:");

// Now, let's change some parameters; let's have each input column output not
// as a square, but as a 4x1 vector.
ctb.BlockWidth(1);
ctb.BlockHeight(4);
ctb.Transform(input, output);

// The output here will be similar, but each maximal input is 4x1:
// -1.0000 -1.0000 -1.0000 -1.0000 -1.0000
// -1.0000 -1.0000 -1.0000  0.1429 -1.0000
// -1.0000 -0.7143 -1.0000  0.4286 -1.0000
// -1.0000 -0.4286 -1.0000  0.7143 -1.0000
// -1.0000 -0.1429 -1.0000  1.0000 -1.0000
// -1.0000 -1.0000 -1.0000 -1.0000 -1.0000
output.print("Output using 4x1 block size:");

Load simple images and reshape into blocks.

// Load some favicons from websites associated with mlpack.
std::vector<std::string> images;
// See the following files:
// - https://datasets.mlpack.org/images/mlpack-favicon.png
// - https://datasets.mlpack.org/images/ensmallen-favicon.png
// - https://datasets.mlpack.org/images/armadillo-favicon.png 
// - https://datasets.mlpack.org/images/bandicoot-favicon.png
images.push_back("mlpack-favicon.png");
images.push_back("ensmallen-favicon.png");
images.push_back("armadillo-favicon.png");
images.push_back("bandicoot-favicon.png");

mlpack::data::ImageInfo info;
info.Channels() = 1; // Force loading in grayscale.

arma::mat matrix;
mlpack::data::Load(images, matrix, info, true);

// Now `matrix` has 4 columns, each of which is an individual image.
// Let's save that as its own image just for visualization.
mlpack::data::ImageInfo outInfo(matrix.n_cols, matrix.n_rows, 1);
mlpack::data::Save("favicons-matrix.png", matrix, outInfo, true);

// Use ColumnsToBlocks to create a 2x2 block matrix holding each image.
mlpack::ColumnsToBlocks ctb(2, 2);
ctb.BufValue(0.0); // Use 0 for the margin value.
ctb.BufSize(2); // Use 2-pixel margins.

arma::mat blocks;
ctb.Transform(matrix, blocks);

mlpack::data::ImageInfo blockOutInfo(blocks.n_cols, blocks.n_rows, 1);
mlpack::data::Save("favicons-blocks.png", blocks, blockOutInfo, true);

The resulting images (before and after using ColumnsToBlocks) are shown below.

Before:

four favicons each as a column in a matrix, unintelligible

After:

four favicons each as a block in a larger image, much better

See also


πŸ”— Distribution utilities

Example:

const double d1 = mlpack::Digamma(0.25);
const double d2 = mlpack::Digamma(1.0);

const double t1 = mlpack::Trigamma(0.25);
const double t2 = mlpack::Trigamma(1.0);

std::cout << "Digamma(0.25):  " << d1 << "." << std::endl;
std::cout << "Digamma(1.0):   " << d2 << "." << std::endl;
std::cout << "Trigamma(0.25): " << t1 << "." << std::endl;
std::cout << "Trigamma(1.0):  " << t2 << "." << std::endl;

πŸ”— RandVector()

Example:

// Generate a random 10-dimensional vector.
arma::vec v;
v.set_size(10);
RandVector(v);
v.print("Random 10-dimensional vector: ");

std::cout << "Random 10-dimensional vector: " << std::endl;
std::cout << v.t();
std::cout << "L2-norm of vector (should be 1): " << arma::norm(v, 2) << "."
    << std::endl;

πŸ”— Logarithmic utilities

mlpack contains a few functions that are useful for working with logarithms, or vectors containing logarithms.




πŸ”— MultiplyCube2Cube()

Example usage:

// Generate two random cubes.
arma::cube x(10, 100, 5, arma::fill::randu); // 5 matrices, each 10x100.
arma::cube y(12, 100, 5, arma::fill::randu); // 5 matrices, each 12x100.

arma::cube z = mlpack::MultiplyCube2Cube(x, y, false, true);

// Output size should be 10x12x5.
std::cout << "Output size: " << z.n_rows << "x" << z.n_cols << "x" << z.n_slices
    << "." << std::endl;

πŸ”— MultiplyMat2Cube()

Example usage:

// Generate random inputs.
arma::mat  x(10, 100,    arma::fill::randu); // Random 10x100 matrix.
arma::cube y(12, 100, 5, arma::fill::randu); // 5 matrices, each 12x100.

arma::cube z = mlpack::MultiplyMat2Cube(x, y, false, true);

// Output size should be 10x12x5.
std::cout << "Output size: " << z.n_rows << "x" << z.n_cols << "x" << z.n_slices
    << "." << std::endl;

πŸ”— MultiplyCube2Mat()

Example usage:

// Generate two random cubes.
arma::cube x(12, 50, 5, arma::fill::randu); // 5 matrices, each 12x50.
arma::mat  y(12, 60,    arma::fill::randu); // Random 12x60 matrix.

arma::cube z = mlpack::MultiplyCube2Mat(x, y, true, false);

// Output size should be 50x60x5.
std::cout << "Output size: " << z.n_rows << "x" << z.n_cols << "x" << z.n_slices
    << "." << std::endl;

πŸ”— Quantile()

Example usage:

// 70% of points from N(0, 1) are less than q1 = 0.524.
double q1 = mlpack::Quantile(0.7);

// 90% of points from N(0, 1) are less than q2 = 1.282.
double q2 = mlpack::Quantile(0.9);

// 50% of points from N(1, 1) are less than q3 = 1.0.
double q3 = mlpack::Quantile(0.5, 1.0); // Quantile of 1.0 for N(1, 1) is 1.0.

// 10% of points from N(1, 0.1) are less than q4 = 0.871.
double q4 = mlpack::Quantile(0.1, 1.0, 0.1);

std::cout << "Quantile(0.7): " << q1 << "." << std::endl;
std::cout << "Quantile(0.9): " << q2 << "." << std::endl;
std::cout << "Quantile(0.5, 1.0): " << q3 << "." << std::endl;
std::cout << "Quantile(0.1, 1.0, 0.1): " << q4 << "." << std::endl;

πŸ”— RNG and random number utilities

On top of the random number generation support that Armadillo provides via randu(), randn(), and randi(), mlpack provides a few additional thread-safe random number generation functions for generating random scalar values.

Examples:

mlpack::RandomSeed(123); // Set a specific random seed.

const double r1 = mlpack::Random();             // In the range [0, 1).
const double r2 = mlpack::Random(3, 4);         // In the range [3, 4).
const double r3 = mlpack::RandBernoulli(0.25);  // P(1) = 0.25.
const int    r4 = mlpack::RandInt(10);          // In the range [0, 10).
const int    r5 = mlpack::RandInt(5, 10);       // In the range [5, 10).
const double r6 = mlpack::RandNormal();         // r6 ~ N(0, 1).
const double r7 = mlpack::RandNormal(2.0, 3.0); // r7 ~ N(2, 3).

std::cout << "Random():            " << r1 << "." << std::endl;
std::cout << "Random(3, 4):        " << r2 << "." << std::endl;
std::cout << "RandBernoulli(0.25): " << r3 << "." << std::endl;
std::cout << "RandInt(10):         " << r4 << "." << std::endl;
std::cout << "RandInt(5, 10):      " << r5 << "." << std::endl;
std::cout << "RandNormal():        " << r6 << "." << std::endl;
std::cout << "RandNormal(2, 3):    " << r7 << "." << std::endl;

πŸ”— RandomBasis()

The RandomBasis() function generates a random d-dimensional orthogonal basis.

Example:

arma::mat basis;

// Generate a 10-dimensional random basis.
mlpack::RandomBasis(basis, 10);

// Each two vectors are orthogonal.
std::cout << "Dot product of basis vectors 2 and 4: "
    << arma::dot(basis.col(2), basis.col(4))
    << " (should be zero or very close!)." << std::endl;

πŸ”— ShuffleData()

Shuffle a column-major dataset and associated labels/responses, optionally with weights. This preserves the connection of each data point to its label (and optionally its weight).

Note: when inputData is a cube (e.g. arma::cube or similar), the columns of the cube will be shuffled.

Example usage:

// See https://datasets.mlpack.org/iris.csv.
arma::mat dataset;
mlpack::data::Load("iris.csv", dataset, true);
// See https://datasets.mlpack.org/iris.labels.csv.
arma::Row<size_t> labels;
mlpack::data::Load("iris.labels.csv", labels, true);

// Now shuffle the points in the iris dataset.
arma::mat shuffledDataset;
arma::Row<size_t> shuffledLabels;
mlpack::ShuffleData(dataset, labels, shuffledDataset, shuffledLabels);

std::cout << "Before shuffling, the first point was: " << std::endl;
std::cout << "  " << dataset.col(0).t();
std::cout << "with label " << labels[0] << "." << std::endl;
std::cout << std::endl;
std::cout << "After shuffling, the first point is: " << std::endl;
std::cout << "  " << shuffledDataset.col(0).t();
std::cout << "with label " << shuffledLabels[0] << "." << std::endl;

// Generate random weights, then shuffle those also.
arma::rowvec weights(dataset.n_cols, arma::fill::randu);
arma::rowvec shuffledWeights;
mlpack::ShuffleData(dataset, labels, weights, shuffledDataset, shuffledLabels,
    shuffledWeights);

std::cout << std::endl << std::endl;
std::cout << "Before shuffling with weights, the first point was: "
    << std::endl;
std::cout << "  " << dataset.col(0).t();
std::cout << "with label " << labels[0] << " and weight " << weights[0] << "."
    << std::endl;
std::cout << std::endl;
std::cout << "After shuffling with weights, the first point is: " << std::endl;
std::cout << "  " << shuffledDataset.col(0).t();
std::cout << "with label " << shuffledLabels[0] << " and weight "
    << shuffledWeights[0] << "." << std::endl;

πŸ”— Distributions

mlpack has support for a number of different distributions, each supporting the same API. These can be used with, for instance, the HMM class.

πŸ”— DiscreteDistribution

DiscreteDistribution represents a multidimensional categorical distribution (or generalized Bernoulli distribution) where integer-valued vectors (e.g. [0, 3, 4]) are associated with specific probabilities in each dimension.

Example: a 3-dimensional DiscreteDistribution will have a specific probability value associated with each integer value in each dimension. So, for the vector [0, 3, 4], P(0) in dimension 0 could be, e.g., 0.3, P(3) in dimension 1 could be, e.g., 0.4, and P(4) in dimension 2 could be, e.g., 0.6. Then, P([0, 3, 4]) would be 0.3 * 0.4 * 0.6 = 0.072.


Constructors


Access and modify properties of distribution


Compute probabilities of points


Sample from the distribution


Fit the distribution to observations


Example usage:

// Create a single-dimension Bernoulli distribution: P([0]) = 0.3, P([1]) = 0.7.
mlpack::DiscreteDistribution bernoulli(2);
bernoulli.Probabilities(0)[0] = 0.3;
bernoulli.Probabilities(0)[1] = 0.7;

const double p1 = bernoulli.Probability(arma::vec("0")); // p1 = 0.3.
const double p2 = bernoulli.Probability(arma::vec("1")); // p2 = 0.7.

// Create a 3-dimensional discrete distribution by specifying the probabilities
// manually.
arma::vec probDim0 = arma::vec("0.1 0.3 0.5 0.1"); // 4 possible values.
arma::vec probDim1 = arma::vec("0.7 0.3");         // 2 possible values.
arma::vec probDim2 = arma::vec("0.4 0.4 0.2");     // 3 possible values.
std::vector<arma::vec> probs { probDim0, probDim1, probDim2 };
mlpack::DiscreteDistribution d(probs);

arma::vec obs("2 0 1");
const double p3 = d.Probability(obs); // p3 = 0.5 * 0.7 * 0.4 = 0.14.

// Estimate a 10-dimensional discrete distribution.
// Each dimension takes values between 0 and 9.
arma::mat observations = arma::randi<arma::mat>(10, 1000,
    arma::distr_param(0, 9));

// Create a distribution with 10 observations in each of the 10 dimensions.
mlpack::DiscreteDistribution d2(
    arma::Col<size_t>("10 10 10 10 10 10 10 10 10 10"));
d2.Train(observations);

// Compute the probabilities of each point.
arma::vec probabilities;
d2.Probability(observations, probabilities);
std::cout << "Average probability: " << arma::mean(probabilities) << "."
    << std::endl;

πŸ”— GaussianDistribution

GaussianDistribution is a standard multivariate Gaussian distribution with parameterized mean and covariance.


Constructors


Access and modify properties of distribution


Compute probabilities of points


Sample from the distribution


Fit the distribution to observations


Example usage:

// Create a Gaussian distribution in 3 dimensions with zero mean and unit
// covariance.
mlpack::GaussianDistribution g(3);

// Compute the probability of the point [0, 0.5, 0.25].
const double p = g.Probability(arma::vec("0 0.5 0.25"));

// Modify the mean in dimension 0.
g.Mean()[0] = 0.5;

// Set a random covariance.
arma::mat newCov(3, 3, arma::fill::randu);
newCov *= newCov.t(); // Ensure covariance is positive semidefinite.
g.Covariance(std::move(newCov)); // Set new covariance.

// Compute the probability of the same point [0, 0.5, 0.25].
const double p2 = g.Probability(arma::vec("0 0.5 0.25"));

// Create a Gaussian distribution that is estimated from random samples in 50
// dimensions.
arma::mat samples(50, 10000, arma::fill::randn); // Normally distributed.

mlpack::GaussianDistribution g2(50);
g2.Train(samples);

// Compute the probability of all of the samples.
arma::vec probabilities;
g2.Probability(samples, probabilities);

std::cout << "Average probability is: " << arma::mean(probabilities) << "."
    << std::endl;

πŸ”— Metrics

mlpack includes a number of distance metrics for its distance-based techniques. These all implement the same API, providing one Evaluate() method, and can be used with a variety of different techniques, including:

Supported metrics:

πŸ”— LMetric

The LMetric template class implements a generalized L-metric (L1-metric, L2-metric, etc.). The class has two template parameters:

LMetric<Power, TakeRoot>

Several convenient typedefs are available:


The static Evaluate() method can be used to compute the distance between two vectors.

Note: The vectors given to Evaluate() can have any type so long as the type implements the Armadillo API (e.g. arma::fvec, arma::sp_fvec, etc.).


Example usage:

// Create two vectors: [0, 1.0, 5.0] and [1.0, 3.0, 5.0].
arma::vec a("0.0 1.0 5.0");
arma::vec b("1.0 3.0 5.0");

const double d1 = mlpack::ManhattanDistance::Evaluate(a, b);        // d1 = 3.0
const double d2 = mlpack::EuclideanDistance::Evaluate(a, b);        // d2 = 2.24
const double d3 = mlpack::SquaredEuclideanDistance::Evaluate(a, b); // d3 = 5.0
const double d4 = mlpack::ChebyshevDistance::Evaluate(a, b);        // d4 = 2.0
const double d5 = mlpack::LMetric<4>::Evaluate(a, b);               // d5 = 2.03
const double d6 = mlpack::LMetric<3, false>::Evaluate(a, b);        // d6 = 9.0

std::cout << "Manhattan distance:         " << d1 << "." << std::endl;
std::cout << "Euclidean distance:         " << d2 << "." << std::endl;
std::cout << "Squared Euclidean distance: " << d3 << "." << std::endl;
std::cout << "Chebyshev distance:         " << d4 << "." << std::endl;
std::cout << "L4-distance:                " << d5 << "." << std::endl;
std::cout << "Cubed L3-distance:          " << d6 << "." << std::endl;

// Compute the distance between two random 10-dimensional vectors in a matrix.
arma::mat m(10, 100, arma::fill::randu);

const double d7 = mlpack::EuclideanDistance::Evaluate(m.col(0), m.col(7));

std::cout << std::endl;
std::cout << "Distance between two random vectors: " << d7 << "." << std::endl;
std::cout << std::endl;

// Compute the distance between two 32-bit precision `float` vectors.
arma::fvec fa("0.0 1.0 5.0");
arma::fvec fb("1.0 3.0 5.0");

const double d8 = mlpack::EuclideanDistance::Evaluate(fa, fb); // d8 = 2.236

std::cout << "Euclidean distance (fvec): " << d8 << "." << std::endl;

πŸ”— Kernels

mlpack includes a number of Mercer kernels for its kernel-based techniques. These all implement the same API, providing one Evaluate() method, and can be used with a variety of different techniques, including:

Supported kernels:

πŸ”— GaussianKernel

The GaussianKernel class implements the standard Gaussian kernel (also called the radial basis function kernel or RBF kernel).

The Gaussian kernel is defined as: k(x1, x2) = exp(-|| x1 - x2 ||^2 / (2 * bw^2)) where bw is the bandwidth parameter of the kernel.


Constructors and properties


Kernel evaluation


Other utilities


Example usage:

// Create a Gaussian kernel with default bandwidth.
mlpack::GaussianKernel g;

// Create a Gaussian kernel with bandwidth 5.0.
mlpack::GaussianKernel g2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = g.Evaluate(x1, x2);
const double k2 = g2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Evaluate the kernel value when the distance between two points is already
// computed.
const double distance = 1.5;
const double k3 = g.Evaluate(distance);

// Change the bandwidth of the kernel to 2.5.
g.Bandwidth(2.5);
const double k4 = g.Evaluate(x1, x2);
std::cout << "Kernel value with bw=2.5: " << k4 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = g.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k5 = g.Evaluate(fx1, fx2);
const double k6 = g2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k5
    << " (bw=2.5), " << k6 << " (bw=5.0)." << std::endl;

πŸ”— CauchyKernel

The CauchyKernel class implements the Cauchy kernel, a kernel function with a longer tail than the Gaussian kernel, defined as: k(x1, x2) = 1 / (1 + (|| x1 - x2 ||^2 / bw^2)) where bw is the bandwidth parameter of the kernel.


Constructors and properties


Kernel evaluation


Example usage:

// Create a Cauchy kernel with default bandwidth.
mlpack::CauchyKernel c;

// Create a Cauchy kernel with bandwidth 5.0.
mlpack::CauchyKernel c2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = c.Evaluate(x1, x2);
const double k2 = c2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Change the bandwidth of the kernel to 2.5.
c.Bandwidth(2.5);
const double k3 = c.Evaluate(x1, x2);
std::cout << "Kernel value with bw=2.5: " << k3 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = c.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k4 = c.Evaluate(fx1, fx2);
const double k5 = c2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k4
    << " (bw=2.5), " << k5 << " (bw=5.0)." << std::endl;

πŸ”— CosineSimilarity

The CosineSimilarity class implements the dot-product cosine similarity, defined as: k(x1, x2) = (x1^T x2) / (|| x1 || * || x2 ||). The value of the kernel is limited to the range [-1, 1]. The cosine similarity is often used in text mining tasks.


Constructor

Note: because the CosineSimilarity kernel has no parameters, it is not necessary to create an object and the Evaluate() function (below) can be called statically.


Kernel evaluation


Example usage:

// Create a cosine similarity kernel.
mlpack::CosineSimilarity c;

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = c.Evaluate(x1, x2);
const double k2 = c.Evaluate(x1, x1);
const double k3 = c.Evaluate(x2, x2);
std::cout << "Cosine similarity values:" << std::endl;
std::cout << "  - k(x1, x2): " << k1 << "." << std::endl;
std::cout << "  - k(x1, x1): " << k2 << "." << std::endl;
std::cout << "  - k(x2, x2): " << k3 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix,
// using the static Evaluate() function.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = mlpack::CosineSimilarity::Evaluate(x1, r.col(i));
std::cout << "Average cosine similarity for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the cosine similarity between two sparse 32-bit floating point
// vectors.
arma::sp_fvec x3, x4;
x3.sprandu(100, 1, 0.2);
x4.sprandu(100, 1, 0.2);
const double k4 = mlpack::CosineSimilarity::Evaluate(x3, x4);
std::cout << "Cosine similarity between two random sparse 32-bit floating "
    << "point vectors: " << k4 << "." << std::endl;

πŸ”— EpanechnikovKernel

The EpanechnikovKernel implements the parabolic or Epanechnikov kernel, defined as the following function: k(x1, x2) = max(0, (3 / 4) * (1 - (|| x1 - x2 ||_2 / bw)^2)), where bw is the bandwidth parameter of the kernel.

The kernel takes the value 0 when || x1 - x2 ||_2 (the Euclidean distance between x1 and x2) is greater than or equal to bw.


Constructors and properties


Kernel evaluation


Other utilities


Example usage:

// Create an Epanechnikov kernel with default bandwidth.
mlpack::EpanechnikovKernel e;

// Create an Epanechnikov kernel with bandwidth 5.0.
mlpack::EpanechnikovKernel e2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = e.Evaluate(x1, x2);
const double k2 = e2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Evaluate the kernel value when the distance between two points is already
// computed.
const double distance = 1.5;
const double k3 = e.Evaluate(distance);

// Change the bandwidth of the kernel to 2.5.
e.Bandwidth(2.5);
const double k4 = e.Evaluate(x1, x2);
std::cout << "Kernel value with bw=2.5: " << k4 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = e.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k5 = e.Evaluate(fx1, fx2);
const double k6 = e2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k5
    << " (bw=2.5), " << k6 << " (bw=5.0)." << std::endl;

πŸ”— HyperbolicTangentKernel

The HyperbolicTangentKernel implements the hyperbolic tangent kernel, which is defined by the following equation: f(x1, x2) = tanh(s * (x1^T x2) + t) where s is the scale parameter and t is the offset parameter.

The hyperbolic tangent kernel is not a positive definite Mercer kernel and thus does not satisfy the theoretical requirements of many kernel methods. See this discussion for more details. In practice, for many kernel methods, it may still provide compelling results despite this theoretical limitation.


Constructors and properties


Kernel evaluation


Example usage:

// Create a hyperbolic tangent kernel with default scale and offset.
mlpack::HyperbolicTangentKernel h;

// Create a hyperbolic tangent kernel with scale 2.0 and offset 1.0.
mlpack::HyperbolicTangentKernel h2(2.0, 1.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = h.Evaluate(x1, x2);
const double k2 = h2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (s=1.0, t=0.0), " << k2
    << " (s=2.0, t=1.0)." << std::endl;

// Change the scale and offset of the kernel.
h.Scale(2.5);
h.Offset(-1.0);
const double k3 = h.Evaluate(x1, x2);
std::cout << "Kernel value with s=2.5, t=-1.0: " << k3 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = h.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k4 = h.Evaluate(fx1, fx2);
const double k5 = h2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k4
    << " (s=2.5, t=-1.0), " << k5 << " (s=2.0, t=1.0)." << std::endl;

πŸ”— LaplacianKernel

The LaplacianKernel class implements the Laplacian kernel, also known as the exponential kernel, defined by the following equation: k(x1, x2) = exp(-|| x1 - x2 || / bw) where bw is the bandwidth parameter.


Constructors and properties


Kernel evaluation


Other utilities


Example usage:

// Create a Laplacian kernel with default bandwidth.
mlpack::LaplacianKernel l;

// Create a Laplacian kernel with bandwidth 5.0.
mlpack::LaplacianKernel l2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = l.Evaluate(x1, x2);
const double k2 = l2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Evaluate the kernel value when the distance between two points is already
// computed.
const double distance = 1.5;
const double k3 = l.Evaluate(distance);

// Change the bandwidth of the kernel to 2.5.
l.Bandwidth(2.5);
const double k4 = l.Evaluate(x1, x2);
std::cout << "Kernel value with bw=2.5: " << k4 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = l.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k5 = l.Evaluate(fx1, fx2);
const double k6 = l2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k5
    << " (bw=2.5), " << k6 << " (bw=5.0)." << std::endl;

πŸ”— LinearKernel

The LinearKernel class implements the simple linear dot product kernel, defined by the following equation: k(x1, x2) = x1^T x2.

The use of the linear kernel for kernel methods generally results in the non-kernelized version of the algorithm; for instance, a kernel support vector machine using the linear kernel amounts to a linear SVM.


Constructor

Note: because the LinearKernel kernel has no parameters, it is not necessary to create an object and the Evaluate() function (below) can be called statically.


Kernel evaluation


Example usage:

// Create a linear kernel.
mlpack::LinearKernel l;

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = l.Evaluate(x1, x2); // Identical to arma::dot(x1, x2).
const double k2 = l.Evaluate(x1, x1);
const double k3 = l.Evaluate(x2, x2);
std::cout << "Linear kernel values:" << std::endl;
std::cout << "  - k(x1, x2): " << k1 << "." << std::endl;
std::cout << "  - k(x1, x1): " << k2 << "." << std::endl;
std::cout << "  - k(x2, x2): " << k3 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix,
// using the static Evaluate() function.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = mlpack::LinearKernel::Evaluate(x1, r.col(i));
std::cout << "Average linear kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the cosine similarity between two sparse 32-bit floating point
// vectors.
arma::sp_fvec x3, x4;
x3.sprandu(100, 1, 0.2);
x4.sprandu(100, 1, 0.2);
const double k4 = mlpack::LinearKernel::Evaluate(x3, x4);
std::cout << "Linear kernel value between two random sparse 32-bit floating "
    << "point vectors: " << k4 << "." << std::endl;

πŸ”— PolynomialKernel

The PolynomialKernel class implements the standard polynomial kernel, which is defined by the following equation: k(x1, x2) = (x1^T x2 + t)^d where d is the degree of the polynomial and t is the offset.

The use of the polynomial kernel has a similar effect to the use of polynomial (interaction) features in standard machine learning methods.


Constructors and properties


Kernel evaluation


Example usage:

// Create a polynomial kernel with default degree (2) and offset (0).
mlpack::PolynomialKernel p;

// Create a polynomial kernel with degree 3.0 and offset -1.0.
mlpack::PolynomialKernel p2(3.0, -1.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = p.Evaluate(x1, x2);
const double k2 = p2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Change the degree of the kernel to 2.5 and the offset to 1.0.
p.Degree(2.5);
p.Offset(1.0);
const double k3 = p.Evaluate(x1, x2);
std::cout << "Kernel value with d=2.5, t=1.0: " << k3 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = p.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k4 = p.Evaluate(fx1, fx2);
const double k5 = p2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k4
    << " (d=2.5, t=1.0), " << k5 << " (d=3.0, t=-1.0)." << std::endl;

πŸ”— PSpectrumStringKernel

The PSpectrumStringKernel class implements the length-p string spectrum kernel, proposed by Leslie, Eskin, and Noble (pdf). The kernel finds the contiguous subsequence match count between two strings.

Due to mlpack’s use of Armadillo, which requires that all matrix data be numeric, this class operates by internally storing all strings, and passing in numeric vectors such as [0 1] that reference string index 1 in dataset index 0. In turn, this means that the data points given to the PSpectrumStringKernel are simply IDs and have no geometric meaning.


Constructors and properties


Kernel evaluation


Example usage:

// Create two example datasets:
//      ["hello", "goodbye", "package"],
//      ["mlpack", "is", "really", "great"]
std::vector<std::vector<std::string>> datasets;
datasets.push_back({ "hello", "goodbye", "package" });
datasets.push_back({ "mlpack", "is", "really", "great" });

// Create a p-spectrum string kernel with a substring length of 2,
// and another with a substring length of 3.
mlpack::PSpectrumStringKernel p(datasets, 2);
mlpack::PSpectrumStringKernel p2(datasets, 3);

// Evaluate the kernel value between "mlpack" and "package".
arma::uvec x1("1 0"); // "mlpack": dataset 1, string 0.
arma::uvec x2("0 2"); // "package": dataset 0, string 2.
const double k1 = p.Evaluate(x1, x2);
const double k2 = p2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (p=2), " << k2 << " (p=3)."
    << std::endl;

πŸ”— SphericalKernel

The SphericalKernel class implements the simple spherical kernel, also known as the uniform kernel, or rectangular window kernel. The value of the SphericalKernel is 1 when the Euclidean distance between two points x1 and x2 is less than the bandwidth bw, and 0 otherwise: k(x1, x2) = 1(|| x1 - x2 || <= bw).


Constructors and properties


Kernel evaluation


Other utilities


Example usage:

// Create a spherical kernel with default bandwidth.
mlpack::SphericalKernel s;

// Create a spherical kernel with bandwidth 5.0.
mlpack::SphericalKernel s2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 2.5");
arma::vec x2("2.5 1.0 0.5");
const double k1 = s.Evaluate(x1, x2);
const double k2 = s2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Evaluate the kernel value when the distance between two points is already
// computed.
const double distance = 0.9;
const double k3 = s.Evaluate(distance);

// Change the bandwidth of the kernel to 3.0.
s.Bandwidth(3.0);
const double k4 = s.Evaluate(x1, x2);
std::cout << "Kernel value with bw=3.0: " << k4 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix, using
// a kernel bandwidth of 2.5.
s.Bandwidth(2.5);
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = s.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 2.5");
arma::fvec fx2("2.5 1.0 0.5");
const double k5 = s.Evaluate(fx1, fx2);
const double k6 = s2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k5
    << " (bw=2.5), " << k6 << " (bw=5.0)." << std::endl;

πŸ”— TriangularKernel

The TriangularKernel class implements the simple triangular kernel, defined by the following equation: k(x1, x2) = max(0, 1 - || x1 - x2 || / bw) where bw is the bandwidth of the kernel.


Constructors and properties


Kernel evaluation


Other utilities


Example usage:

// Create a triangular kernel with default bandwidth.
mlpack::TriangularKernel t;

// Create a triangular kernel with bandwidth 5.0.
mlpack::TriangularKernel t2(5.0);

// Evaluate the kernel value between two 3-dimensional points.
arma::vec x1("0.5 1.0 1.5");
arma::vec x2("1.5 1.0 0.5");
const double k1 = t.Evaluate(x1, x2);
const double k2 = t2.Evaluate(x1, x2);
std::cout << "Kernel values: " << k1 << " (bw=1.0), " << k2 << " (bw=5.0)."
    << std::endl;

// Evaluate the kernel value when the distance between two points is already
// computed.
const double distance = 0.75;
const double k3 = t.Evaluate(distance);

// Change the bandwidth of the kernel to 2.5.
t.Bandwidth(2.5);
const double k4 = t.Evaluate(x1, x2);
std::cout << "Kernel value with bw=2.5: " << k4 << "." << std::endl;

// Evaluate the kernel value between x1 and all points in a random matrix.
arma::mat r(3, 100, arma::fill::randu);
arma::vec kernelValues(100);
for (size_t i = 0; i < r.n_cols; ++i)
  kernelValues[i] = t.Evaluate(x1, r.col(i));
std::cout << "Average kernel value for random points: "
    << arma::mean(kernelValues) << "." << std::endl;

// Compute the kernel value between two 32-bit floating-point vectors.
arma::fvec fx1("0.5 1.0 1.5");
arma::fvec fx2("1.5 1.0 0.5");
const double k5 = t.Evaluate(fx1, fx2);
const double k6 = t2.Evaluate(fx1, fx2);
std::cout << "Kernel values between two floating-point vectors: " << k5
    << " (bw=2.5), " << k6 << " (bw=5.0)." << std::endl;

πŸ”— Implement a custom kernel

mlpack supports custom kernels, so long as they implement an appropriate Evaluate() function.

See The KernelType Policy in mlpack for more information.