Changeset 461 for library/bdm

Show
Ignore:
Timestamp:
07/31/09 13:06:49 (15 years ago)
Author:
vbarta
Message:

mpdf (& its dependencies) reformat: now using shared_ptr, moved virtual method bodies to .cpp

Location:
library/bdm
Files:
9 modified

Legend:

Unmodified
Added
Removed
  • library/bdm/base/bdmbase.cpp

    r436 r461  
    135135        for ( int i = 0;i < N;i++ ) X.set_col ( i, this->sample() ); 
    136136        return X; 
    137 }; 
    138  
     137} 
     138 
     139vec mpdf::samplecond(const vec &cond) { 
     140    condition(cond); 
     141    vec temp = shep->sample(); 
     142    return temp; 
     143} 
     144 
     145mat mpdf::samplecond_m(const vec &cond, int N) { 
     146    condition(cond); 
     147    mat temp(shep->dimension(), N); 
     148    vec smp(shep->dimension()); 
     149    for (int i = 0; i < N; i++) { 
     150        smp = shep->sample(); 
     151        temp.set_col(i, smp); 
     152    } 
     153 
     154    return temp; 
     155} 
     156 
     157double mpdf::evallogcond(const vec &dt, const vec &cond) { 
     158    double tmp; 
     159    condition(cond); 
     160    tmp = shep->evallog(dt); 
     161    // it_assert_debug(std::isfinite(tmp), "Infinite value"); 
     162    return tmp; 
     163} 
     164 
     165vec mpdf::evallogcond_m(const mat &Dt, const vec &cond) { 
     166    condition(cond); 
     167    return shep->evallog_m(Dt); 
     168} 
     169 
     170vec mpdf::evallogcond_m(const Array<vec> &Dt, const vec &cond) { 
     171    condition(cond); 
     172    return shep->evallog_m(Dt); 
     173} 
     174 
     175void mpdf::from_setting(const Setting &set){ 
     176    if (set.exists("rv")) { 
     177        RV *r = UI::build<RV>(set, "rv"); 
     178        set_rv(*r); 
     179        delete r; 
     180    } 
     181 
     182    if (set.exists("rvc")) { 
     183        RV *r = UI::build<RV>(set, "rvc"); 
     184        set_rvc(*r);  
     185        delete r; 
     186    } 
     187} 
    139188 
    140189std::ostream &operator<< ( std::ostream &os, const RV &rv ) { 
  • library/bdm/base/bdmbase.h

    r436 r461  
    1818#include "../itpp_ext.h" 
    1919#include "../bdmroot.h" 
     20#include "../shared_ptr.h" 
    2021#include "user_info.h" 
    21  
    2222 
    2323using namespace libconfig; 
     
    350350  //! random variable in condition 
    351351  RV rvc; 
     352 
     353private: 
    352354  //! pointer to internal epdf 
    353   epdf* ep; 
     355  shared_ptr<epdf> shep; 
     356 
    354357public: 
    355358  //! \name Constructors 
    356359  //! @{ 
    357360 
    358   mpdf() : dimc(0), rvc() {}; 
    359   //! copy constructor does not set pointer \c ep - has to be done in offsprings! 
    360   mpdf(const mpdf &m) : dimc(m.dimc), rvc(m.rvc) {}; 
    361   //!@} 
    362  
    363   //! \name Matematical operations 
    364   //!@{ 
    365  
    366   //! Returns a sample from the density conditioned on \c cond, \f$x \sim epdf(rv|cond)\f$. \param cond is numeric value of \c rv 
    367   virtual vec samplecond(const vec &cond) { 
    368     this->condition(cond); 
    369     vec temp = ep->sample(); 
    370     return temp; 
    371   }; 
    372   //! Returns \param N samples from the density conditioned on \c cond, \f$x \sim epdf(rv|cond)\f$. \param cond is numeric value of \c rv 
    373   virtual mat samplecond_m(const vec &cond, int N) { 
    374     this->condition(cond); 
    375     mat temp(ep->dimension(), N); 
    376     vec smp(ep->dimension()); 
    377     for (int i = 0; i < N; i++) {smp = ep->sample() ; temp.set_col(i, smp);} 
    378     return temp; 
    379   }; 
    380   //! Update \c ep so that it represents this mpdf conditioned on \c rvc = cond 
    381   virtual void condition(const vec &cond) {it_error("Not implemented");}; 
     361    mpdf():dimc(0), rvc() { } 
     362 
     363    mpdf(const mpdf &m):dimc(m.dimc), rvc(m.rvc), shep(m.shep) { } 
     364    //!@} 
     365 
     366    //! \name Matematical operations 
     367    //!@{ 
     368 
     369    //! Returns a sample from the density conditioned on \c cond, \f$x \sim epdf(rv|cond)\f$. \param cond is numeric value of \c rv 
     370    virtual vec samplecond(const vec &cond); 
     371 
     372    //! Returns \param N samples from the density conditioned on \c cond, \f$x \sim epdf(rv|cond)\f$. \param cond is numeric value of \c rv 
     373    virtual mat samplecond_m(const vec &cond, int N); 
     374 
     375    //! Update \c ep so that it represents this mpdf conditioned on \c rvc = cond 
     376    virtual void condition(const vec &cond) {it_error("Not implemented");}; 
    382377 
    383378  //! Shortcut for conditioning and evaluation of the internal epdf. In some cases,  this operation can be implemented efficiently. 
    384   virtual double evallogcond(const vec &dt, const vec &cond) { 
    385     double tmp; 
    386     this->condition(cond); 
    387     tmp = ep->evallog(dt); 
    388  //   it_assert_debug(std::isfinite(tmp), "Infinite value"); 
    389     return tmp; 
    390   }; 
    391  
    392   //! Matrix version of evallogcond 
    393   virtual vec evallogcond_m(const mat &Dt, const vec &cond) {this->condition(cond); return ep->evallog_m(Dt);}; 
    394   //! Array<vec> version of evallogcond 
    395   virtual vec evallogcond_m(const Array<vec> &Dt, const vec &cond) {this->condition(cond); return ep->evallog_m(Dt);}; 
    396  
    397   //! \name Access to attributes 
    398   //! @{ 
    399  
    400   RV _rv() {return ep->_rv();} 
    401   RV _rvc() {it_assert_debug(isnamed(), ""); return rvc;} 
    402   int dimension() {return ep->dimension();} 
    403   int dimensionc() {return dimc;} 
    404   epdf& _epdf() {return *ep;} 
    405   epdf* _e() {return ep;} 
    406   //! Load from structure with elements: 
    407   //!  \code 
    408   //! { rv = {class="RV", names=(...),}; // RV describing meaning of random variable 
    409   //!   rvc= {class="RV", names=(...),}; // RV describing meaning of random variable in condition 
    410   //!   // elements of offsprings 
    411   //! } 
    412   //! \endcode 
    413   //!@} 
    414   void from_setting(const Setting &set){ 
    415           if (set.exists("rv")){ 
    416                   RV* r = UI::build<RV>(set,"rv"); 
    417                   set_rv(*r);  
    418                   delete r; 
    419           } 
    420           if (set.exists("rvc")){ 
    421                   RV* r = UI::build<RV>(set,"rvc"); 
    422                   set_rvc(*r);  
    423                   delete r; 
    424           } 
    425   } 
    426   //!@} 
    427  
    428   //! \name Connection to other objects 
    429   //!@{ 
    430   void set_rvc(const RV &rvc0) {rvc = rvc0;} 
    431   void set_rv(const RV &rv0) {ep->set_rv(rv0);} 
    432   bool isnamed() {return (ep->isnamed()) && (dimc == rvc._dsize());} 
    433   //!@} 
     379    virtual double evallogcond(const vec &dt, const vec &cond); 
     380 
     381    //! Matrix version of evallogcond 
     382    virtual vec evallogcond_m(const mat &Dt, const vec &cond); 
     383 
     384    //! Array<vec> version of evallogcond 
     385    virtual vec evallogcond_m(const Array<vec> &Dt, const vec &cond); 
     386 
     387    //! \name Access to attributes 
     388    //! @{ 
     389 
     390    RV _rv() { return shep->_rv(); } 
     391    RV _rvc() { it_assert_debug(isnamed(), ""); return rvc; } 
     392    int dimension() { return shep->dimension(); } 
     393    int dimensionc() { return dimc; } 
     394 
     395    epdf *e() { return shep.get(); } 
     396 
     397    void set_ep(shared_ptr<epdf> ep) { shep = ep; } 
     398 
     399    //! Load from structure with elements: 
     400    //!  \code 
     401    //! { rv = {class="RV", names=(...),}; // RV describing meaning of random variable 
     402    //!   rvc= {class="RV", names=(...),}; // RV describing meaning of random variable in condition 
     403    //!   // elements of offsprings 
     404    //! } 
     405    //! \endcode 
     406    //!@} 
     407    void from_setting(const Setting &set); 
     408    //!@} 
     409 
     410    //! \name Connection to other objects 
     411    //!@{ 
     412    void set_rvc(const RV &rvc0) { rvc = rvc0; } 
     413    void set_rv(const RV &rv0) { shep->set_rv(rv0); } 
     414    bool isnamed() { return (shep->isnamed()) && (dimc == rvc._dsize()); } 
     415    //!@} 
    434416}; 
    435417 
     
    628610*/ 
    629611class mepdf : public mpdf { 
    630         bool owning_ep; // flag trigers deleting ep by destructor 
    631612public: 
    632613        //!Default constructor 
    633         mepdf(){}; 
    634         mepdf ( epdf* em, bool owning_ep0=false ) :mpdf ( ) {ep= em ;owning_ep=owning_ep0;dimc=0;}; 
    635         mepdf (const epdf* em ) :mpdf ( ) {ep=const_cast<epdf*>( em );}; 
     614        mepdf() { } 
     615 
     616        mepdf(shared_ptr<epdf> em) :mpdf ( ) { set_ep(em); dimc=0; } 
     617 
    636618        void condition ( const vec &cond ) {} 
    637         ~mepdf(){if (owning_ep) delete ep;} 
     619 
    638620  //! Load from structure with elements: 
    639621  //!  \code 
     
    644626  //!@} 
    645627        void from_setting(const Setting &set){ 
    646                 epdf* e = UI::build<epdf>(set,"epdf"); 
    647                 ep=     e;  
    648                 owning_ep=true; 
     628                epdf *e = UI::build<epdf>(set, "epdf"); 
     629                set_ep(e); 
    649630        } 
    650631}; 
  • library/bdm/estim/particles.cpp

    r384 r461  
    1414                //generate new samples from paramater evolution model; 
    1515                _samples ( i ) = par->samplecond ( _samples ( i ) ); 
    16                 lls ( i ) = par->_e()->evallog ( _samples ( i ) ); 
     16                lls ( i ) = par->e()->evallog ( _samples ( i ) ); 
    1717                lls ( i ) *= obs->evallogcond ( dt,_samples ( i ) ); 
    1818 
  • library/bdm/estim/particles.h

    r384 r461  
    216216                //generate new samples from paramater evolution model; 
    217217                _samples ( i ) = par->samplecond ( _samples ( i ) ); 
    218                 llsP ( i ) = par->_e()->evallog ( _samples ( i ) ); 
     218                llsP ( i ) = par->e()->evallog ( _samples ( i ) ); 
    219219                BMs ( i )->condition ( _samples ( i ) ); 
    220220                BMs ( i )->bayes ( dt ); 
  • library/bdm/shared_ptr.h

    r420 r461  
    2828template <typename T> 
    2929class shared_ptr { 
     30    template<class U> friend class shared_ptr; 
     31 
    3032private: 
    3133    T *payload; 
     
    5153    //! If other is empty, constructs an empty shared_ptr; otherwise, 
    5254    //! constructs a shared_ptr that shares ownership with other. 
    53     shared_ptr(const shared_ptr &other): 
     55    shared_ptr(const shared_ptr<T> &other): 
     56        payload(other.payload), 
     57        refCnt(other.refCnt) 
     58    { 
     59        add_ref(); 
     60    } 
     61       
     62    //! If other is empty, constructs an empty shared_ptr; otherwise, 
     63    //! constructs a shared_ptr that shares ownership with other. 
     64    template<typename U> 
     65    shared_ptr(const shared_ptr<U> &other): 
    5466        payload(other.payload), 
    5567        refCnt(other.refCnt) 
  • library/bdm/stat/emix.h

    r395 r461  
    1616#define LOG2  0.69314718055995   
    1717 
     18#include "../shared_ptr.h" 
    1819#include "exp_family.h" 
    1920 
     
    5152                rvc = nom0->_rv().subt ( rv ); 
    5253                dimc = rvc._dsize(); 
    53                 ep = new epdf; 
    54                 ep->set_parameters(rv._dsize()); 
    55                 ep->set_rv(rv); 
     54                set_ep(shared_ptr<epdf>(new epdf)); 
     55                e()->set_parameters(rv._dsize()); 
     56                e()->set_rv(rv); 
    5657                 
    5758                //prepare data structures 
     
    6667        double evallogcond ( const vec &val, const vec &cond ) { 
    6768                double tmp; 
    68                 vec nom_val ( ep->dimension() + dimc ); 
     69                vec nom_val ( e()->dimension() + dimc ); 
    6970                dl.pushup_cond ( nom_val,val,cond ); 
    7071                tmp = exp ( nom->evallog ( nom_val ) - den->evallog ( cond ) ); 
     
    228229        //! Data link for each mpdfs 
    229230        Array<datalink_m2m*> dls; 
     231 
    230232        //! dummy ep 
    231         epdf dummy; 
     233        shared_ptr<epdf> dummy; 
     234 
    232235public: 
    233236        /*!\brief Constructor from list of mFacs, 
    234237        */ 
    235         mprod (){}; 
    236         mprod (Array<mpdf*> mFacs ){set_elements( mFacs );}; 
     238        mprod():dummy(new epdf()) { } 
     239        mprod (Array<mpdf*> mFacs): 
     240            dummy(new epdf()) { 
     241            set_elements(mFacs); 
     242        } 
     243 
    237244        void set_elements(Array<mpdf*> mFacs , bool own=false) { 
    238245                 
     
    241248                epdfs.set_size(mFacs.length()); 
    242249                                 
    243                 ep=&dummy; 
     250                set_ep(dummy); 
    244251                RV rv=getrv ( true ); 
    245252                set_rv ( rv ); 
    246                 dummy.set_parameters ( rv._dsize() ); 
    247                 setrvc ( ep->_rv(),rvc ); 
     253                dummy->set_parameters(rv._dsize()); 
     254                setrvc(e()->_rv(), rvc); 
    248255                // rv and rvc established = > we can link them with mpdfs 
    249256                for ( int i = 0;i < mpdfs.length(); i++ ) { 
     
    253260 
    254261                for ( int i=0; i<mpdfs.length(); i++ ) { 
    255                         epdfs ( i ) =& ( mpdfs ( i )->_epdf() ); 
     262                        epdfs(i) = mpdfs(i)->e(); 
    256263                } 
    257264        }; 
     
    291298        //TODO smarter... 
    292299        vec samplecond ( const vec &cond ) { 
    293                 //! Ugly hack to help to discover if mpfs are not in proper order. Correct solution = check that explicitely. 
    294                 vec smp= std::numeric_limits<double>::infinity() * ones ( ep->dimension() ); 
     300            //! Ugly hack to help to discover if mpfs are not in proper order. Correct solution = check that explicitely. 
     301            vec smp= std::numeric_limits<double>::infinity() * ones ( e()->dimension() ); 
    295302                vec smpi; 
    296303                // Hard assumption here!!! We are going backwards, to assure that samples that are needed from smp are already generated! 
     
    417424        //! Component (epdfs) 
    418425        Array<mpdf*> Coms; 
     426 
    419427        //!Internal epdf 
    420         emix Epdf; 
     428        shared_ptr<emix> iepdf; 
     429 
    421430public: 
    422431        //!Default constructor 
    423         mmix ( ) : mpdf ( ), Epdf () {ep = &Epdf;}; 
     432        mmix():iepdf(new emix()) { 
     433            set_ep(iepdf); 
     434        } 
     435 
    424436        //! Set weights \c w and components \c R 
    425437        void set_parameters ( const vec &w, const Array<mpdf*> &Coms ) { 
     
    427439 
    428440                for ( int i = 0;i < Coms.length();i++ ) { 
    429                         Eps ( i ) = & ( Coms ( i )->_epdf() ); 
    430                 } 
    431                 Epdf.set_parameters ( w, Eps ); 
    432         }; 
     441                        Eps(i) = Coms(i)->e(); 
     442                } 
     443 
     444                iepdf->set_parameters(w, Eps); 
     445        } 
    433446 
    434447        void condition ( const vec &cond ) { 
  • library/bdm/stat/exp_family.cpp

    r404 r461  
    211211} 
    212212 
    213 //MGamma 
    214 void mgamma::set_parameters ( double k0, const vec &beta0 ) { 
    215         k=k0; 
    216         ep = &epdf; 
    217         epdf.set_parameters ( k*ones ( beta0.length()),beta0 ); 
    218         dimc = ep->dimension(); 
    219 }; 
     213void mgamma::set_parameters(double k0, const vec &beta0) { 
     214        k = k0; 
     215        iepdf->set_parameters(k * ones(beta0.length()), beta0); 
     216        dimc = e()->dimension(); 
     217} 
    220218 
    221219ivec eEmp::resample (RESAMPLING_METHOD method) { 
  • library/bdm/stat/exp_family.h

    r450 r461  
    1515 
    1616 
     17#include "../shared_ptr.h" 
    1718#include "../base/bdmbase.h" 
    1819#include "../math/chmat.h" 
     
    530531                protected: 
    531532                        //! Internal epdf that arise by conditioning on \c rvc 
    532                         enorm<sq_T> epdf; 
     533                        shared_ptr<enorm<sq_T> > iepdf; 
    533534                        mat A; 
    534535                        vec mu_const; 
     
    537538                        //! \name Constructors 
    538539                        //!@{ 
    539                         mlnorm ( ) :mEF (),epdf ( ),A ( ),_mu ( epdf._mu() ) {ep =&epdf; }; 
    540                         mlnorm ( const  mat &A, const vec &mu0, const sq_T &R ) :epdf ( ),_mu ( epdf._mu() ) 
    541                         { 
    542                                 ep =&epdf; set_parameters ( A,mu0,R ); 
    543                         }; 
     540                        mlnorm():iepdf(new enorm<sq_T>()), _mu(iepdf->_mu()) { set_ep(iepdf); }; 
     541                        mlnorm (const mat &A, const vec &mu0, const sq_T &R ) :iepdf(new enorm<sq_T>()), _mu(iepdf->_mu()) 
     542                        { 
     543                                set_ep(iepdf); set_parameters ( A,mu0,R ); 
     544                        } 
     545 
    544546                        //! Set \c A and \c R 
    545547                        void set_parameters ( const  mat &A, const vec &mu0, const sq_T &R ); 
     
    553555                        mat& _A() {return A;} 
    554556                        //!access function 
    555                         mat _R() {return epdf._R().to_mat();} 
     557                        mat _R() { return iepdf->_R().to_mat(); } 
    556558 
    557559                        template<class sq_M> 
     
    577579                protected: 
    578580                        //! Internal epdf that arise by conditioning on \c rvc 
    579                         enorm<sq_T> epdf; 
     581                        shared_ptr<enorm<sq_T> > iepdf; 
    580582                        vec &mu; 
    581583                        fnc* g; 
    582584                public: 
    583585                        //!default constructor 
    584                         mgnorm() :mu ( epdf._mu() ) {ep=&epdf;} 
     586                        mgnorm():iepdf(new enorm<sq_T>()), mu(iepdf->_mu()) { set_ep(iepdf); } 
    585587                        //!set mean function 
    586                         void set_parameters ( fnc* g0, const sq_T &R0 ) {g=g0; epdf.set_parameters ( zeros ( g->dimension() ), R0 );} 
     588                        void set_parameters ( fnc* g0, const sq_T &R0 ) {g=g0; iepdf->set_parameters ( zeros ( g->dimension() ), R0 );} 
    587589                        void condition ( const vec &cond ) {mu=g->eval ( cond );}; 
    588590 
     
    662664                public: 
    663665                        mlstudent ( ) :mlnorm<ldmat> (), 
    664                                         Lambda (),      _R ( epdf._R() ) {} 
     666                                        Lambda (),      _R ( iepdf->_R() ) {} 
    665667                        void set_parameters ( const mat &A0, const vec &mu0, const ldmat &R0, const ldmat& Lambda0 ) 
    666668                        { 
     
    668670                                it_assert_debug ( R0.rows() ==A0.rows(),"" ); 
    669671 
    670                                 epdf.set_parameters ( mu0,Lambda ); // 
     672                                iepdf->set_parameters ( mu0,Lambda ); // 
    671673                                A = A0; 
    672674                                mu_const = mu0; 
     
    705707                protected: 
    706708                        //! Internal epdf that arise by conditioning on \c rvc 
    707                         egamma epdf; 
     709                        shared_ptr<egamma> iepdf; 
     710 
    708711                        //! Constant \f$k\f$ 
    709712                        double k; 
    710                         //! cache of epdf.beta 
     713 
     714                        //! cache of iepdf.beta 
    711715                        vec &_beta; 
    712716 
    713717                public: 
    714718                        //! Constructor 
    715                         mgamma ( ) : mEF ( ), epdf (), _beta ( epdf._beta() ) {ep=&epdf;}; 
     719                        mgamma():iepdf(new egamma()), k(0), 
     720                            _beta(iepdf->_beta()) { 
     721                            set_ep(iepdf); 
     722                        } 
     723 
    716724                        //! Set value of \c k 
    717                         void set_parameters ( double k, const vec &beta0 ); 
     725                        void set_parameters(double k, const vec &beta0); 
     726 
    718727                        void condition ( const vec &val ) {_beta=k/val;}; 
    719728                        //! Load from structure with elements: 
     
    749758                protected: 
    750759                        //! Internal epdf that arise by conditioning on \c rvc 
    751                         eigamma epdf; 
     760                        shared_ptr<eigamma> iepdf; 
     761 
    752762                        //! Constant \f$k\f$ 
    753763                        double k; 
    754                         //! cache of epdf.alpha 
     764 
     765                        //! cache of iepdf.alpha 
    755766                        vec &_alpha; 
    756                         //! cache of epdf.beta 
     767 
     768                        //! cache of iepdf.beta 
    757769                        vec &_beta; 
    758770 
     
    760772                        //! \name Constructors 
    761773                        //!@{ 
    762                         migamma ( ) : mEF (), epdf ( ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;}; 
    763                         migamma ( const migamma &m ) : mEF (), epdf ( m.epdf ), _alpha ( epdf._alpha() ), _beta ( epdf._beta() ) {ep=&epdf;}; 
     774                        migamma():iepdf(new eigamma()), 
     775                            k(0), 
     776                            _alpha(iepdf->_alpha()), 
     777                            _beta(iepdf->_beta()) { 
     778                            set_ep(iepdf); 
     779                        } 
     780 
     781                        migamma(const migamma &m):iepdf(m.iepdf), 
     782                            k(0), 
     783                            _alpha(iepdf->_alpha()), 
     784                            _beta(iepdf->_beta()) { 
     785                            set_ep(iepdf); 
     786                        } 
    764787                        //!@} 
    765788 
     
    768791                        { 
    769792                                k=k0; 
    770                                 epdf.set_parameters ( ( 1.0/ ( k*k ) +2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ ); 
     793                                iepdf->set_parameters ( ( 1.0/ ( k*k ) +2.0 ) *ones ( len ) /*alpha*/, ones ( len ) /*beta*/ ); 
    771794                                dimc = dimension(); 
    772795                        }; 
     
    907930        { 
    908931                protected: 
    909                         elognorm eno; 
     932                        shared_ptr<elognorm> eno; 
     933 
    910934                        //! parameter 1/2*sigma^2 
    911935                        double sig2; 
     936 
    912937                        //! access 
    913938                        vec &mu; 
    914939                public: 
    915940                        //! Constructor 
    916                         mlognorm ( ) : eno (), mu ( eno._mu() ) {ep=&eno;}; 
     941                        mlognorm():eno(new elognorm()), 
     942                            sig2(0), 
     943                            mu(eno->_mu()) { 
     944                            set_ep(eno); 
     945                        } 
     946 
    917947                        //! Set value of \c k 
    918948                        void set_parameters ( int size, double k ) 
    919949                        { 
    920950                                sig2 = 0.5*log ( k*k+1 ); 
    921                                 eno.set_parameters ( zeros ( size ),2*sig2*eye ( size ) ); 
     951                                eno->set_parameters ( zeros ( size ),2*sig2*eye ( size ) ); 
    922952 
    923953                                dimc = size; 
     
    10531083        { 
    10541084                protected: 
    1055                         eiWishartCh eiW; 
     1085                        shared_ptr<eiWishartCh> eiW; 
    10561086                        //!square root of \f$ \nu-p-1 \f$ - needed for computation of \f$ \Psi \f$ from conditions 
    10571087                        double sqd; 
     
    10601090                        double l; 
    10611091                        int p; 
    1062                 public: 
     1092 
     1093                public: 
     1094                        rwiWishartCh():eiW(new eiWishartCh()), 
     1095                            sqd(0), l(0), p(0) { 
     1096                            set_ep(eiW); 
     1097                        } 
     1098 
    10631099                        void set_parameters ( int p0, double k, vec ref0, double l0  ) 
    10641100                        { 
     
    10691105                                refl=pow(ref0,1-l); 
    10701106                                 
    1071                                 eiW.set_parameters ( eye ( p ),delta ); 
    1072                                 ep=&eiW; 
    1073                                 dimc=eiW.dimension(); 
     1107                                eiW->set_parameters ( eye ( p ),delta ); 
     1108                                dimc=eiW->dimension(); 
    10741109                        } 
    10751110                        void condition ( const vec &c ) { 
     
    10811116                                } 
    10821117 
    1083                                 eiW._setY ( sqd*z ); 
     1118                                eiW->_setY ( sqd*z ); 
    10841119                        } 
    10851120        }; 
     
    12391274                it_assert_debug ( A0.rows() ==R0.rows(),"" ); 
    12401275 
    1241                 epdf.set_parameters ( zeros ( A0.rows() ),R0 ); 
     1276                iepdf->set_parameters(zeros(A0.rows()), R0); 
    12421277                A = A0; 
    12431278                mu_const = mu0; 
    1244                 dimc=A0.cols(); 
     1279                dimc = A0.cols(); 
    12451280        } 
    12461281 
     
    13321367                os << "A:"<< ml.A<<endl; 
    13331368                os << "mu:"<< ml.mu_const<<endl; 
    1334                 os << "R:" << ml.epdf._R().to_mat() <<endl; 
     1369                os << "R:" << ml.iepdf->_R().to_mat() <<endl; 
    13351370                return os; 
    13361371        }; 
  • library/bdm/stat/merger.cpp

    r423 r461  
    156156                                { 
    157157                                        // no need for conditioning or marginalization 
    158                                         lw_src = mpdfs ( i )->_epdf().evallog_m ( Smp  ); 
     158                                        lw_src = mpdfs ( i )->e()->evallog_m ( Smp  ); 
    159159                                } 
    160160                                else 
     
    202202                                        { 
    203203                                                mpdfs ( i )->condition ( dls ( i )->get_cond ( Smp ( k ) ) ); 
    204                                                 lw_src ( k ) += mpdfs ( i )->_epdf().evallog ( dls ( i )->pushdown ( Smp ( k ) ) ); 
     204                                                lw_src ( k ) += mpdfs ( i )->e()->evallog ( dls ( i )->pushdown ( Smp ( k ) ) ); 
    205205                                        } 
    206206