libDAI
Classes | Typedefs | Enumerations | Functions | Variables
dai Namespace Reference

Namespace for libDAI. More...

Classes

class  BBP
 Implements BBP (Back-Belief-Propagation) [EaG09]. More...
 
class  BBPCostFunction
 Predefined cost functions that can be used with BBP. More...
 
class  BipartiteGraph
 Represents the neighborhood structure of nodes in an undirected, bipartite graph. More...
 
class  BP
 Approximate inference algorithm "(Loopy) Belief Propagation". More...
 
class  BP_dual
 Calculates both types of BP messages and their normalizers from an InfAlg. More...
 
class  CBP
 Class for CBP (Conditioned Belief Propagation) [EaG09]. More...
 
class  ClusterGraph
 A ClusterGraph is a hypergraph with variables as nodes, and "clusters" (sets of variables) as hyperedges. More...
 
class  CobwebGraph
 A CobwebGraph is a special type of region graph used by the GLC algorithm. More...
 
class  CondProbEstimation
 Estimates the parameters of a conditional probability table, using pseudocounts. More...
 
class  DAG
 Represents the neighborhood structure of nodes in a directed cyclic graph. More...
 
class  DAIAlg
 Combines the abstract base class InfAlg with a graphical model (e.g., a FactorGraph or RegionGraph). More...
 
class  DEdge
 Represents a directed edge. More...
 
class  EMAlg
 EMAlg performs Expectation Maximization to learn factor parameters. More...
 
class  Evidence
 Stores a data set consisting of multiple samples, where each sample is the observed joint state of some variables. More...
 
class  ExactInf
 Exact inference algorithm using brute force enumeration (mainly useful for testing purposes) More...
 
class  Exception
 Error handling in libDAI is done by throwing an instance of the Exception class. More...
 
class  FactorGraph
 Represents a factor graph. More...
 
class  FBP
 Approximate inference algorithm "Fractional Belief Propagation" [WiH03]. More...
 
struct  fo_abs
 Function object that takes the absolute value. More...
 
struct  fo_absdiff
 Function object that returns the absolute difference of x and y. More...
 
struct  fo_divides0
 Function object similar to std::divides(), but different in that dividing by zero results in zero. More...
 
struct  fo_exp
 Function object that takes the exponent. More...
 
struct  fo_Hellinger
 Function object useful for calculating the Hellinger distance. More...
 
struct  fo_id
 Function object that returns the value itself. More...
 
struct  fo_inv
 Function object that takes the inverse. More...
 
struct  fo_inv0
 Function object that takes the inverse, except that 1/0 is defined to be 0. More...
 
struct  fo_KL
 Function object useful for calculating the KL distance. More...
 
struct  fo_log
 Function object that takes the logarithm. More...
 
struct  fo_log0
 Function object that takes the logarithm, except that log(0) is defined to be 0. More...
 
struct  fo_max
 Function object that returns the maximum of two values. More...
 
struct  fo_min
 Function object that returns the minimum of two values. More...
 
struct  fo_plog0p
 Function object that returns p*log0(p) More...
 
struct  fo_pow
 Function object that returns x to the power y. More...
 
class  FRegion
 An FRegion is a factor with a counting number. More...
 
class  GraphAL
 Represents the neighborhood structure of nodes in an undirected graph. More...
 
class  GraphEL
 Represents an undirected graph, implemented as a std::set of undirected edges. More...
 
class  greedyVariableElimination
 Helper object for dai::ClusterGraph::VarElim() More...
 
class  HAK
 Approximate inference algorithm: implementation of single-loop ("Generalized Belief Propagation") and double-loop algorithms by Heskes, Albers and Kappen [HAK03]. More...
 
class  hash_map
 hash_map is an alias for std::tr1::unordered_map. More...
 
class  IndexFor
 Tool for looping over the states of several variables. More...
 
class  InfAlg
 InfAlg is an abstract base class, defining the common interface of all inference algorithms in libDAI. More...
 
class  JTree
 Exact inference algorithm using junction tree. More...
 
class  LC
 Approximate inference algorithm "Loop Corrected Belief Propagation" [MoK07]. More...
 
class  MaximizationStep
 A MaximizationStep groups together several parameter estimation tasks (SharedParameters objects) into a single unit. More...
 
class  MF
 Approximate inference algorithm "Mean Field". More...
 
class  MR
 Approximate inference algorithm by Montanari and Rizzo [MoR05]. More...
 
class  multifor
 multifor makes it easy to perform a dynamic number of nested for loops. More...
 
struct  Neighbor
 Describes the neighbor relationship of two nodes in a graph. More...
 
class  ParameterEstimation
 Base class for parameter estimation methods. More...
 
class  Permute
 Tool for calculating permutations of linear indices of multi-dimensional arrays. More...
 
class  PropertySet
 Represents a set of properties, mapping keys (of type PropertyKey) to values (of type PropertyValue) More...
 
class  Region
 A Region is a set of variables with a counting number. More...
 
class  RegionGraph
 A RegionGraph combines a bipartite graph consisting of outer regions (type FRegion) and inner regions (type Region) with a FactorGraph. More...
 
class  RootedTree
 Represents a rooted tree, implemented as a vector of directed edges. More...
 
class  sequentialVariableElimination
 Helper object for dai::ClusterGraph::VarElim() More...
 
class  SharedParameters
 Represents a single factor or set of factors whose parameters should be estimated. More...
 
class  SmallSet
 Represents a set; the implementation is optimized for a small number of elements. More...
 
class  State
 Makes it easy to iterate over all possible joint states of variables within a VarSet. More...
 
class  TFactor
 Represents a (probability) factor. More...
 
class  TProb
 Represents a vector with entries of type T. More...
 
class  TreeEP
 Approximate inference algorithm "Tree Expectation Propagation" [MiQ04]. More...
 
class  TRWBP
 Approximate inference algorithm "Tree-Reweighted Belief Propagation" [WJW03]. More...
 
class  UEdge
 Represents an undirected edge. More...
 
class  Var
 Represents a discrete random variable. More...
 
class  VarSet
 Represents a set of variables. More...
 
