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

Revision 404, 34.5 kB (checked in by smidl, 15 years ago)

Change in epdf: evallog returns -inf for points out of support. Merger is aware of it now.

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