Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:25

0001 #include "RecoJets/JetAlgorithms/interface/CATopJetHelper.h"
0002 
0003 struct GreaterByPtCandPtr {
0004   bool operator()(const edm::Ptr<reco::Candidate>& t1, const edm::Ptr<reco::Candidate>& t2) const {
0005     return t1->pt() > t2->pt();
0006   }
0007 };
0008 
0009 reco::CATopJetProperties CATopJetHelper::operator()(reco::Jet const& ihardJet) const {
0010   reco::CATopJetProperties properties;
0011   // Get subjets
0012   reco::Jet::Constituents subjets = ihardJet.getJetConstituents();
0013   properties.nSubJets = subjets.size();  // number of subjets
0014   properties.topMass = ihardJet.mass();  // jet mass
0015   properties.wMass = 99999.;             // best W mass
0016   properties.minMass = 999999.;          // minimum mass pairing
0017 
0018   // Require at least three subjets in all cases, if not, untagged
0019   if (properties.nSubJets >= 3) {
0020     // Take the highest 3 pt subjets for cuts
0021     sort(subjets.begin(), subjets.end(), GreaterByPtCandPtr());
0022 
0023     // Now look at the subjets that were formed
0024     for (int isub = 0; isub < 2; ++isub) {
0025       // Get this subjet
0026       reco::Jet::Constituent icandJet = subjets[isub];
0027 
0028       // Now look at the "other" subjets than this one, form the minimum invariant mass
0029       // pairing, as well as the "closest" combination to the W mass
0030       for (int jsub = isub + 1; jsub < 3; ++jsub) {
0031         // Get the second subjet
0032         reco::Jet::Constituent jcandJet = subjets[jsub];
0033 
0034         reco::Candidate::LorentzVector wCand = icandJet->p4() + jcandJet->p4();
0035 
0036         // Get the candidate mass
0037         double imw = wCand.mass();
0038 
0039         // Find the combination closest to the W mass
0040         if (fabs(imw - WMass_) < fabs(properties.wMass - WMass_)) {
0041           properties.wMass = imw;
0042         }
0043         // Find the minimum mass pairing.
0044         if (fabs(imw) < properties.minMass) {
0045           properties.minMass = imw;
0046         }
0047       }  // end second loop over subjets
0048     }    // end first loop over subjets
0049   }      // endif 3 subjets
0050 
0051   if (properties.minMass == 999999) {
0052     properties.minMass = -1;
0053   }
0054 
0055   return properties;
0056 }