class  WeightedGraph
 Represents an undirected weighted graph, with weights of type T, implemented as a std::map mapping undirected edges to weights. More...
 

Typedefs

typedef boost::mt19937 _rnd_gen_type
 Type of global random number generator. More...
 
typedef DAIAlg< FactorGraphDAIAlgFG
 Base class for inference algorithms that operate on a FactorGraph. More...
 
typedef DAIAlg< RegionGraphDAIAlgRG
 Base class for inference algorithms that operate on a RegionGraph. More...
 
typedef DAIAlg< CobwebGraphDAIAlgCG
 Base class for GLC that operates on CobwebGraph. More...
 
typedef TFactor< RealFactor
 Represents a factor with values of type dai::Real. More...
 
typedef std::vector< NeighborNeighbors
 Describes the set of neighbors of some node in a graph. More...
 
typedef std::pair< size_t, size_t > Edge
 Represents an edge in a graph: an Edge(i,j) corresponds to the edge between node i and node j. More...
 
typedef TProb< RealProb
 Represents a vector with entries of type dai::Real. More...
 
typedef std::string PropertyKey
 Type of the key of a Property. More...
 
typedef boost::any PropertyValue
 Type of the value of a Property. More...
 
typedef std::pair< PropertyKey, PropertyValueProperty
 A Property is a pair of a key and a corresponding value. More...
 
typedef double Real
 Real number (alias for double, which could be changed to long double if necessary) More...
 
typedef mpz_class BigInt
 Arbitrary precision integer number. More...
 

Enumerations

enum  ProbNormType { NORMPROB, NORMLINF }
 Enumerates different ways of normalizing a probability measure. More...
 
enum  ProbDistType {
  DISTL1, DISTLINF, DISTTV, DISTKL,
  DISTHEL
}
 Enumerates different distance measures between probability measures. More...
 

Functions

std::map< std::string, InfAlg * > & builtinInfAlgs ()
 Returns a map that contains for each built-in inference algorithm its name and a pointer to an object of that type. More...
 
std::set< std::string > builtinInfAlgNames ()
 Returns a set of names of all available inference algorithms. More...
 
InfAlgnewInfAlg (const std::string &name, const FactorGraph &fg, const PropertySet &opts)
 Constructs a new inference algorithm. More...
 
InfAlgnewInfAlgFromString (const std::string &nameOpts, const FactorGraph &fg)
 Constructs a new inference algorithm. More...
 
InfAlgnewInfAlgFromString (const std::string &nameOpts, const FactorGraph &fg, const std::map< std::string, std::string > &aliases)
 Constructs a new inference algorithm. More...
 
std::pair< std::string, PropertySetparseNameProperties (const std::string &s)
 Extracts the name and property set from a string s in the format "name[key1=val1,key2=val2,...]" or "name". More...
 
std::pair< std::string, PropertySetparseNameProperties (const std::string &s, const std::map< std::string, std::string > &aliases)
 Extracts the name and property set from a string s in the format "name[key1=val1,key2=val2,...]" or "name", performing alias substitution. More...
 
std::map< std::string, std::string > readAliasesFile (const std::string &filename)
 Reads aliases from file named filename. More...
 
size_t getFactorEntryForState (const FactorGraph &fg, size_t I, const vector< size_t > &state)
 Returns the entry of the I'th factor corresponding to a global state. More...
 
Real numericBBPTest (const InfAlg &bp, const std::vector< size_t > *state, const PropertySet &bbp_props, const BBPCostFunction &cfn, Real h)
 
vector< size_t > complement (vector< size_t > &xis, size_t n_states)
 Given a sorted vector of states xis and total state count n_states, return a vector of states not in xis. More...
 
Real unSoftMax (Real a, Real b)
 Computes $\frac{\exp(a)}{\exp(a)+\exp(b)}$. More...
 
Real logSumExp (Real a, Real b)
 Computes log of sum of exponents, i.e., $\log\left(\exp(a) + \exp(b)\right)$. More...
 
Real dist (const vector< Factor > &b1, const vector< Factor > &b2, size_t nv)
 Compute sum of pairwise L-infinity distances of the first nv factors in each vector. More...
 
static vector< FactormixBeliefs (Real p, const vector< Factor > &b, const vector< Factor > &c)
 Calculates a vector of mixtures p * b + (1-p) * c. More...
 
std::pair< size_t, size_t > BBPFindClampVar (const InfAlg &in_bp, bool clampingVar, const PropertySet &bbp_props, const BBPCostFunction &cfn, Real *maxVarOut)
 
size_t eliminationCost_MinNeighbors (const ClusterGraph &cl, size_t i)
 Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinNeighbors" criterion. More...
 
size_t eliminationCost_MinWeight (const ClusterGraph &cl, size_t i)
 Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinWeight" criterion. More...
 
size_t eliminationCost_MinFill (const ClusterGraph &cl, size_t i)
 Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinFill" criterion. More...
 
size_t eliminationCost_WeightedMinFill (const ClusterGraph &cl, size_t i)
 Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "WeightedMinFill" criterion. More...
 
FactorGraph createFG (const GraphAL &G, FactorType ft, size_t states, const PropertySet &props)
 Creates a factor graph from a pairwise interactions graph. More...
 
FactorGraph createHOIFG (size_t N, size_t M, size_t k, Real beta)
 Return a random factor graph with higher-order interactions. More...
 
BipartiteGraph createRandomBipartiteGraph (size_t N1, size_t N2, size_t d1, size_t d2)
 Creates a regular random bipartite graph. More...
 
int powmod (int x, int n, int p)
 Returns x**n % p, assuming p is prime. More...
 
size_t order (int x, int p)
 Returns order of x in GF(p) with p prime. More...
 
bool isPrime (size_t n)
 Returns whether n is a prime number. More...
 
BipartiteGraph createSmallLDPCGraph ()
 Constructs a regular LDPC graph with N=6, j=2, K=4, k=3. More...
 
BipartiteGraph createGroupStructuredLDPCGraph (size_t p, size_t j, size_t k)
 Creates group-structured LDPC code. More...
 
