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

Revision 725, 41.4 kB (checked in by smidl, 15 years ago)

Sampling from egiw

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