distributions.h
Classes, functions, and variables in this file
doublelog_rising_factorial (double base, unsigned int exponent)
structsymmetric_dirichlet
boolinit (symmetric_dirichlet< V > & distribution, const symmetric_dirichlet< V > & src)
boolinit (symmetric_dirichlet< V > & distribution, const dirichlet< V > & src)
boolread (symmetric_dirichlet< V > & distribution, FILE * in)
boolwrite (const symmetric_dirichlet< V > & distribution, FILE * out)
boolsample (const symmetric_dirichlet< V > & distribution, Destination & dst)
structdirichlet
boolinit (dirichlet< V > & distribution, const symmetric_dirichlet< V > & src)
boolinit (dirichlet< V > & distribution, const dirichlet< V > & src)
boolread (dirichlet< V > & distribution, FILE * in)
boolwrite (const dirichlet< V > & distribution, FILE * out)
boolsample (const dirichlet< V > & distribution, Destination & dst)
structis_dirichlet
structsparse_categorical
boolinit (sparse_categorical< K, V > & distribution, const sparse_categorical< K, V > & src)
boolsample (const sparse_categorical< K, V > & distribution, K & output)
structdense_categorical
boolinit (dense_categorical< V > & distribution, unsigned int atom_count)
boolinit (dense_categorical< V > & distribution, const dense_categorical< V > & src)
boolread (dense_categorical< V > & distribution, Stream & stream)
boolwrite (const dense_categorical< V > & distribution, Stream & stream)
boolsample (const dense_categorical< V > & distribution, unsigned int & output)
structconstant
structuniform_distribution
structsequence_distribution
boolinit (sequence_distribution< ElementDistribution > & distribution, const sequence_distribution< ElementDistribution > & src)
boolread (sequence_distribution< ElementDistribution > & distribution, Stream & stream)
boolwrite (const sequence_distribution< ElementDistribution > & distribution, Stream & stream)
boolsample (const sequence_distribution< ElementDistribution > & distribution, SequenceType & output)
double log_rising_factorial(
doublebase,
unsigned intexponent)

Returns the natural logarithm of the rising factorial:

\[ \log a^{(n)} = \log \left\{ \prod_{i=0}^{n-1} a^i \right\}, \]
where base is $ a $ and exponent is $ n $.

struct symmetric_dirichlet

template<typename V>

A structure representing a symmetric Dirichlet distribution, which is a regular Dirichlet distribution, where all the elements of the concentration parameter vector pi have the same value.

The following example constructs a Dirichlet distribution with dimension 3, and concentration parameter [0.1, 0.1, 0.1]. It proceeds to generate two samples from this distribution. The expected output is [0.004223, 0.995777, 0.000000] [0.964576, 0.035400, 0.000024].

#include <math/distributions.h>
using namespace core;

template<typename V>
struct custom_vector {
    V* elements;

    custom_vector(unsigned int length) {
        elements = (V*) malloc(sizeof(V) * length);
    }

    ~custom_vector() { free(elements); }

    inline V get(unsigned int index) const {
        return elements[index];
    }

    void set(unsigned int index, const V& value) {
        elements[index] = value;
    }
};

int main() {
    set_seed(100);
    symmetric_dirichlet<double> dir(3, 0.1);

    custom_vector<double> output(3);
    for (unsigned int i = 0; i < 2; i++) {
        sample(dir, output);
        print(output.elements, 3, stdout); print(' ', stdout);
    }
}

V

satisfies is_arithmetic.

Public members
Vpi
unsigned intatom_count
Vlog_prob
Vtotal
symmetric_dirichlet (const symmetric_dirichlet< V > & prior)
symmetric_dirichlet (unsigned int atom_count, const V & prior)
voidensure_atom_count (unsigned int new_atom_count)
Vsum () const
Vmax () const
Vget_at_index (unsigned int index) const
Vget_for_atom (const K & atom) const
const symmetric_dirichlet< V > &get_parameters () const
voidsample (Destination & dst) const
Vlog_probability (unsigned int item) const
static voidfree (symmetric_dirichlet< V > & distribution)
typedefV value_type
V symmetric_dirichlet::pi

The value of each element of the concentration parameter.

unsigned int symmetric_dirichlet::atom_count

The number of dimensions of the Dirichlet distribution.

V symmetric_dirichlet::log_prob

The negative logarithm of symmetric_dirichlet::atom_count.

V symmetric_dirichlet::total

The sum of the elements in the concentration parameter vector. In the case of a symmetric Dirichlet distribution, this simply symmetric_dirichlet::pi * symmetric_dirichlet::atom_count.

symmetric_dirichlet::symmetric_dirichlet(
const symmetric_dirichlet< V > &prior)

Constructs a symmetric Dirichlet distribution by copying the fields from the given symmetric Dirichlet distribution.

