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

Revision 665, 39.9 kB (checked in by smidl, 15 years ago)

Compilation and minor extensions

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