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

Revision 796, 41.6 kB (checked in by smidl, 14 years ago)

UI changes for MixEF_init

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