rivet is hosted by Hepforge, IPPP Durham
CMS_2015_I1370682.cc
Go to the documentation of this file.
00001 #include "Rivet/Analysis.hh"
00002 #include "Rivet/Math/LorentzTrans.hh"
00003 #include "Rivet/Projections/FinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 
00006 namespace Rivet {
00007 
00008 
00009   namespace { //< only visible in this compilation unit
00010 
00011     /// @brief Pseudo top finder
00012     ///
00013     /// Find top quark in the particle level.
00014     /// The definition is based on the agreement at the LHC working group.
00015     class PseudoTop : public FinalState {
00016     public:
00017       /// @name Standard constructors and destructors.
00018       //@{
00019 
00020       /// The default constructor. May specify the minimum and maximum
00021       /// pseudorapidity \f$ \eta \f$ and the min \f$ p_T \f$ (in GeV).
00022       PseudoTop(double lepR = 0.1, double lepMinPt = 20, double lepMaxEta = 2.4,
00023                 double jetR = 0.4, double jetMinPt = 30, double jetMaxEta = 4.7)
00024         : FinalState(-MAXDOUBLE, MAXDOUBLE, 0*GeV),
00025           _lepR(lepR), _lepMinPt(lepMinPt), _lepMaxEta(lepMaxEta),
00026           _jetR(jetR), _jetMinPt(jetMinPt), _jetMaxEta(jetMaxEta)
00027       {
00028         setName("PseudoTop");
00029       }
00030 
00031       enum TTbarMode {CH_NONE=-1, CH_FULLHADRON = 0, CH_SEMILEPTON, CH_FULLLEPTON};
00032       enum DecayMode {CH_HADRON = 0, CH_MUON, CH_ELECTRON};
00033 
00034       TTbarMode mode() const {
00035         if (!_isValid) return CH_NONE;
00036         if (_mode1 == CH_HADRON && _mode2 == CH_HADRON) return CH_FULLHADRON;
00037         else if ( _mode1 != CH_HADRON && _mode2 != CH_HADRON) return CH_FULLLEPTON;
00038         else return CH_SEMILEPTON;
00039       }
00040       DecayMode mode1() const {return _mode1;}
00041       DecayMode mode2() const {return _mode2;}
00042 
00043       /// Clone on the heap.
00044       virtual unique_ptr<Projection> clone() const {
00045         return unique_ptr<Projection>(new PseudoTop(*this));
00046       }
00047 
00048       //@}
00049 
00050     public:
00051       Particle t1() const {return _t1;}
00052       Particle t2() const {return _t2;}
00053       Particle b1() const {return _b1;}
00054       Particle b2() const {return _b2;}
00055       ParticleVector wDecays1() const {return _wDecays1;}
00056       ParticleVector wDecays2() const {return _wDecays2;}
00057       Jets jets() const {return _jets;}
00058       Jets bjets() const {return _bjets;}
00059       Jets ljets() const {return _ljets;}
00060 
00061     protected:
00062       // Apply the projection to the event
00063       void project(const Event& e); // override; ///< @todo Re-enable when C++11 allowed
00064       void cleanup(std::map<double, std::pair<size_t, size_t> >& v, const bool doCrossCleanup=false) const;
00065 
00066     private:
00067       const double _lepR, _lepMinPt, _lepMaxEta;
00068       const double _jetR, _jetMinPt, _jetMaxEta;
00069 
00070       //constexpr ///< @todo Re-enable when C++11 allowed
00071       static double _tMass; // = 172.5*GeV; ///< @todo Re-enable when C++11 allowed
00072       //constexpr ///< @todo Re-enable when C++11 allowed
00073       static double _wMass; // = 80.4*GeV; ///< @todo Re-enable when C++11 allowed
00074 
00075     private:
00076       bool _isValid;
00077       DecayMode _mode1, _mode2;
00078 
00079       Particle _t1, _t2;
00080       Particle _b1, _b2;
00081       ParticleVector _wDecays1, _wDecays2;
00082       Jets _jets, _bjets, _ljets;
00083 
00084     };
00085 
00086     // More implementation below the analysis code
00087 
00088   }
00089 
00090 
00091 
00092   /// Pseudo-top analysis from CMS
00093   class CMS_2015_I1370682 : public Analysis {
00094   public:
00095 
00096     CMS_2015_I1370682()
00097       : Analysis("CMS_2015_I1370682"),
00098         _applyCorrection(true),
00099         _doShapeOnly(false)
00100     {    }
00101 
00102 
00103     void init() {
00104       declare(PseudoTop(0.1, 20, 2.4, 0.5, 30, 2.4), "ttbar");
00105 
00106       // Lepton + Jet channel
00107       _hSL_topPt         = bookHisto1D("d15-x01-y01"); // 1/sigma dsigma/dpt(top)
00108       _hSL_topPtTtbarSys = bookHisto1D("d16-x01-y01"); // 1/sigma dsigma/dpt*(top)
00109       _hSL_topY          = bookHisto1D("d17-x01-y01"); // 1/sigma dsigma/dy(top)
00110       _hSL_ttbarDelPhi   = bookHisto1D("d18-x01-y01"); // 1/sigma dsigma/ddeltaphi(t,tbar)
00111       _hSL_topPtLead     = bookHisto1D("d19-x01-y01"); // 1/sigma dsigma/dpt(t1)
00112       _hSL_topPtSubLead  = bookHisto1D("d20-x01-y01"); // 1/sigma dsigma/dpt(t2)
00113       _hSL_ttbarPt       = bookHisto1D("d21-x01-y01"); // 1/sigma dsigma/dpt(ttbar)
00114       _hSL_ttbarY        = bookHisto1D("d22-x01-y01"); // 1/sigma dsigma/dy(ttbar)
00115       _hSL_ttbarMass     = bookHisto1D("d23-x01-y01"); // 1/sigma dsigma/dm(ttbar)
00116 
00117       // Dilepton channel
00118       _hDL_topPt         = bookHisto1D("d24-x01-y01"); // 1/sigma dsigma/dpt(top)
00119       _hDL_topPtTtbarSys = bookHisto1D("d25-x01-y01"); // 1/sigma dsigma/dpt*(top)
00120       _hDL_topY          = bookHisto1D("d26-x01-y01"); // 1/sigma dsigma/dy(top)
00121       _hDL_ttbarDelPhi   = bookHisto1D("d27-x01-y01"); // 1/sigma dsigma/ddeltaphi(t,tbar)
00122       _hDL_topPtLead     = bookHisto1D("d28-x01-y01"); // 1/sigma dsigma/dpt(t1)
00123       _hDL_topPtSubLead  = bookHisto1D("d29-x01-y01"); // 1/sigma dsigma/dpt(t2)
00124       _hDL_ttbarPt       = bookHisto1D("d30-x01-y01"); // 1/sigma dsigma/dpt(ttbar)
00125       _hDL_ttbarY        = bookHisto1D("d31-x01-y01"); // 1/sigma dsigma/dy(ttbar)
00126       _hDL_ttbarMass     = bookHisto1D("d32-x01-y01"); // 1/sigma dsigma/dm(ttbar)
00127 
00128     }
00129 
00130 
00131     void analyze(const Event& event) {
00132 
00133       // Get the ttbar candidate
00134       const PseudoTop& ttbar = apply<PseudoTop>(event, "ttbar");
00135       if ( ttbar.mode() == PseudoTop::CH_NONE ) vetoEvent;
00136 
00137       const FourMomentum& t1P4 = ttbar.t1().momentum();
00138       const FourMomentum& t2P4 = ttbar.t2().momentum();
00139       const double pt1 = std::max(t1P4.pT(), t2P4.pT());
00140       const double pt2 = std::min(t1P4.pT(), t2P4.pT());
00141       const double dPhi = deltaPhi(t1P4, t2P4);
00142       const FourMomentum ttP4 = t1P4 + t2P4;
00143       const FourMomentum t1P4AtCM = LorentzTransform::mkFrameTransformFromBeta(ttP4.betaVec()).transform(t1P4);
00144 
00145       const double weight = event.weight();
00146 
00147       if ( ttbar.mode() == PseudoTop::CH_SEMILEPTON ) {
00148         const Particle lCand1 = ttbar.wDecays1()[0]; // w1 dau0 is the lepton in the PseudoTop
00149         if (lCand1.pT() < 33*GeV || lCand1.abseta() > 2.1) vetoEvent;
00150         _hSL_topPt->fill(t1P4.pT(), weight);
00151         _hSL_topPt->fill(t2P4.pT(), weight);
00152         _hSL_topPtTtbarSys->fill(t1P4AtCM.pT(), weight);
00153         _hSL_topY->fill(t1P4.rapidity(), weight);
00154         _hSL_topY->fill(t2P4.rapidity(), weight);
00155         _hSL_ttbarDelPhi->fill(dPhi, weight);
00156         _hSL_topPtLead->fill(pt1, weight);
00157         _hSL_topPtSubLead->fill(pt2, weight);
00158         _hSL_ttbarPt->fill(ttP4.pT(), weight);
00159         _hSL_ttbarY->fill(ttP4.rapidity(), weight);
00160         _hSL_ttbarMass->fill(ttP4.mass(), weight);
00161       }
00162       else if ( ttbar.mode() == PseudoTop::CH_FULLLEPTON ) {
00163         const Particle lCand1 = ttbar.wDecays1()[0]; // dau0 are the lepton in the PseudoTop
00164         const Particle lCand2 = ttbar.wDecays2()[0]; // dau0 are the lepton in the PseudoTop
00165         if (lCand1.pT() < 20*GeV || lCand1.abseta() > 2.4) vetoEvent;
00166         if (lCand2.pT() < 20*GeV || lCand2.abseta() > 2.4) vetoEvent;
00167         _hDL_topPt->fill(t1P4.pT(), weight);
00168         _hDL_topPt->fill(t2P4.pT(), weight);
00169         _hDL_topPtTtbarSys->fill(t1P4AtCM.pT(), weight);
00170         _hDL_topY->fill(t1P4.rapidity(), weight);
00171         _hDL_topY->fill(t2P4.rapidity(), weight);
00172         _hDL_ttbarDelPhi->fill(dPhi, weight);
00173         _hDL_topPtLead->fill(pt1, weight);
00174         _hDL_topPtSubLead->fill(pt2, weight);
00175         _hDL_ttbarPt->fill(ttP4.pT(), weight);
00176         _hDL_ttbarY->fill(ttP4.rapidity(), weight);
00177         _hDL_ttbarMass->fill(ttP4.mass(), weight);
00178       }
00179 
00180     }
00181 
00182 
00183     void finalize() {
00184       if ( _applyCorrection ) {
00185         // Correction functions for TOP-12-028 paper, (parton bin height)/(pseudotop bin height)
00186         const double ch15[] = { 5.473609, 4.941048, 4.173346, 3.391191, 2.785644, 2.371346, 2.194161, 2.197167, };
00187         const double ch16[] = { 5.470905, 4.948201, 4.081982, 3.225532, 2.617519, 2.239217, 2.127878, 2.185918, };
00188         const double ch17[] = { 10.003667, 4.546519, 3.828115, 3.601018, 3.522194, 3.524694, 3.600951, 3.808553, 4.531891, 9.995370, };
00189         const double ch18[] = { 4.406683, 4.054041, 3.885393, 4.213646, };
00190         const double ch19[] = { 6.182537, 5.257703, 4.422280, 3.568402, 2.889408, 2.415878, 2.189974, 2.173210, };
00191         const double ch20[] = { 5.199874, 4.693318, 3.902882, 3.143785, 2.607877, 2.280189, 2.204124, 2.260829, };
00192         const double ch21[] = { 6.053523, 3.777506, 3.562251, 3.601356, 3.569347, 3.410472, };
00193         const double ch22[] = { 11.932351, 4.803773, 3.782709, 3.390775, 3.226806, 3.218982, 3.382678, 3.773653, 4.788191, 11.905338, };
00194         const double ch23[] = { 7.145255, 5.637595, 4.049882, 3.025917, 2.326430, 1.773824, 1.235329, };
00195 
00196         const double ch24[] = { 2.268193, 2.372063, 2.323975, 2.034655, 1.736793, };
00197         const double ch25[] = { 2.231852, 2.383086, 2.341894, 2.031318, 1.729672, 1.486993, };
00198         const double ch26[] = { 3.993526, 2.308249, 2.075136, 2.038297, 2.036302, 2.078270, 2.295817, 4.017713, };
00199         const double ch27[] = { 2.205978, 2.175010, 2.215376, 2.473144, };
00200         const double ch28[] = { 2.321077, 2.371895, 2.338871, 2.057821, 1.755382, };
00201         const double ch29[] = { 2.222707, 2.372591, 2.301688, 1.991162, 1.695343, };
00202         const double ch30[] = { 2.599677, 2.026855, 2.138620, 2.229553, };
00203         const double ch31[] = { 5.791779, 2.636219, 2.103642, 1.967198, 1.962168, 2.096514, 2.641189, 5.780828, };
00204         const double ch32[] = { 2.006685, 2.545525, 2.477745, 2.335747, 2.194226, 2.076500, };
00205 
00206         applyCorrection(_hSL_topPt, ch15);
00207         applyCorrection(_hSL_topPtTtbarSys, ch16);
00208         applyCorrection(_hSL_topY, ch17);
00209         applyCorrection(_hSL_ttbarDelPhi, ch18);
00210         applyCorrection(_hSL_topPtLead, ch19);
00211         applyCorrection(_hSL_topPtSubLead, ch20);
00212         applyCorrection(_hSL_ttbarPt, ch21);
00213         applyCorrection(_hSL_ttbarY, ch22);
00214         applyCorrection(_hSL_ttbarMass, ch23);
00215 
00216         applyCorrection(_hDL_topPt, ch24);
00217         applyCorrection(_hDL_topPtTtbarSys, ch25);
00218         applyCorrection(_hDL_topY, ch26);
00219         applyCorrection(_hDL_ttbarDelPhi, ch27);
00220         applyCorrection(_hDL_topPtLead, ch28);
00221         applyCorrection(_hDL_topPtSubLead, ch29);
00222         applyCorrection(_hDL_ttbarPt, ch30);
00223         applyCorrection(_hDL_ttbarY, ch31);
00224         applyCorrection(_hDL_ttbarMass, ch32);
00225       }
00226 
00227       if ( _doShapeOnly ) {
00228         normalize(_hSL_topPt        );
00229         normalize(_hSL_topPtTtbarSys);
00230         normalize(_hSL_topY         );
00231         normalize(_hSL_ttbarDelPhi  );
00232         normalize(_hSL_topPtLead    );
00233         normalize(_hSL_topPtSubLead );
00234         normalize(_hSL_ttbarPt      );
00235         normalize(_hSL_ttbarY       );
00236         normalize(_hSL_ttbarMass    );
00237 
00238         normalize(_hDL_topPt        );
00239         normalize(_hDL_topPtTtbarSys);
00240         normalize(_hDL_topY         );
00241         normalize(_hDL_ttbarDelPhi  );
00242         normalize(_hDL_topPtLead    );
00243         normalize(_hDL_topPtSubLead );
00244         normalize(_hDL_ttbarPt      );
00245         normalize(_hDL_ttbarY       );
00246         normalize(_hDL_ttbarMass    );
00247       }
00248       else {
00249         const double s = 1./sumOfWeights();
00250         scale(_hSL_topPt        , s);
00251         scale(_hSL_topPtTtbarSys, s);
00252         scale(_hSL_topY         , s);
00253         scale(_hSL_ttbarDelPhi  , s);
00254         scale(_hSL_topPtLead    , s);
00255         scale(_hSL_topPtSubLead , s);
00256         scale(_hSL_ttbarPt      , s);
00257         scale(_hSL_ttbarY       , s);
00258         scale(_hSL_ttbarMass    , s);
00259         scale(_hDL_topPt        , s);
00260         scale(_hDL_topPtTtbarSys, s);
00261         scale(_hDL_topY         , s);
00262         scale(_hDL_ttbarDelPhi  , s);
00263         scale(_hDL_topPtLead    , s);
00264         scale(_hDL_topPtSubLead , s);
00265         scale(_hDL_ttbarPt      , s);
00266         scale(_hDL_ttbarY       , s);
00267         scale(_hDL_ttbarMass    , s);
00268       }
00269 
00270     }
00271 
00272 
00273     void applyCorrection(Histo1DPtr h, const double* cf) {
00274       vector<YODA::HistoBin1D>& bins = h->bins();
00275       for (size_t i=0, n=bins.size(); i<n; ++i ) {
00276         const double s = cf[i];
00277         YODA::HistoBin1D& bin = bins[i];
00278         bin.scaleW(s);
00279       }
00280     }
00281 
00282 
00283   private:
00284 
00285     const bool _applyCorrection, _doShapeOnly;
00286     Histo1DPtr _hSL_topPt, _hSL_topPtTtbarSys, _hSL_topY, _hSL_ttbarDelPhi, _hSL_topPtLead,
00287       _hSL_topPtSubLead, _hSL_ttbarPt, _hSL_ttbarY, _hSL_ttbarMass;
00288     Histo1DPtr _hDL_topPt, _hDL_topPtTtbarSys, _hDL_topY, _hDL_ttbarDelPhi, _hDL_topPtLead,
00289       _hDL_topPtSubLead, _hDL_ttbarPt, _hDL_ttbarY, _hDL_ttbarMass;
00290 
00291   };
00292 
00293 
00294 
00295   DECLARE_RIVET_PLUGIN(CMS_2015_I1370682);
00296 
00297 
00298   ///////////////
00299 
00300   // More PseudoTop implementation
00301   namespace {
00302 
00303 
00304     double PseudoTop::_tMass = 172.5*GeV;
00305     double PseudoTop::_wMass = 80.4*GeV;
00306 
00307 
00308     void PseudoTop::cleanup(map<double, pair<size_t, size_t> >& v, const bool doCrossCleanup) const {
00309       vector<map<double, pair<size_t, size_t> >::iterator> toErase;
00310       set<size_t> usedLeg1, usedLeg2;
00311       if ( !doCrossCleanup ) {
00312         /// @todo Reinstate when C++11 allowed: for (auto key = v.begin(); key != v.end(); ++key) {
00313         for (map<double, pair<size_t, size_t> >::iterator key = v.begin(); key != v.end(); ++key) {
00314           const size_t leg1 = key->second.first;
00315           const size_t leg2 = key->second.second;
00316           if (usedLeg1.find(leg1) == usedLeg1.end() and
00317               usedLeg2.find(leg2) == usedLeg2.end()) {
00318             usedLeg1.insert(leg1);
00319             usedLeg2.insert(leg2);
00320           } else {
00321             toErase.push_back(key);
00322           }
00323         }
00324       }
00325       else {
00326         /// @todo Reinstate when C++11 allowed: for (auto key = v.begin(); key != v.end(); ++key) {
00327         for (map<double, pair<size_t, size_t> >::iterator key = v.begin(); key != v.end(); ++key) {
00328           const size_t leg1 = key->second.first;
00329           const size_t leg2 = key->second.second;
00330           if (usedLeg1.find(leg1) == usedLeg1.end() and
00331               usedLeg1.find(leg2) == usedLeg1.end()) {
00332             usedLeg1.insert(leg1);
00333             usedLeg1.insert(leg2);
00334           } else {
00335             toErase.push_back(key);
00336           }
00337         }
00338       }
00339       /// @todo Reinstate when C++11 allowed:  for (auto& key : toErase) v.erase(key);
00340       for (size_t i = 0; i < toErase.size(); ++i) v.erase(toErase[i]);
00341     }
00342 
00343 
00344     void PseudoTop::project(const Event& e) {
00345       // Leptons : do the lepton clustering anti-kt R=0.1 using stable photons and leptons not from hadron decay
00346       // Neutrinos : neutrinos not from hadron decay
00347       // MET : vector sum of all invisible particles in x-y plane
00348       // Jets : anti-kt R=0.4 using all particles excluding neutrinos and particles used in lepton clustering
00349       //        add ghost B hadrons during the jet clustering to identify B jets.
00350 
00351       // W->lv : dressed lepton and neutrino pairs
00352       // W->jj : light flavored dijet
00353       // W candidate : select lv or jj pairs which minimise |mW1-80.4|+|mW2-80.4|
00354       //               lepton-neutrino pair will be selected with higher priority
00355 
00356       // t->Wb : W candidate + b jet
00357       // t candidate : select Wb pairs which minimise |mtop1-172.5|+|mtop2-172.5|
00358 
00359       _isValid = false;
00360       _theParticles.clear();
00361       _wDecays1.clear();
00362       _wDecays2.clear();
00363       _jets.clear();
00364       _bjets.clear();
00365       _ljets.clear();
00366       _mode1 = _mode2 = CH_HADRON;
00367 
00368       // Collect final state particles
00369       Particles pForLep, pForJet;
00370       Particles neutrinos; // Prompt neutrinos
00371       /// @todo Avoid this unsafe jump into HepMC -- all this can be done properly via VisibleFS and HeavyHadrons projections
00372       for (const GenParticle* p : Rivet::particles(e.genEvent())) {
00373         const int status = p->status();
00374         const int pdgId = p->pdg_id();
00375         if (status == 1) {
00376           Particle rp = *p;
00377           if (!PID::isHadron(pdgId) && !rp.fromHadron()) {
00378             // Collect particles not from hadron decay
00379             if (rp.isNeutrino()) {
00380               // Prompt neutrinos are kept in separate collection
00381               neutrinos.push_back(rp);
00382             } else if (pdgId == 22 || rp.isLepton()) {
00383               // Leptons and photons for the dressing
00384               pForLep.push_back(rp);
00385             }
00386           } else if (!rp.isNeutrino()) {
00387             // Use all particles from hadron decay
00388             pForJet.push_back(rp);
00389           }
00390         } else if (PID::isHadron(pdgId) && PID::hasBottom(pdgId)) {
00391           // NOTE: Consider B hadrons with pT > 5GeV - not in CMS proposal
00392           //if ( p->momentum().perp() < 5 ) continue;
00393 
00394           // Do unstable particles, to be used in the ghost B clustering
00395           // Use last B hadrons only
00396           bool isLast = true;
00397           for (GenParticle* pp : Rivet::particles(p->end_vertex(), HepMC::children)) {
00398             if (PID::hasBottom(pp->pdg_id())) {
00399               isLast = false;
00400               break;
00401             }
00402           }
00403           if (!isLast) continue;
00404 
00405           // Rescale momentum by 10^-20
00406           Particle ghost(pdgId, FourMomentum(p->momentum())*1e-20/p->momentum().rho());
00407           pForJet.push_back(ghost);
00408         }
00409       }
00410 
00411       // Start object building from trivial thing - prompt neutrinos
00412       sortByPt(neutrinos);
00413 
00414       // Proceed to lepton dressing
00415       const PseudoJets lep_pjs = mkPseudoJets(pForLep);
00416       const fastjet::JetDefinition lep_jdef(fastjet::antikt_algorithm, _lepR);
00417       const Jets leps_all = mkJets(fastjet::ClusterSequence(lep_pjs, lep_jdef).inclusive_jets());
00418       const Jets leps_sel = sortByPt(filterBy(leps_all, Cuts::pT > _lepMinPt));
00419       // FastJets fjLep(FastJets::ANTIKT, _lepR);
00420       // fjLep.calc(pForLep);
00421 
00422       Jets leptons;
00423       vector<int> leptonsId;
00424       set<int> dressedIdxs;
00425       for (const Jet& lep : leps_sel) {
00426         if (lep.abseta() > _lepMaxEta) continue;
00427         double leadingPt = -1;
00428         int leptonId = 0;
00429         for (const Particle& p : lep.particles()) {
00430           /// @warning Barcodes aren't future-proof in HepMC
00431           dressedIdxs.insert(p.genParticle()->barcode());
00432           if (p.isLepton() && p.pT() > leadingPt) {
00433             leadingPt = p.pT();
00434             leptonId = p.pid();
00435           }
00436         }
00437         if (leptonId == 0) continue;
00438         leptons.push_back(lep);
00439         leptonsId.push_back(leptonId);
00440       }
00441 
00442       // Re-use particles not used in lepton dressing
00443       for (const Particle& rp : pForLep) {
00444         /// @warning Barcodes aren't future-proof in HepMC
00445         const int barcode = rp.genParticle()->barcode();
00446         // Skip if the particle is used in dressing
00447         if (dressedIdxs.find(barcode) != dressedIdxs.end()) continue;
00448         // Put back to be used in jet clustering
00449         pForJet.push_back(rp);
00450       }
00451 
00452       // Then do the jet clustering
00453       const PseudoJets jet_pjs = mkPseudoJets(pForJet);
00454       const fastjet::JetDefinition jet_jdef(fastjet::antikt_algorithm, _jetR);
00455       const Jets jets_all = mkJets(fastjet::ClusterSequence(jet_pjs, jet_jdef).inclusive_jets());
00456       const Jets jets_sel = sortByPt(filterBy(jets_all, Cuts::pT > _jetMinPt));
00457       // FastJets fjJet(FastJets::ANTIKT, _jetR);
00458       //fjJet.useInvisibles(); // NOTE: CMS proposal to remove neutrinos (AB: wouldn't work anyway, since they were excluded from clustering inputs)
00459       // fjJet.calc(pForJet);
00460       for (const Jet& jet : jets_sel) {
00461         if (jet.abseta() > _jetMaxEta) continue;
00462         _jets.push_back(jet);
00463         bool isBJet = false;
00464         for (const Particle& rp : jet.particles()) {
00465           if (PID::hasBottom(rp.pdgId())) {
00466             isBJet = true;
00467             break;
00468           }
00469         }
00470         if ( isBJet ) _bjets.push_back(jet);
00471         else _ljets.push_back(jet);
00472       }
00473 
00474       // Every building blocks are ready. Continue to pseudo-W and pseudo-top combination
00475 
00476       if (_bjets.size() < 2) return; // Ignore single top for now
00477       map<double, pair<size_t, size_t> > wLepCandIdxs;
00478       map<double, pair<size_t, size_t> > wHadCandIdxs;
00479 
00480       // Collect leptonic-decaying W's
00481       for (size_t iLep = 0, nLep = leptons.size(); iLep < nLep; ++iLep) {
00482         const Jet& lep = leptons.at(iLep);
00483         for (size_t iNu = 0, nNu = neutrinos.size(); iNu < nNu; ++iNu) {
00484           const Particle& nu = neutrinos.at(iNu);
00485           const double m = (lep.momentum()+nu.momentum()).mass();
00486           const double dm = std::abs(m-_wMass);
00487           wLepCandIdxs[dm] = make_pair(iLep, iNu);
00488         }
00489       }
00490 
00491       // Continue to hadronic decaying W's
00492       for (size_t i = 0, nLjet = _ljets.size(); i < nLjet; ++i) {
00493         const Jet& ljet1 = _ljets[i];
00494         for (size_t j = i+1; j < nLjet; ++j) {
00495           const Jet& ljet2 = _ljets[j];
00496           const double m = (ljet1.momentum()+ljet2.momentum()).mass();
00497           const double dm = std::abs(m-_wMass);
00498           wHadCandIdxs[dm] = make_pair(i, j);
00499         }
00500       }
00501 
00502       // Cleanup W candidate, choose pairs with minimum dm if they share decay products
00503       cleanup(wLepCandIdxs);
00504       cleanup(wHadCandIdxs, true);
00505       const size_t nWLepCand = wLepCandIdxs.size();
00506       const size_t nWHadCand = wHadCandIdxs.size();
00507 
00508       if (nWLepCand + nWHadCand < 2) return; // We skip single top
00509 
00510       int w1Q = 1, w2Q = -1;
00511       int w1dau1Id = 1, w2dau1Id = -1;
00512       FourMomentum w1dau1LVec, w1dau2LVec;
00513       FourMomentum w2dau1LVec, w2dau2LVec;
00514       if (nWLepCand == 0) { // Full hadronic case
00515         const pair<size_t, size_t>& idPair1 = wHadCandIdxs.begin()->second;
00516         const pair<size_t, size_t>& idPair2 = (++wHadCandIdxs.begin())->second;  ///< @todo Reinstate std::next
00517         const Jet& w1dau1 = _ljets[idPair1.first];
00518         const Jet& w1dau2 = _ljets[idPair1.second];
00519         const Jet& w2dau1 = _ljets[idPair2.first];
00520         const Jet& w2dau2 = _ljets[idPair2.second];
00521         w1dau1LVec = w1dau1.momentum();
00522         w1dau2LVec = w1dau2.momentum();
00523         w2dau1LVec = w2dau1.momentum();
00524         w2dau2LVec = w2dau2.momentum();
00525       } else if (nWLepCand == 1) { // Semi-leptonic case
00526         const pair<size_t, size_t>& idPair1 = wLepCandIdxs.begin()->second;
00527         const pair<size_t, size_t>& idPair2 = wHadCandIdxs.begin()->second;
00528         const Jet& w1dau1 = leptons[idPair1.first];
00529         const Particle& w1dau2 = neutrinos[idPair1.second];
00530         const Jet& w2dau1 = _ljets[idPair2.first];
00531         const Jet& w2dau2 = _ljets[idPair2.second];
00532         w1dau1LVec = w1dau1.momentum();
00533         w1dau2LVec = w1dau2.momentum();
00534         w2dau1LVec = w2dau1.momentum();
00535         w2dau2LVec = w2dau2.momentum();
00536         w1dau1Id = leptonsId[idPair1.first];
00537         w1Q = w1dau1Id > 0 ? -1 : 1;
00538         w2Q = -w1Q;
00539         switch (w1dau1Id) {
00540         case 13: case -13: _mode1 = CH_MUON; break;
00541         case 11: case -11: _mode1 = CH_ELECTRON; break;
00542         }
00543       } else { // Full leptonic case
00544         const pair<size_t, size_t>& idPair1 = wLepCandIdxs.begin()->second;
00545         const pair<size_t, size_t>& idPair2 = (++wLepCandIdxs.begin())->second;  ///< @todo Reinstate std::next
00546         const Jet& w1dau1 = leptons[idPair1.first];
00547         const Particle& w1dau2 = neutrinos[idPair1.second];
00548         const Jet& w2dau1 = leptons[idPair2.first];
00549         const Particle& w2dau2 = neutrinos[idPair2.second];
00550         w1dau1LVec = w1dau1.momentum();
00551         w1dau2LVec = w1dau2.momentum();
00552         w2dau1LVec = w2dau1.momentum();
00553         w2dau2LVec = w2dau2.momentum();
00554         w1dau1Id = leptonsId[idPair1.first];
00555         w2dau1Id = leptonsId[idPair2.first];
00556         w1Q = w1dau1Id > 0 ? -1 : 1;
00557         w2Q = w2dau1Id > 0 ? -1 : 1;
00558         switch (w1dau1Id) {
00559         case 13: case -13: _mode1 = CH_MUON; break;
00560         case 11: case -11: _mode1 = CH_ELECTRON; break;
00561         }
00562         switch (w2dau1Id) {
00563         case 13: case -13: _mode2 = CH_MUON; break;
00564         case 11: case -11: _mode2 = CH_ELECTRON; break;
00565         }
00566       }
00567       const FourMomentum w1LVec = w1dau1LVec+w1dau2LVec;
00568       const FourMomentum w2LVec = w2dau1LVec+w2dau2LVec;
00569 
00570       // Combine b jets
00571       double sumDm = 1e9;
00572       FourMomentum b1LVec, b2LVec;
00573       for (size_t i = 0, n = _bjets.size(); i < n; ++i) {
00574         const Jet& bjet1 = _bjets[i];
00575         const double mtop1 = (w1LVec+bjet1.momentum()).mass();
00576         const double dmtop1 = std::abs(mtop1-_tMass);
00577         for (size_t j=0; j<n; ++j) {
00578           if (i == j) continue;
00579           const Jet& bjet2 = _bjets[j];
00580           const double mtop2 = (w2LVec+bjet2.momentum()).mass();
00581           const double dmtop2 = std::abs(mtop2-_tMass);
00582 
00583           if (sumDm <= dmtop1+dmtop2) continue;
00584 
00585           sumDm = dmtop1+dmtop2;
00586           b1LVec = bjet1.momentum();
00587           b2LVec = bjet2.momentum();
00588         }
00589       }
00590       if (sumDm >= 1e9) return; // Failed to make top, but this should not happen.
00591 
00592       const FourMomentum t1LVec = w1LVec + b1LVec;
00593       const FourMomentum t2LVec = w2LVec + b2LVec;
00594 
00595       // Put all of them into candidate collection
00596       _t1 = Particle(w1Q*6, t1LVec);
00597       _b1 = Particle(w1Q*5, b1LVec);
00598       _wDecays1.push_back(Particle(w1dau1Id, w1dau1LVec));
00599       _wDecays1.push_back(Particle(-w1dau1Id+w1Q, w1dau2LVec));
00600 
00601       _t2 = Particle(w2Q*6, t2LVec);
00602       _b2 = Particle(w2Q*5, b2LVec);
00603       _wDecays2.push_back(Particle(w2dau1Id, w2dau1LVec));
00604       _wDecays2.push_back(Particle(-w2dau1Id+w2Q, w2dau2LVec));
00605 
00606       _isValid = true;
00607     }
00608 
00609   }
00610 
00611 
00612 }