CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/TopQuarkAnalysis/TopPairBSM/src/BoostedTopProducer.cc

Go to the documentation of this file.
00001 #include "BoostedTopProducer.h"
00002 #include "PhysicsTools/CandUtils/interface/AddFourMomenta.h"
00003 
00004 #include <string>
00005 #include <sstream>
00006 #include <iostream>
00007 
00008 using std::string;
00009 using std::cout;
00010 using std::endl;
00011 
00012 //
00013 // constants, enums and typedefs
00014 //
00015 
00016 //
00017 // static data member definitions
00018 //
00019 
00020 //
00021 // constructors and destructor
00022 //
00023 BoostedTopProducer::BoostedTopProducer(const edm::ParameterSet& iConfig) :
00024   eleLabel_   (iConfig.getParameter<edm::InputTag>  ("electronLabel")),
00025   muoLabel_   (iConfig.getParameter<edm::InputTag>  ("muonLabel")),
00026   jetLabel_   (iConfig.getParameter<edm::InputTag>  ("jetLabel")),
00027   metLabel_   (iConfig.getParameter<edm::InputTag>  ("metLabel")),
00028   solLabel_   (iConfig.getParameter<edm::InputTag>  ("solLabel")),
00029   caloIsoCut_ (iConfig.getParameter<double>         ("caloIsoCut") ),
00030   mTop_       (iConfig.getParameter<double>         ("mTop") )
00031 {
00032   //register products
00033   produces<std::vector<reco::CompositeCandidate> > ();
00034 }
00035 
00036 
00037 BoostedTopProducer::~BoostedTopProducer()
00038 {
00039 }
00040 
00041 
00042 //
00043 // member functions
00044 //
00045 
00046 // ------------ method called to produce the data  ------------
00047 void
00048 BoostedTopProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00049 {
00050 
00051   using namespace edm;
00052 
00053   bool debug = false;
00054 
00055   // -----------------------------------------------------
00056   // get the bare PAT objects
00057   // -----------------------------------------------------
00058    edm::Handle<std::vector<pat::Muon> > muonHandle;
00059    iEvent.getByLabel(muoLabel_,muonHandle);
00060    std::vector<pat::Muon> const & muons = *muonHandle;
00061   
00062    edm::Handle<std::vector<pat::Jet> > jetHandle;
00063    iEvent.getByLabel(jetLabel_,jetHandle);
00064    std::vector<pat::Jet> const & jets = *jetHandle;
00065 
00066    edm::Handle<std::vector<pat::Electron> > electronHandle;
00067    iEvent.getByLabel(eleLabel_,electronHandle);
00068    std::vector<pat::Electron> const & electrons = *electronHandle;
00069 
00070    edm::Handle<std::vector<pat::MET> > metHandle;
00071    iEvent.getByLabel(metLabel_,metHandle);
00072    std::vector<pat::MET> const & mets = *metHandle;
00073 
00074    // -----------------------------------------------------
00075    // Event Preselection:
00076    //    <= 1 isolated electron or muon
00077    //    >= 1 electron or muon
00078    //    >= 2 jets
00079    //    >= 1 missing et
00080    //
00081    // To explain:
00082    //    We want to look at leptons within "top jets" in some
00083    //    cases. This means the isolation will kill those events.
00084    //    However, if there IS an isolated lepton, we want only
00085    //    one of them. 
00086    // 
00087    //    So to select the prompt W lepton, the logic is:
00088    //    1. If there is an isolated lepton, accept it as the W lepton.
00089    //    2. Else, take the highest Pt lepton (possibly non-isolated)
00090    // 
00091    // -----------------------------------------------------
00092    bool preselection = true;
00093   
00094    // This will hold the prompt W lepton candidate, and a
00095    // maximum pt decision variable
00096    double maxWLeptonPt = -1;
00097    //reco::Candidate const * Wlepton = 0;
00098 
00099    // ----------------------
00100    // Find isolated muons, and highest pt lepton
00101    // ----------------------
00102    std::vector<pat::Muon>::const_iterator isolatedMuon     = muons.end();
00103    std::vector<pat::Muon>::const_iterator muon = muons.end();
00104    bool nIsolatedMuons = 0;
00105    std::vector<pat::Muon>::const_iterator muonIt = muons.begin(),
00106      muonEnd = muons.end();
00107    for (; muonIt != muonEnd; ++muonIt ) {
00108 
00109      // Find highest pt lepton
00110      double pt = muonIt->pt();
00111      if ( pt > maxWLeptonPt ) {
00112        maxWLeptonPt = pt;
00113        muon = muonIt;
00114      }
00115 
00116      // Find any isolated muons
00117      double caloIso = muonIt->caloIso();
00118      if ( caloIso >= 0 && caloIso < caloIsoCut_ ) {
00119        nIsolatedMuons++;
00120        isolatedMuon = muonIt;
00121      }
00122    }
00123 
00124    // ----------------------
00125    // Find isolated electrons, and highest pt lepton
00126    // ----------------------
00127    std::vector<pat::Electron>::const_iterator isolatedElectron     = electrons.end();
00128    std::vector<pat::Electron>::const_iterator electron = electrons.end();
00129    bool nIsolatedElectrons = 0;
00130    std::vector<pat::Electron>::const_iterator electronIt = electrons.begin(),
00131      electronEnd = electrons.end();
00132    for (; electronIt != electronEnd; ++electronIt ) {
00133 
00134      // Find highest pt lepton
00135      double pt = electronIt->pt();
00136      if ( pt > maxWLeptonPt ) {
00137        maxWLeptonPt = pt;
00138        electron = electronIt;
00139      }
00140 
00141      // Find any isolated electrons
00142      double caloIso = electronIt->caloIso();
00143      if ( caloIso >= 0 && caloIso < caloIsoCut_ ) {
00144        nIsolatedElectrons++;
00145        isolatedElectron = electronIt;
00146      }
00147    }
00148 
00149 
00150    // ----------------------
00151    // Now decide on the "prompt" lepton from the W:
00152    // Choose isolated leptons over all, and if no isolated,
00153    // then take highest pt lepton. 
00154    // ----------------------
00155    bool isMuon = true;
00156    if      ( isolatedMuon     != muonEnd     ) { muon     = isolatedMuon;     isMuon = true; }
00157    else if ( isolatedElectron != electronEnd ) { electron = isolatedElectron; isMuon = false; } 
00158    else {
00159      // Set to the highest pt lepton
00160      if      ( muon != muonEnd && electron == electronEnd ) isMuon = true;
00161      else if ( muon == muonEnd && electron != electronEnd ) isMuon = false;
00162      else if ( muon != muonEnd && electron != electronEnd ) {
00163        isMuon =  muon->pt() > electron->pt();
00164      }
00165    }
00166 
00167    // ----------------------
00168    // Veto events that have more than one isolated lepton
00169    // ----------------------
00170    int nIsolatedLeptons = nIsolatedMuons + nIsolatedElectrons;
00171    if ( nIsolatedLeptons > 1 ) {
00172      preselection = false;
00173    }
00174 
00175    // ----------------------
00176    // Veto events that have no prompt lepton candidates
00177    // ----------------------
00178    if ( muon == muonEnd && electron == electronEnd ) {
00179      preselection = false;
00180    }
00181 
00182    // ----------------------
00183    // Veto events with < 2 jets or no missing et
00184    // ----------------------
00185    if ( jets.size() < 2 ||
00186         mets.size() == 0 ) {
00187      preselection = false;
00188    }
00189 
00190    bool write = false;
00191 
00192 
00193    
00194    // -----------------------------------------------------
00195    //
00196    // CompositeCandidates to store the event solution.
00197    // This will take one of two forms:
00198    //    a) lv jj jj   Full reconstruction.
00199    //       
00200    //   ttbar->
00201    //       (hadt -> (hadW -> hadp + hadq) + hadb) + 
00202    //       (lept -> (lepW -> lepton + neutrino) + lepb)
00203    // 
00204    //    b) lv jj (j)  Partial reconstruction, associate 
00205    //                  at least 1 jet to the lepton 
00206    //                  hemisphere, and at least one jet in 
00207    //                  the opposite hemisphere.
00208    //
00209    //    ttbar->
00210    //        (hadt -> (hadJet1 [+ hadJet2] ) ) +
00211    //        (lept -> (lepW -> lepton + neutrino) + lepJet1 )
00212    //
00213    // There will also be two subcategories of (b) that 
00214    // will correspond to physics cases:
00215    // 
00216    //    b1)           Lepton is isolated: Moderate ttbar mass.
00217    //    b2)           Lepton is nonisolated: High ttbar mass. 
00218    //
00219    // -----------------------------------------------------
00220    reco::CompositeCandidate ttbar("ttbar");
00221    AddFourMomenta addFourMomenta;
00222 
00223 
00224    // Main decisions after preselection
00225    if ( preselection ) {
00226 
00227      if ( debug ) cout << "Preselection is satisfied" << endl;
00228 
00229      if ( debug ) cout << "Jets.size() = " << jets.size() << endl;
00230 
00231      // This will be modified for the z solution, so make a copy
00232      pat::MET              neutrino( mets[0] );
00233 
00234 
00235      // 1. First examine the low mass case with 4 jets and widely separated
00236      //    products. We take out the TtSemiLeptonicEvent from the TQAF and
00237      //    form the ttbar invariant mass.
00238      if ( jets.size() >= 4 ) {
00239 
00240        if ( debug ) cout << "Getting ttbar semileptonic solution" << endl;
00241 
00242        // get the ttbar semileptonic event solution if there are more than 3 jets
00243        edm::Handle< TtSemiLeptonicEvent > eSol;
00244        iEvent.getByLabel(solLabel_, eSol);
00245 
00246        // Have solution, continue
00247        if ( eSol.isValid() ) {
00248          if ( debug ) cout << "Got a nonzero size solution vector" << endl;
00249          // Just set the ttbar solution to the best ttbar solution from
00250          // TtSemiEvtSolutionMaker
00251          ttbar = eSol->eventHypo(TtSemiLeptonicEvent::kMVADisc);
00252          write = true;
00253        }
00254        // No ttbar solution with 4 jets, something is weird, print a warning
00255        else {
00256          edm::LogWarning("DataNotFound") << "BoostedTopProducer: Cannot find TtSemiEvtSolution\n";
00257        }
00258      }
00259      // 2. With 2 or 3 jets, we decide based on the separation between
00260      // the lepton and the closest jet in that hemisphere whether to
00261      // consider it "moderate" or "high" mass. 
00262      else if ( jets.size() == 2 || jets.size() == 3 ) {
00263 
00264        // ------------------------------------------------------------------
00265        // First create a leptonic W candidate
00266        // ------------------------------------------------------------------
00267        reco::CompositeCandidate lepW("lepW");
00268        
00269        if ( isMuon ) {
00270          if ( debug ) cout << "Adding muon as daughter" << endl;
00271          lepW.addDaughter(  *muon,     "muon" );
00272        } else {
00273          if ( debug ) cout << "Adding electron as daughter" << endl;
00274          lepW.addDaughter( *electron, "electron" );
00275        }
00276        if ( debug ) cout << "Adding neutrino as daughter" << endl;
00277        lepW.addDaughter  ( neutrino, "neutrino");
00278        addFourMomenta.set( lepW );
00279 
00280        //bool nuzHasComplex = false;
00281        METzCalculator zcalculator;
00282        
00283        zcalculator.SetMET( neutrino );
00284        if ( isMuon ) 
00285          zcalculator.SetMuon( *muon );
00286        else
00287          zcalculator.SetMuon( *electron ); // This name is misleading, should be setLepton
00288        double neutrinoPz = zcalculator.Calculate(1);// closest to the lepton Pz
00289        //if (zcalculator.IsComplex()) nuzHasComplex = true;
00290        // Set the neutrino pz
00291        neutrino.setPz( neutrinoPz );
00292 
00293        if ( debug ) cout << "Set neutrino pz to " << neutrinoPz << endl;
00294 
00295        // ------------------------------------------------------------------
00296        // Next ensure that there is a jet within the hemisphere of the 
00297        // leptonic W, and one in the opposite hemisphere
00298        // ------------------------------------------------------------------
00299        reco::CompositeCandidate hadt("hadt");
00300        reco::CompositeCandidate lept("lept");
00301        if ( debug ) cout << "Adding lepW as daughter" << endl;
00302        lept.addDaughter( lepW, "lepW" );
00303 
00304        std::string hadName("hadJet");
00305        std::string lepName("lepJet");
00306 
00307        // Get the W momentum
00308        TLorentzVector p4_W (lepW.px(), lepW.py(), lepW.pz(), lepW.energy() );
00309 
00310        // Loop over the jets
00311        std::vector<pat::Jet>::const_iterator jetit = jets.begin(),
00312          jetend = jets.end();
00313        unsigned long ii = 1; // Count by 1 for naming histograms
00314        for ( ; jetit != jetend; ++jetit, ++ii ) {
00315          // Get this jet's momentum
00316          TLorentzVector p4_jet (jetit->px(), jetit->py(), jetit->pz(), jetit->energy() );
00317 
00318          // Calculate psi (like DeltaR, only more invariant under Rapidity)
00319          double psi = Psi( p4_W, p4_jet, mTop_ );
00320 
00321          // Get jets that are in the leptonic hemisphere
00322          if ( psi < TMath::Pi() ) {
00323            // Add this jet to the leptonic top
00324            std::stringstream s;
00325            s << lepName << ii;
00326            if ( debug ) cout << "Adding daughter " << s.str() << endl;
00327            lept.addDaughter( *jetit, s.str() );
00328          }
00329          // Get jets that are in the hadronic hemisphere
00330          if ( psi > TMath::Pi() ) {
00331            // Add this jet to the hadronic top. We don't
00332            // make any W hypotheses in this case, since
00333            // we cannot determine which of the three
00334            // jets are merged.
00335            std::stringstream s;
00336            s << hadName << ii;
00337            if ( debug ) cout << "Adding daughter " << s.str() << endl;
00338            hadt.addDaughter( *jetit, s.str() );
00339            
00340          }
00341        } // end loop over jets
00342 
00343        addFourMomenta.set( lept );
00344        addFourMomenta.set( hadt );
00345 
00346        bool lepWHasJet = lept.numberOfDaughters() >= 2; // W and >= 1 jet
00347        bool hadWHasJet = hadt.numberOfDaughters() >= 1; // >= 1 jet
00348        if ( lepWHasJet && hadWHasJet ) {
00349          if ( debug ) cout << "Adding daughters lept and hadt" << endl;
00350          ttbar.addDaughter( lept, "lept");
00351          ttbar.addDaughter( hadt, "hadt");
00352          addFourMomenta.set( ttbar );
00353          write = true; 
00354        } // end of hadronic jet and leptonic jet
00355 
00356 
00357      } // end if there are 2 or 3 jets 
00358    
00359    } // end if preselection is satisfied
00360 
00361    // Write the solution to the event record   
00362    std::vector<reco::CompositeCandidate> ttbarList;
00363    if ( write ) {
00364      if ( debug ) cout << "Writing out" << endl;
00365      ttbarList.push_back( ttbar );
00366    }
00367    std::auto_ptr<std::vector<reco::CompositeCandidate> > pTtbar ( new std::vector<reco::CompositeCandidate>(ttbarList) );
00368    iEvent.put( pTtbar );
00369 
00370    
00371 }
00372 
00373 // ------------ method called once each job just before starting event loop  ------------
00374 void 
00375 BoostedTopProducer::beginJob(const edm::EventSetup&)
00376 {
00377 }
00378 
00379 // ------------ method called once each job just after ending the event loop  ------------
00380 void 
00381 BoostedTopProducer::endJob() {
00382 }
00383 
00384 double
00385 BoostedTopProducer::Psi(TLorentzVector p1, TLorentzVector p2, double mass) {
00386 
00387         TLorentzVector ptot = p1 + p2;
00388         Double_t theta1 = TMath::ACos( (p1.Vect().Dot(ptot.Vect()))/(p1.P()*ptot.P()) );
00389         Double_t theta2 = TMath::ACos( (p2.Vect().Dot(ptot.Vect()))/(p2.P()*ptot.P()) );
00390         //Double_t sign = 1.;
00391         //if ( (theta1+theta2) > (TMath::Pi()/2) ) sign = -1.;
00392         double th1th2 = theta1 + theta2;
00393         double psi = (p1.P()+p2.P())*TMath::Abs(TMath::Sin(th1th2))/(2.* mass );
00394         if ( th1th2 > (TMath::Pi()/2) )
00395                 psi = (p1.P()+p2.P())*( 1. + TMath::Abs(TMath::Cos(th1th2)))/(2.* mass );
00396         
00397         return psi;
00398 }
00399 
00400 //define this as a plug-in
00401 DEFINE_FWK_MODULE(BoostedTopProducer);