void createParityCheck (Real *result, size_t n, Real eps)
 
 DAI_ENUM (LDPCType, SMALL, GROUP, RANDOM)
 Possible LDPC structures. More...
 
Factor calcMarginal (const InfAlg &obj, const VarSet &vs, bool reInit)
 Calculates the marginal probability distribution for vs using inference algorithm obj. More...
 
vector< FactorcalcPairBeliefs (const InfAlg &obj, const VarSet &vs, bool reInit, bool accurate=false)
 Calculates beliefs for all pairs of variables in vs using inference algorithm obj. More...
 
std::vector< size_t > findMaximum (const InfAlg &obj)
 Calculates the joint state of all variables that has maximum probability, according to the inference algorithm obj. More...
 
Factor createFactorIsing (const Var &x, Real h)
 Returns a binary unnormalized single-variable factor $ \exp(hx) $ where $ x = \pm 1 $. More...
 
Factor createFactorIsing (const Var &x1, const Var &x2, Real J)
 Returns a binary unnormalized pairwise factor $ \exp(J x_1 x_2) $ where $ x_1, x_2 = \pm 1 $. More...
 
Factor createFactorExpGauss (const VarSet &vs, Real beta)
 Returns a random factor on the variables vs with strength beta. More...
 
Factor createFactorPotts (const Var &x1, const Var &x2, Real J)
 Returns a pairwise Potts factor $ \exp( J \delta_{x_1, x_2} ) $. More...
 
Factor createFactorDelta (const Var &v, size_t state)
 Returns a Kronecker delta point mass. More...
 
Factor createFactorDelta (const VarSet &vs, size_t state)
 Returns a Kronecker delta point mass. More...
 
std::ostream & operator<< (std::ostream &os, const FactorGraph &fg)
 Writes a FactorGraph to an output stream. More...
 
std::istream & operator>> (std::istream &is, FactorGraph &fg)
 Reads a FactorGraph from an input stream. More...
 
GraphAL createGraphFull (size_t N)
 Creates a fully-connected graph with N nodes. More...
 
GraphAL createGraphGrid (size_t N1, size_t N2, bool periodic)
 Creates a two-dimensional rectangular grid of N1 by N2 nodes, which can be periodic. More...
 
GraphAL createGraphGrid3D (size_t N1, size_t N2, size_t N3, bool periodic)
 Creates a three-dimensional rectangular grid of N1 by N2 by N3 nodes, which can be periodic. More...
 
GraphAL createGraphLoop (size_t N)
 Creates a graph consisting of a single loop of N nodes. More...
 
GraphAL createGraphTree (size_t N)
 Creates a random tree-structured graph of N nodes. More...
 
GraphAL createGraphRegular (size_t N, size_t d)
 Creates a random regular graph of N nodes with uniform connectivity d. More...
 
template<class T >
TFactor< T > & makePositive (TFactor< T > &f, T epsilon)
 Sets factor entries that lie between 0 and epsilon to epsilon. More...
 
template<class T >
TFactor< T > & makeZero (TFactor< T > &f, T epsilon)
 Sets factor entries that are smaller (in absolute value) than epsilon to 0. More...
 
void ReadUaiAieFactorGraphFile (const char *filename, size_t verbose, std::vector< Var > &vars, std::vector< Factor > &factors, std::vector< Permute > &permutations)
 Reads factor graph (as a pair of a variable vector and factor vector) from a file in the UAI approximate inference challenge format. More...
 
std::vector< std::map< size_t, size_t > > ReadUaiAieEvidenceFile (const char *filename, size_t verbose)
 Reads evidence (a mapping from observed variable labels to the observed values) from a file in the UAI approximate inference challenge format. More...
 
std::pair< size_t, BigIntboundTreewidth (const FactorGraph &fg, greedyVariableElimination::eliminationCostFunction fn, size_t maxStates)
 
std::ostream & operator<< (std::ostream &os, const Property &p)
 Writes a Property object (key-value pair) to an output stream. More...
 
std::ostream & operator<< (std::ostream &os, const PropertySet &ps)
 Writes a PropertySet object to an output stream. More...
 
std::istream & operator>> (std::istream &is, PropertySet &ps)
 Reads a PropertySet object from an input stream, storing values as strings. More...
 
ostream & operator<< (ostream &os, const RegionGraph &rg)
 Send RegionGraph to output stream. More...
 
bool isnan (Real x)
 Returns true if argument is NAN (Not A Number) More...
 
double toc ()
 Returns wall clock time in seconds. More...
 
_rnd_gen_type _rnd_gen (42U)
 Global random number generator. More...
 
boost::uniform_real< Real_uni_dist (0, 1)
 Uniform distribution with values between 0 and 1 (0 inclusive, 1 exclusive). More...
 
boost::variate_generator< _rnd_gen_type &, boost::uniform_real< Real > > _uni_rnd (_rnd_gen, _uni_dist)
 Global uniform random random number. More...
 
boost::variate_generator< _rnd_gen_type &, boost::normal_distribution< Real > > _normal_rnd (_rnd_gen, _normal_dist)
 Global random number generator with standard normal distribution. More...
 
void rnd_seed (size_t seed)
 Sets the random seed. More...
 
Real rnd_uniform ()
 Returns a real number, distributed uniformly on [0,1) More...
 
Real rnd_stdnormal ()
 Returns a real number from a standard-normal distribution. More...
 
int rnd_int (int min, int max)
 Returns a random integer in interval [min, max]. More...
 
std::vector< std::string > tokenizeString (const std::string &s, bool singleDelim, const std::string &delim="\t\n")
 Split a string into tokens delimited by one of the characters in delim. More...
 
size_t calcLinearState (const VarSet &vs, const std::map< Var, size_t > &state)
 Calculates the linear index in the Cartesian product of the variables in vs that corresponds to a particular joint assignment of the variables, specified by state. More...
 
std::map< Var, size_t > calcState (const VarSet &vs, size_t linearState)
 Calculates the joint assignment of the variables in vs corresponding to the linear index linearState. More...
 
mxArray * Factors2mx (const std::vector< Factor > &Ps)
 Convert vector<Factor> structure to a cell vector of CPTAB-like structs. More...
 
