CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/DQMOffline/RecoB/interface/FlavourHistorgrams.h

Go to the documentation of this file.
00001 
00002 #ifndef FlavourHistograms_H
00003 #define FlavourHistograms_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 "TH1F.h"
00012 #include "TCanvas.h"
00013 #include "TROOT.h"
00014 #include "TSystem.h"
00015 #include "TStyle.h"
00016 
00017 #include "DQMOffline/RecoB/interface/Tools.h"
00018 #include "DQMOffline/RecoB/interface/HistoProviderDQM.h"
00019 #include <iostream>
00020 #include <string>
00021 
00022 //class DQMStore;
00023 
00024 //
00025 // class to describe Histo
00026 //
00027 template <class T>
00028 class FlavourHistograms {
00029 
00030 public:
00031 
00032   FlavourHistograms (const std::string& baseNameTitle_ , const std::string& baseNameDescription_ ,
00033                      const int& nBins_ , const double& lowerBound_ , const double& upperBound_ ,
00034                      const bool& statistics_ , const bool& plotLog_ , const bool& plotNormalized_ ,
00035                      const std::string& plotFirst_ , const bool& update, const std::string& folder, const unsigned int& mc) ;
00036 
00037   virtual ~FlavourHistograms () ;
00038 
00039 
00040   // define arrays (if needed)
00041 //   void defineArray ( int * dimension , int max , int indexToPlot ) ;
00042   
00043   // fill entry
00044   // For single variables and arrays (for arrays only a single index can be filled)
00045   void fill ( const int & flavour,  const T & variable) const;
00046   void fill ( const int & flavour,  const T & variable, const T & w) const;
00047 
00048   // For single variables and arrays
00049   void fill ( const int & flavour,  const T * variable) const;
00050 
00051 
00052   void settitle(const char* title) ;
00053   
00054   void plot (TPad * theCanvas = 0) ;
00055 
00056   void epsPlot(const std::string& name);
00057 
00058   // needed for efficiency computations -> this / b
00059   // (void : alternative would be not to overwrite the histos but to return a cloned HistoDescription)
00060   void divide ( const FlavourHistograms<T> & bHD ) const ;
00061 
00062   inline void SetMaximum(const double& max) { theMax = max;}
00063   inline void SetMinimum(const double& min) { theMin = min;}
00064 
00065 
00066   // trivial access functions
00067   inline std::string  baseNameTitle       () const { return theBaseNameTitle       ; }
00068   inline std::string  baseNameDescription () const { return theBaseNameDescription ; }
00069   inline int    nBins               () const { return theNBins               ; }
00070   inline double lowerBound          () const { return theLowerBound          ; } 
00071   inline double upperBound          () const { return theUpperBound          ; }
00072   inline bool   statistics          () const { return theStatistics          ; }
00073   inline bool   plotLog             () const { return thePlotLog             ; }
00074   inline bool   plotNormalized      () const { return thePlotNormalized      ; }
00075   inline std::string  plotFirst           () const { return thePlotFirst           ; }
00076   inline int* arrayDimension () const { return theArrayDimension; }
00077   inline int maxDimension () const {return theMaxDimension; }
00078   inline int indexToPlot () const {return theIndexToPlot; }
00079 
00080   // access to the histos
00081   inline TH1F * histo_all  () const { return theHisto_all->getTH1F()  ; }    
00082   inline TH1F * histo_d    () const { return theHisto_d ->getTH1F()   ; }    
00083   inline TH1F * histo_u    () const { return theHisto_u->getTH1F()    ; }
00084   inline TH1F * histo_s    () const { return theHisto_s->getTH1F()    ; }
00085   inline TH1F * histo_c    () const { return theHisto_c->getTH1F()    ; }
00086   inline TH1F * histo_b    () const { return theHisto_b->getTH1F()    ; }
00087   inline TH1F * histo_g    () const { return theHisto_g->getTH1F()    ; }
00088   inline TH1F * histo_ni   () const { return theHisto_ni->getTH1F()   ; }
00089   inline TH1F * histo_dus  () const { return theHisto_dus->getTH1F()  ; }
00090   inline TH1F * histo_dusg () const { return theHisto_dusg->getTH1F() ; }
00091 
00092   std::vector<TH1F*> getHistoVector() const;
00093   
00094 
00095 protected:
00096 
00097   void fillVariable ( const int & flavour , const T & var , const T & w) const;
00098   
00099   //
00100   // the data members
00101   //
00102   
00103 //   T *   theVariable   ;
00104 
00105   // for arrays
00106   int * theArrayDimension ;
00107   int   theMaxDimension ;
00108   int   theIndexToPlot ; // in case that not the complete array has to be plotted
00109 
00110   std::string  theBaseNameTitle ;
00111   std::string  theBaseNameDescription ;
00112   int    theNBins ;
00113   double theLowerBound ;
00114   double theUpperBound ;
00115   bool   theStatistics ;
00116   bool   thePlotLog    ;
00117   bool   thePlotNormalized ;
00118   std::string  thePlotFirst  ; // one character; gives flavour to plot first: l (light) , c , b
00119   double theMin, theMax;
00120 
00121   // the histos
00122   MonitorElement *theHisto_all  ;    
00123   MonitorElement *theHisto_d    ;    
00124   MonitorElement *theHisto_u    ;
00125   MonitorElement *theHisto_s    ;
00126   MonitorElement *theHisto_c    ;
00127   MonitorElement *theHisto_b    ;
00128   MonitorElement *theHisto_g    ;
00129   MonitorElement *theHisto_ni   ;
00130   MonitorElement *theHisto_dus  ;
00131   MonitorElement *theHisto_dusg ;
00132 
00133   //  DQMStore * dqmStore_; 
00134 
00135 
00136   // the canvas to plot
00137   TCanvas* theCanvas ;
00138   private:
00139   FlavourHistograms(){}
00140 
00141   unsigned int mcPlots_;
00142 
00143 } ;
00144 
00145 template <class T>
00146 FlavourHistograms<T>::FlavourHistograms (const std::string& baseNameTitle_ , const std::string& baseNameDescription_ ,
00147                                          const int& nBins_ , const double& lowerBound_ , const double& upperBound_ ,
00148                                          const bool& statistics_ , const bool& plotLog_ , const bool& plotNormalized_ ,
00149                                          const std::string& plotFirst_, const bool& update, const std::string& folder, const unsigned int& mc) :
00150   // BaseFlavourHistograms () ,
00151   // theVariable ( variable_ ) ,
00152   theMaxDimension(-1), theIndexToPlot(-1), theBaseNameTitle ( baseNameTitle_ ) , theBaseNameDescription ( baseNameDescription_ ) ,
00153   theNBins ( nBins_ ) , theLowerBound ( lowerBound_ ) , theUpperBound ( upperBound_ ) ,
00154   theStatistics ( statistics_ ) , thePlotLog ( plotLog_ ) , thePlotNormalized ( plotNormalized_ ) ,
00155   thePlotFirst ( plotFirst_ ), theMin(-1.), theMax(-1.), mcPlots_(mc)
00156 {
00157   // defaults for array dimensions
00158   theArrayDimension = 0  ;
00159     
00160   // check plot order string 
00161   if ( thePlotFirst == "l" || thePlotFirst == "c" || thePlotFirst == "b" ) {
00162     // OK
00163   }
00164   else {
00165     // not correct: print warning and set default (l)
00166     std::cout << "FlavourHistograms::FlavourHistograms : thePlotFirst was not correct : " << thePlotFirst << std::endl ;
00167     std::cout << "FlavourHistograms::FlavourHistograms : Set it to default value (l)! " << std::endl ;
00168     thePlotFirst = "l" ;
00169   }
00170 
00171   if (!update) {
00172     // book histos
00173     HistoProviderDQM prov("Btag",folder);
00174     if(mcPlots_%2==0) theHisto_all   = (prov.book1D( theBaseNameTitle + "ALL"  , theBaseNameDescription + " all jets"  , theNBins , theLowerBound , theUpperBound )) ; 
00175     else theHisto_all = 0;
00176     if (mcPlots_) {
00177       if(mcPlots_>2) {
00178         theHisto_d     = (prov.book1D ( theBaseNameTitle + "D"    , theBaseNameDescription + " d-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00179         theHisto_u     = (prov.book1D ( theBaseNameTitle + "U"    , theBaseNameDescription + " u-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00180         theHisto_s     = (prov.book1D ( theBaseNameTitle + "S"    , theBaseNameDescription + " s-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00181         theHisto_g     = (prov.book1D ( theBaseNameTitle + "G"    , theBaseNameDescription + " g-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00182         theHisto_dus   = (prov.book1D ( theBaseNameTitle + "DUS"  , theBaseNameDescription + " dus-jets"  , theNBins , theLowerBound , theUpperBound )) ; 
00183       }
00184       else{
00185         theHisto_d = 0;
00186         theHisto_u = 0;
00187         theHisto_s = 0;
00188         theHisto_g = 0;
00189         theHisto_dus = 0;
00190       }
00191       theHisto_c     = (prov.book1D ( theBaseNameTitle + "C"    , theBaseNameDescription + " c-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00192       theHisto_b     = (prov.book1D ( theBaseNameTitle + "B"    , theBaseNameDescription + " b-jets"    , theNBins , theLowerBound , theUpperBound )) ; 
00193       theHisto_ni    = (prov.book1D ( theBaseNameTitle + "NI"   , theBaseNameDescription + " ni-jets"   , theNBins , theLowerBound , theUpperBound )) ; 
00194       theHisto_dusg  = (prov.book1D ( theBaseNameTitle + "DUSG" , theBaseNameDescription + " dusg-jets" , theNBins , theLowerBound , theUpperBound )) ;
00195     }else{
00196       theHisto_d = 0;
00197       theHisto_u = 0;
00198       theHisto_s = 0;
00199       theHisto_c = 0;
00200       theHisto_b = 0;
00201       theHisto_g = 0;
00202       theHisto_ni = 0;
00203       theHisto_dus = 0;
00204       theHisto_dusg = 0;
00205     }
00206 
00207       // statistics if requested
00208     if ( theStatistics ) {
00209       if(theHisto_all) theHisto_all ->getTH1F()->Sumw2() ; 
00210       if (mcPlots_ ) {  
00211         if (mcPlots_>2 ) {
00212           theHisto_d   ->getTH1F()->Sumw2() ; 
00213           theHisto_u   ->getTH1F()->Sumw2() ; 
00214           theHisto_s   ->getTH1F()->Sumw2() ;
00215           theHisto_g   ->getTH1F()->Sumw2() ; 
00216           theHisto_dus ->getTH1F()->Sumw2() ; 
00217         } 
00218         theHisto_c   ->getTH1F()->Sumw2() ; 
00219         theHisto_b   ->getTH1F()->Sumw2() ; 
00220         theHisto_ni  ->getTH1F()->Sumw2() ; 
00221         theHisto_dusg->getTH1F()->Sumw2() ;
00222       }
00223     }
00224   } else {
00225     HistoProviderDQM prov("Btag",folder);
00226     if(theHisto_all) theHisto_all   = prov.access(theBaseNameTitle + "ALL" ) ; 
00227     if (mcPlots_) {  
00228       if (mcPlots_>2 ) {
00229         theHisto_d     = prov.access(theBaseNameTitle + "D"   ) ; 
00230         theHisto_u     = prov.access(theBaseNameTitle + "U"   ) ; 
00231         theHisto_s     = prov.access(theBaseNameTitle + "S"   ) ;
00232         theHisto_g     = prov.access(theBaseNameTitle + "G"   ) ; 
00233         theHisto_dus   = prov.access(theBaseNameTitle + "DUS" ) ; 
00234       } 
00235       theHisto_c     = prov.access(theBaseNameTitle + "C"   ) ; 
00236       theHisto_b     = prov.access(theBaseNameTitle + "B"   ) ; 
00237       theHisto_ni    = prov.access(theBaseNameTitle + "NI"  ) ; 
00238       theHisto_dusg  = prov.access(theBaseNameTitle + "DUSG") ;
00239     }
00240   }
00241 
00242   // defaults for other data members
00243   theCanvas = 0 ;
00244 }
00245 
00246 
00247 template <class T>
00248 FlavourHistograms<T>::~FlavourHistograms () {
00249   // delete the canvas*/
00250   delete theCanvas ;
00251 }
00252 
00253 
00254 // define arrays (if needed)
00255 // template <class T>
00256 // void FlavourHistograms<T>::defineArray ( int * dimension , int max , int indexToPlot ) {
00257 //   // indexToPlot < 0 if all to be plotted
00258 //   theArrayDimension = dimension ;
00259 //   theMaxDimension   = max ;
00260 //   theIndexToPlot    = indexToPlot ;
00261 // }
00262   
00263 // fill entry
00264 template <class T> void
00265 FlavourHistograms<T>::fill ( const int & flavour,  const T & variable) const 
00266 {
00267   // For single variables and arrays (for arrays only a single index can be filled)
00268   fillVariable ( flavour , variable, 1.) ;
00269 }
00270 
00271 template <class T> void
00272 FlavourHistograms<T>::fill ( const int & flavour,  const T & variable, const T & w) const 
00273 {
00274   // For single variables and arrays (for arrays only a single index can be filled)
00275   fillVariable ( flavour , variable, w) ;
00276 }
00277 
00278 template <class T> void
00279 FlavourHistograms<T>::fill ( const int & flavour,  const T * variable) const
00280 {
00281   if ( theArrayDimension == 0 ) {       
00282     // single variable
00283     fillVariable ( flavour , *variable, 1.) ;
00284   } else {
00285     // array      
00286     int iMax = (*theArrayDimension > theMaxDimension) ? theMaxDimension : *theArrayDimension ;
00287     //
00288     for ( int i = 0 ; i != iMax ; ++i ) {
00289       // check if only one index to be plotted (<0: switched off -> plot all)
00290       if ( ( theIndexToPlot < 0 ) || ( i == theIndexToPlot ) ) { 
00291         fillVariable ( flavour , *(variable+i) , 1.) ;
00292       }
00293     }
00294 
00295     // if single index to be filled but not enough entries: fill 0.0 (convention!)
00296     if ( theIndexToPlot >= iMax ) { 
00297       // cout << "==>> The index to be filled is too big -> fill 0.0 : " << theBaseNameTitle << " : " << theIndexToPlot << " >= " << iMax << endl ;
00298       const T& theZero = static_cast<T> ( 0.0 ) ;
00299       fillVariable ( flavour , theZero , 1.) ;
00300     }
00301   }
00302 } 
00303 
00304 
00305 template <class T>
00306 void FlavourHistograms<T>::settitle(const char* title) {
00307  if(theHisto_all) theHisto_all ->setAxisTitle(title) ;
00308     if (mcPlots_) {  
00309       if (mcPlots_>2 ) {
00310         theHisto_d   ->setAxisTitle(title) ;
00311         theHisto_u   ->setAxisTitle(title) ;
00312         theHisto_s   ->setAxisTitle(title) ;
00313         theHisto_g   ->setAxisTitle(title) ;
00314         theHisto_dus ->setAxisTitle(title) ;
00315       }
00316       theHisto_c   ->setAxisTitle(title) ;
00317       theHisto_b   ->setAxisTitle(title) ;
00318       theHisto_ni  ->setAxisTitle(title) ;
00319       theHisto_dusg->setAxisTitle(title) ;
00320     }
00321 }
00322 
00323 
00324 template <class T>
00325 void FlavourHistograms<T>::plot (TPad * theCanvas /* = 0 */) {
00326 
00327 //fixme:
00328   bool btppNI = false;
00329   bool btppColour = true;
00330 
00331   if (theCanvas)
00332     theCanvas->cd();
00333   
00334   RecoBTag::setTDRStyle()->cd();
00335   gPad->UseCurrentStyle();
00336 //   if ( !btppTitle ) gStyle->SetOptTitle ( 0 ) ;
00337 //   
00338 //   // here: plot histograms in a canvas
00339 //   theCanvas = new TCanvas ( "C" + theBaseNameTitle , "C" + theBaseNameDescription ,
00340 //                          btppXCanvas , btppYCanvas ) ;
00341 //   theCanvas->SetFillColor ( 0 ) ;
00342 //   theCanvas->cd  ( 1 ) ;
00343   gPad->SetLogy  ( 0 ) ;
00344   if ( thePlotLog ) gPad->SetLogy ( 1 ) ; 
00345   gPad->SetGridx ( 0 ) ;
00346   gPad->SetGridy ( 0 ) ;
00347   gPad->SetTitle ( 0 ) ;
00348 
00349   MonitorElement * histo[4];
00350   int col[4], lineStyle[4], markerStyle[4];
00351   int lineWidth = 1 ;
00352 
00353   const double markerSize = gPad->GetWh() * gPad->GetHNDC() / 500.;
00354 
00355   // default (l)
00356   histo[0] = theHisto_dusg ;
00357   //CW histo_1 = theHisto_dus ;
00358   histo[1] = theHisto_b ;
00359   histo[2] = theHisto_c ;
00360   histo[3]= 0 ;
00361 
00362   double max = theMax;
00363   if (theMax<=0.) {
00364     max = theHisto_dusg->getTH1F()->GetMaximum();
00365     if (theHisto_b->getTH1F()->GetMaximum() > max) max = theHisto_b->getTH1F()->GetMaximum();
00366     if (theHisto_c->getTH1F()->GetMaximum() > max) max = theHisto_c->getTH1F()->GetMaximum();
00367   }
00368 
00369   if (btppNI) {
00370     histo[3] = theHisto_ni ;
00371     if (theHisto_ni->getTH1F()->GetMaximum() > max) max = theHisto_ni->getTH1F()->GetMaximum();
00372   }
00373 
00374   if ( btppColour ) { // print colours 
00375     col[0] = 4 ;
00376     col[1] = 2 ;
00377     col[2] = 6 ;
00378     col[3] = 3 ;
00379     lineStyle[0] = 1 ;
00380     lineStyle[1] = 1 ;
00381     lineStyle[2] = 1 ;
00382     lineStyle[3] = 1 ;
00383     markerStyle[0] = 20 ;
00384     markerStyle[1] = 21 ;
00385     markerStyle[2] = 22 ;
00386     markerStyle[3] = 23 ;
00387    lineWidth = 1 ;
00388   }
00389   else { // different marker/line styles
00390     col[1] = 1 ;
00391     col[2] = 1 ;
00392     col[3] = 1 ;
00393     col[0] = 1 ;
00394     lineStyle[0] = 2 ;
00395     lineStyle[1] = 1 ;
00396     lineStyle[2] = 3 ;
00397     lineStyle[3] = 4 ;
00398     markerStyle[0] = 20 ;
00399     markerStyle[1] = 21 ;
00400     markerStyle[2] = 22 ;
00401     markerStyle[3] = 23 ;
00402   }
00403 
00404   // if changing order (NI stays always last)
00405   
00406   // c to plot first   
00407   if ( thePlotFirst == "c" ) {
00408     histo[0] = theHisto_c ;
00409     if ( btppColour  ) col[0] = 6 ;
00410     if ( !btppColour ) lineStyle[0] = 3 ;
00411     histo[2] = theHisto_dusg ;
00412     if ( btppColour  ) col[2] = 4 ;
00413     if ( !btppColour ) lineStyle[2] = 2 ;
00414   }
00415 
00416   // b to plot first   
00417   if ( thePlotFirst == "b" ) {
00418     histo[0] = theHisto_b ;
00419     if ( btppColour  ) col[0] = 2 ;
00420     if ( !btppColour ) lineStyle[0] = 1 ;
00421     histo[1] = theHisto_dusg ;
00422     if ( btppColour  ) col[1] = 4 ;
00423     if ( !btppColour ) lineStyle[1] = 2 ;
00424   }
00425 
00426 
00427   histo[0] ->getTH1F()->GetXaxis()->SetTitle ( theBaseNameDescription.c_str() ) ;
00428   histo[0] ->getTH1F()->GetYaxis()->SetTitle ( "Arbitrary Units" ) ;
00429   histo[0] ->getTH1F()->GetYaxis()->SetTitleOffset(1.25) ;
00430 
00431   for (int i=0; i != 4; ++i) {
00432     if (histo[i]== 0 ) continue;
00433     histo[i] ->getTH1F()->SetStats ( false ) ;
00434     histo[i] ->getTH1F()->SetLineStyle ( lineStyle[i] ) ;
00435     histo[i] ->getTH1F()->SetLineWidth ( lineWidth ) ;
00436     histo[i] ->getTH1F()->SetLineColor ( col[i] ) ;
00437     histo[i] ->getTH1F()->SetMarkerStyle ( markerStyle[i] ) ;
00438     histo[i] ->getTH1F()->SetMarkerColor ( col[i] ) ;
00439     histo[i] ->getTH1F()->SetMarkerSize ( markerSize ) ;
00440   }
00441 
00442   if ( thePlotNormalized ) {
00443 //   cout <<histo[0]->GetEntries() <<" "<< histo[1]->GetEntries() <<" "<< histo[2]->GetEntries()<<" "<<histo[3]->GetEntries()<<endl;
00444     if (histo[0]->getTH1F()->GetEntries() != 0) {histo[0]->getTH1F() ->DrawNormalized() ;} else {    histo[0]->getTH1F() ->SetMaximum(1.0);
00445 histo[0] ->getTH1F()->Draw() ;}
00446     if (histo[1]->getTH1F()->GetEntries() != 0) histo[1] ->getTH1F()->DrawNormalized("Same") ;
00447     if (histo[2]->getTH1F()->GetEntries() != 0) histo[2]->getTH1F() ->DrawNormalized("Same") ;
00448     if ((histo[3] != 0) && (histo[3]->getTH1F()->GetEntries() != 0))  histo[3] ->getTH1F()->DrawNormalized("Same") ;
00449   }
00450   else {
00451     histo[0]->getTH1F()->SetMaximum(max*1.05);
00452     if (theMin!=-1.) histo[0]->getTH1F()->SetMinimum(theMin);
00453     histo[0]->getTH1F()->Draw() ;
00454     histo[1]->getTH1F()->Draw("Same") ;
00455     histo[2]->getTH1F()->Draw("Same") ;
00456     if ( histo[3] != 0 ) histo[3]->getTH1F()->Draw("Same") ;
00457   }
00458 
00459 }
00460 
00461 template <class T>
00462 void FlavourHistograms<T>::epsPlot(const std::string& name)
00463 {
00464    TCanvas tc(theBaseNameTitle.c_str() , theBaseNameDescription.c_str());
00465    
00466    plot(&tc);
00467    tc.Print((name + theBaseNameTitle + ".eps").c_str());
00468 }
00469 
00470 
00471 // needed for efficiency computations -> this / b
00472 // (void : alternative would be not to overwrite the histos but to return a cloned HistoDescription)
00473 template <class T>
00474 void FlavourHistograms<T>::divide ( const FlavourHistograms<T> & bHD ) const {
00475   // divide histos using binomial errors
00476   //
00477   // ATTENTION: It's the responsability of the user to make sure that the HistoDescriptions
00478   //            involved in this operation have been constructed with the statistics option switched on!!
00479   //
00480   if(theHisto_all) theHisto_all  ->getTH1F()-> Divide ( theHisto_all->getTH1F()  , bHD.histo_all () , 1.0 , 1.0 , "b" ) ;    
00481     if (mcPlots_) {
00482       if (mcPlots_>2 ) {
00483         theHisto_d    ->getTH1F()-> Divide ( theHisto_d ->getTH1F()   , bHD.histo_d   () , 1.0 , 1.0 , "b" ) ;    
00484         theHisto_u    ->getTH1F()-> Divide ( theHisto_u ->getTH1F()   , bHD.histo_u   () , 1.0 , 1.0 , "b" ) ;
00485         theHisto_s    ->getTH1F()-> Divide ( theHisto_s ->getTH1F()   , bHD.histo_s   () , 1.0 , 1.0 , "b" ) ;
00486         theHisto_g    ->getTH1F()-> Divide ( theHisto_g  ->getTH1F()  , bHD.histo_g   () , 1.0 , 1.0 , "b" ) ;
00487         theHisto_dus  ->getTH1F()-> Divide ( theHisto_dus->getTH1F()  , bHD.histo_dus () , 1.0 , 1.0 , "b" ) ;
00488       }
00489       theHisto_c    ->getTH1F()-> Divide ( theHisto_c ->getTH1F()   , bHD.histo_c   () , 1.0 , 1.0 , "b" ) ;
00490       theHisto_b    ->getTH1F()-> Divide ( theHisto_b ->getTH1F()   , bHD.histo_b   () , 1.0 , 1.0 , "b" ) ;
00491       theHisto_ni   ->getTH1F()-> Divide ( theHisto_ni->getTH1F()   , bHD.histo_ni  () , 1.0 , 1.0 , "b" ) ;
00492       theHisto_dusg ->getTH1F()-> Divide ( theHisto_dusg->getTH1F() , bHD.histo_dusg() , 1.0 , 1.0 , "b" ) ;
00493     }
00494 }
00495   
00496 
00497 template <class T>
00498 void FlavourHistograms<T>::fillVariable ( const int & flavour , const T & var , const T & w) const {
00499 
00500   // all, except for the Jet Multiplicity which is not filled for each jets but for each events  
00501   if( (theBaseNameDescription != "Jet Multiplicity" || flavour == -1) && theHisto_all) theHisto_all->Fill ( var ,w) ;
00502 
00503   // flavour specific
00504   if (!mcPlots_ || (theBaseNameDescription == "Jet Multiplicity" && flavour == -1)) return;
00505 
00506   switch(flavour) {
00507     case 1:
00508       if (mcPlots_>2 ) {
00509         theHisto_d->Fill( var ,w);
00510         if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
00511       }
00512       if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
00513       return;
00514     case 2:
00515       if (mcPlots_>2 ) {
00516         theHisto_u->Fill( var ,w);
00517         if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
00518       }
00519       if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
00520       return;
00521     case 3:
00522       if (mcPlots_>2 ) {
00523         theHisto_s->Fill( var ,w);
00524         if(theBaseNameDescription != "Jet Multiplicity") theHisto_dus->Fill( var ,w);
00525       }
00526       if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
00527       return;
00528     case 4:
00529       theHisto_c->Fill( var ,w);
00530       return;
00531     case 5:
00532       theHisto_b->Fill( var ,w);
00533       return;
00534     case 21:
00535       if (mcPlots_>2 ) theHisto_g->Fill( var ,w);
00536       if(theBaseNameDescription != "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
00537       return;
00538     case 123:
00539       if (mcPlots_>2 && theBaseNameDescription == "Jet Multiplicity") theHisto_dus->Fill( var ,w);
00540       return;
00541     case 12321:
00542       if (theBaseNameDescription == "Jet Multiplicity") theHisto_dusg->Fill( var ,w);
00543       return;
00544     default:
00545       theHisto_ni->Fill( var ,w);
00546       return;
00547   }
00548 }
00549 
00550 template <class T>
00551 std::vector<TH1F*> FlavourHistograms<T>::getHistoVector() const
00552 {
00553   std::vector<TH1F*> histoVector;
00554   if(theHisto_all) histoVector.push_back ( theHisto_all->getTH1F() );
00555     if (mcPlots_) {  
00556       if (mcPlots_>2 ) {
00557         histoVector.push_back ( theHisto_d->getTH1F()   );
00558         histoVector.push_back ( theHisto_u->getTH1F()   );
00559         histoVector.push_back ( theHisto_s->getTH1F()   );
00560         histoVector.push_back ( theHisto_g ->getTH1F()  );
00561         histoVector.push_back ( theHisto_dus->getTH1F() );
00562       }
00563       histoVector.push_back ( theHisto_c->getTH1F()   );
00564       histoVector.push_back ( theHisto_b->getTH1F()   );
00565       histoVector.push_back ( theHisto_ni->getTH1F()  );
00566       histoVector.push_back ( theHisto_dusg->getTH1F());
00567     }
00568   return histoVector;
00569 }
00570 #endif