symmetric_dirichlet::symmetric_dirichlet(
unsigned intatom_count,
const V &prior)

Constructs a symmetric Dirichlet distribution with the given dimension and concentration parameter.

void symmetric_dirichlet::ensure_atom_count(
unsigned intnew_atom_count)

Checks if symmetric_dirichlet::atom_count is smaller than the given new_atom_count. If so, symmetric_dirichlet::atom_count is set to new_atom_count and symmetric_dirichlet::total is updated accordingly.

V symmetric_dirichlet::sum()

Returns the sum of the elements in the concentration parameter vector (symmetric_dirichlet::total).

V symmetric_dirichlet::max()

Returns the maximum of the elements in the concentration parameter vector. In the case of the symmetric Dirichlet distribution, this is simply symmetric_dirichlet::pi.

V symmetric_dirichlet::get_at_index(
unsigned intindex) const

Returns the element at the given index in the concentration parameter vector. In the case of the symmetric Dirichlet distribution, this function always returns symmetric_dirichlet::pi.

template<typename K>
V symmetric_dirichlet::get_for_atom(
const K &atom) const

Returns the element at the given atom in the concentration parameter vector. The atom is a non-empty element drawn from a categorical distribution with a Dirichlet prior. In the case of the symmetric Dirichlet distribution, this function always returns symmetric_dirichlet::pi.

const symmetric_dirichlet< V > & symmetric_dirichlet::get_parameters()

Returns the parameters that may be used to construct this distribution, using either the constructor or init.

template<typename Destination>
void symmetric_dirichlet::sample(
Destination &dst) const
Samples from this symmetric Dirichlet distribution and puts the result in dst.
Destination

a vector type that implements the public member function set(unsigned int, const V&).

V symmetric_dirichlet::log_probability(
unsigned intitem) const

Returns the log probability of a single observation atom, drawn from a Dirichlet-categorical distribution, where the Dirichlet is represented by this object. In the case of the symmetric Dirichlet distribution, this function always returns symmetric_dirichlet::log_prob.

static void symmetric_dirichlet::free(
symmetric_dirichlet< V > &distribution)

Frees the given symmetric Dirichlet distribution. Since symmetric_dirichlet does not allocate additional memory, this function is a no-op.

typedef V symmetric_dirichlet::value_type

The type of the probabilities.

template<typename V>
bool init(
symmetric_dirichlet< V > &distribution,
const symmetric_dirichlet< V > &src)

Initializes the given symmetric_dirichlet distribution by copying its fields from the given symmetric_dirichlet src.

template<typename V>
bool init(
symmetric_dirichlet< V > &distribution,
const dirichlet< V > &src)

Since a symmetric Dirichlet is a special case of a Dirichlet distribution, this function prints an error and exits.

template<typename V>
bool read(
symmetric_dirichlet< V > &distribution,
FILE *in)

Reads a symmetric_dirichlet distribution from in.

template<typename V>
bool write(
const symmetric_dirichlet< V > &distribution,
FILE *out)

Writes the symmetric_dirichlet distribution to out.

template<typename V, typename Destination>
bool sample(
const symmetric_dirichlet< V > &distribution,
Destination &dst)
Samples from the given symmetric Dirichlet distribution and puts the result in dst.
Destination

a vector type that implements the public member function set(unsigned int, const V&).

struct dirichlet

template<typename V>

A structure representing a finite Dirichlet distribution, which is a generalization of the symmetric Dirichlet distribution (see struct symmetric_dirichlet).

The following example constructs a Dirichlet distribution with dimension 3, and concentration parameter [10, 1, 2]. It proceeds to generate two samples from this distribution. The expected output is [0.491806, 0.023958, 0.484236] [0.804916, 0.019964, 0.175120].

#include <math/distributions.h>
using namespace core;

template<typename V>
struct custom_vector {
    V* elements;

    custom_vector(unsigned int length) {
        elements = (V*) malloc(sizeof(V) * length);
    }

    ~custom_vector() { free(elements); }

    inline V get(unsigned int index) const {
        return elements[index];
    }

    void set(unsigned int index, const V& value) {
        elements[index] = value;
    }
};

int main() {
    set_seed(100);
    double alpha[] = {10.0, 1.0, 2.0};
    dirichlet<double> dir(3, alpha);

    custom_vector<double> output(3);
    for (unsigned int i = 0; i < 2; i++) {
        sample(dir, output);
        print(output.elements, 3, stdout); print(' ', stdout);
    }
}

V

satisfies is_arithmetic.

