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

Revision 900, 46.3 kB (checked in by smidl, 14 years ago)

particle bug fixing

  • 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//! Global Uniform_RNG
25extern Uniform_RNG UniRNG;
26//! Global Normal_RNG
27extern Normal_RNG NorRNG;
28//! Global Gamma_RNG
29extern Gamma_RNG GamRNG;
30
31/*!
32* \brief General conjugate exponential family posterior density.
33
34* More?...
35*/
36
37class eEF : public epdf {
38public:
39//      eEF() :epdf() {};
40        //! default constructor
41        eEF () : epdf () {};
42        //! logarithm of the normalizing constant, \f$\mathcal{I}\f$
43        virtual double lognc() const = 0;
44
45        //!Evaluate normalized log-probability
46        virtual double evallog_nn ( const vec &val ) const NOT_IMPLEMENTED(0);
47
48        //!Evaluate normalized log-probability
49        virtual double evallog ( const vec &val ) const {
50                double tmp;
51                tmp = evallog_nn ( val ) - lognc();
52                return tmp;
53        }
54        //!Evaluate normalized log-probability for many samples
55        virtual vec evallog_mat ( const mat &Val ) const {
56                vec x ( Val.cols() );
57                for ( int i = 0; i < Val.cols(); i++ ) {
58                        x ( i ) = evallog_nn ( Val.get_col ( i ) ) ;
59                }
60                return x - lognc();
61        }
62        //!Evaluate normalized log-probability for many samples
63        virtual vec evallog_mat ( const Array<vec> &Val ) const {
64                vec x ( Val.length() );
65                for ( int i = 0; i < Val.length(); i++ ) {
66                        x ( i ) = evallog_nn ( Val ( i ) ) ;
67                }
68                return x - lognc();
69        }
70
71        //!Power of the density, used e.g. to flatten the density
72        virtual void pow ( double p ) NOT_IMPLEMENTED_VOID;
73};
74
75
76//! Estimator for Exponential family
77class BMEF : public BM {
78protected:
79        //! forgetting factor
80        double frg;
81        //! cached value of lognc() in the previous step (used in evaluation of \c ll )
82        double last_lognc;
83public:
84        //! Default constructor (=empty constructor)
85        BMEF ( double frg0 = 1.0 ) : BM (), frg ( frg0 ) {}
86        //! Copy constructor
87        BMEF ( const BMEF &B ) : BM ( B ), frg ( B.frg ), last_lognc ( B.last_lognc ) {}
88        //!get statistics from another model
89        virtual void set_statistics ( const BMEF* BM0 ) NOT_IMPLEMENTED_VOID;
90
91        //! Weighted update of sufficient statistics (Bayes rule)
92        virtual void bayes_weighted ( const vec &data, const vec &cond = empty_vec, const double w = 1.0 ) {};
93        //original Bayes
94        void bayes ( const vec &yt, const vec &cond = empty_vec );
95
96        //!Flatten the posterior according to the given BMEF (of the same type!)
97        virtual void flatten ( const BMEF * B ) NOT_IMPLEMENTED_VOID;
98
99
100        void to_setting ( Setting &set ) const
101        {                       
102                BM::to_setting( set );
103                UI::save(frg, set, "frg");
104                // TODO DOPLNIT? ALE MOMENTALNE CHYBI FROM_SETTING PRO INSPIRACI
105        } 
106        void from_setting( const Setting &set) {
107                BM::from_setting(set);
108                if ( !UI::get ( frg, set, "frg" ) )
109                        frg = 1.0;
110       
111        }
112
113        void validate() {
114                BM::validate(); 
115        }
116
117};
118
119/*! Dirac delta density with predefined transformation
120
121Density of the type:\f[ f(x_t | y_t) = \delta (x_t - g(y_t)) \f]
122where \f$ x_t \f$ is the \c rv, \f$ y_t \f$ is the \c rvc and g is a deterministic transformation of class fn.
123*/
124class mgdirac: public pdf{
125        protected:
126        shared_ptr<fnc> g;
127        public:
128                vec samplecond(const vec &cond) {
129                        bdm_assert_debug(cond.length()==g->dimensionc(),"given cond in not compatible with g");
130                        vec tmp = g->eval(cond);
131                        return tmp;
132                } 
133                double evallogcond ( const vec &yt, const vec &cond ){
134                        return std::numeric_limits< double >::max();
135                }
136                void from_setting(const Setting& set){
137                        pdf::from_setting(set);
138                        g=UI::build<fnc>(set,"g",UI::compulsory);
139                        validate();
140                }
141                void to_setting(Setting &set) const{
142                        pdf::to_setting(set);
143                        UI::save(g.get(), set, "g");
144                }
145                void validate() {
146                        dim = g->dimension();
147                        dimc = g->dimensionc();
148                }
149};
150UIREGISTER(mgdirac);
151
152
153template<class sq_T, template <typename> class TEpdf>
154class mlnorm;
155
156/*!
157* \brief Gaussian density with positive definite (decomposed) covariance matrix.
158
159* More?...
160*/
161template<class sq_T>
162class enorm : public eEF {
163protected:
164        //! mean value
165        vec mu;
166        //! Covariance matrix in decomposed form
167        sq_T R;
168public:
169        //!\name Constructors
170        //!@{
171
172        enorm () : eEF (), mu (), R () {};
173        enorm ( const vec &mu, const sq_T &R ) {
174                set_parameters ( mu, R );
175        }
176        void set_parameters ( const vec &mu, const sq_T &R );
177        /*! Create Normal density
178        \f[ f(rv) = N(\mu, R) \f]
179        from structure
180        \code
181        class = 'enorm<ldmat>', (OR) 'enorm<chmat>', (OR) 'enorm<fsqmat>';
182        mu    = [];                  // mean value
183        R     = [];                  // variance, square matrix of appropriate dimension
184        \endcode
185        */
186        void from_setting ( const Setting &root );
187        void to_setting ( Setting &root ) const ;
188       
189        void validate() {
190                bdm_assert ( mu.length() == R.rows(), "mu and R parameters do not match" );
191                dim = mu.length();
192        }
193        //!@}
194
195        //! \name Mathematical operations
196        //!@{
197
198        //! dupdate in exponential form (not really handy)
199        void dupdate ( mat &v, double nu = 1.0 );
200
201        //! evaluate bhattacharya distance
202        double bhattacharyya(const enorm<sq_T> &e2){
203                bdm_assert(dim == e2.dimension(), "enorms of differnt dimensions");
204                sq_T P=R;
205                P.add(e2._R());
206               
207                double tmp = 0.125*P.invqform(mu - e2._mu()) + 0.5*(P.logdet() - 0.5*(R.logdet() + e2._R().logdet()));
208                return tmp;
209        }
210       
211        vec sample() const;
212
213        double evallog_nn ( const vec &val ) const;
214        double lognc () const;
215        vec mean() const {
216                return mu;
217        }
218        vec variance() const {
219                return diag ( R.to_mat() );
220        }
221        mat covariance() const {
222                return R.to_mat();
223        }
224        //      mlnorm<sq_T>* condition ( const RV &rvn ) const ; <=========== fails to cmpile. Why?
225        shared_ptr<pdf> condition ( const RV &rvn ) const;
226
227        // target not typed to mlnorm<sq_T, enorm<sq_T> > &
228        // because that doesn't compile (perhaps because we
229        // haven't finished defining enorm yet), but the type
230        // is required
231        void condition ( const RV &rvn, pdf &target ) const;
232
233        shared_ptr<epdf> marginal ( const RV &rvn ) const;
234        void marginal ( const RV &rvn, enorm<sq_T> &target ) const;
235        //!@}
236
237        //! \name Access to attributes
238        //!@{
239
240        vec& _mu() {
241                return mu;
242        }
243        const vec& _mu() const {
244                return mu;
245        }
246        void set_mu ( const vec mu0 ) {
247                mu = mu0;
248        }
249        sq_T& _R() {
250                return R;
251        }
252        const sq_T& _R() const {
253                return R;
254        }
255        //!@}
256
257};
258UIREGISTER2 ( enorm, chmat );
259SHAREDPTR2 ( enorm, chmat );
260UIREGISTER2 ( enorm, ldmat );
261SHAREDPTR2 ( enorm, ldmat );
262UIREGISTER2 ( enorm, fsqmat );
263SHAREDPTR2 ( enorm, fsqmat );
264
265typedef enorm<ldmat> egauss;
266UIREGISTER(egauss);
267
268
269//forward declaration
270class mstudent;
271
272/*! distribution of multivariate Student t density
273
274Based on article by Genest and Zidek,
275*/
276template<class sq_T>
277class estudent : public eEF{
278        protected:
279                //! mena value
280                vec mu;
281                //! matrix H
282                sq_T H;
283                //! degrees of freedom
284                double delta;
285        public:
286                double evallog_nn(const vec &val) const{
287                        double tmp = -0.5*H.logdet() - 0.5*(delta + dim) * log(1+ H.invqform(val - mu)/delta);
288                        return tmp;
289                }
290                double lognc() const {
291                        //log(pi) = 1.14472988584940
292                        double tmp = -lgamma(0.5*(delta+dim))+lgamma(0.5*delta) + 0.5*dim*(log(delta) + 1.14472988584940);
293                        return tmp;
294                }
295                void marginal (const RV &rvm, estudent<sq_T> &marg) const {
296                        ivec ind = rvm.findself_ids(rv); // indices of rvm in rv
297                        marg._mu() = mu(ind);
298                        marg._H() = sq_T(H,ind);
299                        marg._delta() = delta;
300                        marg.validate();
301                }
302                shared_ptr<epdf> marginal(const RV &rvm) const {
303                        shared_ptr<estudent<sq_T> > tmp = new estudent<sq_T>;
304                        marginal(rvm, *tmp);
305                        return tmp;
306                }
307                vec sample() const NOT_IMPLEMENTED(vec(0))
308               
309                vec mean() const {return mu;}
310                mat covariance() const {return delta/(delta-2)*H.to_mat();}
311                vec variance() const {return diag(covariance());}
312                        //! \name access
313                //! @{
314                        //! access function
315                        vec& _mu() {return mu;}
316                        //! access function
317                        sq_T& _H() {return H;}
318                        //! access function
319                        double& _delta() {return delta;}
320                        //!@}
321                        //! todo
322                        void from_setting(const Setting &set){
323                                epdf::from_setting(set);
324                                mat H0;
325                                UI::get(H0,set, "H");
326                                H= H0; // conversion!!
327                                UI::get(delta,set,"delta");
328                                UI::get(mu,set,"mu");
329                        }
330                        void to_setting(Setting &set) const{
331                                epdf::to_setting(set);
332                                UI::save(H.to_mat(), set, "H");
333                                UI::save(delta, set, "delta");
334                                UI::save(mu, set, "mu");
335                        }
336                        void validate() {
337                                dim = H.rows();
338                        }
339};
340UIREGISTER2(estudent,fsqmat);
341UIREGISTER2(estudent,ldmat);
342UIREGISTER2(estudent,chmat);
343
344/*!
345* \brief Gauss-inverse-Wishart density stored in LD form
346
347* For \f$p\f$-variate densities, given rv.count() should be \f$p\times\f$ V.rows().
348*
349*/
350class egiw : public eEF {
351        //! \var log_level_enums tri
352        //! TODO DOPLNIT
353        LOG_LEVEL(egiw,tri);
354
355protected:
356        //! Extended information matrix of sufficient statistics
357        ldmat V;
358        //! Number of data records (degrees of freedom) of sufficient statistics
359        double nu;
360        //! Dimension of the output
361        int dimx;
362        //! Dimension of the regressor
363        int nPsi;
364public:
365        //!\name Constructors
366        //!@{
367        egiw() : eEF() {};
368        egiw ( int dimx0, ldmat V0, double nu0 = -1.0 ) : eEF() {
369                set_parameters ( dimx0, V0, nu0 );
370                validate();
371        };
372
373        void set_parameters ( int dimx0, ldmat V0, double nu0 = -1.0 );
374        //!@}
375
376        vec sample() const;
377        mat sample_mat ( int n ) const;
378        vec mean() const;
379        vec variance() const;
380        void sample_mat ( mat &Mi, chmat &Ri ) const;
381
382        void factorize ( mat &M, ldmat &Vz, ldmat &Lam ) const;
383        //! LS estimate of \f$\theta\f$
384        vec est_theta() const;
385
386        //! Covariance of the LS estimate
387        ldmat est_theta_cov() const;
388
389        //! expected values of the linear coefficient and the covariance matrix are written to \c M and \c R , respectively
390        void mean_mat ( mat &M, mat&R ) const;
391        //! In this instance, val= [theta, r]. For multivariate instances, it is stored columnwise val = [theta_1 theta_2 ... r_1 r_2 ]
392        double evallog_nn ( const vec &val ) const;
393        double lognc () const;
394        void pow ( double p ) {
395                V *= p;
396                nu *= p;
397        };
398
399        //! marginal density (only student for now)
400        shared_ptr<epdf> marginal(const RV &rvm) const {
401                bdm_assert(dimx==1, "Not supported");
402                //TODO - this is too trivial!!!
403                ivec ind = rvm.findself_ids(rv);
404                if (min(ind)==0) { //assume it si
405                        shared_ptr<estudent<ldmat> > tmp = new estudent<ldmat>;
406                        mat M;
407                        ldmat Vz;
408                        ldmat Lam;
409                        factorize(M,Vz,Lam);
410                       
411                        tmp->_mu() = M.get_col(0);
412                        ldmat H;
413                        Vz.inv(H);
414                        H *=Lam._D()(0)/nu;
415                        tmp->_H() = H;
416                        tmp->_delta() = nu;
417                        tmp->validate();
418                        return tmp;
419                }
420                return NULL;
421        }
422        //! \name Access attributes
423        //!@{
424
425        ldmat& _V() {
426                return V;
427        }
428        const ldmat& _V() const {
429                return V;
430        }
431        double& _nu()  {
432                return nu;
433        }
434        const double& _nu() const {
435                return nu;
436        }
437        const int & _dimx() const {
438                return dimx;
439        }
440
441        /*! Create Gauss-inverse-Wishart density
442        \f[ f(rv) = GiW(V,\nu) \f]
443        from structure
444        \code
445        class = 'egiw';
446        V     = [];               // square matrix
447        dV    = [];               // vector of diagonal of V (when V not given)
448        nu    = [];               // scalar \nu ((almost) degrees of freedom)
449                                                          // when missing, it will be computed to obtain proper pdf
450        dimx  = [];               // dimension of the wishart part
451        rv = RV({'name'})         // description of RV
452        rvc = RV({'name'})        // description of RV in condition
453        log_level = 'tri';                // set the level of logged details
454        \endcode
455
456        \sa log_level_enums
457        */
458        void from_setting ( const Setting &set ) {
459                epdf::from_setting ( set );
460                UI::get ( dimx, set, "dimx", UI::compulsory );
461                if ( !UI::get ( nu, set, "nu", UI::optional ) ) {
462                        nu = -1;
463                }
464                mat V;
465                if ( !UI::get ( V, set, "V", UI::optional ) ) {
466                        vec dV;
467                        UI::get ( dV, set, "dV", UI::compulsory );
468                        set_parameters ( dimx, ldmat ( dV ), nu );
469                        validate();
470
471                } else {
472                        set_parameters ( dimx, V, nu );
473                        validate();
474                }
475        }
476
477        void to_setting ( Setting& set ) const {
478                epdf::to_setting ( set );
479                UI::save ( dimx, set, "dimx" );
480                UI::save ( V.to_mat(), set, "V" );
481                UI::save ( nu, set, "nu" );
482        };
483
484        void validate() {
485        dim = dimx * ( dimx + nPsi );   
486       
487            // check sizes, rvs etc.
488                // also check if RV are meaningful!!!
489                // meaningful =  rv for theta  and rv for r are split!
490        }
491        void log_register ( bdm::logger& L, const string& prefix );
492
493        void log_write() const;
494        //!@}
495};
496UIREGISTER ( egiw );
497SHAREDPTR ( egiw );
498
499/*! \brief Dirichlet posterior density
500
501Continuous Dirichlet density of \f$n\f$-dimensional variable \f$x\f$
502\f[
503f(x|\beta) = \frac{\Gamma[\gamma]}{\prod_{i=1}^{n}\Gamma(\beta_i)} \prod_{i=1}^{n}x_i^{\beta_i-1}
504\f]
505where \f$\gamma=\sum_i \beta_i\f$.
506*/
507class eDirich: public eEF {
508protected:
509        //!sufficient statistics
510        vec beta;
511public:
512        //!\name Constructors
513        //!@{
514
515        eDirich () : eEF () {};
516        eDirich ( const eDirich &D0 ) : eEF () {
517                set_parameters ( D0.beta );
518                validate();
519        };
520        eDirich ( const vec &beta0 ) {
521                set_parameters ( beta0 );
522                validate();
523        };
524        void set_parameters ( const vec &beta0 ) {
525                beta = beta0;
526                dim = beta.length();
527        }
528        //!@}
529
530        //! using sampling procedure from wikipedia
531        vec sample() const {
532                vec y ( beta.length() );
533                for ( int i = 0; i < beta.length(); i++ ) {
534                        GamRNG.setup ( beta ( i ), 1 );
535#pragma omp critical
536                        y ( i ) = GamRNG();
537                }
538                return y / sum ( y );
539        }
540
541        vec mean() const {
542                return beta / sum ( beta );
543        };
544        vec variance() const {
545                double gamma = sum ( beta );
546                return elem_mult ( beta, ( gamma - beta ) ) / ( gamma*gamma* ( gamma + 1 ) );
547        }
548        //! In this instance, val is ...
549        double evallog_nn ( const vec &val ) const {
550                double tmp;
551                tmp = ( beta - 1 ) * log ( val );
552                return tmp;
553        }
554
555        double lognc () const {
556                double tmp;
557                double gam = sum ( beta );
558                double lgb = 0.0;
559                for ( int i = 0; i < beta.length(); i++ ) {
560                        lgb += lgamma ( beta ( i ) );
561                }
562                tmp = lgb - lgamma ( gam );
563                return tmp;
564        }
565
566        //!access function
567        vec& _beta()  {
568                return beta;
569        }
570        /*! configuration structure
571        \code
572        class = 'eDirich';
573        beta  = [];           //parametr beta
574        \endcode
575        */
576        void from_setting ( const Setting &set ) {
577                epdf::from_setting ( set );
578                UI::get ( beta, set, "beta", UI::compulsory );
579                validate();
580        }
581        void validate() {
582                //check rv
583                dim = beta.length();
584        }
585
586        void to_setting ( Setting &set ) const
587        {                       
588                eEF::to_setting( set );
589                UI::save( beta, set, "beta" );
590        } 
591};
592UIREGISTER ( eDirich );
593
594/*! Random Walk on Dirichlet
595Using simple assignment
596 \f[ \beta = rvc / k + \beta_c \f]
597 hence, mean value = rvc, variance = (k+1)*mean*mean;
598
599 The greater k is, the greater is the variance of the random walk;
600
601 \f$ \beta_c \f$ is used as regularizing element to avoid corner cases, i.e. when one element of rvc is zero.
602 By default is it set to 0.1;
603*/
604
605class mDirich: public pdf_internal<eDirich> {
606protected:
607        //! constant \f$ k \f$ of the random walk
608        double k;
609        //! cache of beta_i
610        vec &_beta;
611        //! stabilizing coefficient \f$ \beta_c \f$
612        vec betac;
613public:
614        mDirich() : pdf_internal<eDirich>(), _beta ( iepdf._beta() ) {};
615        void condition ( const vec &val ) {
616                _beta =  val / k + betac;
617        };
618        /*! Create Dirichlet random walk
619        \f[ f(rv|rvc) = Di(rvc*k) \f]
620        from structure
621        \code
622        class = 'mDirich';
623        k = 1;                      // multiplicative constant k
624        --- optional ---
625        rv = RV({'name'},size)      // description of RV
626        beta0 = [];                 // initial value of beta
627        betac = [];                 // initial value of beta
628        \endcode
629        */
630        void from_setting ( const Setting &set ) {
631                pdf::from_setting ( set ); // reads rv and rvc
632                if ( _rv()._dsize() > 0 ) {
633                        rvc = _rv().copy_t ( -1 );
634                }
635                vec beta0;
636                if ( !UI::get ( beta0, set, "beta0", UI::optional ) ) {
637                        beta0 = ones ( _rv()._dsize() );
638                }
639                if ( !UI::get ( betac, set, "betac", UI::optional ) ) {
640                        betac = 0.1 * ones ( _rv()._dsize() );
641                }
642                _beta = beta0;
643
644                UI::get ( k, set, "k", UI::compulsory );
645                validate();
646        }
647        void validate() {
648                pdf_internal<eDirich>::validate();
649                bdm_assert ( _beta.length() == betac.length(), "beta0 and betac are not compatible" );
650                if ( _rv()._dsize() > 0 ) {
651                        bdm_assert ( ( _rv()._dsize() == dimension() ) , "Size of rv does not match with beta" );
652                }
653                dimc = _beta.length();
654        };
655};
656UIREGISTER ( mDirich );
657
658//! \brief Estimator for Multinomial density
659class multiBM : public BMEF {
660protected:
661        //! Conjugate prior and posterior
662        eDirich est;
663        //! Pointer inside est to sufficient statistics
664        vec &beta;
665public:
666        //!Default constructor
667        multiBM () : BMEF (), est (), beta ( est._beta() ) {
668                if ( beta.length() > 0 ) {
669                        last_lognc = est.lognc();
670                } else {
671                        last_lognc = 0.0;
672                }
673        }
674        //!Copy constructor
675        multiBM ( const multiBM &B ) : BMEF ( B ), est ( B.est ), beta ( est._beta() ) {}
676        //! Sets sufficient statistics to match that of givefrom mB0
677        void set_statistics ( const BM* mB0 ) {
678                const multiBM* mB = dynamic_cast<const multiBM*> ( mB0 );
679                beta = mB->beta;
680        }
681        void bayes ( const vec &yt, const vec &cond = empty_vec );
682
683        double logpred ( const vec &yt ) const;
684
685        void flatten ( const BMEF* B );
686
687        //! return correctly typed posterior (covariant return)
688        const eDirich& posterior() const {
689                return est;
690        };
691        //! constructor function
692        void set_parameters ( const vec &beta0 ) {
693                est.set_parameters ( beta0 );
694                est.validate();
695                if ( evalll ) {
696                        last_lognc = est.lognc();
697                }
698        }
699
700        void to_setting ( Setting &set ) const {
701                BMEF::to_setting ( set );
702                UI::save( &est, set, "prior" );
703        }
704};
705UIREGISTER( multiBM );
706
707/*!
708 \brief Gamma posterior density
709
710 Multivariate Gamma density as product of independent univariate densities.
711 \f[
712 f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
713 \f]
714*/
715
716class egamma : public eEF {
717protected:
718        //! Vector \f$\alpha\f$
719        vec alpha;
720        //! Vector \f$\beta\f$
721        vec beta;
722public :
723        //! \name Constructors
724        //!@{
725        egamma () : eEF (), alpha ( 0 ), beta ( 0 ) {};
726        egamma ( const vec &a, const vec &b ) {
727                set_parameters ( a, b );
728                validate();
729        };
730        void set_parameters ( const vec &a, const vec &b ) {
731                alpha = a, beta = b;
732        };
733        //!@}
734
735        vec sample() const;
736        double evallog ( const vec &val ) const;
737        double lognc () const;
738        //! Returns pointer to internal alpha. Potentially dengerous: use with care!
739        vec& _alpha() {
740                return alpha;
741        }
742        //! Returns pointer to internal beta. Potentially dengerous: use with care!
743        vec& _beta() {
744                return beta;
745        }
746        vec mean() const {
747                return elem_div ( alpha, beta );
748        }
749        vec variance() const {
750                return elem_div ( alpha, elem_mult ( beta, beta ) );
751        }
752
753        /*! Create Gamma density
754        \f[ f(rv|rvc) = \Gamma(\alpha, \beta) \f]
755        from structure
756        \code
757        class = 'egamma';
758        alpha = [...];         // vector of alpha
759        beta = [...];          // vector of beta
760        rv = RV({'name'})      // description of RV
761        \endcode
762        */
763        void from_setting ( const Setting &set ) {
764                epdf::from_setting ( set ); // reads rv
765                UI::get ( alpha, set, "alpha", UI::compulsory );
766                UI::get ( beta, set, "beta", UI::compulsory );
767                validate();
768        }
769        void validate() {
770                bdm_assert ( alpha.length() == beta.length(), "parameters do not match" );
771                dim = alpha.length();
772        }
773};
774UIREGISTER ( egamma );
775SHAREDPTR ( egamma );
776
777/*!
778 \brief Inverse-Gamma posterior density
779
780 Multivariate inverse-Gamma density as product of independent univariate densities.
781 \f[
782 f(x|\alpha,\beta) = \prod f(x_i|\alpha_i,\beta_i)
783 \f]
784
785Vector \f$\beta\f$ has different meaning (in fact it is 1/beta as used in definition of iG)
786
787 Inverse Gamma can be converted to Gamma using \f[
788 x\sim iG(a,b) => 1/x\sim G(a,1/b)
789\f]
790This relation is used in sampling.
791 */
792
793class eigamma : public egamma {
794protected:
795public :
796        //! \name Constructors
797        //! All constructors are inherited
798        //!@{
799        //!@}
800
801        vec sample() const {
802                return 1.0 / egamma::sample();
803        };
804        //! Returns poiter to alpha and beta. Potentially dangerous: use with care!
805        vec mean() const {
806                return elem_div ( beta, alpha - 1 );
807        }
808        vec variance() const {
809                vec mea = mean();
810                return elem_div ( elem_mult ( mea, mea ), alpha - 2 );
811        }
812};
813/*
814//! Weighted mixture of epdfs with external owned components.
815class emix : public epdf {
816protected:
817        int n;
818        vec &w;
819        Array<epdf*> Coms;
820public:
821//! Default constructor
822        emix ( const RV &rv, vec &w0): epdf(rv), n(w0.length()), w(w0), Coms(n) {};
823        void set_parameters( int &i, double wi, epdf* ep){w(i)=wi;Coms(i)=ep;}
824        vec mean(){vec pom; for(int i=0;i<n;i++){pom+=Coms(i)->mean()*w(i);} return pom;};
825};
826*/
827
828//!  Uniform distributed density on a rectangular support
829
830class euni: public epdf {
831protected:
832//! lower bound on support
833        vec low;
834//! upper bound on support
835        vec high;
836//! internal
837        vec distance;
838//! normalizing coefficients
839        double nk;
840//! cache of log( \c nk )
841        double lnk;
842public:
843        //! \name Constructors
844        //!@{
845        euni () : epdf () {}
846        euni ( const vec &low0, const vec &high0 ) {
847                set_parameters ( low0, high0 );
848        }
849        void set_parameters ( const vec &low0, const vec &high0 ) {
850                distance = high0 - low0;
851                low = low0;
852                high = high0;
853                nk = prod ( 1.0 / distance );
854                lnk = log ( nk );
855        }
856        //!@}
857
858        double evallog ( const vec &val ) const  {
859                if ( any ( val < low ) && any ( val > high ) ) {
860                        return -inf;
861                } else return lnk;
862        }
863        vec sample() const {
864                vec smp ( dim );
865#pragma omp critical
866                UniRNG.sample_vector ( dim , smp );
867                return low + elem_mult ( distance, smp );
868        }
869        //! set values of \c low and \c high
870        vec mean() const {
871                return ( high - low ) / 2.0;
872        }
873        vec variance() const {
874                return ( pow ( high, 2 ) + pow ( low, 2 ) + elem_mult ( high, low ) ) / 3.0;
875        }
876        /*! Create Uniform density
877        \f[ f(rv) = U(low,high) \f]
878        from structure
879         \code
880         class = 'euni'
881         high = [...];          // vector of upper bounds
882         low = [...];           // vector of lower bounds
883         rv = RV({'name'});     // description of RV
884         \endcode
885         */
886        void from_setting ( const Setting &set ) {
887                epdf::from_setting ( set ); // reads rv and rvc
888
889                UI::get ( high, set, "high", UI::compulsory );
890                UI::get ( low, set, "low", UI::compulsory );
891                set_parameters ( low, high );
892                validate();
893        }
894        void validate() {
895                bdm_assert ( high.length() == low.length(), "Incompatible high and low vectors" );
896                dim = high.length();
897                bdm_assert ( min ( distance ) > 0.0, "bad support" );
898        }
899};
900UIREGISTER ( euni );
901
902//! Uniform density with conditional mean value
903class mguni : public pdf_internal<euni> {
904        //! function of the mean value
905        shared_ptr<fnc> mean;
906        //! distance from mean to both sides
907        vec delta;
908public:
909        void condition ( const vec &cond ) {
910                vec mea = mean->eval ( cond );
911                iepdf.set_parameters ( mea - delta, mea + delta );
912        }
913        //! load from
914        void from_setting ( const Setting &set ) {
915                pdf::from_setting ( set ); //reads rv and rvc
916                UI::get ( delta, set, "delta", UI::compulsory );
917                mean = UI::build<fnc> ( set, "mean", UI::compulsory );
918                iepdf.set_parameters ( -delta, delta );
919        }
920       
921        void validate(){
922                dimc = mean->dimensionc();
923                iepdf.validate();
924        }
925
926};
927UIREGISTER ( mguni );
928/*!
929 \brief Normal distributed linear function with linear function of mean value;
930
931 Mean value \f$ \mu=A*\mbox{rvc}+\mu_0 \f$.
932*/
933template < class sq_T, template <typename> class TEpdf = enorm >
934class mlnorm : public pdf_internal< TEpdf<sq_T> > {
935protected:
936        //! Internal epdf that arise by conditioning on \c rvc
937        mat A;
938        //! Constant additive term
939        vec mu_const;
940//                      vec& _mu; //cached epdf.mu; !!!!!! WHY NOT?
941public:
942        //! \name Constructors
943        //!@{
944        mlnorm() : pdf_internal< TEpdf<sq_T> >() {};
945        mlnorm ( const mat &A, const vec &mu0, const sq_T &R ) : pdf_internal< TEpdf<sq_T> >() {
946                set_parameters ( A, mu0, R );
947                validate();
948        }
949
950        //! Set \c A and \c R
951        void set_parameters ( const  mat &A0, const vec &mu0, const sq_T &R0 ) {
952                this->iepdf.set_parameters ( zeros ( A0.rows() ), R0 );
953                A = A0;
954                mu_const = mu0;
955        }
956
957        //!@}
958        //! Set value of \c rvc . Result of this operation is stored in \c epdf use function \c _ep to access it.
959        void condition ( const vec &cond ) {
960                this->iepdf._mu() = A * cond + mu_const;
961//R is already assigned;
962        }
963
964        //!access function
965        const vec& _mu_const() const {
966                return mu_const;
967        }
968        //!access function
969        const mat& _A() const {
970                return A;
971        }
972        //!access function
973        mat _R() const {
974                return this->iepdf._R().to_mat();
975        }
976        //!access function
977        sq_T __R() const {
978                return this->iepdf._R();
979        }
980
981        //! Debug stream
982        template<typename sq_M>
983        friend std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_M, enorm> &ml );
984
985        /*! Create Normal density with linear function of mean value
986         \f[ f(rv|rvc) = N(A*rvc+const, R) \f]
987        from structure
988        \code
989        class = 'mlnorm<ldmat>', (OR) 'mlnorm<chmat>', (OR) 'mlnorm<fsqmat>';
990        A     = [];                  // matrix or vector of appropriate dimension
991        const = [];                  // vector of constant term
992        R     = [];                  // square matrix of appropriate dimension
993        \endcode
994        */
995        void from_setting ( const Setting &set ) {
996                pdf::from_setting ( set );
997
998                UI::get ( A, set, "A", UI::compulsory );
999                UI::get ( mu_const, set, "const", UI::compulsory );
1000                mat R0;
1001                UI::get ( R0, set, "R", UI::compulsory );
1002                set_parameters ( A, mu_const, R0 );
1003                validate();
1004        };
1005        void to_setting (Setting &set) const {
1006                pdf::to_setting(set);
1007                UI::save ( A, set, "A");
1008                UI::save ( mu_const, set, "const");
1009                UI::save ( _R(), set, "R");
1010        }
1011        void validate() {
1012                pdf_internal<TEpdf<sq_T> >::validate();
1013                bdm_assert ( A.rows() == mu_const.length(), "mlnorm: A vs. mu mismatch" );
1014                bdm_assert ( A.rows() == _R().rows(), "mlnorm: A vs. R mismatch" );
1015                this->dimc = A.cols();
1016
1017        }
1018};
1019UIREGISTER2 ( mlnorm, ldmat );
1020SHAREDPTR2 ( mlnorm, ldmat );
1021UIREGISTER2 ( mlnorm, fsqmat );
1022SHAREDPTR2 ( mlnorm, fsqmat );
1023UIREGISTER2 ( mlnorm, chmat );
1024SHAREDPTR2 ( mlnorm, chmat );
1025
1026//! pdf with general function for mean value
1027template<class sq_T>
1028class mgnorm : public pdf_internal< enorm< sq_T > > {
1029private:
1030//                      vec &mu; WHY NOT?
1031        shared_ptr<fnc> g;
1032
1033public:
1034        //!default constructor
1035        mgnorm() : pdf_internal<enorm<sq_T> >() { }
1036        //!set mean function
1037        inline void set_parameters ( const shared_ptr<fnc> &g0, const sq_T &R0 );
1038        inline void condition ( const vec &cond );
1039
1040
1041        /*! Create Normal density with given function of mean value
1042        \f[ f(rv|rvc) = N( g(rvc), R) \f]
1043        from structure
1044        \code
1045        class = 'mgnorm';
1046        g.class =  'fnc';      // function for mean value evolution
1047        g._fields_of_fnc = ...;
1048
1049        R = [1, 0;             // covariance matrix
1050                0, 1];
1051                --OR --
1052        dR = [1, 1];           // diagonal of cavariance matrix
1053
1054        rv = RV({'name'})      // description of RV
1055        rvc = RV({'name'})     // description of RV in condition
1056        \endcode
1057        */
1058
1059        void from_setting ( const Setting &set ) {
1060                pdf::from_setting ( set );
1061                shared_ptr<fnc> g = UI::build<fnc> ( set, "g", UI::compulsory );
1062
1063                mat R;
1064                vec dR;
1065                if ( UI::get ( dR, set, "dR" ) )
1066                        R = diag ( dR );
1067                else
1068                        UI::get ( R, set, "R", UI::compulsory );
1069
1070                set_parameters ( g, R );
1071                validate();
1072        }
1073        void validate() {
1074                bdm_assert ( g->dimension() == this->dimension(), "incompatible function" );
1075                this->dim = g->dimension();
1076                this->dimc = g->dimensionc();
1077                this->iepdf.validate();
1078        }
1079};
1080
1081UIREGISTER2 ( mgnorm, chmat );
1082SHAREDPTR2 ( mgnorm, chmat );
1083
1084
1085/*! (Approximate) Student t density with linear function of mean value
1086
1087The internal epdf of this class is of the type of a Gaussian (enorm).
1088However, each conditioning is trying to assure the best possible approximation by taking into account the zeta function. See [] for reference.
1089
1090Perhaps a moment-matching technique?
1091*/
1092class mlstudent : public mlnorm<ldmat, enorm> {
1093protected:
1094        //! Variable \f$ \Lambda \f$ from theory
1095        ldmat Lambda;
1096        //! Reference to variable \f$ R \f$
1097        ldmat &_R;
1098        //! Variable \f$ R_e \f$
1099        ldmat Re;
1100public:
1101        mlstudent () : mlnorm<ldmat, enorm> (),
1102                        Lambda (),      _R ( iepdf._R() ) {}
1103        //! constructor function
1104        void set_parameters ( const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0 ) {
1105                iepdf.set_parameters ( mu0, R0 );// was Lambda, why?
1106                A = A0;
1107                mu_const = mu0;
1108                Re = R0;
1109                Lambda = Lambda0;
1110        }
1111
1112        void condition ( const vec &cond ); 
1113
1114        void validate() {
1115                bdm_assert ( A.rows() == mu_const.length(), "mlstudent: A vs. mu mismatch" );
1116                bdm_assert ( _R.rows() == A.rows(), "mlstudent: A vs. R mismatch" );
1117
1118        }
1119};
1120
1121/*!
1122\brief  Gamma random walk
1123
1124Mean value, \f$\mu\f$, of this density is given by \c rvc .
1125Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
1126This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
1127
1128The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
1129*/
1130class mgamma : public pdf_internal<egamma> {
1131protected:
1132
1133        //! Constant \f$k\f$
1134        double k;
1135
1136        //! cache of iepdf.beta
1137        vec &_beta;
1138
1139public:
1140        //! Constructor
1141        mgamma() : pdf_internal<egamma>(), k ( 0 ),
1142                        _beta ( iepdf._beta() ) {
1143        }
1144
1145        //! Set value of \c k
1146        void set_parameters ( double k, const vec &beta0 );
1147
1148        void condition ( const vec &val ) {
1149                _beta = k / val;
1150        };
1151        /*! Create Gamma density with conditional mean value
1152        \f[ f(rv|rvc) = \Gamma(k, k/rvc) \f]
1153        from structure
1154        \code
1155          class = 'mgamma';
1156          beta = [...];          // vector of initial alpha
1157          k = 1.1;               // multiplicative constant k
1158          rv = RV({'name'})      // description of RV
1159          rvc = RV({'name'})     // description of RV in condition
1160         \endcode
1161        */
1162        void from_setting ( const Setting &set ) {
1163                pdf::from_setting ( set ); // reads rv and rvc
1164                vec betatmp; // ugly but necessary
1165                UI::get ( betatmp, set, "beta", UI::compulsory );
1166                UI::get ( k, set, "k", UI::compulsory );
1167                set_parameters ( k, betatmp );
1168                validate();
1169        }
1170        void validate() {
1171                pdf_internal<egamma>::validate();
1172
1173                dim = _beta.length();
1174                dimc = _beta.length();
1175        }
1176};
1177UIREGISTER ( mgamma );
1178SHAREDPTR ( mgamma );
1179
1180/*!
1181\brief  Inverse-Gamma random walk
1182
1183Mean value, \f$ \mu \f$, of this density is given by \c rvc .
1184Standard deviation of the random walk is proportional to one \f$ k \f$-th the mean.
1185This is achieved by setting \f$ \alpha=\mu/k^2+2 \f$ and \f$ \beta=\mu(\alpha-1)\f$.
1186
1187The standard deviation of the walk is then: \f$ \mu/\sqrt(k)\f$.
1188 */
1189class migamma : public pdf_internal<eigamma> {
1190protected:
1191        //! Constant \f$k\f$
1192        double k;
1193
1194        //! cache of iepdf.alpha
1195        vec &_alpha;
1196
1197        //! cache of iepdf.beta
1198        vec &_beta;
1199
1200public:
1201        //! \name Constructors
1202        //!@{
1203        migamma() : pdf_internal<eigamma>(),
1204                        k ( 0 ),
1205                        _alpha ( iepdf._alpha() ),
1206                        _beta ( iepdf._beta() ) {
1207        }
1208
1209        migamma ( const migamma &m ) : pdf_internal<eigamma>(),
1210                        k ( 0 ),
1211                        _alpha ( iepdf._alpha() ),
1212                        _beta ( iepdf._beta() ) {
1213        }
1214        //!@}
1215
1216        //! Set value of \c k
1217        void set_parameters ( int len, double k0 ) {
1218                k = k0;
1219                iepdf.set_parameters ( ( 1.0 / ( k*k ) + 2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ );               
1220        };
1221
1222        void    validate (){
1223                iepdf.validate();
1224                dimc = dimension();
1225};
1226
1227        void condition ( const vec &val ) {
1228                _beta = elem_mult ( val, ( _alpha - 1.0 ) );
1229        };
1230};
1231
1232
1233/*!
1234\brief  Gamma random walk around a fixed point
1235
1236Mean 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
1237\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
1238
1239Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
1240This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
1241
1242The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
1243*/
1244class mgamma_fix : public mgamma {
1245protected:
1246        //! parameter l
1247        double l;
1248        //! reference vector
1249        vec refl;
1250public:
1251        //! Constructor
1252        mgamma_fix () : mgamma (), refl () {};
1253        //! Set value of \c k
1254        void set_parameters ( double k0 , vec ref0, double l0 ) {
1255                mgamma::set_parameters ( k0, ref0 );
1256                refl = pow ( ref0, 1.0 - l0 );
1257                l = l0; 
1258        };
1259
1260        void    validate (){
1261                mgamma::validate();
1262                dimc = dimension();
1263        };
1264
1265        void condition ( const vec &val ) {
1266                vec mean = elem_mult ( refl, pow ( val, l ) );
1267                _beta = k / mean;
1268        };
1269};
1270
1271
1272/*!
1273\brief  Inverse-Gamma random walk around a fixed point
1274
1275Mean 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
1276\f[ \mu = \mu_{t-1} ^{l} p^{1-l}\f]
1277
1278==== Check == vv =
1279Standard deviation of the random walk is proportional to one \f$k\f$-th the mean.
1280This is achieved by setting \f$\alpha=k\f$ and \f$\beta=k/\mu\f$.
1281
1282The standard deviation of the walk is then: \f$\mu/\sqrt(k)\f$.
1283 */
1284class migamma_ref : public migamma {
1285protected:
1286        //! parameter l
1287        double l;
1288        //! reference vector
1289        vec refl;
1290public:
1291        //! Constructor
1292        migamma_ref () : migamma (), refl () {};
1293       
1294        //! Set value of \c k
1295        void set_parameters ( double k0 , vec ref0, double l0 ) {
1296                migamma::set_parameters ( ref0.length(), k0 );
1297                refl = pow ( ref0, 1.0 - l0 );
1298                l = l0;
1299        };
1300
1301        void validate(){
1302                migamma::validate();
1303                dimc = dimension();
1304        };
1305       
1306        void condition ( const vec &val ) {
1307                vec mean = elem_mult ( refl, pow ( val, l ) );
1308                migamma::condition ( mean );
1309        };
1310
1311
1312        /*! Create inverse-Gamma density with conditional mean value
1313        \f[ f(rv|rvc) = i\Gamma(k, k/(rvc^l \circ ref^{(1-l)}) \f]
1314        from structure
1315        \code
1316        class = 'migamma_ref';
1317        ref = [1e-5; 1e-5; 1e-2 1e-3];            // reference vector
1318        l = 0.999;                                // constant l
1319        k = 0.1;                                  // constant k
1320        rv = RV({'name'})                         // description of RV
1321        rvc = RV({'name'})                        // description of RV in condition
1322        \endcode
1323         */
1324        void from_setting ( const Setting &set );
1325
1326        // TODO dodelat void to_setting( Setting &set ) const;
1327};
1328
1329
1330UIREGISTER ( migamma_ref );
1331SHAREDPTR ( migamma_ref );
1332
1333/*! Log-Normal probability density
1334 only allow diagonal covariances!
1335
1336Density of the form \f$ \log(x)\sim \mathcal{N}(\mu,\sigma^2) \f$ , i.e.
1337\f[
1338x \sim \frac{1}{x\sigma\sqrt{2\pi}}\exp{-\frac{1}{2\sigma^2}(\log(x)-\mu)}
1339\f]
1340
1341Function from_setting loads mu and R in the same way as it does for enorm<>!
1342*/
1343class elognorm: public enorm<ldmat> {
1344public:
1345        vec sample() const {
1346                return exp ( enorm<ldmat>::sample() );
1347        };
1348        vec mean() const {
1349                vec var = enorm<ldmat>::variance();
1350                return exp ( mu - 0.5*var );
1351        };
1352
1353};
1354
1355/*!
1356\brief  Log-Normal random walk
1357
1358Mean value, \f$\mu\f$, is...
1359
1360 */
1361class mlognorm : public pdf_internal<elognorm> {
1362protected:
1363        //! parameter 1/2*sigma^2
1364        double sig2;
1365
1366        //! access
1367        vec &mu;
1368public:
1369        //! Constructor
1370        mlognorm() : pdf_internal<elognorm>(),
1371                        sig2 ( 0 ),
1372                        mu ( iepdf._mu() ) {
1373        }
1374
1375        //! Set value of \c k
1376        void set_parameters ( int size, double k ) {
1377                sig2 = 0.5 * log ( k * k + 1 );
1378                iepdf.set_parameters ( zeros ( size ), 2*sig2*eye ( size ) );
1379        };
1380       
1381        void validate(){
1382                iepdf.validate();
1383                dimc = iepdf.dimension();
1384        }
1385
1386        void condition ( const vec &val ) {
1387                mu = log ( val ) - sig2;//elem_mult ( refl,pow ( val,l ) );
1388        };
1389
1390        /*! Create logNormal random Walk
1391        \f[ f(rv|rvc) = log\mathcal{N}( \log(rvc)-0.5\log(k^2+1), k I) \f]
1392        from structure
1393        \code
1394        class = 'mlognorm';
1395        k   = 0.1;               // "variance" k
1396        mu0 = 0.1;               // Initial value of mean
1397        rv  = RV({'name'})       // description of RV
1398        rvc = RV({'name'})       // description of RV in condition
1399        \endcode
1400        */
1401        void from_setting ( const Setting &set );
1402
1403        // TODO dodelat void to_setting( Setting &set ) const;
1404
1405};
1406
1407UIREGISTER ( mlognorm );
1408SHAREDPTR ( mlognorm );
1409
1410/*! inverse Wishart density defined on Choleski decomposition
1411
1412*/
1413class eWishartCh : public epdf {
1414protected:
1415        //! Upper-Triagle of Choleski decomposition of \f$ \Psi \f$
1416        chmat Y;
1417        //! dimension of matrix  \f$ \Psi \f$
1418        int p;
1419        //! degrees of freedom  \f$ \nu \f$
1420        double delta;
1421public:
1422        //! Set internal structures
1423        void set_parameters ( const mat &Y0, const double delta0 ) {
1424                Y = chmat ( Y0 );
1425                delta = delta0;
1426                p = Y.rows();   
1427        }
1428        //! Set internal structures
1429        void set_parameters ( const chmat &Y0, const double delta0 ) {
1430                Y = Y0;
1431                delta = delta0;
1432                p = Y.rows();
1433                }
1434       
1435        virtual void validate (){
1436                dim = p * p;
1437        }
1438
1439        //! Sample matrix argument
1440        mat sample_mat() const {
1441                mat X = zeros ( p, p );
1442
1443                //sample diagonal
1444                for ( int i = 0; i < p; i++ ) {
1445                        GamRNG.setup ( 0.5* ( delta - i ) , 0.5 );   // no +1 !! index if from 0
1446#pragma omp critical
1447                        X ( i, i ) = sqrt ( GamRNG() );
1448                }
1449                //do the rest
1450                for ( int i = 0; i < p; i++ ) {
1451                        for ( int j = i + 1; j < p; j++ ) {
1452#pragma omp critical
1453                                X ( i, j ) = NorRNG.sample();
1454                        }
1455                }
1456                return X*Y._Ch();// return upper triangular part of the decomposition
1457        }
1458
1459        vec sample () const {
1460                return vec ( sample_mat()._data(), p*p );
1461        }
1462
1463        virtual vec mean() const NOT_IMPLEMENTED(0);
1464
1465        //! return expected variance (not covariance!)
1466        virtual vec variance() const NOT_IMPLEMENTED(0);
1467
1468        virtual double evallog ( const vec &val ) const NOT_IMPLEMENTED(0);
1469
1470        //! fast access function y0 will be copied into Y.Ch.
1471        void setY ( const mat &Ch0 ) {
1472                copy_vector ( dim, Ch0._data(), Y._Ch()._data() );
1473        }
1474
1475        //! fast access function y0 will be copied into Y.Ch.
1476        void _setY ( const vec &ch0 ) {
1477                copy_vector ( dim, ch0._data(), Y._Ch()._data() );
1478        }
1479
1480        //! access function
1481        const chmat& getY() const {
1482                return Y;
1483        }
1484};
1485
1486//! Inverse Wishart on Choleski decomposition
1487/*! Being computed by conversion from `standard' Wishart
1488*/
1489class eiWishartCh: public epdf {
1490protected:
1491        //! Internal instance of Wishart density
1492        eWishartCh W;
1493        //! size of Ch
1494        int p;
1495        //! parameter delta
1496        double delta;
1497public:
1498        //! constructor function
1499        void set_parameters ( const mat &Y0, const double delta0 ) {
1500                delta = delta0;
1501                W.set_parameters ( inv ( Y0 ), delta0 );
1502                p = Y0.rows();
1503        }
1504       
1505        virtual void    validate (){
1506                W.validate();
1507                dim = W.dimension();
1508        }
1509       
1510       
1511        vec sample() const {
1512                mat iCh;
1513                iCh = inv ( W.sample_mat() );
1514                return vec ( iCh._data(), dim );
1515        }
1516        //! access function
1517        void _setY ( const vec &y0 ) {
1518                mat Ch ( p, p );
1519                mat iCh ( p, p );
1520                copy_vector ( dim, y0._data(), Ch._data() );
1521
1522                iCh = inv ( Ch );
1523                W.setY ( iCh );
1524        }
1525
1526        virtual double evallog ( const vec &val ) const {
1527                chmat X ( p );
1528                const chmat& Y = W.getY();
1529
1530                copy_vector ( p*p, val._data(), X._Ch()._data() );
1531                chmat iX ( p );
1532                X.inv ( iX );
1533                // compute
1534//                              \frac{ |\Psi|^{m/2}|X|^{-(m+p+1)/2}e^{-tr(\Psi X^{-1})/2} }{ 2^{mp/2}\Gamma_p(m/2)},
1535                mat M = Y.to_mat() * iX.to_mat();
1536
1537                double log1 = 0.5 * p * ( 2 * Y.logdet() ) - 0.5 * ( delta + p + 1 ) * ( 2 * X.logdet() ) - 0.5 * trace ( M );
1538                //Fixme! Multivariate gamma omitted!! it is ok for sampling, but not otherwise!!
1539
1540                /*                              if (0) {
1541                                                        mat XX=X.to_mat();
1542                                                        mat YY=Y.to_mat();
1543
1544                                                        double log2 = 0.5*p*log(det(YY))-0.5*(delta+p+1)*log(det(XX))-0.5*trace(YY*inv(XX));
1545                                                        cout << log1 << "," << log2 << endl;
1546                                                }*/
1547                return log1;
1548        };
1549
1550        virtual vec mean() const NOT_IMPLEMENTED(0);
1551
1552        //! return expected variance (not covariance!)
1553        virtual vec variance() const NOT_IMPLEMENTED(0);
1554};
1555
1556//! Random Walk on inverse Wishart
1557class rwiWishartCh : public pdf_internal<eiWishartCh> {
1558protected:
1559        //!square root of \f$ \nu-p-1 \f$ - needed for computation of \f$ \Psi \f$ from conditions
1560        double sqd;
1561        //!reference point for diagonal
1562        vec refl;
1563        //! power of the reference
1564        double l;
1565        //! dimension
1566        int p;
1567
1568public:
1569        rwiWishartCh() : sqd ( 0 ), l ( 0 ), p ( 0 ) {}
1570        //! constructor function
1571        void set_parameters ( int p0, double k, vec ref0, double l0 ) {
1572                p = p0;
1573                double delta = 2 / ( k * k ) + p + 3;
1574                sqd = sqrt ( delta - p - 1 );
1575                l = l0;
1576                refl = pow ( ref0, 1 - l );
1577                iepdf.set_parameters ( eye ( p ), delta );
1578        };
1579       
1580        void validate(){
1581                iepdf.validate();
1582                dimc = iepdf.dimension();
1583        }
1584       
1585        void condition ( const vec &c ) {
1586                vec z = c;
1587                int ri = 0;
1588                for ( int i = 0; i < p*p; i += ( p + 1 ) ) {//trace diagonal element
1589                        z ( i ) = pow ( z ( i ), l ) * refl ( ri );
1590                        ri++;
1591                }
1592
1593                iepdf._setY ( sqd*z );
1594        }
1595};
1596
1597//! Switch between various resampling methods.
1598enum RESAMPLING_METHOD { MULTINOMIAL = 0, STRATIFIED = 1, SYSTEMATIC = 3 };
1599
1600//! Shortcut for multinomial.sample(int n). Various simplifications are allowed see RESAMPLING_METHOD
1601void resample(const vec &w, ivec &ind, RESAMPLING_METHOD=SYSTEMATIC);
1602
1603/*!
1604\brief Weighted empirical density
1605
1606Used e.g. in particle filters.
1607*/
1608class eEmp: public epdf {
1609protected :
1610        //! Number of particles
1611        int n;
1612        //! Sample weights \f$w\f$
1613        vec w;
1614        //! Samples \f$x^{(i)}, i=1..n\f$
1615        Array<vec> samples;
1616public:
1617        //! \name Constructors
1618        //!@{
1619        eEmp () : epdf (), w (), samples () {};
1620        //! copy constructor
1621        eEmp ( const eEmp &e ) : epdf ( e ), w ( e.w ), samples ( e.samples ) {};
1622        //!@}
1623
1624        //! Set samples and weights
1625        void set_statistics ( const vec &w0, const epdf &pdf0 );
1626        //! Set samples and weights
1627        void set_statistics ( const epdf &pdf0 , int n ) {
1628                set_statistics ( ones ( n ) / n, pdf0 );
1629        };
1630        //! Set sample
1631        void set_samples ( const epdf* pdf0 );
1632        //! Set sample
1633        void set_parameters ( int n0, bool copy = true ) {
1634                n = n0;
1635                w.set_size ( n0, copy );
1636                samples.set_size ( n0, copy );
1637        };
1638        //! Set samples
1639        void set_parameters ( const Array<vec> &Av ) {
1640                n = Av.size();
1641                w = 1 / n * ones ( n );
1642                samples = Av;
1643        };
1644        virtual void    validate (){
1645          bdm_assert (samples.length()==w.length(),"samples and weigths are of different lengths");
1646          n = w.length();
1647          if (n>0)
1648                pdf::dim = samples ( 0 ).length();
1649        }
1650        //! Potentially dangerous, use with care.
1651        vec& _w()  {
1652                return w;
1653        };
1654        //! Potentially dangerous, use with care.
1655        const vec& _w() const {
1656                return w;
1657        };
1658        //! access function
1659        Array<vec>& _samples() {
1660                return samples;
1661        };
1662        //! access function
1663        const vec& _sample ( int i ) const {
1664                return samples ( i );
1665        };
1666        //! access function
1667        const Array<vec>& _samples() const {
1668                return samples;
1669        };
1670        //! 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.
1671        void resample ( RESAMPLING_METHOD method = SYSTEMATIC );
1672
1673        //! inherited operation : NOT implemented
1674        vec sample() const NOT_IMPLEMENTED(0);
1675
1676        //! inherited operation : NOT implemented
1677        double evallog ( const vec &val ) const NOT_IMPLEMENTED(0);
1678
1679        vec mean() const {
1680                vec pom = zeros ( dim );
1681                for ( int i = 0; i < n; i++ ) {
1682                        pom += samples ( i ) * w ( i );
1683                }
1684                return pom;
1685        }
1686        vec variance() const {
1687                vec pom = zeros ( dim );
1688                for ( int i = 0; i < n; i++ ) {
1689                        pom += pow ( samples ( i ), 2 ) * w ( i );
1690                }
1691                return pom - pow ( mean(), 2 );
1692        }
1693        //! For this class, qbounds are minimum and maximum value of the population!
1694        void qbounds ( vec &lb, vec &ub, double perc = 0.95 ) const;
1695
1696        void to_setting ( Setting &set ) const {
1697                epdf::to_setting( set );
1698                UI::save ( samples, set, "samples" );
1699                UI::save ( w, set, "w" );
1700        }
1701
1702        void from_setting ( const Setting &set ) {
1703                epdf::from_setting( set );
1704               
1705                UI::get( samples, set, "samples", UI::compulsory );
1706                UI::get ( w, set, "w", UI::compulsory );
1707                validate();
1708        }
1709
1710};
1711UIREGISTER(eEmp);
1712
1713
1714////////////////////////
1715
1716template<class sq_T>
1717void enorm<sq_T>::set_parameters ( const vec &mu0, const sq_T &R0 ) {
1718//Fixme test dimensions of mu0 and R0;
1719        mu = mu0;
1720        R = R0;
1721        validate();
1722};
1723
1724template<class sq_T>
1725void enorm<sq_T>::from_setting ( const Setting &set ) {
1726        epdf::from_setting ( set ); //reads rv
1727
1728        UI::get ( mu, set, "mu", UI::compulsory );
1729        mat Rtmp;// necessary for conversion
1730        UI::get ( Rtmp, set, "R", UI::compulsory );
1731        R = Rtmp; // conversion
1732        validate();
1733}
1734
1735template<class sq_T>
1736void enorm<sq_T>::to_setting ( Setting &set ) const {
1737        epdf::to_setting ( set ); //reads rv   
1738        UI::save ( mu, set, "mu");
1739        UI::save ( R.to_mat(), set, "R");
1740}
1741
1742template<class sq_T>
1743void enorm<sq_T>::dupdate ( mat &v, double nu ) {
1744        //
1745};
1746
1747// template<class sq_T>
1748// void enorm<sq_T>::tupdate ( double phi, mat &vbar, double nubar ) {
1749//      //
1750// };
1751
1752template<class sq_T>
1753vec enorm<sq_T>::sample() const {
1754        vec x ( dim );
1755#pragma omp critical
1756        NorRNG.sample_vector ( dim, x );
1757        vec smp = R.sqrt_mult ( x );
1758
1759        smp += mu;
1760        return smp;
1761};
1762
1763// template<class sq_T>
1764// double enorm<sq_T>::eval ( const vec &val ) const {
1765//      double pdfl,e;
1766//      pdfl = evallog ( val );
1767//      e = exp ( pdfl );
1768//      return e;
1769// };
1770
1771template<class sq_T>
1772double enorm<sq_T>::evallog_nn ( const vec &val ) const {
1773        // 1.83787706640935 = log(2pi)
1774        double tmp = -0.5 * ( R.invqform ( mu - val ) );// - lognc();
1775        return  tmp;
1776};
1777
1778template<class sq_T>
1779inline double enorm<sq_T>::lognc () const {
1780        // 1.83787706640935 = log(2pi)
1781        double tmp = 0.5 * ( R.cols() * 1.83787706640935 + R.logdet() );
1782        return tmp;
1783};
1784
1785
1786// template<class sq_T>
1787// vec mlnorm<sq_T>::samplecond (const  vec &cond, double &lik ) {
1788//      this->condition ( cond );
1789//      vec smp = epdf.sample();
1790//      lik = epdf.eval ( smp );
1791//      return smp;
1792// }
1793
1794// template<class sq_T>
1795// mat mlnorm<sq_T>::samplecond (const vec &cond, vec &lik, int n ) {
1796//      int i;
1797//      int dim = rv.count();
1798//      mat Smp ( dim,n );
1799//      vec smp ( dim );
1800//      this->condition ( cond );
1801//
1802//      for ( i=0; i<n; i++ ) {
1803//              smp = epdf.sample();
1804//              lik ( i ) = epdf.eval ( smp );
1805//              Smp.set_col ( i ,smp );
1806//      }
1807//
1808//      return Smp;
1809// }
1810
1811
1812template<class sq_T>
1813shared_ptr<epdf> enorm<sq_T>::marginal ( const RV &rvn ) const {
1814        enorm<sq_T> *tmp = new enorm<sq_T> ();
1815        shared_ptr<epdf> narrow ( tmp );
1816        marginal ( rvn, *tmp );
1817        return narrow;
1818}
1819
1820template<class sq_T>
1821void enorm<sq_T>::marginal ( const RV &rvn, enorm<sq_T> &target ) const {
1822        bdm_assert ( isnamed(), "rv description is not assigned" );
1823        ivec irvn = rvn.dataind ( rv );
1824
1825        sq_T Rn ( R, irvn );  // select rows and columns of R
1826
1827        target.set_rv ( rvn );
1828        target.set_parameters ( mu ( irvn ), Rn );
1829}
1830
1831template<class sq_T>
1832shared_ptr<pdf> enorm<sq_T>::condition ( const RV &rvn ) const {
1833        mlnorm<sq_T> *tmp = new mlnorm<sq_T> ();
1834        shared_ptr<pdf> narrow ( tmp );
1835        condition ( rvn, *tmp );
1836        return narrow;
1837}
1838
1839template<class sq_T>
1840void enorm<sq_T>::condition ( const RV &rvn, pdf &target ) const {
1841        typedef mlnorm<sq_T> TMlnorm;
1842
1843        bdm_assert ( isnamed(), "rvs are not assigned" );
1844        TMlnorm &uptarget = dynamic_cast<TMlnorm &> ( target );
1845
1846        RV rvc = rv.subt ( rvn );
1847        bdm_assert ( ( rvc._dsize() + rvn._dsize() == rv._dsize() ), "wrong rvn" );
1848        //Permutation vector of the new R
1849        ivec irvn = rvn.dataind ( rv );
1850        ivec irvc = rvc.dataind ( rv );
1851        ivec perm = concat ( irvn , irvc );
1852        sq_T Rn ( R, perm );
1853
1854        //fixme - could this be done in general for all sq_T?
1855        mat S = Rn.to_mat();
1856        //fixme
1857        int n = rvn._dsize() - 1;
1858        int end = R.rows() - 1;
1859        mat S11 = S.get ( 0, n, 0, n );
1860        mat S12 = S.get ( 0, n , rvn._dsize(), end );
1861        mat S22 = S.get ( rvn._dsize(), end, rvn._dsize(), end );
1862
1863        vec mu1 = mu ( irvn );
1864        vec mu2 = mu ( irvc );
1865        mat A = S12 * inv ( S22 );
1866        sq_T R_n ( S11 - A *S12.T() );
1867
1868        uptarget.set_rv ( rvn );
1869        uptarget.set_rvc ( rvc );
1870        uptarget.set_parameters ( A, mu1 - A*mu2, R_n );
1871        uptarget.validate();
1872}
1873
1874/*! Dirac delta function distribution */
1875class dirac: public epdf{
1876        protected: 
1877                vec point;
1878        public:
1879                double evallog (const vec &dt) const {return -inf;}
1880                vec mean () const {return point;}
1881                vec variance () const {return zeros(point.length());}
1882                void qbounds ( vec &lb, vec &ub, double percentage = 0.95 ) const { lb = point; ub = point;}
1883                //! access
1884                const vec& _point() {return point;}
1885                void set_point(const vec& p){point =p; dim=p.length();}
1886                vec sample() const {return point;}
1887        };
1888
1889////
1890///////
1891template<class sq_T>
1892void mgnorm<sq_T >::set_parameters ( const shared_ptr<fnc> &g0, const sq_T &R0 ) {
1893        g = g0;
1894        this->iepdf.set_parameters ( zeros ( g->dimension() ), R0 );
1895}
1896
1897template<class sq_T>
1898void mgnorm<sq_T >::condition ( const vec &cond ) {
1899        this->iepdf._mu() = g->eval ( cond );
1900};
1901
1902//! \todo unify this stuff with to_string()
1903template<class sq_T>
1904std::ostream &operator<< ( std::ostream &os,  mlnorm<sq_T> &ml ) {
1905        os << "A:" << ml.A << endl;
1906        os << "mu:" << ml.mu_const << endl;
1907        os << "R:" << ml._R() << endl;
1908        return os;
1909};
1910
1911}
1912#endif //EF_H
Note: See TracBrowser for help on using the browser.