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

Revision 488, 33.5 kB (checked in by smidl, 15 years ago)

changes in mpdf -> compile OK, broken tests!

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