Public members
V *pi
Vpi_sum
unsigned intatom_count
dirichlet (const symmetric_dirichlet< V > & prior)
dirichlet (const dirichlet< V > & prior)
dirichlet (unsigned int atom_count, const V & prior)
dirichlet (unsigned int atom_count, const V * prior)
voidensure_atom_count (unsigned int new_atom_count)
Vsum () const
Vget_at_index (unsigned int index) const
Vget_for_atom (unsigned int atom) const
const dirichlet< V > &get_parameters () const
voidsample (Destination & dst) const
Vlog_probability (unsigned int atom) const
static voidfree (dirichlet< V > & distribution)
typedefV value_type
V * dirichlet::pi

The concentration parameter.

V dirichlet::pi_sum

The sum of all elements in the concentration parameter dirichlet::pi.

unsigned int dirichlet::atom_count

The number of dimensions of this Dirichlet distribution.

dirichlet::dirichlet(
const symmetric_dirichlet< V > &prior)

Constructs a Dirichlet distribution by copying the fields from the given symmetric Dirichlet distribution.

dirichlet::dirichlet(
const dirichlet< V > &prior)

Constructs a Dirichlet distribution by copying the fields from the given Dirichlet distribution.

dirichlet::dirichlet(
unsigned intatom_count,
const V &prior)

Constructs a Dirichlet distribution with the given dimension atom_count and symmetric concentration parameter prior.

dirichlet::dirichlet(
unsigned intatom_count,
const V *prior)

Constructs a Dirichlet distribution with the given dimension atom_count and concentration parameter vector prior.

void dirichlet::ensure_atom_count(
unsigned intnew_atom_count)

Checks if dirichlet::atom_count is smaller than the given new_atom_count. If so, dirichlet::atom_count is set to new_atom_count.

V dirichlet::sum()

Returns the sum of the elements in the concentration parameter vector (dirichlet::pi_sum).

V dirichlet::get_at_index(
unsigned intindex) const

Returns the element at the given index in the concentration parameter vector. This function does not perform any bounds checking.

V dirichlet::get_for_atom(
unsigned intatom) const

Returns the element at the given atom in the concentration parameter vector. The atom is a non-zero unsigned integer drawn from a categorical distribution with a Dirichlet prior. Thus, the atom n corresponds to the index n - 1. This function does not perform any bounds checking.

const dirichlet< V > & dirichlet::get_parameters()

Returns the parameters that may be used to construct this distribution, using either the constructor or init.

template<typename Destination>
void dirichlet::sample(
Destination &dst) const
Samples from this Dirichlet distribution and puts the result in dst.
Destination

a vector type that implements the public member function set(unsigned int, const V&).

V dirichlet::log_probability(
unsigned intatom) const

Returns the log probability of a single observation atom, drawn from a Dirichlet-categorical distribution, where the Dirichlet is represented by this object.

static void dirichlet::free(
dirichlet< V > &distribution)

Frees dirichlet::pi in the given Dirichlet distribution.

typedef V dirichlet::value_type

The type of the probabilities.

template<typename V>
bool init(
dirichlet< V > &distribution,
const symmetric_dirichlet< V > &src)

Initializes the given Dirichlet distribution by copying its fields from the given symmetric_dirichlet distribution src.

template<typename V>
bool init(
dirichlet< V > &distribution,
const dirichlet< V > &src)

Initializes the given Dirichlet distribution by copying its fields from the given Dirichlet distribution src.

template<typename V>
bool read(
dirichlet< V > &distribution,
FILE *in)

Reads a dirichlet distribution from in.

template<typename V>
bool write(
const dirichlet< V > &distribution,
FILE *out)

Writes the given dirichlet distribution to out.

template<typename V, typename Destination>
bool sample(
const dirichlet< V > &distribution,
Destination &dst)
Samples from the given Dirichlet distribution and puts the result in dst.
Destination

a vector type that implements the public member function set(unsigned int, const V&).

struct is_dirichlet

template<typename T>

A type trait that is true_type if and only if T is either symmetric_dirichlet or dirichlet.

struct sparse_categorical

template<typename K, typename V>

This struct represents a categorical distribution, where the probabilities are represented as a sum of uniform and a non-uniform component. The non-uniform component is stored as a core::hash_map from atoms to probabilities. The dense_categorical struct instead stores all probabilities contiguously as a single native array. dense_categorical should be used if the dimension is small or if the probabilities cannot be easily represented sparsely as a sum of uniform and non-uniform components. Unlike dense_categorical, the observations do not necessarily have type unsigned int, and can have generic type K.

The following is an example where a sparse_categorical distribution is constructed over the domain {'a', 'b', 'c', 'd', 'e'}. The probability of each event is specified and 10 samples are drawn. The expected output is d, b, d, d, b, c, a, b, e, e,.

#include <math/distributions.h>
using namespace core;

int main() {
    set_seed(100);
    sparse_categorical<char, double> categorical(5);
    categorical.set('a', 0.1);
    categorical.set('b', 0.4);
    categorical.set('c', 0.1);
    categorical.set('d', 0.2);
    categorical.set('e', 0.2);

    for (unsigned int i = 0; i < 10; i++) {
        char c;
        sample(categorical, c);
        printf("%c, ", c);
    }
}