vector< Factormx2Factors (const mxArray *psi, long verbose)
 Convert cell vector of CPTAB-like structs to vector<Factor> More...
 
Factor mx2Factor (const mxArray *psi)
 Convert CPTAB-like struct to Factor. More...
 
 DAI_ENUM (BBPCostFunctionBase, CFN_GIBBS_B, CFN_GIBBS_B2, CFN_GIBBS_EXP, CFN_GIBBS_B_FACTOR, CFN_GIBBS_B2_FACTOR, CFN_GIBBS_EXP_FACTOR, CFN_VAR_ENT, CFN_FACTOR_ENT, CFN_BETHE_ENT)
 Enumeration of several cost functions that can be used with BBP. More...
 
 DAI_ENUM (FactorType, ISINGGAUSS, ISINGUNIFORM, EXPGAUSS, POTTS)
 Possible factor types. More...
 
size_t BigInt_size_t (const BigInt &N)
 Safe down-cast of big integer to size_t. More...
 
Real log (Real x)
 Returns logarithm of x. More...
 
Real log0 (Real x)
 Returns logarithm of x, or 0 if x == 0. More...
 
Real exp (Real x)
 Returns exponent of x. More...
 
Real pow (Real x, Real y)
 Returns x to the power y. More...
 
template<class T >
abs (const T &t)
 Returns absolute value of t. More...
 
int rnd (int n)
 Returns a random integer in the half-open interval [0, n) More...
 
template<class T >
std::string toString (const T &x)
 Converts a variable of type T to a std::string by using a boost::lexical_cast. More...
 
template<class T >
fromString (const std::string &x)
 Converts a variable of type std::string to T by using a boost::lexical_cast. More...
 
template<class T >
std::ostream & operator<< (std::ostream &os, const std::vector< T > &x)
 Writes a std::vector<> to a std::ostream. More...
 
template<class T >
std::ostream & operator<< (std::ostream &os, const std::set< T > &x)
 Writes a std::set<> to a std::ostream. More...
 
template<class T1 , class T2 >
std::ostream & operator<< (std::ostream &os, const std::map< T1, T2 > &x)
 Writes a std::map<> to a std::ostream. More...
 
template<class T1 , class T2 >
std::ostream & operator<< (std::ostream &os, const std::pair< T1, T2 > &x)
 Writes a std::pair<> to a std::ostream. More...
 
template<class T >
std::vector< T > concat (const std::vector< T > &u, const std::vector< T > &v)
 Concatenates two vectors. More...
 
template<typename T >
RootedTree MinSpanningTree (const WeightedGraph< T > &G, bool usePrim)
 Constructs a minimum spanning tree from the (non-negatively) weighted graph G. More...
 
template<typename T >
RootedTree MaxSpanningTree (const WeightedGraph< T > &G, bool usePrim)
 Constructs a minimum spanning tree from the (non-negatively) weighted graph G. More...
 

Variables

static _builtinInfAlgs allBuiltinInfAlgs
 
const char * FULL_TYPE = "FULL"
 Predefined names of various factor graph types. More...
 
const char * DREG_TYPE = "DREG"
 
const char * LOOP_TYPE = "LOOP"
 
const char * TREE_TYPE = "TREE"
 
const char * GRID_TYPE = "GRID"
 
const char * GRID3D_TYPE = "GRID3D"
 
const char * HOI_TYPE = "HOI"
 
const char * LDPC_TYPE = "LDPC"
 
boost::normal_distribution< Real_normal_dist
 Normal distribution with mean 0 and standard deviation 1. More...
 

Detailed Description

Namespace for libDAI.

Typedef Documentation

typedef boost::mt19937 dai::_rnd_gen_type

Type of global random number generator.

Base class for inference algorithms that operate on a FactorGraph.

Base class for inference algorithms that operate on a RegionGraph.

Base class for GLC that operates on CobwebGraph.

Represents a factor with values of type dai::Real.

Examples:
example_sprinkler.cpp, and uai2010-aie-solver.cpp.
typedef std::vector<Neighbor> dai::Neighbors

Describes the set of neighbors of some node in a graph.

typedef std::pair<size_t,size_t> dai::Edge

Represents an edge in a graph: an Edge(i,j) corresponds to the edge between node i and node j.

Note
If the edge is interpreted as a directed edge, then it points from i to j.
If the edge is part of a bipartite graph, i is understood to correspond to a node of type 1, and j to a node of type 2.
Examples:
example_bipgraph.cpp.
typedef TProb<Real> dai::Prob

Represents a vector with entries of type dai::Real.

typedef std::string dai::PropertyKey

Type of the key of a Property.

typedef boost::any dai::PropertyValue

Type of the value of a Property.

A Property is a pair of a key and a corresponding value.

typedef double dai::Real

Real number (alias for double, which could be changed to long double if necessary)

typedef mpz_class dai::BigInt

Arbitrary precision integer number.

Enumeration Type Documentation

Enumerates different ways of normalizing a probability measure.

  • NORMPROB means that the sum of all entries should be 1;
  • NORMLINF means that the maximum absolute value of all entries should be 1.

Enumerates different distance measures between probability measures.

  • DISTL1 is the $\ell_1$ distance (sum of absolute values of pointwise difference);
  • DISTLINF is the $\ell_\infty$ distance (maximum absolute value of pointwise difference);
  • DISTTV is the total variation distance (half of the $\ell_1$ distance);
  • DISTKL is the Kullback-Leibler distance ( $\sum_i p_i (\log p_i - \log q_i)$).
  • DISTHEL is the Hellinger distance ( $\frac{1}{2}\sum_i (\sqrt{p_i}-\sqrt{q_i})^2$).

Function Documentation

std::map< std::string, InfAlg * > & dai::builtinInfAlgs ( )

Returns a map that contains for each built-in inference algorithm its name and a pointer to an object of that type.

This functionality is obsolete and will be removed in future versions of libDAI

std::set< std::string > dai::builtinInfAlgNames ( )

Returns a set of names of all available inference algorithms.

Examples:
example.cpp.
InfAlg * dai::newInfAlg ( const std::string &  name,
const FactorGraph fg,
const PropertySet opts 
)

