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

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

replace assert_debug by assert in classes interacting with users (from_setting, set_parameters, validate)

  • 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 (mu.length() == R.rows(), "mu and R parameters do not match");
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 (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 (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                        this->iepdf.set_parameters (zeros (A0.rows()), R0);
562                        A = A0;
563                        mu_const = mu0;
564                        this->dimc = A0.cols();
565                }
566                //!@}
567                //! Set value of \c rvc . Result of this operation is stored in \c epdf use function \c _ep to access it.
568                void condition (const vec &cond) {
569                        this->iepdf._mu() = A * cond + mu_const;
570//R is already assigned;
571                }
572
573                //!access function
574                const vec& _mu_const() const {return mu_const;}
575                //!access function
576                const mat& _A() const {return A;}
577                //!access function
578                mat _R() const { return this->iepdf._R().to_mat(); }
579
580                //! Debug stream
581                template<typename sq_M>
582                friend std::ostream &operator<< (std::ostream &os,  mlnorm<sq_M, enorm> &ml);
583
584                void from_setting (const Setting &set) {
585                        mpdf::from_setting (set);
586
587                        UI::get (A, set, "A", UI::compulsory);
588                        UI::get (mu_const, set, "const", UI::compulsory);
589                        mat R0;
590                        UI::get (R0, set, "R", UI::compulsory);
591                        set_parameters (A, mu_const, R0);
592                        validate();
593                };
594                void validate() {
595                        bdm_assert (A.rows() == mu_const.length(), "mlnorm: A vs. mu mismatch");
596                        bdm_assert (A.rows() == _R().rows(), "mlnorm: A vs. R mismatch");
597                       
598                }
599};
600UIREGISTER2 (mlnorm,ldmat);
601SHAREDPTR2 ( mlnorm, ldmat );
602UIREGISTER2 (mlnorm,fsqmat);
603SHAREDPTR2 ( mlnorm, fsqmat );
604UIREGISTER2 (mlnorm, chmat);
605SHAREDPTR2 ( mlnorm, chmat );
606
607//! Mpdf with general function for mean value
608template<class sq_T>
609class mgnorm : public mpdf_internal< enorm< sq_T > >
610{
611        private:
612//                      vec &mu; WHY NOT?
613                shared_ptr<fnc> g;
614
615        public:
616                //!default constructor
617                mgnorm() : mpdf_internal<enorm<sq_T> >() { }
618                //!set mean function
619                inline void set_parameters (const shared_ptr<fnc> &g0, const sq_T &R0);
620                inline void condition (const vec &cond);
621
622
623                /*! UI for mgnorm
624
625                The mgnorm is constructed from a structure with fields:
626                \code
627                system = {
628                        type = "mgnorm";
629                        // function for mean value evolution
630                        g = {type="fnc"; ... }
631
632                        // variance
633                        R = [1, 0,
634                                 0, 1];
635                        // --OR --
636                        dR = [1, 1];
637
638                        // == OPTIONAL ==
639
640                        // description of y variables
641                        y = {type="rv"; names=["y", "u"];};
642                        // description of u variable
643                        u = {type="rv"; names=[];}
644                };
645                \endcode
646
647                Result if
648                */
649
650                void from_setting (const Setting &set) {
651                        shared_ptr<fnc> g = UI::build<fnc> (set, "g", UI::compulsory);
652
653                        mat R;
654                        vec dR;
655                        if (UI::get (dR, set, "dR"))
656                                R = diag (dR);
657                        else
658                                UI::get (R, set, "R", UI::compulsory);
659
660                        set_parameters (g, R);
661                }
662};
663
664UIREGISTER2 (mgnorm, chmat);
665SHAREDPTR2 ( mgnorm, chmat );
666
667
668/*! (Approximate) Student t density with linear function of mean value
669
670The internal epdf of this class is of the type of a Gaussian (enorm).
671However, each conditioning is trying to assure the best possible approximation by taking into account the zeta function. See [] for reference.
672
673Perhaps a moment-matching technique?
674*/
675class mlstudent : public mlnorm<ldmat, enorm>
676{
677        protected:
678                //! Variable \f$ \Lambda \f$ from theory
679                ldmat Lambda;
680                //! Reference to variable \f$ R \f$
681                ldmat &_R;
682                //! Variable \f$ R_e \f$
683                ldmat Re;
684        public:
685                mlstudent () : mlnorm<ldmat, enorm> (),
686                                Lambda (),      _R (iepdf._R()) {}
687                //! constructor function
688                void set_parameters (const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0) {
689                        iepdf.set_parameters (mu0, R0);// was Lambda, why?
690                        A = A0;
691                        mu_const = mu0;
692                        Re = R0;
693                        Lambda = Lambda0;
694                }
695                void condition (const vec &cond) {
696                        iepdf._mu() = A * cond + mu_const;
697                        double zeta;
698                        //ugly hack!
699                        if ( (cond.length() + 1) == Lambda.rows()) {
700                                zeta = Lambda.invqform (concat (cond, vec_1 (1.0)));
701                        } else {
702                                zeta = Lambda.invqform (cond);
703                        }
704                        _R = Re;
705                        _R *= (1 + zeta);// / ( nu ); << nu is in Re!!!!!!
706                };
707
708                void validate() {
709                        bdm_assert (A.rows() == mu_const.length(), "mlstudent: A vs. mu mismatch");
710                        bdm_assert (_R.rows() == A.rows(), "mlstudent: A vs. R mismatch");
711                       
712                }
713};
714/*!
715\brief  Gamma random walk
716
717Mean value, \f$\mu\f$, of this density is given by \c rvc .
718Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
719This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
720
721The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
722*/
723class mgamma : public mpdf_internal<egamma>
724{
725        protected:
726
727                //! Constant \f$k\f$
728                double k;
729
730                //! cache of iepdf.beta
731                vec &_beta;
732
733        public:
734                //! Constructor
735                mgamma() : mpdf_internal<egamma>(), k (0),
736                                _beta (iepdf._beta()) {
737                }
738
739                //! Set value of \c k
740                void set_parameters (double k, const vec &beta0);
741
742                void condition (const vec &val) {_beta = k / val;};
743                //! Load from structure with elements:
744                //!  \code
745                //! { alpha = [...];         // vector of alpha
746                //!   k = 1.1;               // multiplicative constant k
747                //!   rv = {class="RV",...}  // description of RV
748                //!   rvc = {class="RV",...} // description of RV in condition
749                //! }
750                //! \endcode
751                //!@}
752                void from_setting (const Setting &set) {
753                        mpdf::from_setting (set); // reads rv and rvc
754                        vec betatmp; // ugly but necessary
755                        UI::get (betatmp, set, "beta", UI::compulsory);
756                        UI::get (k, set, "k", UI::compulsory);
757                        set_parameters (k, betatmp);
758                }
759};
760UIREGISTER (mgamma);
761SHAREDPTR (mgamma);
762
763/*!
764\brief  Inverse-Gamma random walk
765
766Mean value, \f$ \mu \f$, of this density is given by \c rvc .
767Standard deviation of the random walk is proportional to one \f$ k \f$-th the mean.
768This is achieved by setting \f$ \alpha=\mu/k^2+2 \f$ and \f$ \beta=\mu(\alpha-1)\f$.
769
770The standard deviation of the walk is then: \f$ \mu/\sqrt(k)\f$.
771 */
772class migamma : public mpdf_internal<eigamma>
773{
774        protected:
775                //! Constant \f$k\f$
776                double k;
777
778                //! cache of iepdf.alpha
779                vec &_alpha;
780
781                //! cache of iepdf.beta
782                vec &_beta;
783
784        public:
785                //! \name Constructors
786                //!@{
787                migamma() : mpdf_internal<eigamma>(),
788                                k (0),
789                                _alpha (iepdf._alpha()),
790                                _beta (iepdf._beta()) {
791                }
792
793                migamma (const migamma &m) : mpdf_internal<eigamma>(),
794                                k (0),
795                                _alpha (iepdf._alpha()),
796                                _beta (iepdf._beta()) {
797                }
798                //!@}
799
800                //! Set value of \c k
801                void set_parameters (int len, double k0) {
802                        k = k0;
803                        iepdf.set_parameters ( (1.0 / (k*k) + 2.0) *ones (len) /*alpha*/, ones (len) /*beta*/);
804                        dimc = dimension();
805                };
806                void condition (const vec &val) {
807                        _beta = elem_mult (val, (_alpha - 1.0));
808                };
809};
810
811
812/*!
813\brief  Gamma random walk around a fixed point
814
815Mean 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
816\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
817
818Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
819This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
820
821The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
822*/
823class mgamma_fix : public mgamma
824{
825        protected:
826                //! parameter l
827                double l;
828                //! reference vector
829                vec refl;
830        public:
831                //! Constructor
832                mgamma_fix () : mgamma (), refl () {};
833                //! Set value of \c k
834                void set_parameters (double k0 , vec ref0, double l0) {
835                        mgamma::set_parameters (k0, ref0);
836                        refl = pow (ref0, 1.0 - l0);l = l0;
837                        dimc = dimension();
838                };
839
840                void condition (const vec &val) {vec mean = elem_mult (refl, pow (val, l)); _beta = k / mean;};
841};
842
843
844/*!
845\brief  Inverse-Gamma random walk around a fixed point
846
847Mean 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
848\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
849
850==== Check == vv =
851Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
852This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
853
854The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
855 */
856class migamma_ref : public migamma
857{
858        protected:
859                //! parameter l
860                double l;
861                //! reference vector
862                vec refl;
863        public:
864                //! Constructor
865                migamma_ref () : migamma (), refl () {};
866                //! Set value of \c k
867                void set_parameters (double k0 , vec ref0, double l0) {
868                        migamma::set_parameters (ref0.length(), k0);
869                        refl = pow (ref0, 1.0 - l0);
870                        l = l0;
871                        dimc = dimension();
872                };
873
874                void condition (const vec &val) {
875                        vec mean = elem_mult (refl, pow (val, l));
876                        migamma::condition (mean);
877                };
878
879                /*! UI for migamma_ref
880
881                The migamma_ref is constructed from a structure with fields:
882                \code
883                system = {
884                        type = "migamma_ref";
885                        ref = [1e-5; 1e-5; 1e-2 1e-3];            // reference vector
886                        l = 0.999;                                // constant l
887                        k = 0.1;                                  // constant k
888                       
889                        // == OPTIONAL ==
890                        // description of y variables
891                        y = {type="rv"; names=["y", "u"];};
892                        // description of u variable
893                        u = {type="rv"; names=[];}
894                };
895                \endcode
896
897                Result if
898                 */
899                void from_setting (const Setting &set);
900
901                // TODO dodelat void to_setting( Setting &set ) const;
902};
903
904
905UIREGISTER (migamma_ref);
906SHAREDPTR (migamma_ref);
907
908/*! Log-Normal probability density
909 only allow diagonal covariances!
910
911Density of the form \f$ \log(x)\sim \mathcal{N}(\mu,\sigma^2) \f$ , i.e.
912\f[
913x \sim \frac{1}{x\sigma\sqrt{2\pi}}\exp{-\frac{1}{2\sigma^2}(\log(x)-\mu)}
914\f]
915
916*/
917class elognorm: public enorm<ldmat>
918{
919        public:
920                vec sample() const {return exp (enorm<ldmat>::sample());};
921                vec mean() const {vec var = enorm<ldmat>::variance();return exp (mu - 0.5*var);};
922
923};
924
925/*!
926\brief  Log-Normal random walk
927
928Mean value, \f$\mu\f$, is...
929
930==== Check == vv =
931Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
932This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
933
934The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
935 */
936class mlognorm : public mpdf_internal<elognorm>
937{
938        protected:
939                //! parameter 1/2*sigma^2
940                double sig2;
941
942                //! access
943                vec &mu;
944        public:
945                //! Constructor
946                mlognorm() : mpdf_internal<elognorm>(),
947                                sig2 (0),
948                                mu (iepdf._mu()) {
949                }
950
951                //! Set value of \c k
952                void set_parameters (int size, double k) {
953                        sig2 = 0.5 * log (k * k + 1);
954                        iepdf.set_parameters (zeros (size), 2*sig2*eye (size));
955
956                        dimc = size;
957                };
958
959                void condition (const vec &val) {
960                        mu = log (val) - sig2;//elem_mult ( refl,pow ( val,l ) );
961                };
962
963                /*! UI for mlognorm
964
965                The mlognorm is constructed from a structure with fields:
966                \code
967                system = {
968                        type = "mlognorm";
969                        k = 0.1;                                  // constant k
970                        mu0 = [1., 1.];
971                       
972                        // == OPTIONAL ==
973                        // description of y variables
974                        y = {type="rv"; names=["y", "u"];};
975                        // description of u variable
976                        u = {type="rv"; names=[];}
977                };
978                \endcode
979
980                 */
981                void from_setting (const Setting &set);
982
983                // TODO dodelat void to_setting( Setting &set ) const;
984
985};
986
987UIREGISTER (mlognorm);
988SHAREDPTR (mlognorm);
989
990/*! inverse Wishart density defined on Choleski decomposition
991
992*/
993class eWishartCh : public epdf
994{
995        protected:
996                //! Upper-Triagle of Choleski decomposition of \f$ \Psi \f$
997                chmat Y;
998                //! dimension of matrix  \f$ \Psi \f$
999                int p;
1000                //! degrees of freedom  \f$ \nu \f$
1001                double delta;
1002        public:
1003                //! Set internal structures
1004                void set_parameters (const mat &Y0, const double delta0) {Y = chmat (Y0);delta = delta0; p = Y.rows(); dim = p * p; }
1005                //! Sample matrix argument
1006                mat sample_mat() const {
1007                        mat X = zeros (p, p);
1008
1009                        //sample diagonal
1010                        for (int i = 0;i < p;i++) {
1011                                GamRNG.setup (0.5* (delta - i) , 0.5);   // no +1 !! index if from 0
1012#pragma omp critical
1013                                X (i, i) = sqrt (GamRNG());
1014                        }
1015                        //do the rest
1016                        for (int i = 0;i < p;i++) {
1017                                for (int j = i + 1;j < p;j++) {
1018#pragma omp critical
1019                                        X (i, j) = NorRNG.sample();
1020                                }
1021                        }
1022                        return X*Y._Ch();// return upper triangular part of the decomposition
1023                }
1024                vec sample () const {
1025                        return vec (sample_mat()._data(), p*p);
1026                }
1027                //! fast access function y0 will be copied into Y.Ch.
1028                void setY (const mat &Ch0) {copy_vector (dim, Ch0._data(), Y._Ch()._data());}
1029                //! fast access function y0 will be copied into Y.Ch.
1030                void _setY (const vec &ch0) {copy_vector (dim, ch0._data(), Y._Ch()._data()); }
1031                //! access function
1032                const chmat& getY() const {return Y;}
1033};
1034
1035//! Inverse Wishart on Choleski decomposition
1036/*! Being computed by conversion from `standard' Wishart
1037*/
1038class eiWishartCh: public epdf
1039{
1040        protected:
1041                //! Internal instance of Wishart density
1042                eWishartCh W;
1043                //! size of Ch
1044                int p;
1045                //! parameter delta
1046                double delta;
1047        public:
1048                //! constructor function
1049                void set_parameters (const mat &Y0, const double delta0) {
1050                        delta = delta0;
1051                        W.set_parameters (inv (Y0), delta0);
1052                        dim = W.dimension(); p = Y0.rows();
1053                }
1054                vec sample() const {mat iCh; iCh = inv (W.sample_mat()); return vec (iCh._data(), dim);}
1055                //! access function
1056                void _setY (const vec &y0) {
1057                        mat Ch (p, p);
1058                        mat iCh (p, p);
1059                        copy_vector (dim, y0._data(), Ch._data());
1060
1061                        iCh = inv (Ch);
1062                        W.setY (iCh);
1063                }
1064                virtual double evallog (const vec &val) const {
1065                        chmat X (p);
1066                        const chmat& Y = W.getY();
1067
1068                        copy_vector (p*p, val._data(), X._Ch()._data());
1069                        chmat iX (p);X.inv (iX);
1070                        // compute
1071//                              \frac{ |\Psi|^{m/2}|X|^{-(m+p+1)/2}e^{-tr(\Psi X^{-1})/2} }{ 2^{mp/2}\Gamma_p(m/2)},
1072                        mat M = Y.to_mat() * iX.to_mat();
1073
1074                        double log1 = 0.5 * p * (2 * Y.logdet()) - 0.5 * (delta + p + 1) * (2 * X.logdet()) - 0.5 * trace (M);
1075                        //Fixme! Multivariate gamma omitted!! it is ok for sampling, but not otherwise!!
1076
1077                        /*                              if (0) {
1078                                                                mat XX=X.to_mat();
1079                                                                mat YY=Y.to_mat();
1080
1081                                                                double log2 = 0.5*p*log(det(YY))-0.5*(delta+p+1)*log(det(XX))-0.5*trace(YY*inv(XX));
1082                                                                cout << log1 << "," << log2 << endl;
1083                                                        }*/
1084                        return log1;
1085                };
1086
1087};
1088
1089//! Random Walk on inverse Wishart
1090class rwiWishartCh : public mpdf_internal<eiWishartCh>
1091{
1092        protected:
1093                //!square root of \f$ \nu-p-1 \f$ - needed for computation of \f$ \Psi \f$ from conditions
1094                double sqd;
1095                //!reference point for diagonal
1096                vec refl;
1097                //! power of the reference
1098                double l;
1099                //! dimension
1100                int p;
1101
1102        public:
1103                rwiWishartCh() : sqd (0), l (0), p (0) {}
1104                //! constructor function
1105                void set_parameters (int p0, double k, vec ref0, double l0) {
1106                        p = p0;
1107                        double delta = 2 / (k * k) + p + 3;
1108                        sqd = sqrt (delta - p - 1);
1109                        l = l0;
1110                        refl = pow (ref0, 1 - l);
1111
1112                        iepdf.set_parameters (eye (p), delta);
1113                        dimc = iepdf.dimension();
1114                }
1115                void condition (const vec &c) {
1116                        vec z = c;
1117                        int ri = 0;
1118                        for (int i = 0;i < p*p;i += (p + 1)) {//trace diagonal element
1119                                z (i) = pow (z (i), l) * refl (ri);
1120                                ri++;
1121                        }
1122
1123                        iepdf._setY (sqd*z);
1124                }
1125};
1126
1127//! Switch between various resampling methods.
1128enum RESAMPLING_METHOD { MULTINOMIAL = 0, STRATIFIED = 1, SYSTEMATIC = 3 };
1129/*!
1130\brief Weighted empirical density
1131
1132Used e.g. in particle filters.
1133*/
1134class eEmp: public epdf
1135{
1136        protected :
1137                //! Number of particles
1138                int n;
1139                //! Sample weights \f$w\f$
1140                vec w;
1141                //! Samples \f$x^{(i)}, i=1..n\f$
1142                Array<vec> samples;
1143        public:
1144                //! \name Constructors
1145                //!@{
1146                eEmp () : epdf (), w (), samples () {};
1147                //! copy constructor
1148                eEmp (const eEmp &e) : epdf (e), w (e.w), samples (e.samples) {};
1149                //!@}
1150
1151                //! Set samples and weights
1152                void set_statistics (const vec &w0, const epdf &pdf0);
1153                //! Set samples and weights
1154                void set_statistics (const epdf &pdf0 , int n) {set_statistics (ones (n) / n, pdf0);};
1155                //! Set sample
1156                void set_samples (const epdf* pdf0);
1157                //! Set sample
1158                void set_parameters (int n0, bool copy = true) {n = n0; w.set_size (n0, copy);samples.set_size (n0, copy);};
1159                //! Set samples
1160                void set_parameters (const Array<vec> &Av) {
1161                        bdm_assert(Av.size()>0,"Empty samples"); 
1162                        n = Av.size(); 
1163                        epdf::set_parameters(Av(0).length());
1164                        w=1/n*ones(n);
1165                        samples=Av;
1166                };
1167                //! Potentially dangerous, use with care.
1168                vec& _w()  {return w;};
1169                //! Potentially dangerous, use with care.
1170                const vec& _w() const {return w;};
1171                //! access function
1172                Array<vec>& _samples() {return samples;};
1173                //! access function
1174                const Array<vec>& _samples() const {return samples;};
1175                //! 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.
1176                ivec resample (RESAMPLING_METHOD method = SYSTEMATIC);
1177
1178                //! inherited operation : NOT implemented
1179                vec sample() const {
1180                        bdm_error ("Not implemented");
1181                        return vec();
1182                }
1183
1184                //! inherited operation : NOT implemented
1185                double evallog (const vec &val) const {
1186                        bdm_error ("Not implemented");
1187                        return 0.0;
1188                }
1189
1190                vec mean() const {
1191                        vec pom = zeros (dim);
1192                        for (int i = 0;i < n;i++) {pom += samples (i) * w (i);}
1193                        return pom;
1194                }
1195                vec variance() const {
1196                        vec pom = zeros (dim);
1197                        for (int i = 0;i < n;i++) {pom += pow (samples (i), 2) * w (i);}
1198                        return pom -pow (mean(), 2);
1199                }
1200                //! For this class, qbounds are minimum and maximum value of the population!
1201                void qbounds (vec &lb, vec &ub, double perc = 0.95) const {
1202                        // lb in inf so than it will be pushed below;
1203                        lb.set_size (dim);
1204                        ub.set_size (dim);
1205                        lb = std::numeric_limits<double>::infinity();
1206                        ub = -std::numeric_limits<double>::infinity();
1207                        int j;
1208                        for (int i = 0;i < n;i++) {
1209                                for (j = 0;j < dim; j++) {
1210                                        if (samples (i) (j) < lb (j)) {lb (j) = samples (i) (j);}
1211                                        if (samples (i) (j) > ub (j)) {ub (j) = samples (i) (j);}
1212                                }
1213                        }
1214                }
1215};
1216
1217
1218////////////////////////
1219
1220template<class sq_T>
1221void enorm<sq_T>::set_parameters (const vec &mu0, const sq_T &R0)
1222{
1223//Fixme test dimensions of mu0 and R0;
1224        mu = mu0;
1225        R = R0;
1226        validate();
1227};
1228
1229template<class sq_T>
1230void enorm<sq_T>::from_setting (const Setting &set)
1231{
1232        epdf::from_setting (set); //reads rv
1233
1234        UI::get (mu, set, "mu", UI::compulsory);
1235        mat Rtmp;// necessary for conversion
1236        UI::get (Rtmp, set, "R", UI::compulsory);
1237        R = Rtmp; // conversion
1238        validate();
1239}
1240
1241template<class sq_T>
1242void enorm<sq_T>::dupdate (mat &v, double nu)
1243{
1244        //
1245};
1246
1247// template<class sq_T>
1248// void enorm<sq_T>::tupdate ( double phi, mat &vbar, double nubar ) {
1249//      //
1250// };
1251
1252template<class sq_T>
1253vec enorm<sq_T>::sample() const
1254{
1255        vec x (dim);
1256#pragma omp critical
1257        NorRNG.sample_vector (dim, x);
1258        vec smp = R.sqrt_mult (x);
1259
1260        smp += mu;
1261        return smp;
1262};
1263
1264// template<class sq_T>
1265// double enorm<sq_T>::eval ( const vec &val ) const {
1266//      double pdfl,e;
1267//      pdfl = evallog ( val );
1268//      e = exp ( pdfl );
1269//      return e;
1270// };
1271
1272template<class sq_T>
1273double enorm<sq_T>::evallog_nn (const vec &val) const
1274{
1275        // 1.83787706640935 = log(2pi)
1276        double tmp = -0.5 * (R.invqform (mu - val));// - lognc();
1277        return  tmp;
1278};
1279
1280template<class sq_T>
1281inline double enorm<sq_T>::lognc () const
1282{
1283        // 1.83787706640935 = log(2pi)
1284        double tmp = 0.5 * (R.cols() * 1.83787706640935 + R.logdet());
1285        return tmp;
1286};
1287
1288
1289// template<class sq_T>
1290// vec mlnorm<sq_T>::samplecond (const  vec &cond, double &lik ) {
1291//      this->condition ( cond );
1292//      vec smp = epdf.sample();
1293//      lik = epdf.eval ( smp );
1294//      return smp;
1295// }
1296
1297// template<class sq_T>
1298// mat mlnorm<sq_T>::samplecond (const vec &cond, vec &lik, int n ) {
1299//      int i;
1300//      int dim = rv.count();
1301//      mat Smp ( dim,n );
1302//      vec smp ( dim );
1303//      this->condition ( cond );
1304//
1305//      for ( i=0; i<n; i++ ) {
1306//              smp = epdf.sample();
1307//              lik ( i ) = epdf.eval ( smp );
1308//              Smp.set_col ( i ,smp );
1309//      }
1310//
1311//      return Smp;
1312// }
1313
1314
1315template<class sq_T>
1316shared_ptr<epdf> enorm<sq_T>::marginal ( const RV &rvn ) const
1317{
1318        enorm<sq_T> *tmp = new enorm<sq_T> ();
1319        shared_ptr<epdf> narrow(tmp);
1320        marginal ( rvn, *tmp );
1321        return narrow;
1322}
1323
1324template<class sq_T>
1325void enorm<sq_T>::marginal ( const RV &rvn, enorm<sq_T> &target ) const
1326{
1327        bdm_assert (isnamed(), "rv description is not assigned");
1328        ivec irvn = rvn.dataind (rv);
1329
1330        sq_T Rn (R, irvn);  // select rows and columns of R
1331
1332        target.set_rv ( rvn );
1333        target.set_parameters (mu (irvn), Rn);
1334}
1335
1336template<class sq_T>
1337shared_ptr<mpdf> enorm<sq_T>::condition ( const RV &rvn ) const
1338{
1339        mlnorm<sq_T> *tmp = new mlnorm<sq_T> ();
1340        shared_ptr<mpdf> narrow(tmp);
1341        condition ( rvn, *tmp );
1342        return narrow;
1343}
1344
1345template<class sq_T>
1346void enorm<sq_T>::condition ( const RV &rvn, mpdf &target ) const
1347{
1348        typedef mlnorm<sq_T> TMlnorm;
1349
1350        bdm_assert (isnamed(), "rvs are not assigned");
1351        TMlnorm &uptarget = dynamic_cast<TMlnorm &>(target);
1352
1353        RV rvc = rv.subt (rvn);
1354        bdm_assert ( (rvc._dsize() + rvn._dsize() == rv._dsize()), "wrong rvn");
1355        //Permutation vector of the new R
1356        ivec irvn = rvn.dataind (rv);
1357        ivec irvc = rvc.dataind (rv);
1358        ivec perm = concat (irvn , irvc);
1359        sq_T Rn (R, perm);
1360
1361        //fixme - could this be done in general for all sq_T?
1362        mat S = Rn.to_mat();
1363        //fixme
1364        int n = rvn._dsize() - 1;
1365        int end = R.rows() - 1;
1366        mat S11 = S.get (0, n, 0, n);
1367        mat S12 = S.get (0, n , rvn._dsize(), end);
1368        mat S22 = S.get (rvn._dsize(), end, rvn._dsize(), end);
1369
1370        vec mu1 = mu (irvn);
1371        vec mu2 = mu (irvc);
1372        mat A = S12 * inv (S22);
1373        sq_T R_n (S11 - A *S12.T());
1374
1375        uptarget.set_rv (rvn);
1376        uptarget.set_rvc (rvc);
1377        uptarget.set_parameters (A, mu1 - A*mu2, R_n);
1378}
1379
1380////
1381///////
1382template<class sq_T>
1383void mgnorm<sq_T >::set_parameters (const shared_ptr<fnc> &g0, const sq_T &R0) {
1384        g = g0;
1385        this->iepdf.set_parameters (zeros (g->dimension()), R0);
1386}
1387
1388template<class sq_T>
1389void mgnorm<sq_T >::condition (const vec &cond) {this->iepdf._mu() = g->eval (cond);};
1390
1391//! \todo unify this stuff with to_string()
1392template<class sq_T>
1393std::ostream &operator<< (std::ostream &os,  mlnorm<sq_T> &ml)
1394{
1395        os << "A:" << ml.A << endl;
1396        os << "mu:" << ml.mu_const << endl;
1397        os << "R:" << ml._R() << endl;
1398        return os;
1399};
1400
1401}
1402#endif //EF_H
Note: See TracBrowser for help on using the browser.