Kthe generic type of the observations. K must satisfy either:
  1. is_fundamental,

  2. is_enum,

  3. is_pointer,

  4. implements the public static method unsigned int hash(const T&), the public static method void is_empty(const T&), implements the operators ==, satisfies CopyAssignable, and core::is_moveable. NOTE: The first argument to the == operator may be empty.

V

satisfies is_arithmetic.

Public members
hash_map< K, pair< V, V > >probabilities
unsigned intatom_count
Vprob
Vdense_prob
Vlog_prob
sparse_categorical (unsigned int atom_count)
sparse_categorical (const sparse_categorical< K, V > & src)
boolset (const K & key, const V & probability)
Vprobability (const K & observation) const
Vlog_probability (const K & observation) const
static Vconditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, const array_multiset< K > & items, const array_multiset< K > & conditioned)
static voidfree (sparse_categorical< K, V > & distribution)
typedefV value_type
hash_map< K, pair< V, V > > sparse_categorical::probabilities

A hash_map that encodes the non-uniform component of the categorical distribution. It maps from atoms to pairs, where the first entry in the pair contains the probability and the second entry contains the log probability.

unsigned int sparse_categorical::atom_count

The number of dimensions of the categorical distribution.

V sparse_categorical::prob

Stores the probability of every atom in the uniform component of the categorical distribution (i.e. every atom that is not a key in sparse_categorical::probabilities).

V sparse_categorical::dense_prob

Stores the total probability mass in the non-uniform component of the categorical distribution (i.e. the sum of the probabilities in sparse_categorical::probabilities).

V sparse_categorical::log_prob

The natural logarithm of sparse_categorical::prob.

sparse_categorical::sparse_categorical(
unsigned intatom_count)

Initializes this categorical distribution with the given dimension atom_count, setting all probabilities to zero.

sparse_categorical::sparse_categorical(
const sparse_categorical< K, V > &src)

Initializes this categorical distribution by copying its fields from the given sparse_categorical distribution src.

bool sparse_categorical::set(
const K &key,
const V &probability)

Sets the probability of the given observation key. This function will make the key part of the non-uniform component of the categorical distribution.

V sparse_categorical::probability(
const K &observation) const

Returns the probability of the given observation key.

V sparse_categorical::log_probability(
const K &observation) const

Returns the log probability of the given observation key.

template<typename PriorDistribution>
static V sparse_categorical::conditional(
const PriorDistribution &prior,
const K &item,
const array_multiset< K > &conditioned)
Computes the conditional probability of observing the given item drawn from this constant distribution, conditioned on a collection of observations conditioned. This function assumes all the elements in conditioned are identical, and that item and conditioned have non-zero probability according to prior.
Returns

true if item is equivalent to the first element in conditioned.

false otherwise.

PriorDistribution

the type of the prior distribution. This function does not use prior.

template<typename PriorDistribution>
static V sparse_categorical::log_conditional(
const PriorDistribution &prior,
const K &item,
const array_multiset< K > &conditioned)
Returns the log probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V sparse_categorical::log_conditional(
const PriorDistribution &prior,
const array_multiset< K > &items,
const array_multiset< K > &conditioned)
Returns the probability of observing the given collection of items, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

static void sparse_categorical::free(
sparse_categorical< K, V > &distribution)

Frees the given sparse_categorical distribution by releasing the memory resources associated with sparse_categorical::probabilities, along with all of its elements.

typedef V sparse_categorical::value_type

The type of the probabilities.

template<typename K, typename V>
bool init(
sparse_categorical< K, V > &distribution,
const sparse_categorical< K, V > &src)

Initializes the given sparse_categorical distribution by copying its fields from the given sparse_categorical distribution src.

template<typename K, typename V>
bool sample(
const sparse_categorical< K, V > &distribution,
K &output)
Draws a sample from the given sparse_categorical distribution and stores it in output. It is possible that the observation will be sampled from the uniform component of the sparse categorical distribution, in which case this function will print an error and return false.
K

satisfies core::is_copyable.

struct dense_categorical

template<typename V>

This struct represents a categorical distribution, where the probabilities are stored contiguously in a native array of type V. The sparse_categorical struct instead stores only the non-uniform probabilities in a core::hash_map. sparse_categorical should be used if the dimension is large and the probabilities can be represented sparsely as a sum of uniform and non-uniform components. NOTE: This distribution assumes the observations have type unsigned int and take values in {1, ..., dense_categorical::atom_count}.

