root/bdm/stat/libEF.h @ 280

Revision 280, 23.9 kB (checked in by smidl, 15 years ago)

progress...

  • Property svn:eol-style set to native
Line 
1/*!
2  \file
3  \brief Probability distributions for Exponential Family models.
4  \author Vaclav Smidl.
5
6  -----------------------------------
7  BDM++ - C++ library for Bayesian Decision Making under Uncertainty
8
9  Using IT++ for numerical operations
10  -----------------------------------
11*/
12
13#ifndef EF_H
14#define EF_H
15
16
17#include "libBM.h"
18#include "../math/chmat.h"
19//#include <std>
20
21namespace bdm {
22
23
24//! Global Uniform_RNG
25extern Uniform_RNG UniRNG;
26//! Global Normal_RNG
27extern Normal_RNG NorRNG;
28//! Global Gamma_RNG
29extern Gamma_RNG GamRNG;
30
31/*!
32* \brief General conjugate exponential family posterior density.
33
34* More?...
35*/
36
37class eEF : public epdf {
38public:
39//      eEF() :epdf() {};
40        //! default constructor
41        eEF ( ) :epdf ( ) {};
42        //! logarithm of the normalizing constant, \f$\mathcal{I}\f$
43        virtual double lognc() const =0;
44        //!TODO decide if it is really needed
45        virtual void dupdate ( mat &v ) {it_error ( "Not implemented" );};
46        //!Evaluate normalized log-probability
47        virtual double evallog_nn ( const vec &val ) const{it_error ( "Not implemented" );return 0.0;};
48        //!Evaluate normalized log-probability
49        virtual double evallog ( const vec &val ) const {double tmp;tmp= evallog_nn ( val )-lognc();it_assert_debug ( std::isfinite ( tmp ),"Infinite value" ); return tmp;}
50        //!Evaluate normalized log-probability for many samples
51        virtual vec evallog ( const mat &Val ) const {
52                vec x ( Val.cols() );
53                for ( int i=0;i<Val.cols();i++ ) {x ( i ) =evallog_nn ( Val.get_col ( i ) ) ;}
54                return x-lognc();
55        }
56        //!Power of the density, used e.g. to flatten the density
57        virtual void pow ( double p ) {it_error ( "Not implemented" );};
58};
59
60/*!
61* \brief Exponential family model.
62
63* More?...
64*/
65
66class mEF : public mpdf {
67
68public:
69        //! Default constructor
70        mEF ( ) :mpdf ( ) {};
71};
72
73//! Estimator for Exponential family
74class BMEF : public BM {
75protected:
76        //! forgetting factor
77        double frg;
78        //! cached value of lognc() in the previous step (used in evaluation of \c ll )
79        double last_lognc;
80public:
81        //! Default constructor (=empty constructor)
82        BMEF ( double frg0=1.0 ) :BM ( ), frg ( frg0 ) {}
83        //! Copy constructor
84        BMEF ( const BMEF &B ) :BM ( B ), frg ( B.frg ), last_lognc ( B.last_lognc ) {}
85        //!get statistics from another model
86        virtual void set_statistics ( const BMEF* BM0 ) {it_error ( "Not implemented" );};
87        //! Weighted update of sufficient statistics (Bayes rule)
88        virtual void bayes ( const vec &data, const double w ) {};
89        //original Bayes
90        void bayes ( const vec &dt );
91        //!Flatten the posterior according to the given BMEF (of the same type!)
92        virtual void flatten ( const BMEF * B ) {it_error ( "Not implemented" );}
93        //!Flatten the posterior as if to keep nu0 data
94//      virtual void flatten ( double nu0 ) {it_error ( "Not implemented" );}
95
96        BMEF* _copy_ ( bool changerv=false ) {it_error ( "function _copy_ not implemented for this BM" ); return NULL;};
97};
98
99template<class sq_T>
100class mlnorm;
101
102/*!
103* \brief Gaussian density with positive definite (decomposed) covariance matrix.
104
105* More?...
106*/
107template<class sq_T>
108class enorm : public eEF {
109protected:
110        //! mean value
111        vec mu;
112        //! Covariance matrix in decomposed form
113        sq_T R;
114public:
115        //!\name Constructors
116        //!@{
117
118        enorm ( ) :eEF ( ), mu ( ),R ( ) {};
119        enorm ( const vec &mu,const sq_T &R ) {set_parameters ( mu,R );}
120        void set_parameters ( const vec &mu,const sq_T &R );
121        //!@}
122
123        //! \name Mathematical operations
124        //!@{
125
126        //! dupdate in exponential form (not really handy)
127        void dupdate ( mat &v,double nu=1.0 );
128
129        vec sample() const;
130        mat sample ( int N ) const;
131        double evallog_nn ( const vec &val ) const;
132        double lognc () const;
133        vec mean() const {return mu;}
134        vec variance() const {return diag ( R.to_mat() );}
135//      mlnorm<sq_T>* condition ( const RV &rvn ) const ;
136        mpdf* condition ( const RV &rvn ) const ;
137//      enorm<sq_T>* marginal ( const RV &rv ) const;
138        epdf* marginal ( const RV &rv ) const;
139        //!@}
140
141        //! \name Access to attributes
142        //!@{
143
144        vec& _mu() {return mu;}
145        void set_mu ( const vec mu0 ) { mu=mu0;}
146        sq_T& _R() {return R;}
147        const sq_T& _R() const {return R;}
148        //!@}
149
150};
151
152/*!
153* \brief Gauss-inverse-Wishart density stored in LD form
154
155* For \f$p\f$-variate densities, given rv.count() should be \f$p\times\f$ V.rows().
156*
157*/
158class egiw : public eEF {
159protected:
160        //! Extended information matrix of sufficient statistics
161        ldmat V;
162        //! Number of data records (degrees of freedom) of sufficient statistics
163        double nu;
164        //! Dimension of the output
165        int dimx;
166        //! Dimension of the regressor
167        int nPsi;
168public:
169        //!\name Constructors
170        //!@{
171        egiw() :eEF() {};
172        egiw ( int dimx0, ldmat V0, double nu0=-1.0 ) :eEF() {set_parameters ( dimx0,V0, nu0 );};
173
174        void set_parameters ( int dimx0, ldmat V0, double nu0=-1.0 ) {
175                dimx=dimx0;
176                nPsi = V0.rows()-dimx;
177                dim = dimx* ( dimx+nPsi ); // size(R) + size(Theta)
178
179                V=V0;
180                if ( nu0<0 ) {
181                        nu = 0.1 +nPsi +2*dimx +2; // +2 assures finite expected value of R
182                        // terms before that are sufficient for finite normalization
183                }
184                else {
185                        nu=nu0;
186                }
187        }
188        //!@}
189
190        vec sample() const;
191        vec mean() const;
192        vec variance() const;
193        void mean_mat ( mat &M, mat&R ) const;
194        //! In this instance, val= [theta, r]. For multivariate instances, it is stored columnwise val = [theta_1 theta_2 ... r_1 r_2 ]
195        double evallog_nn ( const vec &val ) const;
196        double lognc () const;
197        void pow ( double p ) {V*=p;nu*=p;};
198
199        //! \name Access attributes
200        //!@{
201
202        ldmat& _V() {return V;}
203        const ldmat& _V() const {return V;}
204        double& _nu()  {return nu;}
205        const double& _nu() const {return nu;}
206        //!@}
207};
208
209/*! \brief Dirichlet posterior density
210
211Continuous Dirichlet density of \f$n\f$-dimensional variable \f$x\f$
212\f[
213f(x|\beta) = \frac{\Gamma[\gamma]}{\prod_{i=1}^{n}\Gamma(\beta_i)} \prod_{i=1}^{n}x_i^{\beta_i-1}
214\f]
215where \f$\gamma=\sum_i \beta_i\f$.
216*/
217class eDirich: public eEF {
218protected:
219        //!sufficient statistics
220        vec beta;
221        //!speedup variable
222        double gamma;
223public:
224        //!\name Constructors
225        //!@{
226
227        eDirich () : eEF ( ) {};
228        eDirich ( const eDirich &D0 ) : eEF () {set_parameters ( D0.beta );};
229        eDirich ( const vec &beta0 ) {set_parameters ( beta0 );};
230        void set_parameters ( const vec &beta0 ) {
231                beta= beta0;
232                dim = beta.length();
233                gamma = sum ( beta );
234        }
235        //!@}
236
237        vec sample() const {it_error ( "Not implemented" );return vec_1 ( 0.0 );};
238        vec mean() const {return beta/gamma;};
239        vec variance() const {return elem_mult ( beta, ( beta+1 ) ) / ( gamma* ( gamma+1 ) );}
240        //! In this instance, val is ...
241        double evallog_nn ( const vec &val ) const {
242                double tmp; tmp= ( beta-1 ) *log ( val );               it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
243                return tmp;
244        };
245        double lognc () const {
246                double tmp;
247                double gam=sum ( beta );
248                double lgb=0.0;
249                for ( int i=0;i<beta.length();i++ ) {lgb+=lgamma ( beta ( i ) );}
250                tmp= lgb-lgamma ( gam );
251                it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
252                return tmp;
253        };
254        //!access function
255        vec& _beta()  {return beta;}
256        //!Set internal parameters
257};
258
259//! \brief Estimator for Multinomial density
260class multiBM : public BMEF {
261protected:
262        //! Conjugate prior and posterior
263        eDirich est;
264        //! Pointer inside est to sufficient statistics
265        vec &beta;
266public:
267        //!Default constructor
268        multiBM ( ) : BMEF ( ),est ( ),beta ( est._beta() ) {if ( beta.length() >0 ) {last_lognc=est.lognc();}else{last_lognc=0.0;}}
269        //!Copy constructor
270        multiBM ( const multiBM &B ) : BMEF ( B ),est ( B.est ),beta ( est._beta() ) {}
271        //! Sets sufficient statistics to match that of givefrom mB0
272        void set_statistics ( const BM* mB0 ) {const multiBM* mB=dynamic_cast<const multiBM*> ( mB0 ); beta=mB->beta;}
273        void bayes ( const vec &dt ) {
274                if ( frg<1.0 ) {beta*=frg;last_lognc=est.lognc();}
275                beta+=dt;
276                if ( evalll ) {ll=est.lognc()-last_lognc;}
277        }
278        double logpred ( const vec &dt ) const {
279                eDirich pred ( est );
280                vec &beta = pred._beta();
281
282                double lll;
283                if ( frg<1.0 )
284                        {beta*=frg;lll=pred.lognc();}
285                else
286                        if ( evalll ) {lll=last_lognc;}
287                        else{lll=pred.lognc();}
288
289                beta+=dt;
290                return pred.lognc()-lll;
291        }
292        void flatten ( const BMEF* B ) {
293                const multiBM* E=dynamic_cast<const multiBM*> ( B );
294                // sum(beta) should be equal to sum(B.beta)
295                const vec &Eb=E->beta;//const_cast<multiBM*> ( E )->_beta();
296                beta*= ( sum ( Eb ) /sum ( beta ) );
297                if ( evalll ) {last_lognc=est.lognc();}
298        }
299        const epdf& posterior() const {return est;};
300        const eDirich* _e() const {return &est;};
301        void set_parameters ( const vec &beta0 ) {
302                est.set_parameters ( beta0 );
303                if ( evalll ) {last_lognc=est.lognc();}
304        }
305};
306
307/*!
308 \brief Gamma posterior density
309
310 Multivariate Gamma density as product of independent univariate densities.
311 \f[
312 f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
313 \f]
314*/
315
316class egamma : public eEF {
317protected:
318        //! Vector \f$\alpha\f$
319        vec alpha;
320        //! Vector \f$\beta\f$
321        vec beta;
322public :
323        //! \name Constructors
324        //!@{
325        egamma ( ) :eEF ( ), alpha(), beta() {};
326        egamma ( const vec &a, const vec &b ) {set_parameters ( a, b );};
327        void set_parameters ( const vec &a, const vec &b ) {alpha=a,beta=b;dim = alpha.length();};
328        //!@}
329
330        vec sample() const;
331        //! TODO: is it used anywhere?
332//      mat sample ( int N ) const;
333        double evallog ( const vec &val ) const;
334        double lognc () const;
335        //! Returns poiter to alpha and beta. Potentially dengerous: use with care!
336        vec& _alpha() {return alpha;}
337        vec& _beta() {return beta;}
338        vec mean() const {return elem_div ( alpha,beta );}
339        vec variance() const {return elem_div ( alpha,elem_mult ( beta,beta ) ); }
340};
341
342/*!
343 \brief Inverse-Gamma posterior density
344
345 Multivariate inverse-Gamma density as product of independent univariate densities.
346 \f[
347 f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
348 \f]
349
350 Inverse Gamma can be converted to Gamma using \f[
351 x\sim iG(a,b) => 1/x\sim G(a,1/b)
352\f]
353This relation is used in sampling.
354 */
355
356class eigamma : public eEF {
357protected:
358        //!internal egamma
359        egamma eg;
360        //! Vector \f$\alpha\f$
361        vec &alpha;
362        //! Vector \f$\beta\f$ (in fact it is 1/beta as used in definition of iG)
363        vec &beta;
364public :
365        //! \name Constructors
366        //!@{
367        eigamma ( ) :eEF ( ), eg(),alpha ( eg._alpha() ), beta ( eg._beta() ) {};
368        eigamma ( const vec &a, const vec &b ) :eEF ( ), eg(),alpha ( eg._alpha() ), beta ( eg._beta() ) {eg.set_parameters ( a,b );};
369        void set_parameters ( const vec &a, const vec &b ) {eg.set_parameters ( a,b );};
370        //!@}
371
372        vec sample() const {return 1.0/eg.sample();};
373        //! TODO: is it used anywhere?
374//      mat sample ( int N ) const;
375        double evallog ( const vec &val ) const {return eg.evallog ( val );};
376        double lognc () const {return eg.lognc();};
377        //! Returns poiter to alpha and beta. Potentially dangerous: use with care!
378        vec mean() const {return elem_div ( beta,alpha-1 );}
379        vec variance() const {vec mea=mean(); return elem_div ( elem_mult ( mea,mea ),alpha-2 );}
380        vec& _alpha() {return alpha;}
381        vec& _beta() {return beta;}
382};
383/*
384//! Weighted mixture of epdfs with external owned components.
385class emix : public epdf {
386protected:
387        int n;
388        vec &w;
389        Array<epdf*> Coms;
390public:
391//! Default constructor
392        emix ( const RV &rv, vec &w0): epdf(rv), n(w0.length()), w(w0), Coms(n) {};
393        void set_parameters( int &i, double wi, epdf* ep){w(i)=wi;Coms(i)=ep;}
394        vec mean(){vec pom; for(int i=0;i<n;i++){pom+=Coms(i)->mean()*w(i);} return pom;};
395        vec sample() {it_error ( "Not implemented" );return 0;}
396};
397*/
398
399//!  Uniform distributed density on a rectangular support
400
401class euni: public epdf {
402protected:
403//! lower bound on support
404        vec low;
405//! upper bound on support
406        vec high;
407//! internal
408        vec distance;
409//! normalizing coefficients
410        double nk;
411//! cache of log( \c nk )
412        double lnk;
413public:
414        //! \name Constructors
415        //!@{
416        euni ( ) :epdf ( ) {}
417        euni ( const vec &low0, const vec &high0 ) {set_parameters ( low0,high0 );}
418        void set_parameters ( const vec &low0, const vec &high0 ) {
419                distance = high0-low0;
420                it_assert_debug ( min ( distance ) >0.0,"bad support" );
421                low = low0;
422                high = high0;
423                nk = prod ( 1.0/distance );
424                lnk = log ( nk );
425                dim = low.length();
426        }
427        //!@}
428
429        double eval ( const vec &val ) const  {return nk;}
430        double evallog ( const vec &val ) const  {return lnk;}
431        vec sample() const {
432                vec smp ( dim );
433#pragma omp critical
434                UniRNG.sample_vector ( dim ,smp );
435                return low+elem_mult ( distance,smp );
436        }
437        //! set values of \c low and \c high
438        vec mean() const {return ( high-low ) /2.0;}
439        vec variance() const {return ( pow ( high,2 ) +pow ( low,2 ) +elem_mult ( high,low ) ) /3.0;}
440};
441
442
443/*!
444 \brief Normal distributed linear function with linear function of mean value;
445
446 Mean value \f$mu=A*rvc+mu_0\f$.
447*/
448template<class sq_T>
449class mlnorm : public mEF {
450protected:
451        //! Internal epdf that arise by conditioning on \c rvc
452        enorm<sq_T> epdf;
453        mat A;
454        vec mu_const;
455        vec& _mu; //cached epdf.mu;
456public:
457        //! \name Constructors
458        //!@{
459        mlnorm ( ) :mEF (),epdf ( ),A ( ),_mu ( epdf._mu() ) {ep =&epdf; };
460        mlnorm ( const  mat &A, const vec &mu0, const sq_T &R ) :epdf ( ),_mu ( epdf._mu() ) {
461                ep =&epdf; set_parameters ( A,mu0,R );
462        };
463        //! Set \c A and \c R
464        void set_parameters ( const  mat &A, const vec &mu0, const sq_T &R );
465        //!@}
466        //! Set value of \c rvc . Result of this operation is stored in \c epdf use function \c _ep to access it.
467        void condition ( const vec &cond );
468
469        //!access function
470        vec& _mu_const() {return mu_const;}
471        //!access function
472        mat& _A() {return A;}
473        //!access function
474        mat _R() {return epdf._R().to_mat();}
475
476        template<class sq_M>
477        friend std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_M> &ml );
478};
479
480//! Mpdf with general function for mean value
481template<class sq_T>
482class mgnorm : public mEF {
483protected:
484        //! Internal epdf that arise by conditioning on \c rvc
485        enorm<sq_T> epdf;
486        vec &mu;
487        fnc* g;
488public:
489        //!default constructor
490        mgnorm() :mu ( epdf._mu() ) {ep=&epdf;}
491        //!set mean function
492        void set_parameters ( fnc* g0, const sq_T &R0 ) {g=g0; epdf.set_parameters ( zeros ( g->_dimy() ), R0 );}
493        void condition ( const vec &cond ) {mu=g->eval ( cond );};
494};
495
496/*! (Approximate) Student t density with linear function of mean value
497
498The internal epdf of this class is of the type of a Gaussian (enorm).
499However, each conditioning is trying to assure the best possible approximation by taking into account the zeta function. See [] for reference.
500
501Perhaps a moment-matching technique?
502*/
503class mlstudent : public mlnorm<ldmat> {
504protected:
505        ldmat Lambda;
506        ldmat &_R;
507        ldmat Re;
508public:
509        mlstudent ( ) :mlnorm<ldmat> (),
510                        Lambda (),      _R ( epdf._R() ) {}
511        void set_parameters ( const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0 ) {
512                it_assert_debug ( A0.rows() ==mu0.length(),"" );
513                it_assert_debug ( R0.rows() ==A0.rows(),"" );
514
515                epdf.set_parameters ( mu0,Lambda ); //
516                A = A0;
517                mu_const = mu0;
518                Re=R0;
519                Lambda = Lambda0;
520        }
521        void condition ( const vec &cond ) {
522                _mu = A*cond + mu_const;
523                double zeta;
524                //ugly hack!
525                if ( ( cond.length() +1 ) ==Lambda.rows() ) {
526                        zeta = Lambda.invqform ( concat ( cond, vec_1 ( 1.0 ) ) );
527                }
528                else {
529                        zeta = Lambda.invqform ( cond );
530                }
531                _R = Re;
532                _R*= ( 1+zeta );// / ( nu ); << nu is in Re!!!!!!
533        };
534
535};
536/*!
537\brief  Gamma random walk
538
539Mean value, \f$\mu\f$, of this density is given by \c rvc .
540Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
541This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
542
543The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
544*/
545class mgamma : public mEF {
546protected:
547        //! Internal epdf that arise by conditioning on \c rvc
548        egamma epdf;
549        //! Constant \f$k\f$
550        double k;
551        //! cache of epdf.beta
552        vec &_beta;
553
554public:
555        //! Constructor
556        mgamma ( ) : mEF ( ), epdf (), _beta ( epdf._beta() ) {ep=&epdf;};
557        //! Set value of \c k
558        void set_parameters ( double k, const vec &beta0 );
559        void condition ( const vec &val ) {_beta=k/val;};
560};
561
562/*!
563\brief  Inverse-Gamma random walk
564
565Mean value, \f$ \mu \f$, of this density is given by \c rvc .
566Standard deviation of the random walk is proportional to one \f$ k \f$-th the mean.
567This is achieved by setting \f$ \alpha=\mu/k^2+2 \f$ and \f$ \beta=\mu(\alpha-1)\f$.
568
569The standard deviation of the walk is then: \f$ \mu/\sqrt(k)\f$.
570 */
571class migamma : public mEF {
572protected:
573        //! Internal epdf that arise by conditioning on \c rvc
574        eigamma epdf;
575        //! Constant \f$k\f$
576        double k;
577        //! cache of epdf.alpha
578        vec &_alpha;
579        //! cache of epdf.beta
580        vec &_beta;
581
582public:
583        //! \name Constructors
584        //!@{
585        migamma ( ) : mEF (), epdf ( ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;};
586        migamma ( const migamma &m ) : mEF (), epdf ( m.epdf ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;};
587        //!@}
588
589        //! Set value of \c k
590        void set_parameters ( int len, double k0 ) {
591                k=k0;
592                epdf.set_parameters ( ( 1.0/ ( k*k ) +2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ );
593                dimc = dimension();
594        };
595        void condition ( const vec &val ) {
596                _beta=elem_mult ( val, ( _alpha-1.0 ) );
597        };
598};
599
600/*!
601\brief  Gamma random walk around a fixed point
602
603Mean value, \f$\mu\f$, of this density is given by a geometric combination of \c rvc and given fixed point, \f$p\f$. \f$l\f$ is the coefficient of the geometric combimation
604\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
605
606Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
607This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
608
609The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
610*/
611class mgamma_fix : public mgamma {
612protected:
613        //! parameter l
614        double l;
615        //! reference vector
616        vec refl;
617public:
618        //! Constructor
619        mgamma_fix ( ) : mgamma ( ),refl () {};
620        //! Set value of \c k
621        void set_parameters ( double k0 , vec ref0, double l0 ) {
622                mgamma::set_parameters ( k0, ref0 );
623                refl=pow ( ref0,1.0-l0 );l=l0;
624                dimc=dimension();
625        };
626
627        void condition ( const vec &val ) {vec mean=elem_mult ( refl,pow ( val,l ) ); _beta=k/mean;};
628};
629
630
631/*!
632\brief  Inverse-Gamma random walk around a fixed point
633
634Mean value, \f$\mu\f$, of this density is given by a geometric combination of \c rvc and given fixed point, \f$p\f$. \f$l\f$ is the coefficient of the geometric combimation
635\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
636
637==== Check == vv =
638Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
639This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
640
641The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
642 */
643class migamma_fix : public migamma {
644protected:
645        //! parameter l
646        double l;
647        //! reference vector
648        vec refl;
649public:
650        //! Constructor
651        migamma_fix ( ) : migamma (),refl ( ) {};
652        //! Set value of \c k
653        void set_parameters ( double k0 , vec ref0, double l0 ) {
654                migamma::set_parameters ( ref0.length(), k0 );
655                refl=pow ( ref0,1.0-l0 );
656                l=l0;
657                dimc = dimension();
658        };
659
660        void condition ( const vec &val ) {
661                vec mean=elem_mult ( refl,pow ( val,l ) );
662                migamma::condition ( mean );
663        };
664};
665//! Switch between various resampling methods.
666enum RESAMPLING_METHOD { MULTINOMIAL = 0, STRATIFIED = 1, SYSTEMATIC = 3 };
667/*!
668\brief Weighted empirical density
669
670Used e.g. in particle filters.
671*/
672class eEmp: public epdf {
673protected :
674        //! Number of particles
675        int n;
676        //! Sample weights \f$w\f$
677        vec w;
678        //! Samples \f$x^{(i)}, i=1..n\f$
679        Array<vec> samples;
680public:
681        //! \name Constructors
682        //!@{
683        eEmp ( ) :epdf ( ),w ( ),samples ( ) {};
684        eEmp ( const eEmp &e ) : epdf ( e ), w ( e.w ), samples ( e.samples ) {};
685        //!@}
686
687        //! Set samples and weights
688        void set_parameters ( const vec &w0, const epdf* pdf0 );
689        //! Set sample
690        void set_samples ( const epdf* pdf0 );
691        //! Set sample
692        void set_n ( int n0, bool copy=true ) {w.set_size ( n0,copy );samples.set_size ( n0,copy );};
693        //! Potentially dangerous, use with care.
694        vec& _w()  {return w;};
695        //! Potentially dangerous, use with care.
696        const vec& _w() const {return w;};
697        //! access function
698        Array<vec>& _samples() {return samples;};
699        //! access function
700        const Array<vec>& _samples() const {return samples;};
701        //! Function performs resampling, i.e. removal of low-weight samples and duplication of high-weight samples such that the new samples represent the same density.
702        ivec resample ( RESAMPLING_METHOD method = SYSTEMATIC );
703        //! inherited operation : NOT implemneted
704        vec sample() const {it_error ( "Not implemented" );return 0;}
705        //! inherited operation : NOT implemneted
706        double evallog ( const vec &val ) const {it_error ( "Not implemented" );return 0.0;}
707        vec mean() const {
708                vec pom=zeros ( dim );
709                for ( int i=0;i<n;i++ ) {pom+=samples ( i ) *w ( i );}
710                return pom;
711        }
712        vec variance() const {
713                vec pom=zeros ( dim );
714                for ( int i=0;i<n;i++ ) {pom+=pow ( samples ( i ),2 ) *w ( i );}
715                return pom-pow ( mean(),2 );
716        }
717};
718
719
720////////////////////////
721
722template<class sq_T>
723void enorm<sq_T>::set_parameters ( const vec &mu0, const sq_T &R0 ) {
724//Fixme test dimensions of mu0 and R0;
725        mu = mu0;
726        R = R0;
727        dim = mu0.length();
728};
729
730template<class sq_T>
731void enorm<sq_T>::dupdate ( mat &v, double nu ) {
732        //
733};
734
735// template<class sq_T>
736// void enorm<sq_T>::tupdate ( double phi, mat &vbar, double nubar ) {
737//      //
738// };
739
740template<class sq_T>
741vec enorm<sq_T>::sample() const {
742        vec x ( dim );
743#pragma omp critical
744        NorRNG.sample_vector ( dim,x );
745        vec smp = R.sqrt_mult ( x );
746
747        smp += mu;
748        return smp;
749};
750
751template<class sq_T>
752mat enorm<sq_T>::sample ( int N ) const {
753        mat X ( dim,N );
754        vec x ( dim );
755        vec pom;
756        int i;
757
758        for ( i=0;i<N;i++ ) {
759#pragma omp critical
760                NorRNG.sample_vector ( dim,x );
761                pom = R.sqrt_mult ( x );
762                pom +=mu;
763                X.set_col ( i, pom );
764        }
765
766        return X;
767};
768
769// template<class sq_T>
770// double enorm<sq_T>::eval ( const vec &val ) const {
771//      double pdfl,e;
772//      pdfl = evallog ( val );
773//      e = exp ( pdfl );
774//      return e;
775// };
776
777template<class sq_T>
778double enorm<sq_T>::evallog_nn ( const vec &val ) const {
779        // 1.83787706640935 = log(2pi)
780        double tmp=-0.5* ( R.invqform ( mu-val ) );// - lognc();
781        return  tmp;
782};
783
784template<class sq_T>
785inline double enorm<sq_T>::lognc () const {
786        // 1.83787706640935 = log(2pi)
787        double tmp=0.5* ( R.cols() * 1.83787706640935 +R.logdet() );
788        return tmp;
789};
790
791template<class sq_T>
792void mlnorm<sq_T>::set_parameters ( const mat &A0, const vec &mu0, const sq_T &R0 ) {
793        it_assert_debug ( A0.rows() ==mu0.length(),"" );
794        it_assert_debug ( A0.rows() ==R0.rows(),"" );
795
796        epdf.set_parameters ( zeros ( A0.rows() ),R0 );
797        A = A0;
798        mu_const = mu0;
799        dimc=A0.cols();
800}
801
802// template<class sq_T>
803// vec mlnorm<sq_T>::samplecond (const  vec &cond, double &lik ) {
804//      this->condition ( cond );
805//      vec smp = epdf.sample();
806//      lik = epdf.eval ( smp );
807//      return smp;
808// }
809
810// template<class sq_T>
811// mat mlnorm<sq_T>::samplecond (const vec &cond, vec &lik, int n ) {
812//      int i;
813//      int dim = rv.count();
814//      mat Smp ( dim,n );
815//      vec smp ( dim );
816//      this->condition ( cond );
817//
818//      for ( i=0; i<n; i++ ) {
819//              smp = epdf.sample();
820//              lik ( i ) = epdf.eval ( smp );
821//              Smp.set_col ( i ,smp );
822//      }
823//
824//      return Smp;
825// }
826
827template<class sq_T>
828void mlnorm<sq_T>::condition ( const vec &cond ) {
829        _mu = A*cond + mu_const;
830//R is already assigned;
831}
832
833template<class sq_T>
834epdf* enorm<sq_T>::marginal ( const RV &rvn ) const {
835        it_assert_debug ( isnamed(), "rv description is not assigned" );
836        ivec irvn = rvn.dataind ( rv );
837
838        sq_T Rn ( R,irvn ); //select rows and columns of R
839
840        enorm<sq_T>* tmp = new enorm<sq_T>;
841        tmp->set_rv ( rvn );
842        tmp->set_parameters ( mu ( irvn ), Rn );
843        return tmp;
844}
845
846template<class sq_T>
847mpdf* enorm<sq_T>::condition ( const RV &rvn ) const {
848
849        it_assert_debug ( isnamed(),"rvs are not assigned" );
850
851        RV rvc = rv.subt ( rvn );
852        it_assert_debug ( ( rvc._dsize() +rvn._dsize() ==rv._dsize() ),"wrong rvn" );
853        //Permutation vector of the new R
854        ivec irvn = rvn.dataind ( rv );
855        ivec irvc = rvc.dataind ( rv );
856        ivec perm=concat ( irvn , irvc );
857        sq_T Rn ( R,perm );
858
859        //fixme - could this be done in general for all sq_T?
860        mat S=Rn.to_mat();
861        //fixme
862        int n=rvn._dsize()-1;
863        int end=R.rows()-1;
864        mat S11 = S.get ( 0,n, 0, n );
865        mat S12 = S.get ( 0, n , rvn._dsize(), end );
866        mat S22 = S.get ( rvn._dsize(), end, rvn._dsize(), end );
867
868        vec mu1 = mu ( irvn );
869        vec mu2 = mu ( irvc );
870        mat A=S12*inv ( S22 );
871        sq_T R_n ( S11 - A *S12.T() );
872
873        mlnorm<sq_T>* tmp=new mlnorm<sq_T> ( );
874        tmp->set_rv ( rvn ); tmp->set_rvc ( rvc );
875        tmp->set_parameters ( A,mu1-A*mu2,R_n );
876        return tmp;
877}
878
879///////////
880
881template<class sq_T>
882std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_T> &ml ) {
883        os << "A:"<< ml.A<<endl;
884        os << "mu:"<< ml.mu_const<<endl;
885        os << "R:" << ml.epdf._R().to_mat() <<endl;
886        return os;
887};
888
889}
890#endif //EF_H
Note: See TracBrowser for help on using the browser.