Returns the natural logarithm of the rising factorial:
base
is $ a $ and exponent
is $ n $.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 | |
---|---|
V | pi |
unsigned int | atom_count |
V | log_prob |
V | total |
symmetric_dirichlet (const symmetric_dirichlet< V > & prior) | |
symmetric_dirichlet (unsigned int atom_count, const V & prior) | |
void | ensure_atom_count (unsigned int new_atom_count) |
V | sum () const |
V | max () const |
V | get_at_index (unsigned int index) const |
V | get_for_atom (const K & atom) const |
const symmetric_dirichlet< V > & | get_parameters () const |
void | sample (Destination & dst) const |
V | log_probability (unsigned int item) const |
static void | free (symmetric_dirichlet< V > & distribution) |
typedef | V value_type |
The value of each element of the concentration parameter.
The number of dimensions of the Dirichlet distribution.
The negative logarithm of symmetric_dirichlet::atom_count
.
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
.
const symmetric_dirichlet< V > & | prior | ) |
Constructs a symmetric Dirichlet distribution by copying the fields from the given symmetric Dirichlet distribution.
Constructs a symmetric Dirichlet distribution with the given dimension and concentration parameter.
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.
Returns the sum of the elements in the concentration parameter vector (symmetric_dirichlet::total).
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.
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.
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.
Returns the parameters that may be used to construct this distribution, using either the constructor or init.
dst
.Destination | a vector type that implements the public member function |
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.
Frees the given symmetric Dirichlet distribution. Since symmetric_dirichlet does not allocate additional memory, this function is a no-op.
The type of the probabilities.
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
.
Since a symmetric Dirichlet is a special case of a Dirichlet distribution, this function prints an error and exits.
Reads a symmetric_dirichlet distribution
from in
.
Writes the symmetric_dirichlet distribution
to out
.
const symmetric_dirichlet< V > & | distribution, | |
Destination & | dst | ) |
distribution
and puts the result in dst
.Destination | a vector type that implements the public member function |
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 |
V | pi_sum |
unsigned int | atom_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) | |
void | ensure_atom_count (unsigned int new_atom_count) |
V | sum () const |
V | get_at_index (unsigned int index) const |
V | get_for_atom (unsigned int atom) const |
const dirichlet< V > & | get_parameters () const |
void | sample (Destination & dst) const |
V | log_probability (unsigned int atom) const |
static void | free (dirichlet< V > & distribution) |
typedef | V value_type |
The concentration parameter.
The sum of all elements in the concentration parameter dirichlet::pi.
The number of dimensions of this Dirichlet distribution.
const symmetric_dirichlet< V > & | prior | ) |
Constructs a Dirichlet distribution by copying the fields from the given symmetric Dirichlet distribution.
Constructs a Dirichlet distribution by copying the fields from the given Dirichlet distribution.
Constructs a Dirichlet distribution with the given dimension atom_count
and symmetric concentration parameter prior
.
Constructs a Dirichlet distribution with the given dimension atom_count
and concentration parameter vector prior
.
Checks if dirichlet::atom_count is smaller than the given new_atom_count
. If so, dirichlet::atom_count is set to new_atom_count
.
Returns the sum of the elements in the concentration parameter vector (dirichlet::pi_sum).
Returns the element at the given index
in the concentration parameter vector. This function does not perform any bounds checking.
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.
Returns the parameters that may be used to construct this distribution, using either the constructor or init.
dst
.Destination | a vector type that implements the public member function |
Returns the log probability of a single observation atom
, drawn from a Dirichlet-categorical distribution, where the Dirichlet is represented by this object.
Frees dirichlet::pi in the given Dirichlet distribution.
The type of the probabilities.
Initializes the given Dirichlet distribution
by copying its fields from the given symmetric_dirichlet distribution src
.
Initializes the given Dirichlet distribution
by copying its fields from the given Dirichlet distribution src
.
Reads a dirichlet distribution
from in
.
Writes the given dirichlet distribution
to out
.
const dirichlet< V > & | distribution, | |
Destination & | dst | ) |
distribution
and puts the result in dst
.Destination | a vector type that implements the public member function |
A type trait that is true_type if and only if T
is either symmetric_dirichlet or dirichlet.
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); } }
K | the generic type of the observations. K must satisfy either:
|
V | satisfies is_arithmetic. |
Public members | |
---|---|
hash_map< K, pair< V, V > > | probabilities |
unsigned int | atom_count |
V | prob |
V | dense_prob |
V | log_prob |
sparse_categorical (unsigned int atom_count) | |
sparse_categorical (const sparse_categorical< K, V > & src) | |
bool | set (const K & key, const V & probability) |
V | probability (const K & observation) const |
V | log_probability (const K & observation) const |
static V | conditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned) |
static V | log_conditional (const PriorDistribution & prior, const K & item, const array_multiset< K > & conditioned) |
static V | log_conditional (const PriorDistribution & prior, const array_multiset< K > & items, const array_multiset< K > & conditioned) |
static void | free (sparse_categorical< K, V > & distribution) |
typedef | V value_type |
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.
The number of dimensions of the categorical distribution.
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).
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).
The natural logarithm of sparse_categorical::prob.
Initializes this categorical distribution with the given dimension atom_count
, setting all probabilities to zero.
const sparse_categorical< K, V > & | src | ) |
Initializes this categorical distribution by copying its fields from the given sparse_categorical distribution src
.
Sets the probability
of the given observation key
. This function will make the key part of the non-uniform component of the categorical distribution.
Returns the probability of the given observation key
.
Returns the log probability of the given observation key
.
const PriorDistribution & | prior, | |
const K & | item, | |
const array_multiset< K > & | conditioned | ) |
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
.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 |
const PriorDistribution & | prior, | |
const K & | item, | |
const array_multiset< K > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< K > & | items, | |
const array_multiset< K > & | conditioned | ) |
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 |
Frees the given sparse_categorical distribution
by releasing the memory resources associated with sparse_categorical::probabilities, along with all of its elements.
The type of the probabilities.
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
.
const sparse_categorical< K, V > & | distribution, | |
K & | output | ) |
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. |
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. |
The native array containing the probabilities of each event.
The number of dimensions of the categorical distribution.
Initializes this categorical distribution with the given dimension atom_count
, setting all probabilities to zero.
const dense_categorical & | src | ) |
Initializes this categorical distribution by copying the fields from the given dense_categorical distribution src
.
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.
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.
It is assumed that the probabilities in dense_categorical::phi sum to 1, and so this function always returns 1.0.
Returns the probability of observing item
drawn from this distribution.
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.
Returns the log probability of observing item
drawn from this distribution.
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.
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.
const PriorDistribution & | prior, | |
unsigned int | item | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | item | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< unsigned int > & | items | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | item, | |
const array_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | item, | |
const array_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | item, | |
const hash_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | holdout, | |
const array_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
unsigned int | holdout, | |
const hash_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< unsigned int > & | items, | |
const array_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< unsigned int > & | items, | |
const hash_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< unsigned int > & | holdout, | |
const array_multiset< unsigned int > & | conditioned | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< unsigned int > & | holdout, | |
const hash_multiset< unsigned int > & | conditioned | ) |
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 |
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.
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
.
Frees the given categorical distribution, by releasing the resources associated with the native array dense_categorical::phi.
The type of the probabilities.
Initializes the given categorical distribution
with the given dimension atom_count
, setting all probabilities to zero.
dense_categorical< V > & | distribution, | |
const dense_categorical< V > & | src | ) |
Initializes the given categorical distribution
, copying its fields from the given dense_categorical distribution src
.
Reads a dense_categorical distribution
from stream
.
Stream | satisfies is_readable. |
const dense_categorical< V > & | distribution, | |
Stream & | stream | ) |
Writes the given dense_categorical distribution
to stream
.
Stream | satisfies is_writeable. |
Draws a sample from the given dense_categorical distribution
and writes it to output
.
K | the type of the observation, which must implement the operator |
const PriorDistribution & | prior, | |
const K & | item, | |
const array_multiset< K > & | conditioned | ) |
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
.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 |
const PriorDistribution & | prior, | |
const array_multiset< K > & | items, | |
const array_multiset< K > & | conditioned | ) |
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
.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 |
const PriorDistribution & | prior, | |
const K & | item, | |
const array_multiset< K > & | conditioned | ) |
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
.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 |
const PriorDistribution & | prior, | |
const array_multiset< K > & | items, | |
const array_multiset< K > & | conditioned | ) |
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
.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 |
const PriorDistribution & | prior, | |
const K & | item | ) |
item
, drawn from a constant distribution, which is itself drawn from the given prior
distribution.PriorDistribution | a distribution type that contains the typedef |
const PriorDistribution & | prior, | |
const array_multiset< K > & | items | ) |
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 |
const PriorDistribution & | prior, | |
const K & | item | ) |
item
, drawn from a constant distribution, which is itself drawn from the given prior
distribution.PriorDistribution | a distribution type that contains the typedef |
const PriorDistribution & | prior, | |
const array_multiset< K > & | items | ) |
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 |
const PriorDistribution & | prior, | |
const array_multiset< K > & | observations, | |
K & | sample | ) |
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 |
This struct represents a discrete uniform distribution, where probabilities have type V
.
Public members | |
---|---|
V | prob |
V | log_prob |
uniform_distribution (unsigned int count) | |
V | probability (const K & observation) const |
V | log_probability (const K & observation) const |
typedef | V value_type |
The probability of each event.
The log probability of each event.
Constructs a uniform distribution with the given dimension/number of distinct events count
.
Returns the probability of the given observation
, which is equivalently uniform_distribution::prob.
Returns the log probability of the given observation
, which is equivalently uniform_distribution::log_prob.
The type of the probabilities.
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 |
Public members | |
---|---|
V | end_probability |
V | log_end_probability |
V | log_not_end_probability |
ElementDistribution | element_distribution |
sequence_distribution (ElementDistribution & element_distribution, V end_probability) | |
V | probability (const SequenceType & sequence) const |
V | log_probability (const SequenceType & sequence) const |
static bool | copy (const sequence_distribution< ElementDistribution > & src, sequence_distribution< ElementDistribution > & dst) |
static void | free (sequence_distribution< ElementDistribution > & distribution) |
typedef | ElementDistribution::value_type V |
After the first event in the sequence, this is the probability that no further events are drawn, and the sequence ends.
The natural logarithm of sequence_distribution::end_probability.
The natural logarithm of (1 - sequence_distribution::end_probability).
The distribution from which each element in the sequence is drawn independently and identically.
ElementDistribution & | element_distribution, | |
V | end_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&. |
const SequenceType & | sequence | ) const |
sequence
under this sequence distribution.SequenceType | a sequence type that contains the integral field |
const SequenceType & | sequence | ) const |
sequence
under this sequence distribution.SequenceType | a sequence type that contains the integral field |
const sequence_distribution< ElementDistribution > & | src, | |
sequence_distribution< ElementDistribution > & | dst | ) |
sequence_distribution< ElementDistribution > & | distribution | ) |
distribution
.ElementDistribution | the type of the distribution of each element in the sequence. The function |
The type of the probabilities.
sequence_distribution< ElementDistribution > & | distribution, | |
const sequence_distribution< ElementDistribution > & | src | ) |
distribution
with the given sequence_distribution src
.ElementDistribution | a distribution type for which the function |
sequence_distribution< ElementDistribution > & | distribution, | |
Stream & | stream | ) |
Reads the given sequence_distribution distribution
from stream
.
Stream | satisfies is_readable. |
const sequence_distribution< ElementDistribution > & | distribution, | |
Stream & | stream | ) |
Writes the given sequence_distribution distribution
to stream
.
Stream | satisfies is_readable. |
const sequence_distribution< ElementDistribution > & | distribution, | |
SequenceType & | output | ) |
distribution
and stores the result in output
.SequenceType | a sequence type for which the function |