In the following example, we define a categorical distribution over the ASCII characters (encoded as positive integers from 1 through 256). ASCII characters (encoded as positive integers from 1 through 256). Positive probability is only assigned to the characters {'a', 'b', 'c', 'd', 'e'}, and we draw 10 samples from the distribution. The expected output of the example is a, c, b, b, d, e, b, e, b, b,.

#include <math/distributions.h>
using namespace core;

template<typename V>
inline void set_probability(
        dense_categorical<V>& categorical,
        unsigned int atom, const V& probability)
{
    categorical.phi[atom - 1] = probability;
}

int main() {
    set_seed(100);
    dense_categorical<double> categorical(256);
    set_probability(categorical, 'a', 0.1);
    set_probability(categorical, 'b', 0.4);
    set_probability(categorical, 'c', 0.1);
    set_probability(categorical, 'd', 0.2);
    set_probability(categorical, 'e', 0.2);

    for (unsigned int i = 0; i < 10; i++) {
        unsigned int output;
        sample(categorical, output);
        printf("%c, ", output);
    }
}

V

satisfies is_arithmetic.

Public members
V *phi
unsigned intatom_count
dense_categorical (unsigned int atom_count)
dense_categorical (const dense_categorical & src)
Vget (unsigned int index) const
Vget_for_atom (unsigned int atom) const
Vsum () const
Vprobability (unsigned int item) const
Vprobability (const array_multiset< unsigned int > & items) const
Vlog_probability (unsigned int item) const
Vlog_probability (const array_multiset< unsigned int > & items) const
unsigned intget_parameters () const
static Vprobability (const PriorDistribution & prior, unsigned int item)
static Vlog_probability (const PriorDistribution & prior, unsigned int item)
static Vlog_probability (const PriorDistribution & prior, const array_multiset< unsigned int > & items)
static Vconditional (const PriorDistribution & prior, unsigned int item, const array_multiset< unsigned int > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, unsigned int item, const array_multiset< unsigned int > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, unsigned int item, const hash_multiset< unsigned int > & conditioned)
static Vlog_conditional_without (const PriorDistribution & prior, unsigned int holdout, const array_multiset< unsigned int > & conditioned)
static Vlog_conditional_without (const PriorDistribution & prior, unsigned int holdout, const hash_multiset< unsigned int > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, const array_multiset< unsigned int > & items, const array_multiset< unsigned int > & conditioned)
static Vlog_conditional (const PriorDistribution & prior, const array_multiset< unsigned int > & items, const hash_multiset< unsigned int > & conditioned)
static Vlog_conditional_without (const PriorDistribution & prior, const array_multiset< unsigned int > & holdout, const array_multiset< unsigned int > & conditioned)
static Vlog_conditional_without (const PriorDistribution & prior, const array_multiset< unsigned int > & holdout, const hash_multiset< unsigned int > & conditioned)
static voidmove (const dense_categorical< V > & src, dense_categorical< V > & dst)
static boolcopy (const dense_categorical< V > & src, dense_categorical< V > & dst)
static voidfree (dense_categorical< V > & distribution)
typedefV value_type
V * dense_categorical::phi

The native array containing the probabilities of each event.

unsigned int dense_categorical::atom_count

The number of dimensions of the categorical distribution.

dense_categorical::dense_categorical(
unsigned intatom_count)

Initializes this categorical distribution with the given dimension atom_count, setting all probabilities to zero.

dense_categorical::dense_categorical(
const dense_categorical &src)

Initializes this categorical distribution by copying the fields from the given dense_categorical distribution src.

V dense_categorical::get(
unsigned intindex) const

Returns the element at the given index of the probability vector dense_categorical::phi. This is equivalent to the probability of observing index + 1 drawn from this distribution.

V dense_categorical::get_for_atom(
unsigned intatom) const

Returns the element at the index atom - 1 of the probability vector dense_categorical::phi. This is equivalent to the probability of observing atom drawn from this distribution.

V dense_categorical::sum()

It is assumed that the probabilities in dense_categorical::phi sum to 1, and so this function always returns 1.0.

V dense_categorical::probability(
unsigned intitem) const

Returns the probability of observing item drawn from this distribution.

V dense_categorical::probability(
const array_multiset< unsigned int > &items) const

Returns the joint probability of observing the given collection of items, drawn independently and identically from this distribution.

V dense_categorical::log_probability(
unsigned intitem) const

Returns the log probability of observing item drawn from this distribution.

V dense_categorical::log_probability(
const array_multiset< unsigned int > &items) const

Returns the joint log probability of observing the given collection of items, drawn independently and identically from this distribution.

unsigned int dense_categorical::get_parameters()

Returns the parameters that may be used to construct this distribution, using either the constructor or init. This particular function simply returns dense_categorical::atom_count.