Constructs a new inference algorithm.

Parameters
nameThe name of the inference algorithm.
fgThe FactorGraph that the algorithm should be applied to.
optsA PropertySet specifying the options for the algorithm.
Returns
Returns a pointer to the new InfAlg object; it is the responsibility of the caller to delete it later.
Exceptions
UNKNOWN_DAI_ALGORITHMif the requested name is not known/compiled in.
InfAlg * dai::newInfAlgFromString ( const std::string &  nameOpts,
const FactorGraph fg 
)

Constructs a new inference algorithm.

Parameters
nameOptsThe name and options of the inference algorithm (should be in the format "name[key1=val1,key2=val2,...,keyn=valn]").
fgThe FactorGraph that the algorithm should be applied to.
Returns
Returns a pointer to the new InfAlg object; it is the responsibility of the caller to delete it later.
Exceptions
UNKNOWN_DAI_ALGORITHMif the requested name is not known/compiled in.
Examples:
example_imagesegmentation.cpp.
InfAlg * dai::newInfAlgFromString ( const std::string &  nameOpts,
const FactorGraph fg,
const std::map< std::string, std::string > &  aliases 
)

Constructs a new inference algorithm.

Parameters
nameOptsThe name and options of the inference algorithm (should be in the format "name[key1=val1,key2=val2,...,keyn=valn]").
fgThe FactorGraph that the algorithm should be applied to.
aliasesMaps names to strings in the format "name[key1=val1,key2=val2,...,keyn=valn]"; if not empty, alias substitution will be performed when parsing nameOpts by invoking parseNameProperties(const std::string &,const std::map<std::string,std::string> &)
See also
newInfAlgFromString(const std::string &, const FactorGraph &)
std::pair< std::string, PropertySet > dai::parseNameProperties ( const std::string &  s)

Extracts the name and property set from a string s in the format "name[key1=val1,key2=val2,...]" or "name".

std::pair< std::string, PropertySet > dai::parseNameProperties ( const std::string &  s,
const std::map< std::string, std::string > &  aliases 
)

Extracts the name and property set from a string s in the format "name[key1=val1,key2=val2,...]" or "name", performing alias substitution.

Alias substitution is performed as follows: as long as name appears as a key in aliases, it is substituted by its value. Properties in s override those of the alias (in case of recursion, the "outer" properties override those of the "inner" aliases).

std::map< std::string, std::string > dai::readAliasesFile ( const std::string &  filename)

Reads aliases from file named filename.

Parameters
filenameName of the alias file
Returns
A map that maps aliases to the strings they should be substituted with.
See also
Aliases file format
size_t dai::getFactorEntryForState ( const FactorGraph fg,
size_t  I,
const vector< size_t > &  state 
)

Returns the entry of the I'th factor corresponding to a global state.

vector<size_t> dai::complement ( vector< size_t > &  xis,
size_t  n_states 
)

Given a sorted vector of states xis and total state count n_states, return a vector of states not in xis.

Real dai::unSoftMax ( Real  a,
Real  b 
)

Computes $\frac{\exp(a)}{\exp(a)+\exp(b)}$.

Real dai::logSumExp ( Real  a,
Real  b 
)

Computes log of sum of exponents, i.e., $\log\left(\exp(a) + \exp(b)\right)$.

Real dai::dist ( const vector< Factor > &  b1,
const vector< Factor > &  b2,
size_t  nv 
)

Compute sum of pairwise L-infinity distances of the first nv factors in each vector.

static vector<Factor> dai::mixBeliefs ( Real  p,
const vector< Factor > &  b,
const vector< Factor > &  c 
)
static

Calculates a vector of mixtures p * b + (1-p) * c.

size_t dai::eliminationCost_MinNeighbors ( const ClusterGraph cl,
size_t  i 
)

Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinNeighbors" criterion.

The cost is measured as "number of neigboring nodes in the current adjacency graph", where the adjacency graph has the variables as its nodes and connects nodes i1 and i2 iff i1 and i2 occur together in some common cluster.

size_t dai::eliminationCost_MinWeight ( const ClusterGraph cl,
size_t  i 
)

Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinWeight" criterion.

The cost is measured as "product of weights of neighboring nodes in the current adjacency graph", where the adjacency graph has the variables as its nodes and connects nodes i1 and i2 iff i1 and i2 occur together in some common cluster. The weight of a node is the number of states of the corresponding variable.

size_t dai::eliminationCost_MinFill ( const ClusterGraph cl,
size_t  i 
)

Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "MinFill" criterion.

The cost is measured as "number of added edges in the adjacency graph", where the adjacency graph has the variables as its nodes and connects nodes i1 and i2 iff i1 and i2 occur together in some common cluster.

Examples:
example.cpp.
size_t dai::eliminationCost_WeightedMinFill ( const ClusterGraph cl,
size_t  i 
)

Calculates cost of eliminating the i 'th variable from cluster graph cl according to the "WeightedMinFill" criterion.

The cost is measured as "total weight of added edges in the adjacency graph", where the adjacency graph has the variables as its nodes and connects nodes i1 and i2 iff i1 and i2 occur together in some common cluster. The weight of an edge is the product of the number of states of the variables corresponding with its nodes.

FactorGraph dai::createFG ( const GraphAL G,
FactorType  ft,
size_t  states,
const PropertySet props 
)

Creates a factor graph from a pairwise interactions graph.

Parameters
GGraph describing interactions between variables
ftType of factors to use for interactions
statesNumber of states of the variables
propsAdditional properties for generating the interactions
FactorGraph dai::createHOIFG ( size_t  N,
size_t  M,
size_t  k,
Real  beta 
)

Return a random factor graph with higher-order interactions.

Parameters
Nnumber of variables
Mnumber of factors
knumber of variables that each factor depends on
betastandard-deviation of Gaussian log-factor entries
BipartiteGraph dai::createRandomBipartiteGraph ( size_t  N1,
size_t  N2,
size_t  d1,
size_t  d2 
)

Creates a regular random bipartite graph.

