root/bdm/stat/libEF.h @ 330

Revision 330, 29.9 kB (checked in by dedecius, 15 years ago)

Added methods egiw::est_theta() for LS estimate of theta and egiw::est_theta_cov() for covariance of the LS estimate.

  • 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
25//! Global Uniform_RNG
26        extern Uniform_RNG UniRNG;
27//! Global Normal_RNG
28        extern Normal_RNG NorRNG;
29//! Global Gamma_RNG
30        extern Gamma_RNG GamRNG;
31
32        /*!
33        * \brief General conjugate exponential family posterior density.
34
35        * More?...
36        */
37
38        class eEF : public epdf
39        {
40                public:
41//      eEF() :epdf() {};
42                        //! default constructor
43                        eEF ( ) :epdf ( ) {};
44                        //! logarithm of the normalizing constant, \f$\mathcal{I}\f$
45                        virtual double lognc() const =0;
46                        //!TODO decide if it is really needed
47                        virtual void dupdate ( mat &v ) {it_error ( "Not implemented" );};
48                        //!Evaluate normalized log-probability
49                        virtual double evallog_nn ( const vec &val ) const{it_error ( "Not implemented" );return 0.0;};
50                        //!Evaluate normalized log-probability
51                        virtual double evallog ( const vec &val ) const {double tmp;tmp= evallog_nn ( val )-lognc();it_assert_debug ( std::isfinite ( tmp ),"Infinite value" ); return tmp;}
52                        //!Evaluate normalized log-probability for many samples
53                        virtual vec evallog ( const mat &Val ) const
54                        {
55                                vec x ( Val.cols() );
56                                for ( int i=0;i<Val.cols();i++ ) {x ( i ) =evallog_nn ( Val.get_col ( i ) ) ;}
57                                return x-lognc();
58                        }
59                        //!Power of the density, used e.g. to flatten the density
60                        virtual void pow ( double p ) {it_error ( "Not implemented" );};
61        };
62
63        /*!
64        * \brief Exponential family model.
65
66        * More?...
67        */
68
69        class mEF : public mpdf
70        {
71
72                public:
73                        //! Default constructor
74                        mEF ( ) :mpdf ( ) {};
75        };
76
77//! Estimator for Exponential family
78        class BMEF : public BM
79        {
80                protected:
81                        //! forgetting factor
82                        double frg;
83                        //! cached value of lognc() in the previous step (used in evaluation of \c ll )
84                        double last_lognc;
85                public:
86                        //! Default constructor (=empty constructor)
87                        BMEF ( double frg0=1.0 ) :BM ( ), frg ( frg0 ) {}
88                        //! Copy constructor
89                        BMEF ( const BMEF &B ) :BM ( B ), frg ( B.frg ), last_lognc ( B.last_lognc ) {}
90                        //!get statistics from another model
91                        virtual void set_statistics ( const BMEF* BM0 ) {it_error ( "Not implemented" );};
92                        //! Weighted update of sufficient statistics (Bayes rule)
93                        virtual void bayes ( const vec &data, const double w ) {};
94                        //original Bayes
95                        void bayes ( const vec &dt );
96                        //!Flatten the posterior according to the given BMEF (of the same type!)
97                        virtual void flatten ( const BMEF * B ) {it_error ( "Not implemented" );}
98                        //!Flatten the posterior as if to keep nu0 data
99//      virtual void flatten ( double nu0 ) {it_error ( "Not implemented" );}
100
101                        BMEF* _copy_ () const {it_error ( "function _copy_ not implemented for this BM" ); return NULL;};
102        };
103
104        template<class sq_T>
105        class mlnorm;
106
107        /*!
108        * \brief Gaussian density with positive definite (decomposed) covariance matrix.
109
110        * More?...
111        */
112        template<class sq_T>
113        class enorm : public eEF
114        {
115                protected:
116                        //! mean value
117                        vec mu;
118                        //! Covariance matrix in decomposed form
119                        sq_T R;
120                public:
121                        //!\name Constructors
122                        //!@{
123
124                        enorm ( ) :eEF ( ), mu ( ),R ( ) {};
125                        enorm ( const vec &mu,const sq_T &R ) {set_parameters ( mu,R );}
126                        void set_parameters ( const vec &mu,const sq_T &R );
127                        //!@}
128
129                        //! \name Mathematical operations
130                        //!@{
131
132                        //! dupdate in exponential form (not really handy)
133                        void dupdate ( mat &v,double nu=1.0 );
134
135                        vec sample() const;
136                        mat sample ( int N ) const;
137                        double evallog_nn ( const vec &val ) const;
138                        double lognc () const;
139                        vec mean() const {return mu;}
140                        vec variance() const {return diag ( R.to_mat() );}
141//      mlnorm<sq_T>* condition ( const RV &rvn ) const ; <=========== fails to cmpile. Why?
142                        mpdf* condition ( const RV &rvn ) const ;
143        enorm<sq_T>* marginal ( const RV &rv ) const;
144//                      epdf* marginal ( const RV &rv ) const;
145                        //!@}
146
147                        //! \name Access to attributes
148                        //!@{
149
150                        vec& _mu() {return mu;}
151                        void set_mu ( const vec mu0 ) { mu=mu0;}
152                        sq_T& _R() {return R;}
153                        const sq_T& _R() const {return R;}
154                        //!@}
155
156        };
157
158        /*!
159        * \brief Gauss-inverse-Wishart density stored in LD form
160
161        * For \f$p\f$-variate densities, given rv.count() should be \f$p\times\f$ V.rows().
162        *
163        */
164        class egiw : public eEF
165        {
166                protected:
167                        //! Extended information matrix of sufficient statistics
168                        ldmat V;
169                        //! Number of data records (degrees of freedom) of sufficient statistics
170                        double nu;
171                        //! Dimension of the output
172                        int dimx;
173                        //! Dimension of the regressor
174                        int nPsi;
175                public:
176                        //!\name Constructors
177                        //!@{
178                        egiw() :eEF() {};
179                        egiw ( int dimx0, ldmat V0, double nu0=-1.0 ) :eEF() {set_parameters ( dimx0,V0, nu0 );};
180
181                        void set_parameters ( int dimx0, ldmat V0, double nu0=-1.0 )
182                        {
183                                dimx=dimx0;
184                                nPsi = V0.rows()-dimx;
185                                dim = dimx* ( dimx+nPsi ); // size(R) + size(Theta)
186
187                                V=V0;
188                                if ( nu0<0 )
189                                {
190                                        nu = 0.1 +nPsi +2*dimx +2; // +2 assures finite expected value of R
191                                        // terms before that are sufficient for finite normalization
192                                }
193                                else
194                                {
195                                        nu=nu0;
196                                }
197                        }
198                        //!@}
199
200                        vec sample() const;
201                        vec mean() const;
202                        vec variance() const;
203
204                        //! LS estimate of \f$\theta\f$
205                        vec est_theta() const;
206
207                        //! Covariance of the LS estimate
208                        ldmat est_theta_cov() const;
209
210                        void mean_mat ( mat &M, mat&R ) const;
211                        //! In this instance, val= [theta, r]. For multivariate instances, it is stored columnwise val = [theta_1 theta_2 ... r_1 r_2 ]
212                        double evallog_nn ( const vec &val ) const;
213                        double lognc () const;
214                        void pow ( double p ) {V*=p;nu*=p;};
215
216                        //! \name Access attributes
217                        //!@{
218
219                        ldmat& _V() {return V;}
220                        const ldmat& _V() const {return V;}
221                        double& _nu()  {return nu;}
222                        const double& _nu() const {return nu;}
223                        //!@}
224        };
225
226        /*! \brief Dirichlet posterior density
227
228        Continuous Dirichlet density of \f$n\f$-dimensional variable \f$x\f$
229        \f[
230        f(x|\beta) = \frac{\Gamma[\gamma]}{\prod_{i=1}^{n}\Gamma(\beta_i)} \prod_{i=1}^{n}x_i^{\beta_i-1}
231        \f]
232        where \f$\gamma=\sum_i \beta_i\f$.
233        */
234        class eDirich: public eEF
235        {
236                protected:
237                        //!sufficient statistics
238                        vec beta;
239                public:
240                        //!\name Constructors
241                        //!@{
242
243                        eDirich () : eEF ( ) {};
244                        eDirich ( const eDirich &D0 ) : eEF () {set_parameters ( D0.beta );};
245                        eDirich ( const vec &beta0 ) {set_parameters ( beta0 );};
246                        void set_parameters ( const vec &beta0 )
247                        {
248                                beta= beta0;
249                                dim = beta.length();
250                        }
251                        //!@}
252
253                        vec sample() const {it_error ( "Not implemented" );return vec_1 ( 0.0 );};
254                        vec mean() const {return beta/sum(beta);};
255                        vec variance() const {double gamma =sum(beta); return elem_mult ( beta, ( beta+1 ) ) / ( gamma* ( gamma+1 ) );}
256                        //! In this instance, val is ...
257                        double evallog_nn ( const vec &val ) const
258                        {
259                                double tmp; tmp= ( beta-1 ) *log ( val );               it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
260                                return tmp;
261                        };
262                        double lognc () const
263                        {
264                                double tmp;
265                                double gam=sum ( beta );
266                                double lgb=0.0;
267                                for ( int i=0;i<beta.length();i++ ) {lgb+=lgamma ( beta ( i ) );}
268                                tmp= lgb-lgamma ( gam );
269                                it_assert_debug ( std::isfinite ( tmp ),"Infinite value" );
270                                return tmp;
271                        };
272                        //!access function
273                        vec& _beta()  {return beta;}
274                        //!Set internal parameters
275        };
276
277//! \brief Estimator for Multinomial density
278        class multiBM : public BMEF
279        {
280                protected:
281                        //! Conjugate prior and posterior
282                        eDirich est;
283                        //! Pointer inside est to sufficient statistics
284                        vec &beta;
285                public:
286                        //!Default constructor
287                        multiBM ( ) : BMEF ( ),est ( ),beta ( est._beta() )
288                        {
289                                if ( beta.length() >0 ) {last_lognc=est.lognc();}
290                                else{last_lognc=0.0;}
291                        }
292                        //!Copy constructor
293                        multiBM ( const multiBM &B ) : BMEF ( B ),est ( B.est ),beta ( est._beta() ) {}
294                        //! Sets sufficient statistics to match that of givefrom mB0
295                        void set_statistics ( const BM* mB0 ) {const multiBM* mB=dynamic_cast<const multiBM*> ( mB0 ); beta=mB->beta;}
296                        void bayes ( const vec &dt )
297                        {
298                                if ( frg<1.0 ) {beta*=frg;last_lognc=est.lognc();}
299                                beta+=dt;
300                                if ( evalll ) {ll=est.lognc()-last_lognc;}
301                        }
302                        double logpred ( const vec &dt ) const
303                        {
304                                eDirich pred ( est );
305                                vec &beta = pred._beta();
306
307                                double lll;
308                                if ( frg<1.0 )
309                                        {beta*=frg;lll=pred.lognc();}
310                                else
311                                        if ( evalll ) {lll=last_lognc;}
312                                        else{lll=pred.lognc();}
313
314                                beta+=dt;
315                                return pred.lognc()-lll;
316                        }
317                        void flatten ( const BMEF* B )
318                        {
319                                const multiBM* E=dynamic_cast<const multiBM*> ( B );
320                                // sum(beta) should be equal to sum(B.beta)
321                                const vec &Eb=E->beta;//const_cast<multiBM*> ( E )->_beta();
322                                beta*= ( sum ( Eb ) /sum ( beta ) );
323                                if ( evalll ) {last_lognc=est.lognc();}
324                        }
325                        const epdf& posterior() const {return est;};
326                        const eDirich* _e() const {return &est;};
327                        void set_parameters ( const vec &beta0 )
328                        {
329                                est.set_parameters ( beta0 );
330                                if ( evalll ) {last_lognc=est.lognc();}
331                        }
332        };
333
334        /*!
335         \brief Gamma posterior density
336
337         Multivariate Gamma density as product of independent univariate densities.
338         \f[
339         f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
340         \f]
341        */
342
343        class egamma : public eEF
344        {
345                protected:
346                        //! Vector \f$\alpha\f$
347                        vec alpha;
348                        //! Vector \f$\beta\f$
349                        vec beta;
350                public :
351                        //! \name Constructors
352                        //!@{
353                        egamma ( ) :eEF ( ), alpha ( 0 ), beta ( 0 ) {};
354                        egamma ( const vec &a, const vec &b ) {set_parameters ( a, b );};
355                        void set_parameters ( const vec &a, const vec &b ) {alpha=a,beta=b;dim = alpha.length();};
356                        //!@}
357
358                        vec sample() const;
359                        //! TODO: is it used anywhere?
360//      mat sample ( int N ) const;
361                        double evallog ( const vec &val ) const;
362                        double lognc () const;
363                        //! Returns poiter to alpha and beta. Potentially dengerous: use with care!
364                        vec& _alpha() {return alpha;}
365                        vec& _beta() {return beta;}
366                        vec mean() const {return elem_div ( alpha,beta );}
367                        vec variance() const {return elem_div ( alpha,elem_mult ( beta,beta ) ); }
368        };
369
370        /*!
371         \brief Inverse-Gamma posterior density
372
373         Multivariate inverse-Gamma density as product of independent univariate densities.
374         \f[
375         f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
376         \f]
377
378        Vector \f$\beta\f$ has different meaning (in fact it is 1/beta as used in definition of iG)
379
380         Inverse Gamma can be converted to Gamma using \f[
381         x\sim iG(a,b) => 1/x\sim G(a,1/b)
382        \f]
383        This relation is used in sampling.
384         */
385
386        class eigamma : public egamma
387        {
388                protected:
389                public :
390                        //! \name Constructors
391                        //! All constructors are inherited
392                        //!@{
393                        //!@}
394
395                        vec sample() const {return 1.0/egamma::sample();};
396                        //! Returns poiter to alpha and beta. Potentially dangerous: use with care!
397                        vec mean() const {return elem_div ( beta,alpha-1 );}
398                        vec variance() const {vec mea=mean(); return elem_div ( elem_mult ( mea,mea ),alpha-2 );}
399        };
400        /*
401        //! Weighted mixture of epdfs with external owned components.
402        class emix : public epdf {
403        protected:
404                int n;
405                vec &w;
406                Array<epdf*> Coms;
407        public:
408        //! Default constructor
409                emix ( const RV &rv, vec &w0): epdf(rv), n(w0.length()), w(w0), Coms(n) {};
410                void set_parameters( int &i, double wi, epdf* ep){w(i)=wi;Coms(i)=ep;}
411                vec mean(){vec pom; for(int i=0;i<n;i++){pom+=Coms(i)->mean()*w(i);} return pom;};
412                vec sample() {it_error ( "Not implemented" );return 0;}
413        };
414        */
415
416//!  Uniform distributed density on a rectangular support
417
418        class euni: public epdf
419        {
420                protected:
421//! lower bound on support
422                        vec low;
423//! upper bound on support
424                        vec high;
425//! internal
426                        vec distance;
427//! normalizing coefficients
428                        double nk;
429//! cache of log( \c nk )
430                        double lnk;
431                public:
432                        //! \name Constructors
433                        //!@{
434                        euni ( ) :epdf ( ) {}
435                        euni ( const vec &low0, const vec &high0 ) {set_parameters ( low0,high0 );}
436                        void set_parameters ( const vec &low0, const vec &high0 )
437                        {
438                                distance = high0-low0;
439                                it_assert_debug ( min ( distance ) >0.0,"bad support" );
440                                low = low0;
441                                high = high0;
442                                nk = prod ( 1.0/distance );
443                                lnk = log ( nk );
444                                dim = low.length();
445                        }
446                        //!@}
447
448                        double eval ( const vec &val ) const  {return nk;}
449                        double evallog ( const vec &val ) const  {return lnk;}
450                        vec sample() const
451                        {
452                                vec smp ( dim );
453#pragma omp critical
454                                UniRNG.sample_vector ( dim ,smp );
455                                return low+elem_mult ( distance,smp );
456                        }
457                        //! set values of \c low and \c high
458                        vec mean() const {return ( high-low ) /2.0;}
459                        vec variance() const {return ( pow ( high,2 ) +pow ( low,2 ) +elem_mult ( high,low ) ) /3.0;}
460        };
461
462
463        /*!
464         \brief Normal distributed linear function with linear function of mean value;
465
466         Mean value \f$mu=A*rvc+mu_0\f$.
467        */
468        template<class sq_T>
469        class mlnorm : public mEF
470        {
471                protected:
472                        //! Internal epdf that arise by conditioning on \c rvc
473                        enorm<sq_T> epdf;
474                        mat A;
475                        vec mu_const;
476                        vec& _mu; //cached epdf.mu;
477                public:
478                        //! \name Constructors
479                        //!@{
480                        mlnorm ( ) :mEF (),epdf ( ),A ( ),_mu ( epdf._mu() ) {ep =&epdf; };
481                        mlnorm ( const  mat &A, const vec &mu0, const sq_T &R ) :epdf ( ),_mu ( epdf._mu() )
482                        {
483                                ep =&epdf; set_parameters ( A,mu0,R );
484                        };
485                        //! Set \c A and \c R
486                        void set_parameters ( const  mat &A, const vec &mu0, const sq_T &R );
487                        //!@}
488                        //! Set value of \c rvc . Result of this operation is stored in \c epdf use function \c _ep to access it.
489                        void condition ( const vec &cond );
490
491                        //!access function
492                        vec& _mu_const() {return mu_const;}
493                        //!access function
494                        mat& _A() {return A;}
495                        //!access function
496                        mat _R() {return epdf._R().to_mat();}
497
498                        template<class sq_M>
499                        friend std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_M> &ml );
500        };
501
502//! Mpdf with general function for mean value
503        template<class sq_T>
504        class mgnorm : public mEF
505        {
506                protected:
507                        //! Internal epdf that arise by conditioning on \c rvc
508                        enorm<sq_T> epdf;
509                        vec &mu;
510                        fnc* g;
511                public:
512                        //!default constructor
513                        mgnorm() :mu ( epdf._mu() ) {ep=&epdf;}
514                        //!set mean function
515                        void set_parameters ( fnc* g0, const sq_T &R0 ) {g=g0; epdf.set_parameters ( zeros ( g->dimension() ), R0 );}
516                        void condition ( const vec &cond ) {mu=g->eval ( cond );};
517        };
518
519        /*! (Approximate) Student t density with linear function of mean value
520
521        The internal epdf of this class is of the type of a Gaussian (enorm).
522        However, each conditioning is trying to assure the best possible approximation by taking into account the zeta function. See [] for reference.
523
524        Perhaps a moment-matching technique?
525        */
526        class mlstudent : public mlnorm<ldmat>
527        {
528                protected:
529                        ldmat Lambda;
530                        ldmat &_R;
531                        ldmat Re;
532                public:
533                        mlstudent ( ) :mlnorm<ldmat> (),
534                                        Lambda (),      _R ( epdf._R() ) {}
535                        void set_parameters ( const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0 )
536                        {
537                                it_assert_debug ( A0.rows() ==mu0.length(),"" );
538                                it_assert_debug ( R0.rows() ==A0.rows(),"" );
539
540                                epdf.set_parameters ( mu0,Lambda ); //
541                                A = A0;
542                                mu_const = mu0;
543                                Re=R0;
544                                Lambda = Lambda0;
545                        }
546                        void condition ( const vec &cond )
547                        {
548                                _mu = A*cond + mu_const;
549                                double zeta;
550                                //ugly hack!
551                                if ( ( cond.length() +1 ) ==Lambda.rows() )
552                                {
553                                        zeta = Lambda.invqform ( concat ( cond, vec_1 ( 1.0 ) ) );
554                                }
555                                else
556                                {
557                                        zeta = Lambda.invqform ( cond );
558                                }
559                                _R = Re;
560                                _R*= ( 1+zeta );// / ( nu ); << nu is in Re!!!!!!
561                        };
562
563        };
564        /*!
565        \brief  Gamma random walk
566
567        Mean value, \f$\mu\f$, of this density is given by \c rvc .
568        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
569        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
570
571        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
572        */
573        class mgamma : public mEF
574        {
575                protected:
576                        //! Internal epdf that arise by conditioning on \c rvc
577                        egamma epdf;
578                        //! Constant \f$k\f$
579                        double k;
580                        //! cache of epdf.beta
581                        vec &_beta;
582
583                public:
584                        //! Constructor
585                        mgamma ( ) : mEF ( ), epdf (), _beta ( epdf._beta() ) {ep=&epdf;};
586                        //! Set value of \c k
587                        void set_parameters ( double k, const vec &beta0 );
588                        void condition ( const vec &val ) {_beta=k/val;};
589        };
590
591        /*!
592        \brief  Inverse-Gamma random walk
593
594        Mean value, \f$ \mu \f$, of this density is given by \c rvc .
595        Standard deviation of the random walk is proportional to one \f$ k \f$-th the mean.
596        This is achieved by setting \f$ \alpha=\mu/k^2+2 \f$ and \f$ \beta=\mu(\alpha-1)\f$.
597
598        The standard deviation of the walk is then: \f$ \mu/\sqrt(k)\f$.
599         */
600        class migamma : public mEF
601        {
602                protected:
603                        //! Internal epdf that arise by conditioning on \c rvc
604                        eigamma epdf;
605                        //! Constant \f$k\f$
606                        double k;
607                        //! cache of epdf.alpha
608                        vec &_alpha;
609                        //! cache of epdf.beta
610                        vec &_beta;
611
612                public:
613                        //! \name Constructors
614                        //!@{
615                        migamma ( ) : mEF (), epdf ( ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;};
616                        migamma ( const migamma &m ) : mEF (), epdf ( m.epdf ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;};
617                        //!@}
618
619                        //! Set value of \c k
620                        void set_parameters ( int len, double k0 )
621                        {
622                                k=k0;
623                                epdf.set_parameters ( ( 1.0/ ( k*k ) +2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ );
624                                dimc = dimension();
625                        };
626                        void condition ( const vec &val )
627                        {
628                                _beta=elem_mult ( val, ( _alpha-1.0 ) );
629                        };
630        };
631
632        /*!
633        \brief  Gamma random walk around a fixed point
634
635        Mean 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
636        \f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
637
638        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
639        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
640
641        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
642        */
643        class mgamma_fix : public mgamma
644        {
645                protected:
646                        //! parameter l
647                        double l;
648                        //! reference vector
649                        vec refl;
650                public:
651                        //! Constructor
652                        mgamma_fix ( ) : mgamma ( ),refl () {};
653                        //! Set value of \c k
654                        void set_parameters ( double k0 , vec ref0, double l0 )
655                        {
656                                mgamma::set_parameters ( k0, ref0 );
657                                refl=pow ( ref0,1.0-l0 );l=l0;
658                                dimc=dimension();
659                        };
660
661                        void condition ( const vec &val ) {vec mean=elem_mult ( refl,pow ( val,l ) ); _beta=k/mean;};
662        };
663
664
665        /*!
666        \brief  Inverse-Gamma random walk around a fixed point
667
668        Mean 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
669        \f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
670
671        ==== Check == vv =
672        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
673        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
674
675        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
676         */
677        class migamma_ref : public migamma
678        {
679                protected:
680                        //! parameter l
681                        double l;
682                        //! reference vector
683                        vec refl;
684                public:
685                        //! Constructor
686                        migamma_ref ( ) : migamma (),refl ( ) {};
687                        //! Set value of \c k
688                        void set_parameters ( double k0 , vec ref0, double l0 )
689                        {
690                                migamma::set_parameters ( ref0.length(), k0 );
691                                refl=pow ( ref0,1.0-l0 );
692                                l=l0;
693                                dimc = dimension();
694                        };
695
696                        void condition ( const vec &val )
697                        {
698                                vec mean=elem_mult ( refl,pow ( val,l ) );
699                                migamma::condition ( mean );
700                        };
701        };
702
703        /*! Log-Normal probability density
704         only allow diagonal covariances!
705
706        Density of the form \f$ \log(x)\sim \mathcal{N}(\mu,\sigma^2), i.e.
707        \f[
708        x \sim \frac{1}{x\sigma\sqrt{2\pi}}\exp{-\frac{1}{2\sigma^2}(\log(x)-\mu)}
709        \f]
710
711        */
712        class elognorm: public enorm<ldmat>
713        {
714                public:
715                        vec sample() const {return exp ( enorm<ldmat>::sample() );};
716                        vec mean() const {vec var=enorm<ldmat>::variance();return exp ( mu - 0.5*var );};
717
718        };
719
720        /*!
721        \brief  Log-Normal random walk
722
723        Mean value, \f$\mu\f$, is...
724
725        ==== Check == vv =
726        Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
727        This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
728
729        The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
730         */
731        class mlognorm : public mpdf
732        {
733                protected:
734                        elognorm eno;
735                        //! parameter 1/2*sigma^2
736                        double sig2;
737                        //! access
738                        vec &mu;
739                public:
740                        //! Constructor
741                        mlognorm ( ) : eno (), mu ( eno._mu() ) {ep=&eno;};
742                        //! Set value of \c k
743                        void set_parameters ( int size, double k )
744                        {
745                                sig2 = 0.5*log ( k*k+1 );
746                                eno.set_parameters ( zeros ( size ),2*sig2*eye ( size ) );
747
748                                dimc = size;
749                        };
750
751                        void condition ( const vec &val )
752                        {
753                                mu=log ( val )-sig2;//elem_mult ( refl,pow ( val,l ) );
754                        };
755        };
756
757        /*! inverse Wishart density defined on Choleski decomposition
758
759        */
760        class eWishartCh : public epdf
761        {
762                protected:
763                        //! Upper-Triagle of Choleski decomposition of \f$ \Psi \f$
764                        chmat Y;
765                        //! dimension of matrix  \f$ \Psi \f$
766                        int p;
767                        //! degrees of freedom  \f$ \nu \f$
768                        double delta;
769                public:
770                        void set_parameters ( const mat &Y0, const double delta0 ) {Y=chmat ( Y0 );delta=delta0; p=Y.rows(); dim = p*p; }
771                        mat sample_mat() const
772                        {
773                                mat X=zeros ( p,p );
774
775                                //sample diagonal
776                                for ( int i=0;i<p;i++ )
777                                {
778                                        GamRNG.setup ( 0.5* ( delta-i ) , 0.5 ); // no +1 !! index if from 0
779#pragma omp critical
780                                        X ( i,i ) =sqrt ( GamRNG() );
781                                }
782                                //do the rest
783                                for ( int i=0;i<p;i++ )
784                                {
785                                        for ( int j=i+1;j<p;j++ )
786                                        {
787#pragma omp critical
788                                                X ( i,j ) =NorRNG.sample();
789                                        }
790                                }
791                                return X*Y._Ch();// return upper triangular part of the decomposition
792                        }
793                        vec sample () const
794                        {
795                                return vec ( sample_mat()._data(),p*p );
796                        }
797                        //! fast access function y0 will be copied into Y.Ch.
798                        void setY ( const mat &Ch0 ) {copy_vector ( dim,Ch0._data(), Y._Ch()._data() );}
799                        //! fast access function y0 will be copied into Y.Ch.
800                        void _setY ( const vec &ch0 ) {copy_vector ( dim, ch0._data(), Y._Ch()._data() ); }
801                        //! access function
802                        const chmat& getY()const {return Y;}
803        };
804
805        class eiWishartCh: public epdf
806        {
807                protected:
808                        eWishartCh W;
809                        int p;
810                        double delta;
811                public:
812                        void set_parameters ( const mat &Y0, const double delta0) {
813                                delta = delta0;
814                                W.set_parameters ( inv ( Y0 ),delta0 ); 
815                                dim = W.dimension(); p=Y0.rows();
816                        }
817                        vec sample() const {mat iCh; iCh=inv ( W.sample_mat() ); return vec ( iCh._data(),dim );}
818                        void _setY ( const vec &y0 )
819                        {
820                                mat Ch ( p,p );
821                                mat iCh ( p,p );
822                                copy_vector ( dim, y0._data(), Ch._data() );
823                               
824                                iCh=inv ( Ch );
825                                W.setY ( iCh );
826                        }
827                        virtual double evallog ( const vec &val ) const {
828                                chmat X(p);
829                                const chmat& Y=W.getY();
830                                 
831                                copy_vector(p*p,val._data(),X._Ch()._data());
832                                chmat iX(p);X.inv(iX);
833                                // compute 
834//                              \frac{ |\Psi|^{m/2}|X|^{-(m+p+1)/2}e^{-tr(\Psi X^{-1})/2} }{ 2^{mp/2}\Gamma_p(m/2)},
835                                mat M=Y.to_mat()*iX.to_mat();
836                               
837                                double log1 = 0.5*p*(2*Y.logdet())-0.5*(delta+p+1)*(2*X.logdet())-0.5*trace(M); 
838                                //Fixme! Multivariate gamma omitted!! it is ok for sampling, but not otherwise!!
839                               
840/*                              if (0) {
841                                        mat XX=X.to_mat();
842                                        mat YY=Y.to_mat();
843                                       
844                                        double log2 = 0.5*p*log(det(YY))-0.5*(delta+p+1)*log(det(XX))-0.5*trace(YY*inv(XX));
845                                        cout << log1 << "," << log2 << endl;
846                                }*/
847                                return log1;                           
848                        };
849                       
850        };
851
852        class rwiWishartCh : public mpdf
853        {
854                protected:
855                        eiWishartCh eiW;
856                        //!square root of \f$ \nu-p-1 \f$ - needed for computation of \f$ \Psi \f$ from conditions
857                        double sqd;
858                        //reference point for diagonal
859                        vec refl;
860                        double l;
861                        int p;
862                public:
863                        void set_parameters ( int p0, double k, vec ref0, double l0  )
864                        {
865                                p=p0;
866                                double delta = 2/(k*k)+p+3;
867                                sqd=sqrt ( delta-p-1 );
868                                l=l0;
869                                refl=pow(ref0,1-l);
870                               
871                                eiW.set_parameters ( eye ( p ),delta );
872                                ep=&eiW;
873                                dimc=eiW.dimension();
874                        }
875                        void condition ( const vec &c ) {
876                                vec z=c;
877                                int ri=0;
878                                for(int i=0;i<p*p;i+=(p+1)){//trace diagonal element
879                                        z(i) = pow(z(i),l)*refl(ri);
880                                        ri++;
881                                }
882
883                                eiW._setY ( sqd*z );
884                        }
885        };
886
887//! Switch between various resampling methods.
888        enum RESAMPLING_METHOD { MULTINOMIAL = 0, STRATIFIED = 1, SYSTEMATIC = 3 };
889        /*!
890        \brief Weighted empirical density
891
892        Used e.g. in particle filters.
893        */
894        class eEmp: public epdf
895        {
896                protected :
897                        //! Number of particles
898                        int n;
899                        //! Sample weights \f$w\f$
900                        vec w;
901                        //! Samples \f$x^{(i)}, i=1..n\f$
902                        Array<vec> samples;
903                public:
904                        //! \name Constructors
905                        //!@{
906                        eEmp ( ) :epdf ( ),w ( ),samples ( ) {};
907                        eEmp ( const eEmp &e ) : epdf ( e ), w ( e.w ), samples ( e.samples ) {};
908                        //!@}
909
910                        //! Set samples and weights
911                        void set_statistics ( const vec &w0, const epdf* pdf0 );
912                        //! Set samples and weights
913                        void set_statistics ( const epdf* pdf0 , int n ) {set_statistics ( ones ( n ) /n,pdf0 );};
914                        //! Set sample
915                        void set_samples ( const epdf* pdf0 );
916                        //! Set sample
917                        void set_parameters ( int n0, bool copy=true ) {n=n0; w.set_size ( n0,copy );samples.set_size ( n0,copy );};
918                        //! Potentially dangerous, use with care.
919                        vec& _w()  {return w;};
920                        //! Potentially dangerous, use with care.
921                        const vec& _w() const {return w;};
922                        //! access function
923                        Array<vec>& _samples() {return samples;};
924                        //! access function
925                        const Array<vec>& _samples() const {return samples;};
926                        //! 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.
927                        ivec resample ( RESAMPLING_METHOD method=SYSTEMATIC );
928                        //! inherited operation : NOT implemneted
929                        vec sample() const {it_error ( "Not implemented" );return 0;}
930                        //! inherited operation : NOT implemneted
931                        double evallog ( const vec &val ) const {it_error ( "Not implemented" );return 0.0;}
932                        vec mean() const
933                        {
934                                vec pom=zeros ( dim );
935                                for ( int i=0;i<n;i++ ) {pom+=samples ( i ) *w ( i );}
936                                return pom;
937                        }
938                        vec variance() const
939                        {
940                                vec pom=zeros ( dim );
941                                for ( int i=0;i<n;i++ ) {pom+=pow ( samples ( i ),2 ) *w ( i );}
942                                return pom-pow ( mean(),2 );
943                        }
944                        //! For this class, qbounds are minimum and maximum value of the population!
945                        void qbounds ( vec &lb, vec &ub, double perc=0.95 ) const
946                        {
947                                // lb in inf so than it will be pushed below;
948                                lb.set_size ( dim );
949                                ub.set_size ( dim );
950                                lb = std::numeric_limits<double>::infinity();
951                                ub = -std::numeric_limits<double>::infinity();
952                                int j;
953                                for ( int i=0;i<n;i++ )
954                                {
955                                        for ( j=0;j<dim; j++ )
956                                        {
957                                                if ( samples ( i ) ( j ) <lb ( j ) ) {lb ( j ) =samples ( i ) ( j );}
958                                                if ( samples ( i ) ( j ) >ub ( j ) ) {ub ( j ) =samples ( i ) ( j );}
959                                        }
960                                }
961                        }
962        };
963
964
965////////////////////////
966
967        template<class sq_T>
968        void enorm<sq_T>::set_parameters ( const vec &mu0, const sq_T &R0 )
969        {
970//Fixme test dimensions of mu0 and R0;
971                mu = mu0;
972                R = R0;
973                dim = mu0.length();
974        };
975
976        template<class sq_T>
977        void enorm<sq_T>::dupdate ( mat &v, double nu )
978        {
979                //
980        };
981
982// template<class sq_T>
983// void enorm<sq_T>::tupdate ( double phi, mat &vbar, double nubar ) {
984//      //
985// };
986
987        template<class sq_T>
988        vec enorm<sq_T>::sample() const
989        {
990                vec x ( dim );
991#pragma omp critical
992                NorRNG.sample_vector ( dim,x );
993                vec smp = R.sqrt_mult ( x );
994
995                smp += mu;
996                return smp;
997        };
998
999        template<class sq_T>
1000        mat enorm<sq_T>::sample ( int N ) const
1001        {
1002                mat X ( dim,N );
1003                vec x ( dim );
1004                vec pom;
1005                int i;
1006
1007                for ( i=0;i<N;i++ )
1008                {
1009#pragma omp critical
1010                        NorRNG.sample_vector ( dim,x );
1011                        pom = R.sqrt_mult ( x );
1012                        pom +=mu;
1013                        X.set_col ( i, pom );
1014                }
1015
1016                return X;
1017        };
1018
1019// template<class sq_T>
1020// double enorm<sq_T>::eval ( const vec &val ) const {
1021//      double pdfl,e;
1022//      pdfl = evallog ( val );
1023//      e = exp ( pdfl );
1024//      return e;
1025// };
1026
1027        template<class sq_T>
1028        double enorm<sq_T>::evallog_nn ( const vec &val ) const
1029        {
1030                // 1.83787706640935 = log(2pi)
1031                double tmp=-0.5* ( R.invqform ( mu-val ) );// - lognc();
1032                return  tmp;
1033        };
1034
1035        template<class sq_T>
1036        inline double enorm<sq_T>::lognc () const
1037        {
1038                // 1.83787706640935 = log(2pi)
1039                double tmp=0.5* ( R.cols() * 1.83787706640935 +R.logdet() );
1040                return tmp;
1041        };
1042
1043        template<class sq_T>
1044        void mlnorm<sq_T>::set_parameters ( const mat &A0, const vec &mu0, const sq_T &R0 )
1045        {
1046                it_assert_debug ( A0.rows() ==mu0.length(),"" );
1047                it_assert_debug ( A0.rows() ==R0.rows(),"" );
1048
1049                epdf.set_parameters ( zeros ( A0.rows() ),R0 );
1050                A = A0;
1051                mu_const = mu0;
1052                dimc=A0.cols();
1053        }
1054
1055// template<class sq_T>
1056// vec mlnorm<sq_T>::samplecond (const  vec &cond, double &lik ) {
1057//      this->condition ( cond );
1058//      vec smp = epdf.sample();
1059//      lik = epdf.eval ( smp );
1060//      return smp;
1061// }
1062
1063// template<class sq_T>
1064// mat mlnorm<sq_T>::samplecond (const vec &cond, vec &lik, int n ) {
1065//      int i;
1066//      int dim = rv.count();
1067//      mat Smp ( dim,n );
1068//      vec smp ( dim );
1069//      this->condition ( cond );
1070//
1071//      for ( i=0; i<n; i++ ) {
1072//              smp = epdf.sample();
1073//              lik ( i ) = epdf.eval ( smp );
1074//              Smp.set_col ( i ,smp );
1075//      }
1076//
1077//      return Smp;
1078// }
1079
1080        template<class sq_T>
1081        void mlnorm<sq_T>::condition ( const vec &cond )
1082        {
1083                _mu = A*cond + mu_const;
1084//R is already assigned;
1085        }
1086
1087        template<class sq_T>
1088        enorm<sq_T>* enorm<sq_T>::marginal ( const RV &rvn ) const
1089        {
1090                it_assert_debug ( isnamed(), "rv description is not assigned" );
1091                ivec irvn = rvn.dataind ( rv );
1092
1093                sq_T Rn ( R,irvn ); //select rows and columns of R
1094
1095                enorm<sq_T>* tmp = new enorm<sq_T>;
1096                tmp->set_rv ( rvn );
1097                tmp->set_parameters ( mu ( irvn ), Rn );
1098                return tmp;
1099        }
1100
1101        template<class sq_T>
1102        mpdf* enorm<sq_T>::condition ( const RV &rvn ) const
1103        {
1104
1105                it_assert_debug ( isnamed(),"rvs are not assigned" );
1106
1107                RV rvc = rv.subt ( rvn );
1108                it_assert_debug ( ( rvc._dsize() +rvn._dsize() ==rv._dsize() ),"wrong rvn" );
1109                //Permutation vector of the new R
1110                ivec irvn = rvn.dataind ( rv );
1111                ivec irvc = rvc.dataind ( rv );
1112                ivec perm=concat ( irvn , irvc );
1113                sq_T Rn ( R,perm );
1114
1115                //fixme - could this be done in general for all sq_T?
1116                mat S=Rn.to_mat();
1117                //fixme
1118                int n=rvn._dsize()-1;
1119                int end=R.rows()-1;
1120                mat S11 = S.get ( 0,n, 0, n );
1121                mat S12 = S.get ( 0, n , rvn._dsize(), end );
1122                mat S22 = S.get ( rvn._dsize(), end, rvn._dsize(), end );
1123
1124                vec mu1 = mu ( irvn );
1125                vec mu2 = mu ( irvc );
1126                mat A=S12*inv ( S22 );
1127                sq_T R_n ( S11 - A *S12.T() );
1128
1129                mlnorm<sq_T>* tmp=new mlnorm<sq_T> ( );
1130                tmp->set_rv ( rvn ); tmp->set_rvc ( rvc );
1131                tmp->set_parameters ( A,mu1-A*mu2,R_n );
1132                return tmp;
1133        }
1134
1135///////////
1136
1137        template<class sq_T>
1138        std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_T> &ml )
1139        {
1140                os << "A:"<< ml.A<<endl;
1141                os << "mu:"<< ml.mu_const<<endl;
1142                os << "R:" << ml.epdf._R().to_mat() <<endl;
1143                return os;
1144        };
1145
1146}
1147#endif //EF_H
Note: See TracBrowser for help on using the browser.