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

Revision 613, 35.5 kB (checked in by smidl, 15 years ago)

documentation

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