Parameters
N1= number of nodes of type 1
d1= size of neighborhoods of nodes of type 1
N2= number of nodes of type 2
d2= size of neighborhoods of nodes of type 2
Note
asserts that N1 * d1 == N2 * d2
int dai::powmod ( int  x,
int  n,
int  p 
)

Returns x**n % p, assuming p is prime.

size_t dai::order ( int  x,
int  p 
)

Returns order of x in GF(p) with p prime.

bool dai::isPrime ( size_t  n)

Returns whether n is a prime number.

BipartiteGraph dai::createSmallLDPCGraph ( )

Constructs a regular LDPC graph with N=6, j=2, K=4, k=3.

BipartiteGraph dai::createGroupStructuredLDPCGraph ( size_t  p,
size_t  j,
size_t  k 
)

Creates group-structured LDPC code.

Use construction described in "A Class of Group-Structured LDPC Codes" by R. M. Tanner, D. Sridhara and T. Fuja Proceedings of ICSTA, 2001

Example parameters: (p,j,k) = (31,3,5) (p,j,k) = (37,3,4) (p,j,k) = (7,2,4) (p,j,k) = (29,2,4)

j and k must be divisors of p-1

dai::DAI_ENUM ( LDPCType  ,
SMALL  ,
GROUP  ,
RANDOM   
)

Possible LDPC structures.

Factor dai::calcMarginal ( const InfAlg obj,
const VarSet vs,
bool  reInit 
)

Calculates the marginal probability distribution for vs using inference algorithm obj.

calcMarginal() works by clamping all variables in vs and calculating the partition sum for each clamped state. Therefore, it can be used in combination with any inference algorithm that can calculate/approximate partition sums.

Parameters
objinstance of inference algorithm to be used
vsvariables for which the marginal should be calculated
reInitshould be set to true if at least one of the possible clamped states would be invalid (leading to a factor graph with zero partition sum).
std::vector< Factor > dai::calcPairBeliefs ( const InfAlg obj,
const VarSet vs,
bool  reInit,
bool  accurate = false 
)

Calculates beliefs for all pairs of variables in vs using inference algorithm obj.

calcPairBeliefs() works by

  • clamping single variables in vs and calculating the partition sum and the single variable beliefs for each clamped state, if accurate == false;
  • clamping pairs of variables in vs and calculating the partition sum for each clamped state, if accurate == true.

Therefore, it can be used in combination with any inference algorithm that can calculate/approximate partition sums (and single variable beliefs, if accurate == true).

Parameters
objinstance of inference algorithm to be used
vsvariables for which the pair beliefs should be calculated
reInitshould be set to true if at least one of the possible clamped states would be invalid (leading to a factor graph with zero partition sum).
accurateif true, uses a slower but more accurate approximation algorithm
std::vector< size_t > dai::findMaximum ( const InfAlg obj)

Calculates the joint state of all variables that has maximum probability, according to the inference algorithm obj.

Note
Before this method is called, obj.run() should have been called.
Factor dai::createFactorIsing ( const Var x,
Real  h 
)

Returns a binary unnormalized single-variable factor $ \exp(hx) $ where $ x = \pm 1 $.

Parameters
xVariable (should be binary)
hField strength
Examples:
example_imagesegmentation.cpp.
Factor dai::createFactorIsing ( const Var x1,
const Var x2,
Real  J 
)

Returns a binary unnormalized pairwise factor $ \exp(J x_1 x_2) $ where $ x_1, x_2 = \pm 1 $.

Parameters
x1First variable (should be binary)
x2Second variable (should be binary)
JCoupling strength
Factor dai::createFactorExpGauss ( const VarSet vs,
Real  beta 
)

Returns a random factor on the variables vs with strength beta.

Each entry are set by drawing a normally distributed random with mean 0 and standard-deviation beta, and taking its exponent.

Parameters
vsVariables
betaFactor strength (inverse temperature)
Factor dai::createFactorPotts ( const Var x1,
const Var x2,
Real  J 
)

Returns a pairwise Potts factor $ \exp( J \delta_{x_1, x_2} ) $.

Parameters
x1First variable
x2Second variable (should have the same number of states as x1)
JFactor strength
Factor dai::createFactorDelta ( const Var v,
size_t  state 
)

Returns a Kronecker delta point mass.

Parameters
vVariable
stateThe state of v that should get value 1
Factor dai::createFactorDelta ( const VarSet vs,
size_t  state 
)

Returns a Kronecker delta point mass.

Parameters
vsSet of variables
stateThe state of vs that should get value 1
std::ostream& dai::operator<< ( std::ostream &  os,
const FactorGraph fg 
)

Writes a FactorGraph to an output stream.

Writes a factor graph to an output stream.

See also
Factor graph (.fg) file format
std::istream& dai::operator>> ( std::istream &  is,
FactorGraph fg 
)

Reads a FactorGraph from an input stream.

Reads a factor graph from an input stream.

See also
Factor graph (.fg) file format
Exceptions
INVALID_FACTORGRAPH_FILEif the input stream is not valid
GraphAL dai::createGraphFull ( size_t  N)

Creates a fully-connected graph with N nodes.

GraphAL dai::createGraphGrid ( size_t  N1,
size_t  N2,
bool  periodic 
)

Creates a two-dimensional rectangular grid of N1 by N2 nodes, which can be periodic.

GraphAL dai::createGraphGrid3D ( size_t  N1,
size_t  N2,
size_t  N3,
bool  periodic 
)

Creates a three-dimensional rectangular grid of N1 by N2 by N3 nodes, which can be periodic.

GraphAL dai::createGraphLoop ( size_t  N)

Creates a graph consisting of a single loop of N nodes.

GraphAL dai::createGraphTree ( size_t  N)

Creates a random tree-structured graph of N nodes.

GraphAL dai::createGraphRegular ( size_t  N,
size_t  d 
)

Creates a random regular graph of N nodes with uniform connectivity d.

Algorithm 1 in [StW99]. Draws a random graph of size N and uniform degree d from an almost uniform probability distribution over these graphs (which becomes uniform in the limit that d is small and N goes to infinity).

template<class T >
TFactor<T>& dai::makePositive ( TFactor< T > &  f,
epsilon 
)

Sets factor entries that lie between 0 and epsilon to epsilon.

