CMS 3D CMS Logo

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