CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/HiggsAnalysis/CombinedLimit/interface/VerticalInterpHistPdf.h

Go to the documentation of this file.
00001 #ifndef ROO_VERTICA_INTERP_HIST
00002 #define ROO_VERTICA_INTERP_HIST
00003 
00007 #include "RooAbsPdf.h"
00008 #include "RooRealProxy.h"
00009 #include "RooListProxy.h"
00010 #include "TH1.h"
00011 #include "../interface/SimpleCacheSentry.h"
00012 #include "../interface/FastTemplate.h"
00013 #include <cmath>
00014 
00015 #define FastVerticalInterpHistPdf_CopyConstructor
00016 #define FastVerticalInterpHistPdf_Serializable
00017 
00018 class FastVerticalInterpHistPdf;
00019 
00020 class VerticalInterpHistPdf : public RooAbsPdf {
00021 public:
00022 
00023   VerticalInterpHistPdf() ;
00024   VerticalInterpHistPdf(const char *name, const char *title, const RooRealVar &x, const RooArgList& funcList, const RooArgList& coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1) ;
00025   VerticalInterpHistPdf(const VerticalInterpHistPdf& other, const char* name=0) ;
00026   virtual TObject* clone(const char* newname) const { return new VerticalInterpHistPdf(*this,newname) ; }
00027   virtual ~VerticalInterpHistPdf() ;
00028 
00029   Bool_t selfNormalized() const { return kTRUE; }
00030 
00031   Double_t evaluate() const ;
00032 
00033   const RooArgList& funcList() const { return _funcList ; }
00034   const RooArgList& coefList() const { return _coefList ; }
00035 
00036   const RooRealProxy &x() const { return _x; }
00037   friend class FastVerticalInterpHistPdf;
00038 protected:
00039   RooRealProxy   _x;
00040   RooListProxy _funcList ;   //  List of component FUNCs
00041   RooListProxy _coefList ;  //  List of coefficients
00042   Double_t     _smoothRegion;
00043   Int_t        _smoothAlgo;
00044   TIterator* _funcIter ;     
00045   TIterator* _coefIter ;    
00046 
00047   // TH1 containing the histogram of this pdf
00048   mutable SimpleCacheSentry _sentry; // !not to be serialized
00049   mutable TH1  *_cacheTotal;     
00050   // For additive morphing, histograms of fNominal, fUp and fDown
00051   // For multiplicative morphing, histograms of fNominal, log(fUp/fNominal), -log(fDown/fNominal);
00052   mutable TH1  **_cacheSingle; 
00053   mutable int  *_cacheSingleGood; 
00054 private:
00055 
00056   ClassDef(VerticalInterpHistPdf,1) // 
00057 
00058 protected:
00059   void setupCaches() const ;
00060   void syncTotal() const ;
00061   void syncComponent(int which) const ;
00062   // return a smooth function that is equal to +/-1 for |x| >= smoothRegion_ and it's null in zero
00063   inline double smoothStepFunc(double x) const { 
00064     if (fabs(x) >= _smoothRegion) return x > 0 ? +1 : -1;
00065     double xnorm = x/_smoothRegion, xnorm2 = xnorm*xnorm;
00066     return 0.125 * xnorm * (xnorm2 * (3.*xnorm2 - 10.) + 15);
00067   }
00068 };
00069 
00070 class FastVerticalInterpHistPdfBase : public RooAbsPdf {
00071 public:
00072 
00073   FastVerticalInterpHistPdfBase() ;
00074   FastVerticalInterpHistPdfBase(const char *name, const char *title, const RooArgSet &obs, const RooArgList& funcList, const RooArgList& coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1) ;
00075   FastVerticalInterpHistPdfBase(const FastVerticalInterpHistPdfBase& other, const char* name=0) ;
00076   virtual TObject* clone(const char* newname) const = 0; 
00077   virtual ~FastVerticalInterpHistPdfBase() ;
00078 
00079   Bool_t selfNormalized() const { return kTRUE; }
00080 
00081   Double_t evaluate() const = 0;
00082 
00083   const RooArgList& funcList() const { return _funcList ; }
00084   const RooArgList& coefList() const { return _coefList ; }
00085 
00087   struct Morph { FastTemplate sum; FastTemplate diff; };
00088 protected:
00089   RooRealProxy   _x;
00090   RooListProxy _funcList ;   //  List of component FUNCs
00091   RooListProxy _coefList ;  //  List of coefficients
00092   Double_t     _smoothRegion;
00093   Int_t        _smoothAlgo;
00094   TIterator* _funcIter ;     
00095   TIterator* _coefIter ;    
00096 
00097   // TH1 containing the histogram of this pdf
00098   mutable bool              _init; 
00099   mutable SimpleCacheSentry _sentry; 
00100 
00101   // For additive morphing, histograms of (fUp-f0)+(fDown-f0) and (fUp-f0)-(fDown-f0)
00102   // For multiplicative morphing, log(fUp/f0)+log(fDown/f0),  log(fUp/f0)-log(fDown/f0)
00103   mutable std::vector<Morph> _morphs;  
00104   mutable std::vector<RooAbsReal *> _morphParams; 
00105 
00106   // Prepare morphing data for a triplet of templates
00107   void syncMorph(Morph &out, const FastTemplate &nominal, FastTemplate &lo, FastTemplate &hi) const;
00108 
00109   // Do the vertical morphing from nominal value and morphs into cache. 
00110   // Do not normalize yet, as that depends on the dimension of the template
00111   void syncTotal(FastTemplate &cache, const FastTemplate &cacheNominal, const FastTemplate &cacheNominalLog) const ;
00112 
00113   // return a smooth function that is equal to +/-1 for |x| >= smoothRegion_ and it's null in zero
00114   inline double smoothStepFunc(double x) const { 
00115     if (fabs(x) >= _smoothRegion) return x > 0 ? +1 : -1;
00116     double xnorm = x/_smoothRegion, xnorm2 = xnorm*xnorm;
00117     return 0.125 * xnorm * (xnorm2 * (3.*xnorm2 - 10.) + 15);
00118   }
00119 
00120 private:
00121   ClassDef(FastVerticalInterpHistPdfBase,2) // 
00122 };
00123 
00124 
00125 class FastVerticalInterpHistPdf : public FastVerticalInterpHistPdfBase {
00126 public:
00127 
00128   FastVerticalInterpHistPdf() : FastVerticalInterpHistPdfBase() {}
00129   FastVerticalInterpHistPdf(const char *name, const char *title, const RooRealVar &x, const RooArgList& funcList, const RooArgList& coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1) :
00130     FastVerticalInterpHistPdfBase(name,title,RooArgSet(x),funcList,coefList,smoothRegion,smoothAlgo),
00131     _x("x","Independent variable",this,const_cast<RooRealVar&>(x)),
00132     _cache(), _cacheNominal(), _cacheNominalLog()  {}
00133   FastVerticalInterpHistPdf(const FastVerticalInterpHistPdf& other, const char* name=0) :
00134     FastVerticalInterpHistPdfBase(other, name),
00135     _x("x",this,other._x),
00136 #ifdef FastVerticalInterpHistPdf_CopyConstructor
00137     _cache(other._cache), _cacheNominal(other._cacheNominal), _cacheNominalLog(other._cacheNominalLog)  {}
00138 #else
00139     _cache(), _cacheNominal(), _cacheNominalLog()  {}
00140 #endif
00141   virtual TObject* clone(const char* newname) const { return new FastVerticalInterpHistPdf(*this,newname) ; }
00142   virtual ~FastVerticalInterpHistPdf() {}
00143 
00144   Double_t evaluate() const ;
00145 
00146   Bool_t hasCache()     const { return _cache.size() > 0; }
00147   Bool_t isCacheReady() const { return _cache.size() > 0 && _init; }
00148 protected:
00149   RooRealProxy   _x;
00150 
00152   mutable FastHisto _cache; 
00153 
00154   mutable FastHisto _cacheNominal; 
00155   mutable FastHisto _cacheNominalLog; 
00156 
00157   void setupCaches() const ;
00158   void syncTotal() const ;
00159   void syncComponent(int which) const ;
00160   void syncNominal() const ;
00161   void syncComponents(int dimension) const ;
00162 
00163 private:
00164   ClassDef(FastVerticalInterpHistPdf,1) // 
00165 };
00166 
00167 class FastVerticalInterpHistPdf2D : public FastVerticalInterpHistPdfBase {
00168 public:
00169 
00170   FastVerticalInterpHistPdf2D() : FastVerticalInterpHistPdfBase() {}
00173   FastVerticalInterpHistPdf2D(const char *name, const char *title, const RooRealVar &x, const RooRealVar &y, bool conditional, const RooArgList& funcList, const RooArgList& coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1) :
00174     FastVerticalInterpHistPdfBase(name,title,RooArgSet(x,y),funcList,coefList,smoothRegion,smoothAlgo),
00175     _x("x","Independent variable",this,const_cast<RooRealVar&>(x)),
00176     _y("y","Independent variable",this,const_cast<RooRealVar&>(y)),
00177     _conditional(conditional),
00178     _cache(), _cacheNominal(), _cacheNominalLog()  {}
00179   FastVerticalInterpHistPdf2D(const FastVerticalInterpHistPdf2D& other, const char* name=0) :
00180     FastVerticalInterpHistPdfBase(other, name),
00181     _x("x",this,other._x), _y("y",this,other._y), _conditional(other._conditional),
00182 #ifdef FastVerticalInterpHistPdf_CopyConstructor
00183     _cache(other._cache), _cacheNominal(other._cacheNominal), _cacheNominalLog(other._cacheNominalLog)  {}
00184 #else
00185     _cache(), _cacheNominal(), _cacheNominalLog()  {}
00186 #endif
00187   virtual TObject* clone(const char* newname) const { return new FastVerticalInterpHistPdf2D(*this,newname) ; }
00188   virtual ~FastVerticalInterpHistPdf2D() {}
00189 
00190   Double_t evaluate() const ;
00191 
00192   Bool_t selfNormalized() const { return kTRUE; }
00193   Bool_t conditional() const { return _conditional; }
00194 
00195   Bool_t hasCache()     const { return _cache.size() > 0; }
00196   Bool_t isCacheReady() const { return _cache.size() > 0 && _init; }
00197 protected:
00198   RooRealProxy _x, _y;
00199   bool _conditional;
00200 
00202   mutable FastHisto2D _cache; 
00203 
00204   mutable FastHisto2D _cacheNominal; 
00205   mutable FastHisto2D _cacheNominalLog; 
00206 
00207   void setupCaches() const ;
00208   void syncTotal() const ;
00209   void syncComponent(int which) const ;
00210   void syncNominal() const ;
00211   void syncComponents(int dimension) const ;
00212 
00213 private:
00214   ClassDef(FastVerticalInterpHistPdf2D,1) // 
00215 };
00216 
00217 
00218 
00219 
00220 
00221 #endif