template<class T >
TFactor<T>& dai::makeZero ( TFactor< T > &  f,
epsilon 
)

Sets factor entries that are smaller (in absolute value) than epsilon to 0.

void dai::ReadUaiAieFactorGraphFile ( const char *  filename,
size_t  verbose,
std::vector< Var > &  vars,
std::vector< Factor > &  factors,
std::vector< Permute > &  permutations 
)

Reads factor graph (as a pair of a variable vector and factor vector) from a file in the UAI approximate inference challenge format.

Parameters
[in]filenameThe filename (usually ends with ".uai")
[in]verboseThe amount of output sent to cout
[out]varsArray of variables read from the file
[out]factorsArray of factors read from the file
[out]permutationsArray of permutations, which permute between libDAI canonical ordering and ordering specified by the file
See also
http://www.cs.huji.ac.il/project/UAI10 and http://graphmod.ics.uci.edu/uai08
std::vector< std::map< size_t, size_t > > dai::ReadUaiAieEvidenceFile ( const char *  filename,
size_t  verbose 
)

Reads evidence (a mapping from observed variable labels to the observed values) from a file in the UAI approximate inference challenge format.

Parameters
[in]filenameThe filename (usually ends with ".uai.evid")
[in]verboseThe amount of output sent to cout
See also
http://www.cs.huji.ac.il/project/UAI10 and http://graphmod.ics.uci.edu/uai08
std::ostream & dai::operator<< ( std::ostream &  os,
const Property p 
)

Writes a Property object (key-value pair) to an output stream.

Note
Not all value types are automatically supported; if a type is unknown, an UNKNOWN_PROPERTY_TYPE exception is thrown. Adding support for a new type can be done by changing this function body.
std::ostream& dai::operator<< ( std::ostream &  os,
const PropertySet ps 
)

Writes a PropertySet object to an output stream.

It uses the format "[key1=val1,key2=val2,...,keyn=valn]".

Note
Only a subset of all possible types is supported (see the implementation of this function). Adding support for more types has to be done by hand.
Exceptions
UNKNOWN_PROPERTY_TYPEif the type of a property value is not supported.
std::istream& dai::operator>> ( std::istream &  is,
PropertySet ps 
)

Reads a PropertySet object from an input stream, storing values as strings.

Reads a PropertySet object from an input stream.

It expects a string in the format "[key1=val1,key2=val2,...,keyn=valn]". Values are stored as strings.

Exceptions
MALFORMED_PROPERTYif the string is not in the expected format
ostream& dai::operator<< ( std::ostream &  os,
const RegionGraph rg 
)

Send RegionGraph to output stream.

Writes a region graph to an output stream.

bool dai::isnan ( Real  x)

Returns true if argument is NAN (Not A Number)

Examples:
uai2010-aie-solver.cpp.
double dai::toc ( )

Returns wall clock time in seconds.

_rnd_gen_type dai::_rnd_gen ( 42U  )

Global random number generator.

boost::uniform_real<Real> dai::_uni_dist ( ,
 
)

Uniform distribution with values between 0 and 1 (0 inclusive, 1 exclusive).

boost::variate_generator<_rnd_gen_type&, boost::uniform_real<Real> > dai::_uni_rnd ( _rnd_gen  ,
_uni_dist   
)

Global uniform random random number.

boost::variate_generator<_rnd_gen_type&, boost::normal_distribution<Real> > dai::_normal_rnd ( _rnd_gen  ,
_normal_dist   
)

Global random number generator with standard normal distribution.

void dai::rnd_seed ( size_t  seed)

Sets the random seed.

Real dai::rnd_uniform ( )

Returns a real number, distributed uniformly on [0,1)

Real dai::rnd_stdnormal ( )

Returns a real number from a standard-normal distribution.

int dai::rnd_int ( int  min,
int  max 
)

Returns a random integer in interval [min, max].

std::vector< std::string > dai::tokenizeString ( const std::string &  s,
bool  singleDelim,
const std::string &  delim = "\t\n" 
)

Split a string into tokens delimited by one of the characters in delim.

Parameters
sthe string to be split into tokens
singleDelimif true, any single delimiter forms a boundary between two tokens; if false, a maximal group of consecutive delimiters forms a boundary between two tokens
delimdelimiter characters
size_t dai::calcLinearState ( const VarSet vs,
const std::map< Var, size_t > &  state 
)

Calculates the linear index in the Cartesian product of the variables in vs that corresponds to a particular joint assignment of the variables, specified by state.

Parameters
vsSet of variables for which the linear state should be calculated;
stateSpecifies the states of some variables.
Returns
The linear index in the Cartesian product of the variables in vs corresponding with the joint assignment specified by state, where variables for which no state is specified are assumed to be in state 0.

The linear index is calculated as follows. The variables in vs are ordered according to their label (in ascending order); say vs corresponds with the set $\{x_{l(0)},x_{l(1)},\dots,x_{l(n-1)}\}$ with $l(0) < l(1) < \dots < l(n-1)$, where variable $x_l$ has label l. Denote by $S_l$ the number of possible values ("states") of variable $x_l$. The argument state corresponds with a mapping s that assigns to each variable $x_l$ a state $s(x_l) \in \{0,1,\dots,S_l-1\}$, where $s(x_l)=0$ if $x_l$ is not specified in state. The linear index $S$ corresponding with state is now calculated by:

\begin{eqnarray*} S &:=& \sum_{i=0}^{n-1} s(x_{l(i)}) \prod_{j=0}^{i-1} S_{l(j)} \\ &= & s(x_{l(0)}) + s(x_{l(1)}) S_{l(0)} + s(x_{l(2)}) S_{l(0)} S_{l(1)} + \dots + s(x_{l(n-1)}) S_{l(0)} \cdots S_{l(n-2)}. \end{eqnarray*}

Note
If vs corresponds with $\{x_l\}_{l\in L}$, and state specifies a state for each variable $x_l$ for $l\in L$, calcLinearState() induces a mapping $\sigma : \prod_{l\in L} X_l \to \{0,1,\dots,\prod_{l\in L} S_l-1\}$ that maps a joint state to a linear index; this is the inverse of the mapping $\sigma^{-1}$ induced by calcState().
See also
calcState()
std::map< Var, size_t > dai::calcState ( const VarSet vs,
size_t  linearState 
)

