CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/DQMOffline/RecoB/interface/FlavourHistorgrams2D.h

Go to the documentation of this file.
00001 
00002 #ifndef FlavourHistograms2D_H
00003 #define FlavourHistograms2D_H
00004 
00005 #include "DQMServices/Core/interface/DQMStore.h"
00006 #include "DQMServices/Core/interface/MonitorElement.h"
00007 #include "FWCore/ServiceRegistry/interface/Service.h" 
00008 
00009 // #include "BTagPlotPrintC.h"
00010 
00011 #include "TH2F.h"
00012 #include "TProfile.h"
00013 #include "TCanvas.h"
00014 #include "TROOT.h"
00015 #include "TSystem.h"
00016 #include "TStyle.h"
00017 
00018 #include "DQMOffline/RecoB/interface/Tools.h"
00019 #include "DQMOffline/RecoB/interface/HistoProviderDQM.h"
00020 
00021 #include <iostream>
00022 #include <vector>
00023 #include <string>
00024 //class DQMStore;
00025 
00026 //
00027 // class to describe Histo
00028 //
00029 template <class T, class G>
00030 class FlavourHistograms2D {
00031 
00032 public:
00033 
00034   FlavourHistograms2D (TString baseNameTitle_ , TString baseNameDescription_ ,
00035                       int nBinsX_ , double lowerBoundX_ , double upperBoundX_ ,
00036                       int nBinsY_ , double lowerBoundY_ , double upperBoundY_ ,
00037                       bool statistics_ , bool update, std::string folder, bool mc,
00038                       bool createProfile) ;
00039 
00040   virtual ~FlavourHistograms2D () ;
00041 
00042 
00043   // define arrays (if needed)
00044 //   void defineArray ( int * dimension , int max , int indexToPlot ) ;
00045   
00046   // fill entry
00047   // For single variables and arrays (for arrays only a single index can be filled)
00048   void fill ( const int & flavour,  const T & variableX, const G & variableY) const;
00049 
00050   // For single variables and arrays
00051   void fill ( const int & flavour,  const T * variableX, const G * variableY) const;
00052 
00053 
00054   void settitle(const char* titleX, const char* titleY) ;
00055   
00056   // needed for efficiency computations -> this / b
00057   // (void : alternative would be not to overwrite the histos but to return a cloned HistoDescription)
00058   void divide ( const FlavourHistograms2D<T, G> & bHD ) const ;
00059 
00060   inline void SetMaximum(const double& max) { theMax = max;}
00061   inline void SetMinimum(const double& min) { theMin = min;}
00062 
00063 
00064   // trivial access functions
00065   inline std::string  baseNameTitle       () const { return theBaseNameTitle       ; }
00066   inline std::string  baseNameDescription () const { return theBaseNameDescription ; }
00067   inline int    nBinsX              () const { return theNBinsX              ; }
00068   inline int    nBinsY              () const { return theNBinsY              ; }
00069   inline double lowerBoundX         () const { return theLowerBoundX         ; } 
00070   inline double upperBoundX         () const { return theUpperBoundX         ; }
00071   inline double lowerBoundY         () const { return theLowerBoundY         ; } 
00072   inline double upperBoundY         () const { return theUpperBoundY         ; }
00073   inline bool   statistics          () const { return theStatistics          ; }
00074 
00075   // access to the histos
00076   inline TH2F * histo_all  () const { return theHisto_all->getTH2F()  ; }    
00077   inline TH2F * histo_d    () const { return theHisto_d ->getTH2F()   ; }    
00078   inline TH2F * histo_u    () const { return theHisto_u->getTH2F()    ; }
00079   inline TH2F * histo_s    () const { return theHisto_s->getTH2F()    ; }
00080   inline TH2F * histo_c    () const { return theHisto_c->getTH2F()    ; }
00081   inline TH2F * histo_b    () const { return theHisto_b->getTH2F()    ; }
00082   inline TH2F * histo_g    () const { return theHisto_g->getTH2F()    ; }
00083   inline TH2F * histo_ni   () const { return theHisto_ni->getTH2F()   ; }
00084   inline TH2F * histo_dus  () const { return theHisto_dus->getTH2F()  ; }
00085   inline TH2F * histo_dusg () const { return theHisto_dusg->getTH2F() ; }
00086 
00087   TProfile * profile_all  () const { return theProfile_all->getTProfile()  ; }    
00088   TProfile * profile_d    () const { return theProfile_d ->getTProfile()   ; }    
00089   TProfile * profile_u    () const { return theProfile_u->getTProfile()    ; }
00090   TProfile * profile_s    () const { return theProfile_s->getTProfile()    ; }
00091   TProfile * profile_c    () const { return theProfile_c->getTProfile()    ; }
00092   TProfile * profile_b    () const { return theProfile_b->getTProfile()    ; }
00093   TProfile * profile_g    () const { return theProfile_g->getTProfile()    ; }
00094   TProfile * profile_ni   () const { return theProfile_ni->getTProfile()   ; }
00095   TProfile * profile_dus  () const { return theProfile_dus->getTProfile()  ; }
00096   TProfile * profile_dusg () const { return theProfile_dusg->getTProfile() ; }
00097 
00098   std::vector<TH2F*> getHistoVector() const;
00099 
00100   std::vector<TProfile*> getProfileVector() const;
00101 
00102   
00103 
00104 protected:
00105 
00106   void fillVariable ( const int & flavour , const T & varX , const G & varY ) const;
00107   
00108   //
00109   // the data members
00110   //
00111   
00112 //   T *   theVariable   ;
00113 
00114   // for arrays
00115   int * theArrayDimension ;
00116   int   theMaxDimension ;
00117   int   theIndexToPlot ; // in case that not the complete array has to be plotted
00118 
00119   std::string  theBaseNameTitle ;
00120   std::string  theBaseNameDescription ;
00121   int    theNBinsX ;
00122   int    theNBinsY ;
00123   double theLowerBoundX ;
00124   double theUpperBoundX ;
00125   double theLowerBoundY ;
00126   double theUpperBoundY ;
00127   bool   theStatistics ;
00128   double theMin, theMax;
00129 
00130   // the histos
00131   MonitorElement *theHisto_all  ;    
00132   MonitorElement *theHisto_d    ;    
00133   MonitorElement *theHisto_u    ;
00134   MonitorElement *theHisto_s    ;
00135   MonitorElement *theHisto_c    ;
00136   MonitorElement *theHisto_b    ;
00137   MonitorElement *theHisto_g    ;
00138   MonitorElement *theHisto_ni   ;
00139   MonitorElement *theHisto_dus  ;
00140   MonitorElement *theHisto_dusg ;
00141 
00142   // the profiles
00143   MonitorElement *theProfile_all  ;    
00144   MonitorElement *theProfile_d    ;    
00145   MonitorElement *theProfile_u    ;
00146   MonitorElement *theProfile_s    ;
00147   MonitorElement *theProfile_c    ;
00148   MonitorElement *theProfile_b    ;
00149   MonitorElement *theProfile_g    ;
00150   MonitorElement *theProfile_ni   ;
00151   MonitorElement *theProfile_dus  ;
00152   MonitorElement *theProfile_dusg ;
00153 
00154   //  DQMStore * dqmStore_; 
00155 
00156 
00157   // the canvas to plot
00158   private:
00159   FlavourHistograms2D(){}
00160 
00161   bool mcPlots_;
00162   bool createProfile_;
00163 
00164 } ;
00165 
00166 
00167 
00168 template <class T, class G>
00169 FlavourHistograms2D<T, G>::FlavourHistograms2D (TString baseNameTitle_ , TString baseNameDescription_ ,
00170                                           int nBinsX_ , double lowerBoundX_ , double upperBoundX_ ,
00171                                           int nBinsY_ , double lowerBoundY_ , double upperBoundY_ ,
00172                                           bool statistics_ ,
00173                                           bool update, std::string folder, bool mc, bool createProfile) :
00174   // BaseFlavourHistograms2D () ,
00175   // theVariable ( variable_ ) ,
00176   theMaxDimension(-1), theIndexToPlot(-1), theBaseNameTitle ( baseNameTitle_ ) , theBaseNameDescription ( baseNameDescription_ ) ,
00177   theNBinsX ( nBinsX_ ) , theNBinsY (nBinsY_), 
00178   theLowerBoundX ( lowerBoundX_ ) , theUpperBoundX ( upperBoundX_ ) ,
00179   theLowerBoundY ( lowerBoundY_ ) , theUpperBoundY ( upperBoundY_ ) ,
00180   theStatistics ( statistics_ ) , theMin(-1.), theMax(-1.), mcPlots_(mc), createProfile_(createProfile)
00181 {
00182   // defaults for array dimensions
00183   theArrayDimension = 0  ;
00184     
00185   if (!update) {
00186     // book histos
00187     HistoProviderDQM prov("Btag",folder);
00188     theHisto_all   = (prov.book2D( theBaseNameTitle + "ALL"  , theBaseNameDescription + " all jets"  , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY )) ; 
00189     if (mcPlots_) {  
00190       theHisto_d     = (prov.book2D ( theBaseNameTitle + "D"    , theBaseNameDescription + " d-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY )) ; 
00191       theHisto_u     = (prov.book2D ( theBaseNameTitle + "U"    , theBaseNameDescription + " u-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00192       theHisto_s     = (prov.book2D ( theBaseNameTitle + "S"    , theBaseNameDescription + " s-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00193       theHisto_c     = (prov.book2D ( theBaseNameTitle + "C"    , theBaseNameDescription + " c-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00194       theHisto_b     = (prov.book2D ( theBaseNameTitle + "B"    , theBaseNameDescription + " b-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00195       theHisto_g     = (prov.book2D ( theBaseNameTitle + "G"    , theBaseNameDescription + " g-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00196       theHisto_ni    = (prov.book2D ( theBaseNameTitle + "NI"   , theBaseNameDescription + " ni-jets"   , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00197       theHisto_dus   = (prov.book2D ( theBaseNameTitle + "DUS"  , theBaseNameDescription + " dus-jets"  , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00198       theHisto_dusg  = (prov.book2D ( theBaseNameTitle + "DUSG" , theBaseNameDescription + " dusg-jets" , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ;
00199     }else{
00200       theHisto_d = 0;
00201       theHisto_u = 0;
00202       theHisto_s = 0;
00203       theHisto_c = 0;
00204       theHisto_b = 0;
00205       theHisto_g = 0;
00206       theHisto_ni = 0;
00207       theHisto_dus = 0;
00208       theHisto_dusg = 0;
00209     }
00210 
00211     if (createProfile_) {
00212       theProfile_all = (prov.bookProfile( theBaseNameTitle + "_Profile_ALL" , theBaseNameDescription + " all jets" , theNBinsX, theLowerBoundX, theUpperBoundX, theNBinsY, theLowerBoundY, theUpperBoundY));
00213       if (mcPlots_) {
00214         theProfile_d     = (prov.bookProfile ( theBaseNameTitle + "_Profile_D"    , theBaseNameDescription + " d-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY )) ; 
00215         theProfile_u     = (prov.bookProfile ( theBaseNameTitle + "_Profile_U"    , theBaseNameDescription + " u-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00216         theProfile_s     = (prov.bookProfile ( theBaseNameTitle + "_Profile_S"    , theBaseNameDescription + " s-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00217         theProfile_c     = (prov.bookProfile ( theBaseNameTitle + "_Profile_C"    , theBaseNameDescription + " c-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00218         theProfile_b     = (prov.bookProfile ( theBaseNameTitle + "_Profile_B"    , theBaseNameDescription + " b-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00219         theProfile_g     = (prov.bookProfile ( theBaseNameTitle + "_Profile_G"    , theBaseNameDescription + " g-jets"    , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00220         theProfile_ni    = (prov.bookProfile ( theBaseNameTitle + "_Profile_NI"   , theBaseNameDescription + " ni-jets"   , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00221         theProfile_dus   = (prov.bookProfile ( theBaseNameTitle + "_Profile_DUS"  , theBaseNameDescription + " dus-jets"  , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ; 
00222         theProfile_dusg  = (prov.bookProfile ( theBaseNameTitle + "_Profile_DUSG" , theBaseNameDescription + " dusg-jets" , theNBinsX , theLowerBoundX , theUpperBoundX , theNBinsY, theLowerBoundY, theUpperBoundY)) ;
00223       } else{
00224           theProfile_d = 0;
00225           theProfile_u = 0;
00226           theProfile_s = 0;
00227           theProfile_c = 0;
00228           theProfile_b = 0;
00229           theProfile_g = 0;
00230           theProfile_ni = 0;
00231           theProfile_dus = 0;
00232           theProfile_dusg = 0;
00233       } 
00234     }  else {
00235           theProfile_all = 0;
00236           theProfile_d = 0;
00237           theProfile_u = 0;
00238           theProfile_s = 0;
00239           theProfile_c = 0;
00240           theProfile_b = 0;
00241           theProfile_g = 0;
00242           theProfile_ni = 0;
00243           theProfile_dus = 0;
00244           theProfile_dusg = 0;
00245     }
00246       // statistics if requested
00247     if ( theStatistics ) {
00248       theHisto_all ->getTH2F()->Sumw2() ; 
00249       if(createProfile)
00250         theProfile_all ->getTProfile()->Sumw2() ; 
00251       if (mcPlots_) {  
00252         theHisto_d   ->getTH2F()->Sumw2() ; 
00253         theHisto_u   ->getTH2F()->Sumw2() ; 
00254         theHisto_s   ->getTH2F()->Sumw2() ; 
00255         theHisto_c   ->getTH2F()->Sumw2() ; 
00256         theHisto_b   ->getTH2F()->Sumw2() ; 
00257         theHisto_g   ->getTH2F()->Sumw2() ; 
00258         theHisto_ni  ->getTH2F()->Sumw2() ; 
00259         theHisto_dus ->getTH2F()->Sumw2() ; 
00260         theHisto_dusg->getTH2F()->Sumw2() ;
00261 
00262         if(createProfile) {
00263           theProfile_d   ->getTProfile()->Sumw2() ; 
00264           theProfile_u   ->getTProfile()->Sumw2() ; 
00265           theProfile_s   ->getTProfile()->Sumw2() ; 
00266           theProfile_c   ->getTProfile()->Sumw2() ; 
00267           theProfile_b   ->getTProfile()->Sumw2() ; 
00268           theProfile_g   ->getTProfile()->Sumw2() ; 
00269           theProfile_ni  ->getTProfile()->Sumw2() ; 
00270           theProfile_dus ->getTProfile()->Sumw2() ; 
00271           theProfile_dusg->getTProfile()->Sumw2() ;
00272         }
00273       }
00274     }
00275   } else {
00276     HistoProviderDQM prov("Btag",folder);
00277     theHisto_all   = prov.access(theBaseNameTitle + "ALL" ) ; 
00278     if (mcPlots_) {  
00279       
00280       theHisto_d     = prov.access(theBaseNameTitle + "D"   ) ; 
00281       theHisto_u     = prov.access(theBaseNameTitle + "U"   ) ; 
00282       theHisto_s     = prov.access(theBaseNameTitle + "S"   ) ; 
00283       theHisto_c     = prov.access(theBaseNameTitle + "C"   ) ; 
00284       theHisto_b     =prov.access(theBaseNameTitle + "B"   ) ; 
00285       theHisto_g     =prov.access(theBaseNameTitle + "G"   ) ; 
00286       theHisto_ni    =prov.access(theBaseNameTitle + "NI"  ) ; 
00287       theHisto_dus   =prov.access(theBaseNameTitle + "DUS" ) ; 
00288       theHisto_dusg  =prov.access(theBaseNameTitle + "DUSG") ;
00289     }
00290 
00291     if(createProfile_) {
00292       theProfile_all = prov.access(theBaseNameTitle + "_Profile_ALL");
00293       if(mcPlots_) {
00294         theProfile_d     = prov.access(theBaseNameTitle + "_Profile_D"   ) ; 
00295         theProfile_u     = prov.access(theBaseNameTitle + "_Profile_U"   ) ; 
00296         theProfile_s     = prov.access(theBaseNameTitle + "_Profile_S"   ) ; 
00297         theProfile_c     = prov.access(theBaseNameTitle + "_Profile_C"   ) ; 
00298         theProfile_b     =prov.access(theBaseNameTitle + "_Profile_B"   ) ; 
00299         theProfile_g     =prov.access(theBaseNameTitle + "_Profile_G"   ) ; 
00300         theProfile_ni    =prov.access(theBaseNameTitle + "_Profile_NI"  ) ; 
00301         theProfile_dus   =prov.access(theBaseNameTitle + "_Profile_DUS" ) ; 
00302         theProfile_dusg  =prov.access(theBaseNameTitle + "_Profile_DUSG") ;
00303       }
00304     }
00305   }
00306 }
00307 
00308 
00309 template <class T, class G>
00310 FlavourHistograms2D<T, G>::~FlavourHistograms2D () {}
00311 
00312 
00313 // define arrays (if needed)
00314 // template <class T, class G>
00315 // void FlavourHistograms2D<T, G>::defineArray ( int * dimension , int max , int indexToPlot ) {
00316 //   // indexToPlot < 0 if all to be plotted
00317 //   theArrayDimension = dimension ;
00318 //   theMaxDimension   = max ;
00319 //   theIndexToPlot    = indexToPlot ;
00320 // }
00321   
00322 // fill entry
00323 template <class T, class G> void
00324 FlavourHistograms2D<T, G>::fill ( const int & flavour,  const T & variableX, const G & variableY) const 
00325 {
00326   // For single variables and arrays (for arrays only a single index can be filled)
00327   fillVariable ( flavour , variableX , variableY ) ;
00328 }
00329 
00330 template <class T, class G> void
00331 FlavourHistograms2D<T, G>::fill ( const int & flavour,  const T * variableX, const G * variableY) const
00332 {
00333   if ( theArrayDimension == 0 ) {       
00334     // single variable
00335     fillVariable ( flavour , *variableX, *variableY ) ;
00336   } else {
00337     // array      
00338     int iMax = *theArrayDimension ;
00339     if ( *theArrayDimension > theMaxDimension ) iMax = theMaxDimension ;
00340     //
00341     for ( int i = 0 ; i != iMax ; ++i ) {
00342       // check if only one index to be plotted (<0: switched off -> plot all)
00343       if ( ( theIndexToPlot < 0 ) || ( i == theIndexToPlot ) ) { 
00344         fillVariable ( flavour , *(variableX+i) , *(variableY+i) ) ;
00345       }
00346     }
00347 
00348     // if single index to be filled but not enough entries: fill 0.0 (convention!)
00349     if ( theIndexToPlot >= iMax ) { 
00350       // cout << "==>> The index to be filled is too big -> fill 0.0 : " << theBaseNameTitle << " : " << theIndexToPlot << " >= " << iMax << endl ;
00351       const T& theZeroT = static_cast<T> ( 0.0) ;
00352       const G& theZeroG = static_cast<T> ( 0.0 );
00353       fillVariable ( flavour , theZeroT , theZeroG ) ;
00354     }
00355   }
00356 } 
00357 
00358 
00359 template <class T, class G>
00360 void FlavourHistograms2D<T, G>::settitle(const char* titleX, const char* titleY) {
00361  if(theHisto_all) theHisto_all ->setAxisTitle(titleX) ;
00362  if(theHisto_all) theHisto_all ->setAxisTitle(titleY, 2) ;
00363     if (mcPlots_) {  
00364       
00365       if(theHisto_d)  theHisto_d   ->setAxisTitle(titleX) ;
00366       if(theHisto_u)  theHisto_u   ->setAxisTitle(titleX) ;
00367       if(theHisto_s)  theHisto_s   ->setAxisTitle(titleX) ;
00368       if(theHisto_c)  theHisto_c   ->setAxisTitle(titleX) ;
00369       if(theHisto_b)  theHisto_b   ->setAxisTitle(titleX) ;
00370       if(theHisto_g)  theHisto_g   ->setAxisTitle(titleX) ;
00371       if(theHisto_ni) theHisto_ni  ->setAxisTitle(titleX) ;
00372       if(theHisto_dus) theHisto_dus ->setAxisTitle(titleX) ;
00373       if(theHisto_dusg)theHisto_dusg->setAxisTitle(titleX) ;
00374       if(theHisto_d)  theHisto_d   ->setAxisTitle(titleY, 2) ;
00375       if(theHisto_u)  theHisto_u   ->setAxisTitle(titleY, 2) ;
00376       if(theHisto_s)  theHisto_s   ->setAxisTitle(titleY, 2) ;
00377       if(theHisto_c)  theHisto_c   ->setAxisTitle(titleY, 2) ;
00378       if(theHisto_b)  theHisto_b   ->setAxisTitle(titleY, 2) ;
00379       if(theHisto_g)  theHisto_g   ->setAxisTitle(titleY, 2) ;
00380       if(theHisto_ni) theHisto_ni  ->setAxisTitle(titleY, 2) ;
00381       if(theHisto_dus) theHisto_dus ->setAxisTitle(titleY, 2) ;
00382       if(theHisto_dusg)theHisto_dusg->setAxisTitle(titleY, 2) ;
00383     }
00384 
00385   if(createProfile_) {
00386     if(theProfile_all) theProfile_all ->setAxisTitle(titleX) ;
00387     if(theProfile_all) theProfile_all ->setAxisTitle(titleY, 2) ;
00388     if (mcPlots_ == true) {  
00389       
00390       if(theProfile_d)  theProfile_d   ->setAxisTitle(titleX) ;
00391       if(theProfile_u)  theProfile_u   ->setAxisTitle(titleX) ;
00392       if(theProfile_s)  theProfile_s   ->setAxisTitle(titleX) ;
00393       if(theProfile_c)  theProfile_c   ->setAxisTitle(titleX) ;
00394       if(theProfile_b)  theProfile_b   ->setAxisTitle(titleX) ;
00395       if(theProfile_g)  theProfile_g   ->setAxisTitle(titleX) ;
00396       if(theProfile_ni) theProfile_ni  ->setAxisTitle(titleX) ;
00397       if(theProfile_dus) theProfile_dus ->setAxisTitle(titleX) ;
00398       if(theProfile_dusg)theProfile_dusg->setAxisTitle(titleX) ;
00399       if(theProfile_d)  theProfile_d   ->setAxisTitle(titleY, 2) ;
00400       if(theProfile_u)  theProfile_u   ->setAxisTitle(titleY, 2) ;
00401       if(theProfile_s)  theProfile_s   ->setAxisTitle(titleY, 2) ;
00402       if(theProfile_c)  theProfile_c   ->setAxisTitle(titleY, 2) ;
00403       if(theProfile_b)  theProfile_b   ->setAxisTitle(titleY, 2) ;
00404       if(theProfile_g)  theProfile_g   ->setAxisTitle(titleY, 2) ;
00405       if(theProfile_ni) theProfile_ni  ->setAxisTitle(titleY, 2) ;
00406       if(theProfile_dus) theProfile_dus ->setAxisTitle(titleY, 2) ;
00407       if(theProfile_dusg)theProfile_dusg->setAxisTitle(titleY, 2) ;
00408     }
00409   }
00410 }
00411 
00412 // needed for efficiency computations -> this / b
00413 // (void : alternative would be not to overwrite the histos but to return a cloned HistoDescription)
00414 template <class T, class G>
00415 void FlavourHistograms2D<T, G>::divide ( const FlavourHistograms2D<T, G> & bHD ) const {
00416   // divide histos using binomial errors
00417   //
00418   // ATTENTION: It's the responsability of the user to make sure that the HistoDescriptions
00419   //            involved in this operation have been constructed with the statistics option switched on!!
00420   //
00421   theHisto_all  ->getTH2F()-> Divide ( theHisto_all->getTH2F()  , bHD.histo_all () , 1.0 , 1.0 , "b" ) ;    
00422     if (mcPlots_) {  
00423       theHisto_d    ->getTH2F()-> Divide ( theHisto_d ->getTH2F()   , bHD.histo_d   () , 1.0 , 1.0 , "b" ) ;    
00424       theHisto_u    ->getTH2F()-> Divide ( theHisto_u ->getTH2F()   , bHD.histo_u   () , 1.0 , 1.0 , "b" ) ;
00425       theHisto_s    ->getTH2F()-> Divide ( theHisto_s ->getTH2F()   , bHD.histo_s   () , 1.0 , 1.0 , "b" ) ;
00426       theHisto_c    ->getTH2F()-> Divide ( theHisto_c ->getTH2F()   , bHD.histo_c   () , 1.0 , 1.0 , "b" ) ;
00427       theHisto_b    ->getTH2F()-> Divide ( theHisto_b ->getTH2F()   , bHD.histo_b   () , 1.0 , 1.0 , "b" ) ;
00428       theHisto_g    ->getTH2F()-> Divide ( theHisto_g  ->getTH2F()  , bHD.histo_g   () , 1.0 , 1.0 , "b" ) ;
00429       theHisto_ni   ->getTH2F()-> Divide ( theHisto_ni->getTH2F()   , bHD.histo_ni  () , 1.0 , 1.0 , "b" ) ;
00430       theHisto_dus  ->getTH2F()-> Divide ( theHisto_dus->getTH2F()  , bHD.histo_dus () , 1.0 , 1.0 , "b" ) ;
00431       theHisto_dusg ->getTH2F()-> Divide ( theHisto_dusg->getTH2F() , bHD.histo_dusg() , 1.0 , 1.0 , "b" ) ;
00432     }
00433 }
00434   
00435 
00436 template <class T, class G>
00437 void FlavourHistograms2D<T, G>::fillVariable ( const int & flavour , const T & varX , const G & varY ) const {
00438   // all
00439   theHisto_all                ->Fill ( varX, varY ) ;
00440   if(createProfile_)
00441     theProfile_all->Fill( varX, varY );
00442   // flavour specific
00443   if (!mcPlots_) return;
00444 
00445   switch( flavour ) {
00446     case 1:
00447       theHisto_d->Fill( varX, varY );
00448       theHisto_dus->Fill( varX, varY );
00449       theHisto_dusg->Fill( varX, varY );
00450       if(createProfile_) {
00451         theProfile_d->Fill(varX, varY);
00452         theProfile_dus->Fill(varX, varY);
00453         theProfile_dusg->Fill(varX, varY);
00454       }
00455       return;
00456     case 2:
00457       theHisto_u->Fill( varX, varY );
00458       theHisto_dus->Fill( varX, varY );
00459       theHisto_dusg->Fill( varX, varY );
00460       if(createProfile_) {
00461         theProfile_u->Fill(varX, varY);
00462         theProfile_dus->Fill(varX, varY);
00463         theProfile_dusg->Fill(varX, varY);
00464       }
00465       return;
00466     case 3:
00467       theHisto_s->Fill( varX, varY );
00468       theHisto_dus->Fill( varX, varY );
00469       theHisto_dusg->Fill( varX, varY );
00470       if(createProfile_) {
00471         theProfile_s->Fill(varX, varY);
00472         theProfile_dus->Fill(varX, varY);
00473         theProfile_dusg->Fill(varX, varY);
00474       }
00475       return;
00476     case 4:
00477       theHisto_c->Fill( varX, varY );
00478       if(createProfile_) theProfile_c->Fill(varX, varY);
00479       return;
00480     case 5:
00481       theHisto_b->Fill( varX, varY );
00482       if(createProfile_) theProfile_b->Fill(varX, varY);
00483       return;
00484     case 21:
00485       theHisto_g->Fill( varX, varY );
00486       theHisto_dusg->Fill( varX, varY );
00487       if(createProfile_) {
00488         theProfile_g->Fill(varX, varY);
00489         theProfile_dusg->Fill(varX, varY);
00490       }
00491       return;
00492     default:
00493       theHisto_ni->Fill( varX, varY );
00494       if(createProfile_) theProfile_ni->Fill(varX, varY);
00495       return;
00496   }
00497 }
00498 
00499 template <class T, class G>
00500 std::vector<TH2F*> FlavourHistograms2D<T, G>::getHistoVector() const
00501 {
00502   std::vector<TH2F*> histoVector;
00503   histoVector.push_back ( theHisto_all->getTH2F() );
00504     if (mcPlots_) {  
00505       histoVector.push_back ( theHisto_d->getTH2F()   );
00506       histoVector.push_back ( theHisto_u->getTH2F()   );
00507       histoVector.push_back ( theHisto_s->getTH2F()   );
00508       histoVector.push_back ( theHisto_c->getTH2F()   );
00509       histoVector.push_back ( theHisto_b->getTH2F()   );
00510       histoVector.push_back ( theHisto_g ->getTH2F()  );
00511       histoVector.push_back ( theHisto_ni->getTH2F()  );
00512       histoVector.push_back ( theHisto_dus->getTH2F() );
00513       histoVector.push_back ( theHisto_dusg->getTH2F());
00514     }
00515   return histoVector;
00516 }
00517 
00518 template <class T, class G>
00519 std::vector<TProfile*> FlavourHistograms2D<T, G>::getProfileVector() const
00520 {
00521   std::vector<TProfile*> profileVector;
00522   if(createProfile_) {
00523     profileVector.push_back ( theProfile_all->getTProfile() );
00524       if (mcPlots_) {  
00525         profileVector.push_back ( theProfile_d->getTProfile()   );
00526         profileVector.push_back ( theProfile_u->getTProfile()   );
00527         profileVector.push_back ( theProfile_s->getTProfile()   );
00528         profileVector.push_back ( theProfile_c->getTProfile()   );
00529         profileVector.push_back ( theProfile_b->getTProfile()   );
00530         profileVector.push_back ( theProfile_g ->getTProfile()  );
00531         profileVector.push_back ( theProfile_ni->getTProfile()  );
00532         profileVector.push_back ( theProfile_dus->getTProfile() );
00533         profileVector.push_back ( theProfile_dusg->getTProfile());
00534       }
00535   }
00536   return profileVector;
00537 }
00538 
00539 #endif