CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/TopQuarkAnalysis/TopPairBSM/interface/JetCombinatorics.h

Go to the documentation of this file.
00001 #ifndef JetCombinatorics_h
00002 #define JetCombinatorics_h
00003 
00015 #include "TLorentzVector.h"
00016 #include "TString.h"
00017 #include "TH1F.h"
00018 #include "TFile.h"
00019 #include "TMath.h"
00020 #include <map>
00021 #include <vector>
00022 #include <iostream>
00023 
00024 class Combo {
00025 
00026 
00027   public:
00028 
00029         Combo() {
00030                 
00031           MW = 84.2;//79.8;
00032           Mtop_h = 180.7;//175.;
00033           Mtop_l = 174.9;
00034           sigmaHadW = 10.5;//2.*7.6;
00035           sigmaHadt = 19.2;//2.*12.5;
00036           sigmaLept = 24.2;//2.*15.6;
00037 
00038           SumEt_ = 0.;
00039           usebtag_ = false;
00040           useMtop_ = true;
00041 
00042           useFlv_ = false;
00043           Wp_flv_ = Wq_flv_ = Hadb_flv_ = Lepb_flv_ = 1.;
00044         }
00045         ~Combo(){};
00046 
00047         void SetWp(TLorentzVector Wp) { Wp_ = Wp; }
00048         void SetWq(TLorentzVector Wq) { Wq_ = Wq; }
00049         void SetHadb(TLorentzVector Hadb) { Hadb_ = Hadb; }
00050         void SetLepW(TLorentzVector LepW) { LepW_ = LepW; }
00051         void SetLepb(TLorentzVector Lepb) { Lepb_ = Lepb; }
00052         // flavor corrections
00053         void ApplyFlavorCorrections(bool option=true){ useFlv_ = option;}
00054         void SetFlvCorrWp( double corr ) { Wp_flv_ = corr; }
00055         void SetFlvCorrWq( double corr ) { Wq_flv_ = corr; }
00056         void SetFlvCorrHadb( double corr ) { Hadb_flv_ = corr; }
00057         void SetFlvCorrLepb( double corr ) { Lepb_flv_ = corr; }
00058         // b tagging
00059         void SetWp_disc(double disc) { Wp_disc_ = disc;}
00060         void SetWq_disc(double disc) { Wq_disc_= disc;}
00061         void SetHadb_disc(double disc) { Hadb_disc_= disc;}
00062         void SetLepb_disc(double disc) { Lepb_disc_= disc;}
00063         void SetbDiscPdf(TString filename) { 
00064           pdffile_ = TFile::Open(filename);
00065           hdisc_b_ = (TH1F*) gDirectory->Get("hdiscNorm_b");
00066           hdisc_cl_ = (TH1F*) gDirectory->Get("hdiscNorm_cl");
00067         }
00068         void SetSigmas(int type=0) {
00069 
00070           // type == 0 take defaults
00071           if (type==1) {
00072             // JES +10%
00073             MW = 87.2;
00074             Mtop_h = 193.2;
00075             Mtop_l = 179.0;
00076             sigmaHadW = 13.0;
00077             sigmaHadt = 22.8;
00078             sigmaLept = 26.3;
00079           }
00080           if (type==-1) {
00081             // JES -10%
00082             MW = 81.6;
00083             Mtop_h = 169.3;
00084             Mtop_l = 171.4;
00085             sigmaHadW =8.9;
00086             sigmaHadt =17.9;
00087             sigmaLept =22.6;
00088           }
00089 
00090         }
00091         void Usebtagging(bool option = true) { usebtag_ = option;}
00092         void SetMinMassLepW( double mass ) { minMassLepW_ = mass; }
00093         void SetMaxMassLepW( double mass ) { maxMassLepW_ = mass; }
00094         void SetMinMassHadW( double mass ) { minMassHadW_ = mass; }
00095         void SetMaxMassHadW( double mass ) { maxMassHadW_ = mass; }
00096         void SetMinMassLepTop( double mass ) { minMassLepTop_ = mass; }
00097         void SetMaxMassLepTop( double mass ) { maxMassLepTop_ = mass; }
00098         void UseMtopConstraint(bool option=true) { useMtop_ = option; }
00099                 
00100         void analyze() {
00101 
00102                 if ( useFlv_ ) {
00103                         Wp_ = Wp_flv_ * Wp_;
00104                         Wq_ = Wq_flv_ * Wq_;
00105                         Hadb_ = Hadb_flv_ * Hadb_;
00106                         Lepb_ = Lepb_flv_ * Lepb_;
00107                 }
00108 
00109                 HadW_ = Wp_ + Wq_;
00110                 HadTop_ = HadW_ + Hadb_;
00111                 LepTop_ = LepW_ + Lepb_;
00112                 TopPair_ = HadTop_ + LepTop_;
00113 
00114                 //double sigmaHadW = 10.5;//2.*7.6;
00115                 //double sigmaHadt = 19.2;//2.*12.5;
00116                 //double sigmaLept = 24.2;//2.*15.6;
00117                 
00118                 double chiHadW = (HadW_.M() - MW)/sigmaHadW;
00119                 double chiHadt = (HadTop_.M() - Mtop_h)/sigmaHadt;
00120                 double chiLept = (LepTop_.M() - Mtop_l)/sigmaLept;
00121 
00122                 if ( useMtop_ ) {
00123                         chi2_ = chiHadW*chiHadW + chiHadt*chiHadt + chiLept*chiLept;
00124                         Ndof_ = 3;
00125                 } else {
00126                         chi2_ = chiHadW*chiHadW + (HadTop_.M() - LepTop_.M())*(HadTop_.M() - LepTop_.M())/(sigmaHadt*sigmaHadt+sigmaLept*sigmaLept);
00127                         Ndof_ = 2;
00128                 }
00129                 
00130                 SumEt_ = HadTop_.Pt();
00131 
00132                 if ( usebtag_ ) {
00133 
00134                         double gauss_norm = (2.)*TMath::Log(sigmaHadW*TMath::Sqrt(2*TMath::Pi())) +
00135                                 (2.)*TMath::Log(sigmaHadt*TMath::Sqrt(2*TMath::Pi())) + (2.)*TMath::Log(sigmaLept*TMath::Sqrt(2*TMath::Pi()));
00136 
00137                         double LR_Wp; double LR_Wq;
00138                         double LR_Hadb; double LR_Lepb;
00139 
00140                         double LR_den = 0;
00141                         LR_den = ( getPdfValue("cl", Wp_disc_) + getPdfValue("b", Wp_disc_));
00142                         if (LR_den == 0 ) LR_Wp = 1e-5;
00143                         else LR_Wp = getPdfValue( "cl", Wp_disc_ )/ LR_den;
00144 
00145                         LR_den = ( getPdfValue("cl", Wq_disc_) + getPdfValue("b", Wq_disc_));
00146                         if (LR_den == 0 ) LR_Wq = 1e-5;
00147                         else LR_Wq = getPdfValue( "cl", Wq_disc_ )/ LR_den;
00148 
00149                         LR_den = ( getPdfValue("cl", Hadb_disc_) + getPdfValue("b", Hadb_disc_));
00150                         if (LR_den == 0 ) LR_Hadb = 1e-5;
00151                         else LR_Hadb = getPdfValue( "b", Hadb_disc_ )/ LR_den;
00152 
00153                         LR_den = ( getPdfValue("cl", Lepb_disc_) + getPdfValue("b", Lepb_disc_));
00154                         if (LR_den == 0 ) LR_Lepb = 1e-5;
00155                         else LR_Lepb = getPdfValue( "b", Lepb_disc_ )/ LR_den;
00156 
00157                         double btag_norm = (-0.25-TMath::Log(4)/2);
00158                         double btag_N2LL = btag_norm*4.*( LR_Wp * TMath::Log(LR_Wp/4) + LR_Wq*TMath::Log(LR_Wq/4) + LR_Hadb*TMath::Log(LR_Hadb/4) + LR_Lepb*TMath::Log(LR_Lepb/4) );
00159                   
00160                         chi2_ += btag_N2LL + gauss_norm;
00161                         Ndof_ += 3;
00162                         pdffile_->Close();
00163                 }
00164         }
00165 
00166         TLorentzVector GetWp() { return Wp_; }
00167         TLorentzVector GetWq() { return Wq_; }
00168         TLorentzVector GetHadW() { return HadW_; }
00169         TLorentzVector GetLepW() { return LepW_; }
00170         TLorentzVector GetHadb() { return Hadb_; }
00171         TLorentzVector GetLepb() { return Lepb_; }
00172         TLorentzVector GetHadTop() { return HadTop_; }
00173         TLorentzVector GetLepTop() { return LepTop_; }
00174         TLorentzVector GetTopPair() { return TopPair_; }
00175         double GetChi2() { return chi2_; }
00176         double GetNdof() { return Ndof_; }
00177         double GetSumEt() { return SumEt_; }
00178         int GetIdHadb() { return IdHadb_;}
00179         int GetIdWp() { return IdWp_; }
00180         int GetIdWq() { return IdWq_; }
00181         int GetIdLepb() { return IdLepb_;}
00182         void SetIdHadb(int id) { IdHadb_ = id;}
00183         void SetIdWp(int id) { IdWp_ = id; }
00184         void SetIdWq(int id) { IdWq_ = id; }
00185         void SetIdLepb(int id) { IdLepb_ = id;}
00186         void Print() {
00187           std::cout << " jet Wp  : px = " << Wp_.Px() << " py = " <<  Wp_.Py() << " pz = " << Wp_.Pz() << " e = " << Wp_.E() << std::endl;
00188           std::cout << " jet Wq  : px = " << Wq_.Px() << " py = " <<  Wq_.Py() << " pz = " << Wq_.Pz() << " e = "<< Wq_.E() << std::endl;
00189           std::cout << " jet Hadb: px = " << Hadb_.Px() << " py = " <<  Hadb_.Py() <<" pz = " << Hadb_.Pz() <<" e = "<< Hadb_.E() << std::endl;
00190           std::cout << " jet Lepb: px = " << Lepb_.Px() << " py = " <<  Lepb_.Py() <<" pz = " << Lepb_.Pz() <<" e = "<< Lepb_.E() << std::endl;
00191           std::cout << " chi-squared = " << chi2_ << " sumEt = " << SumEt_ << std::endl;
00192         }
00193         double getPdfValue(std::string flavor, double disc) {
00194           double pdf= 0;
00195           TH1F *hpdf;
00196           if ( flavor == "b" ) hpdf = hdisc_b_;
00197           else hpdf = hdisc_cl_;
00198           int bin = hpdf->GetXaxis()->FindBin( disc );
00199           pdf = hpdf->GetBinContent( bin );
00200           if ( disc < -10 || disc >50 ) return 0;
00201           //if ( pdf == 0 ) return 1.e-7;
00202           return pdf;
00203         }
00204         
00205   private:
00206         
00207         TLorentzVector Wp_;
00208         TLorentzVector Wq_;
00209         TLorentzVector HadW_;
00210         TLorentzVector Hadb_;
00211         TLorentzVector HadTop_;
00212         TLorentzVector LepW_;
00213         TLorentzVector Lepb_;   
00214         TLorentzVector LepTop_;
00215         TLorentzVector TopPair_;
00216         
00217         bool usebtag_;
00218         bool useMtop_;
00219         double Wp_disc_;
00220         double Wq_disc_;
00221         double Hadb_disc_;
00222         double Lepb_disc_;
00223         TFile *pdffile_;
00224         TH1F *hdisc_b_;
00225         TH1F *hdisc_cl_;
00226 
00227         double Wp_flv_, Wq_flv_, Hadb_flv_, Lepb_flv_;
00228         bool useFlv_;
00229         double chi2_;
00230         double Ndof_;
00231         double SumEt_;
00232         double minMassLepW_;
00233         double maxMassLepW_;
00234         double minMassHadW_;
00235         double maxMassHadW_;
00236         
00237         double minMassLepTop_;
00238         double maxMassLepTop_;
00239 
00240         double MW;
00241         double Mtop_h;
00242         double Mtop_l;
00243         double sigmaHadW;
00244         double sigmaHadt;
00245         double sigmaLept;
00246 
00247 
00248         int IdHadb_;
00249         int IdWp_;
00250         int IdWq_;
00251         int IdLepb_;
00252         
00253 };
00254 
00255 struct minChi2
00256 {
00257   bool operator()(Combo s1, Combo s2) const
00258   {
00259     return s1.GetChi2() <= s2.GetChi2();
00260   }
00261 };
00262 
00263 struct maxSumEt
00264 {
00265   bool operator()(Combo s1, Combo s2) const
00266   {
00267     return s1.GetSumEt() >= s2.GetSumEt();
00268   }
00269 };
00270 
00271 
00272 class JetCombinatorics {
00273 
00274   public:
00275 
00276         JetCombinatorics();
00277         ~JetCombinatorics();
00278 
00279         void Verbose() {
00280           verbosef = true;
00281         }
00282 
00283         std::map< int, std::string > Combinatorics(int k, int max = 6);
00284         std::map< int, std::string > NestedCombinatorics();
00285 
00286         void FourJetsCombinations(std::vector<TLorentzVector> jets, std::vector<double> bdiscriminators );
00287         void SetFlavorCorrections(std::vector<double > vector ) { flavorCorrections_ = vector; }
00288         void SetMaxNJets(int n) { maxNJets_ = n; }
00289         Combo GetCombination(int n=0);
00290         Combo GetCombinationSumEt(int n=0);
00291         int GetNumberOfCombos() { return ( (int)allCombos_.size() ); } 
00292         //void SetCandidate( std::vector< TLorentzVector > JetCandidates );
00293 
00294         void SetSigmas(int type = 0) {
00295           SigmasTypef = type;
00296         }
00297         void SetLeptonicW( TLorentzVector LepW ) { theLepW_ = LepW; }
00298 
00299         void SetMinMassLepW( double mass ) { minMassLepW_ = mass; }
00300         void SetMaxMassLepW( double mass ) { maxMassLepW_ = mass; }
00301         void SetMinMassHadW( double mass ) { minMassHadW_ = mass; }
00302         void SetMaxMassHadW( double mass ) { maxMassHadW_ = mass; }
00303         void SetMinMassLepTop( double mass ) { minMassLepTop_ = mass; }
00304         void SetMaxMassLepTop( double mass ) { maxMassLepTop_ = mass; }
00305 
00306         void UsebTagging( bool option = true ) { UsebTagging_ = option; }
00307         void ApplyFlavorCorrection( bool option = true ) { UseFlv_ = option; }
00308         void UseMtopConstraint( bool option = true) { UseMtop_ = option; }
00309         void SetbTagPdf( TString name ) { bTagPdffilename_ = name; }
00310         void Clear();
00311 
00312         std::vector< TLorentzVector > TwoCombos();
00313         std::vector< TLorentzVector > ThreeCombos();
00314 
00315         void RemoveDuplicates( bool option) { removeDuplicates_ = option; }
00316 
00317         std::vector< TLorentzVector > GetComposites();
00318         void AnalyzeCombos();
00319 
00320 
00321   private:
00322 
00323         //int kcombos_;
00324         //int maxcombos_;
00325         int SigmasTypef;
00326         bool verbosef;
00327         std::map< int, std::string > Template4jCombos_;
00328         std::map< int, std::string > Template5jCombos_;
00329         std::map< int, std::string > Template6jCombos_;
00330         std::map< int, std::string > Template7jCombos_;
00331 
00332         std::vector< double > flavorCorrections_;
00333         int maxNJets_;
00334         bool UsebTagging_;
00335         bool UseMtop_;
00336         TString bTagPdffilename_;
00337         bool UseFlv_;
00338         
00339         TLorentzVector theLepW_;
00340 
00341         double minMassLepW_;
00342         double maxMassLepW_;
00343         double minMassHadW_;
00344         double maxMassHadW_;
00345         double minMassLepTop_;
00346         double maxMassLepTop_;
00347         
00348         std::map< Combo, int, minChi2 > allCombos_;
00349         std::map< Combo, int, maxSumEt > allCombosSumEt_;
00350 
00351         Double_t minPhi_;
00352         double chi2_;
00353         int ndf_;
00354         bool removeDuplicates_;
00355         
00356         std::vector< TLorentzVector > cand1_;
00357         std::vector< TLorentzVector > cand2_;
00358         std::vector< TLorentzVector > cand3_;
00359 
00360         //int nLists_;
00361         
00362         //std::vector< TLorentzVector > composites_;
00363         
00364 };
00365 
00366 #endif