Calculates the joint assignment of the variables in vs corresponding to the linear index linearState.

Parameters
vsSet of variables to which linearState refers
linearStateshould be smaller than vs.nrStates().
Returns
A mapping $s$ that maps each Var $x_l$ in vs to its state $s(x_l)$, as specified by linearState.

The variables in vs are ordered according to their label (in ascending order); say vs corresponds with the set $\{x_{l(0)},x_{l(1)},\dots,x_{l(n-1)}\}$ with $l(0) < l(1) < \dots < l(n-1)$, where variable $x_l$ has label l. Denote by $S_l$ the number of possible values ("states") of variable $x_l$ with label l. The mapping $s$ returned by this function is defined as:

\begin{eqnarray*} s(x_{l(i)}) = \left\lfloor\frac{S \mbox { mod } \prod_{j=0}^{i} S_{l(j)}}{\prod_{j=0}^{i-1} S_{l(j)}}\right\rfloor \qquad \mbox{for all $i=0,\dots,n-1$}. \end{eqnarray*}

where $S$ denotes the value of linearState.

Note
If vs corresponds with $\{x_l\}_{l\in L}$, calcState() induces a mapping $\sigma^{-1} : \{0,1,\dots,\prod_{l\in L} S_l-1\} \to \prod_{l\in L} X_l$ that maps a linear index to a joint state; this is the inverse of the mapping $\sigma$ induced by calcLinearState().
See also
calcLinearState()
mxArray * dai::Factors2mx ( const vector< Factor > &  Ps)

Convert vector<Factor> structure to a cell vector of CPTAB-like structs.

std::vector< Factor > dai::mx2Factors ( const mxArray *  psi,
long  verbose 
)

Convert cell vector of CPTAB-like structs to vector<Factor>

Factor dai::mx2Factor ( const mxArray *  psi)

Convert CPTAB-like struct to Factor.

dai::DAI_ENUM ( BBPCostFunctionBase  ,
CFN_GIBBS_B  ,
CFN_GIBBS_B2  ,
CFN_GIBBS_EXP  ,
CFN_GIBBS_B_FACTOR  ,
CFN_GIBBS_B2_FACTOR  ,
CFN_GIBBS_EXP_FACTOR  ,
CFN_VAR_ENT  ,
CFN_FACTOR_ENT  ,
CFN_BETHE_ENT   
)

Enumeration of several cost functions that can be used with BBP.

Note
This class is meant as a base class for BBPCostFunction, which provides additional functionality.
dai::DAI_ENUM ( FactorType  ,
ISINGGAUSS  ,
ISINGUNIFORM  ,
EXPGAUSS  ,
POTTS   
)

Possible factor types.

size_t dai::BigInt_size_t ( const BigInt N)
inline

Safe down-cast of big integer to size_t.

Real dai::log ( Real  x)
inline

Returns logarithm of x.

Examples:
uai2010-aie-solver.cpp.
Real dai::log0 ( Real  x)
inline

Returns logarithm of x, or 0 if x == 0.

Real dai::exp ( Real  x)
inline

Returns exponent of x.

Examples:
uai2010-aie-solver.cpp.
Real dai::pow ( Real  x,
Real  y 
)
inline

Returns x to the power y.

We use the convention that division by zero yields zero; for powers, this means that if x == 0.0 and y < 0.0, we return 0.0 instead of generating an error.

template<class T >
T dai::abs ( const T &  t)
inline

Returns absolute value of t.

int dai::rnd ( int  n)
inline

Returns a random integer in the half-open interval [0, n)

template<class T >
std::string dai::toString ( const T &  x)

Converts a variable of type T to a std::string by using a boost::lexical_cast.

Examples:
uai2010-aie-solver.cpp.
template<class T >
T dai::fromString ( const std::string &  x)

Converts a variable of type std::string to T by using a boost::lexical_cast.

template<class T >
std::ostream& dai::operator<< ( std::ostream &  os,
const std::vector< T > &  x 
)

Writes a std::vector<> to a std::ostream.

template<class T >
std::ostream& dai::operator<< ( std::ostream &  os,
const std::set< T > &  x 
)

Writes a std::set<> to a std::ostream.

template<class T1 , class T2 >
std::ostream& dai::operator<< ( std::ostream &  os,
const std::map< T1, T2 > &  x 
)

Writes a std::map<> to a std::ostream.

template<class T1 , class T2 >
std::ostream& dai::operator<< ( std::ostream &  os,
const std::pair< T1, T2 > &  x 
)

Writes a std::pair<> to a std::ostream.

template<class T >
std::vector<T> dai::concat ( const std::vector< T > &  u,
const std::vector< T > &  v 
)

Concatenates two vectors.

template<typename T >
RootedTree dai::MinSpanningTree ( const WeightedGraph< T > &  G,
bool  usePrim 
)

Constructs a minimum spanning tree from the (non-negatively) weighted graph G.

Parameters
GWeighted graph that should have non-negative weights.
usePrimIf true, use Prim's algorithm (complexity O(E log(V))), otherwise, use Kruskal's algorithm (complexity O(E log(E))).
Note
Uses implementation from Boost Graph Library.
The vertices of G must be in the range [0,N) where N is the number of vertices of G.
template<typename T >
RootedTree dai::MaxSpanningTree ( const WeightedGraph< T > &  G,
bool  usePrim 
)

Constructs a minimum spanning tree from the (non-negatively) weighted graph G.

Parameters
GWeighted graph that should have non-negative weights.
usePrimIf true, use Prim's algorithm (complexity O(E log(V))), otherwise, use Kruskal's algorithm (complexity O(E log(E))).
Note
Uses implementation from Boost Graph Library.
The vertices of G must be in the range [0,N) where N is the number of vertices of G.

Variable Documentation

const char* dai::FULL_TYPE = "FULL"

Predefined names of various factor graph types.

boost::normal_distribution<Real> dai::_normal_dist

Normal distribution with mean 0 and standard deviation 1.