template<typename PriorDistribution>
static V dense_categorical::probability(
const PriorDistribution &prior,
unsigned intitem)
Returns the probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_probability(
const PriorDistribution &prior,
unsigned intitem)
Returns the log probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_probability(
const PriorDistribution &prior,
const array_multiset< unsigned int > &items)
Returns the probability of observing the given collection of items, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::conditional(
const PriorDistribution &prior,
unsigned intitem,
const array_multiset< unsigned int > &conditioned)
Returns the probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional(
const PriorDistribution &prior,
unsigned intitem,
const array_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional(
const PriorDistribution &prior,
unsigned intitem,
const hash_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given item, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional_without(
const PriorDistribution &prior,
unsigned intholdout,
const array_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given item holdout, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned \ { holdout } (the set of conditioned observations with holdout being subtracted). It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional_without(
const PriorDistribution &prior,
unsigned intholdout,
const hash_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given item holdout, drawn from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned \ { holdout } (the set of conditioned observations with holdout being subtracted). It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional(
const PriorDistribution &prior,
const array_multiset< unsigned int > &items,
const array_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given collection of items, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional(
const PriorDistribution &prior,
const array_multiset< unsigned int > &items,
const hash_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given collection of items, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned. It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional_without(
const PriorDistribution &prior,
const array_multiset< unsigned int > &holdout,
const array_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given collection of items holdout, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned \ holdout (the set of conditioned observations with holdout being subtracted). It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

template<typename PriorDistribution>
static V dense_categorical::log_conditional_without(
const PriorDistribution &prior,
const array_multiset< unsigned int > &holdout,
const hash_multiset< unsigned int > &conditioned)
Returns the log probability of observing the given collection of items holdout, each drawn independently and identically from a categorical distribution, which is itself drawn from the given prior distribution, conditioned on the set of observations conditioned \ holdout (the set of conditioned observations with holdout being subtracted). It is assumed the given prior is a Dirichlet.
PriorDistribution

a distribution type with public member functions V get_for_atom(unsigned int) and V sum().

static void dense_categorical::move(
const dense_categorical< V > &src,
dense_categorical< V > &dst)

Moves the dense_categorical distribution from src into dst. This function simply copies the pointers, and does not initialize a new probability array.

static bool dense_categorical::copy(
const dense_categorical< V > &src,
dense_categorical< V > &dst)

Copies the dense_categorical distribution from src into dst. This function initializes a new probability array in dst and copies into it the contents from the array in src.

static void dense_categorical::free(
dense_categorical< V > &distribution)

Frees the given categorical distribution, by releasing the resources associated with the native array dense_categorical::phi.

typedef V dense_categorical::value_type

The type of the probabilities.

template<typename V>
bool init(
dense_categorical< V > &distribution,
unsigned intatom_count)

Initializes the given categorical distribution with the given dimension atom_count, setting all probabilities to zero.

template<typename V>
bool init(
dense_categorical< V > &distribution,
const dense_categorical< V > &src)

Initializes the given categorical distribution, copying its fields from the given dense_categorical distribution src.

template<typename V, typename Stream>
bool read(
dense_categorical< V > &distribution,
Stream &stream)

Reads a dense_categorical distribution from stream.

Stream

satisfies is_readable.

template<typename V, typename Stream>
bool write(
const dense_categorical< V > &distribution,
Stream &stream)

Writes the given dense_categorical distribution to stream.

Stream

satisfies is_writeable.

template<typename V>
bool sample(
const dense_categorical< V > &distribution,
unsigned int &output)

Draws a sample from the given dense_categorical distribution and writes it to output.

struct constant

template<typename K>
A struct that represents the degenerate/constant distribution.
K

the type of the observation, which must implement the operator ==.

Public members
static boolconditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned)
static boolconditional (const PriorDistribution & prior, const array_multiset< K > & items, const array_multiset< K > & conditioned)
static doublelog_conditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned)
static doublelog_conditional (const PriorDistribution & prior, const array_multiset< K > & items, const array_multiset< K > & conditioned)
static PriorDistribution::value_typeprobability (const PriorDistribution & prior, const K & item)
static PriorDistribution::value_typeprobability (const PriorDistribution & prior, const array_multiset< K > & items)
static PriorDistribution::value_typelog_probability (const PriorDistribution & prior, const K & item)
static PriorDistribution::value_typelog_probability (const PriorDistribution & prior, const array_multiset< K > & items)
static boolsample (const PriorDistribution & prior, const array_multiset< K > & observations, K & sample)
template<typename PriorDistribution>
static bool constant::conditional(
const PriorDistribution &prior,
const K &item,
const array_multiset< K > &conditioned)
Computes the conditional probability of observing the given item drawn from this constant distribution, conditioned on a collection of observations conditioned. This function assumes all the elements in conditioned are identical, and that item and conditioned have non-zero probability according to prior.
Returns

true if item is equivalent to the first element in conditioned.

false otherwise.

PriorDistribution

the type of the prior distribution. This function does not use prior.

template<typename PriorDistribution>
static bool constant::conditional(
const PriorDistribution &prior,
const array_multiset< K > &items,
const array_multiset< K > &conditioned)
Computes the conditional probability of observing the given collection of items, each drawn independently and identically from this constant distribution, conditioned on a collection of observations conditioned. This function assumes all the elements in items are identical, all the elements in conditioned are identical, and that items and conditioned have non-zero probability according to prior.
Returns

true if the first element in items is equivalent to the first element in conditioned.

false otherwise.

PriorDistribution

the type of the prior distribution. This function does not use prior.

template<typename PriorDistribution>
static double constant::log_conditional(
const PriorDistribution &prior,
const K &item,
const array_multiset< K > &conditioned)
Computes the conditional log probability of observing the given item drawn from this constant distribution, conditioned on a collection of observations conditioned. This function assumes all the elements in conditioned are identical, and that item and conditioned have non-zero probability according to prior.
Returns

0 if item is equivalent to the first element in conditioned.

-inf otherwise.

PriorDistribution

the type of the prior distribution. This function does not use prior.

template<typename PriorDistribution>
static double constant::log_conditional(
const PriorDistribution &prior,
const array_multiset< K > &items,
const array_multiset< K > &conditioned)
Computes the conditional probability of observing the given collection of items, each drawn independently and identically from this constant distribution, conditioned on a collection of observations conditioned. This function assumes all the elements in items are identical, all the elements in conditioned are identical, and that items and conditioned have non-zero probability according to prior.
Returns

0 if the first element in items is equivalent to the first element in conditioned.

-inf otherwise.

PriorDistribution

the type of the prior distribution. This function does not use prior.

template<typename PriorDistribution>
static PriorDistribution::value_type constant::probability(
const PriorDistribution &prior,
const K &item)
Returns the probability of observing the given item, drawn from a constant distribution, which is itself drawn from the given prior distribution.
PriorDistribution

a distribution type that contains the typedef value_type implements the public member function value_type probability(const K&).

template<typename PriorDistribution>
static PriorDistribution::value_type constant::probability(
const PriorDistribution &prior,
const array_multiset< K > &items)
Returns the probability of observing the given collection of items, each drawn independently and identically from a constant distribution, which is itself drawn from the given prior distribution.
PriorDistribution

a distribution type that contains the typedef value_type implements the public member function value_type probability(const array_multiset<K>&).

template<typename PriorDistribution>
static PriorDistribution::value_type constant::log_probability(
const PriorDistribution &prior,
const K &item)
Returns the log probability of observing the given item, drawn from a constant distribution, which is itself drawn from the given prior distribution.
PriorDistribution

a distribution type that contains the typedef value_type implements the public member function value_type log_probability(const K&).

template<typename PriorDistribution>
static PriorDistribution::value_type constant::log_probability(
const PriorDistribution &prior,
const array_multiset< K > &items)
Returns the log probability of observing the given collection of items, each drawn independently and identically from a constant distribution, which is itself drawn from the given prior distribution.
PriorDistribution

a distribution type that contains the typedef value_type implements the public member function value_type log_probability(const array_multiset<K>&).

template<typename PriorDistribution>
static bool constant::sample(
const PriorDistribution &prior,
const array_multiset< K > &observations,
K &sample)
Samples an observation from a constant distribution, which is itself drawn from a prior distribution, conditioned on the fact that the given set of observations was sampled from the constant distribution. The sample is written to sample. This function assumes the elements in observations are identical, and have non-zero probability according to prior.
K

satisfies core::is_copyable.

PriorDistribution

the type of the prior distribution. This function does not use prior.

struct uniform_distribution

template<typename V>

This struct represents a discrete uniform distribution, where probabilities have type V.

Public members
Vprob
Vlog_prob
uniform_distribution (unsigned int count)
Vprobability (const K & observation) const
Vlog_probability (const K & observation) const
typedefV value_type
V uniform_distribution::prob

The probability of each event.

V uniform_distribution::log_prob

The log probability of each event.

uniform_distribution::uniform_distribution(
unsigned intcount)

Constructs a uniform distribution with the given dimension/number of distinct events count.

template<typename K>
V uniform_distribution::probability(
const K &observation) const

Returns the probability of the given observation, which is equivalently uniform_distribution::prob.

template<typename K>
V uniform_distribution::log_probability(
const K &observation) const

Returns the log probability of the given observation, which is equivalently uniform_distribution::log_prob.

typedef V uniform_distribution::value_type

The type of the probabilities.

struct sequence_distribution

template<typename ElementDistribution>

This struct represents a distribution over sequences of events, where each event is drawn independently and identically from ElementDistribution, until a special "stop event" is drawn with probability sequence_distribution::end_probability. We assume the first event is not the stop event.

The following example constructs a sequence_distribution where the ElementDistribution is a sparse_categorical distribution over the characters {'a', 'b', 'c', 'd', 'e'}. Five samples are drawn from the sequence distribution. The expected output is 'ddbc', 'eecdecdd', 'bbbb', 'dd', 'babababc',.

#include <math/distributions.h>
using namespace core;

int main() {
    set_seed(100);
    sparse_categorical<char, double> categorical(5);
    categorical.set('a', 0.1);
    categorical.set('b', 0.4);
    categorical.set('c', 0.1);
    categorical.set('d', 0.2);
    categorical.set('e', 0.2);

    sequence_distribution<sparse_categorical<char, double>> seq(categorical, 0.2);

    for (unsigned int i = 0; i < 5; i++) {
        string s;
        sample(seq, s);
        print("'", stdout);
        print(s, stdout);
        print("', ", stdout);
    }
}

ElementDistribution

a distribution type that defines a public typedef value_type that indicates the type of the probabilities. Some operations in this class also require the public member functions value_type probability(const T&), value_type log_probability(const T&).

Public members
Vend_probability
Vlog_end_probability
Vlog_not_end_probability
ElementDistributionelement_distribution
sequence_distribution (ElementDistribution & element_distribution, V end_probability)
Vprobability (const SequenceType & sequence) const
Vlog_probability (const SequenceType & sequence) const
static boolcopy (const sequence_distribution< ElementDistribution > & src, sequence_distribution< ElementDistribution > & dst)
static voidfree (sequence_distribution< ElementDistribution > & distribution)
typedefElementDistribution::value_type V
V sequence_distribution::end_probability

After the first event in the sequence, this is the probability that no further events are drawn, and the sequence ends.

V sequence_distribution::log_end_probability

The natural logarithm of sequence_distribution::end_probability.

V sequence_distribution::log_not_end_probability

The natural logarithm of (1 - sequence_distribution::end_probability).

ElementDistribution sequence_distribution::element_distribution

The distribution from which each element in the sequence is drawn independently and identically.

sequence_distribution::sequence_distribution(
ElementDistribution &element_distribution,
Vend_probability)

Constructs a sequence_distribution with the given sequence_distribution::element_distribution and sequence_distribution::end_probability.

ElementDistribution

a distribution type that can be constructed using a single parameter with type ElementDistribution&.

template<typename SequenceType>
V sequence_distribution::probability(
const SequenceType &sequence) const
Returns the probability of the given observation sequence under this sequence distribution.
SequenceType

a sequence type that contains the integral field length that indicates the number of elements in the sequence, and implements the operator [] that returns an element at a particular index of the sequence.

template<typename SequenceType>
V sequence_distribution::log_probability(
const SequenceType &sequence) const
Returns the log probability of the given observation sequence under this sequence distribution.
SequenceType

a sequence type that contains the integral field length that indicates the number of elements in the sequence, and implements the operator [] that returns an element at a particular index of the sequence.

Copies the sequence_distribution src into dst.
ElementDistribution

satisfies core::is_copyable.

static void sequence_distribution::free(
sequence_distribution< ElementDistribution > &distribution)
Frees the given sequence distribution.
ElementDistribution

the type of the distribution of each element in the sequence. The function core::free is called with an argument of type ElementDistribution&.

typedef ElementDistribution::value_type sequence_distribution::V

The type of the probabilities.

template<typename ElementDistribution>
bool init(
sequence_distribution< ElementDistribution > &distribution,
const sequence_distribution< ElementDistribution > &src)
Initializes the given sequence_distribution distribution with the given sequence_distribution src.
ElementDistribution

a distribution type for which the function bool init(ElementDistribution&, const ElementDistribution) is implemented.

template<typename ElementDistribution, typename Stream>
bool read(
sequence_distribution< ElementDistribution > &distribution,
Stream &stream)

Reads the given sequence_distribution distribution from stream.

Stream

satisfies is_readable.

template<typename ElementDistribution, typename Stream>
bool write(
const sequence_distribution< ElementDistribution > &distribution,
Stream &stream)

Writes the given sequence_distribution distribution to stream.

Stream

satisfies is_readable.

template<typename ElementDistribution, typename SequenceType>
bool sample(
const sequence_distribution< ElementDistribution > &distribution,
SequenceType &output)
Samples from the given sequence_distribution distribution and stores the result in output.
SequenceType

a sequence type for which the function bool init(SequenceType&, unsigned int) is implemented, which initializes the sequence with the given length. The operator [] must also be implemented, which returns a reference to the specified index in the sequence with type T, so that the function bool sample(const ElementDistribution&, T&) can be called. Finally, the function core::free may be called on the elements of the sequence, if the sampling fails.