root/library/bdm/stat/exp_family.h @ 395

Revision 395, 34.4 kB (checked in by smidl, 15 years ago)

merging works for merger_mx

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