Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:22:00

0001 #include "RecoTauTag/RecoTau/interface/TauTagTools.h"
0002 
0003 using namespace reco;
0004 
0005 namespace tautagtools {
0006 
0007   TrackRefVector filteredTracks(const TrackRefVector& theInitialTracks,
0008                                 double tkminPt,
0009                                 int tkminPixelHitsn,
0010                                 int tkminTrackerHitsn,
0011                                 double tkmaxipt,
0012                                 double tkmaxChi2,
0013                                 const Vertex& pv) {
0014     TrackRefVector filteredTracks;
0015     for (TrackRefVector::const_iterator iTk = theInitialTracks.begin(); iTk != theInitialTracks.end(); iTk++) {
0016       if ((**iTk).pt() >= tkminPt && (**iTk).normalizedChi2() <= tkmaxChi2 &&
0017           fabs((**iTk).dxy(pv.position())) <= tkmaxipt && (**iTk).numberOfValidHits() >= tkminTrackerHitsn &&
0018           (**iTk).hitPattern().numberOfValidPixelHits() >= tkminPixelHitsn)
0019         filteredTracks.push_back(*iTk);
0020     }
0021     return filteredTracks;
0022   }
0023   TrackRefVector filteredTracks(const TrackRefVector& theInitialTracks,
0024                                 double tkminPt,
0025                                 int tkminPixelHitsn,
0026                                 int tkminTrackerHitsn,
0027                                 double tkmaxipt,
0028                                 double tkmaxChi2,
0029                                 double tktorefpointmaxDZ,
0030                                 const Vertex& pv,
0031                                 double refpoint_Z) {
0032     TrackRefVector filteredTracks;
0033     for (TrackRefVector::const_iterator iTk = theInitialTracks.begin(); iTk != theInitialTracks.end(); iTk++) {
0034       if (pv.isFake())
0035         tktorefpointmaxDZ = 30.;
0036       if ((**iTk).pt() >= tkminPt && (**iTk).normalizedChi2() <= tkmaxChi2 &&
0037           fabs((**iTk).dxy(pv.position())) <= tkmaxipt && (**iTk).numberOfValidHits() >= tkminTrackerHitsn &&
0038           (**iTk).hitPattern().numberOfValidPixelHits() >= tkminPixelHitsn &&
0039           fabs((**iTk).dz(pv.position())) <= tktorefpointmaxDZ)
0040         filteredTracks.push_back(*iTk);
0041     }
0042     return filteredTracks;
0043   }
0044 
0045   std::vector<reco::CandidatePtr> filteredPFChargedHadrCands(const std::vector<reco::CandidatePtr>& theInitialPFCands,
0046                                                              double ChargedHadrCand_tkminPt,
0047                                                              int ChargedHadrCand_tkminPixelHitsn,
0048                                                              int ChargedHadrCand_tkminTrackerHitsn,
0049                                                              double ChargedHadrCand_tkmaxipt,
0050                                                              double ChargedHadrCand_tkmaxChi2,
0051                                                              const Vertex& pv) {
0052     std::vector<reco::CandidatePtr> filteredPFChargedHadrCands;
0053     for (std::vector<reco::CandidatePtr>::const_iterator iPFCand = theInitialPFCands.begin();
0054          iPFCand != theInitialPFCands.end();
0055          iPFCand++) {
0056       if (std::abs((*iPFCand)->pdgId()) == 211 || std::abs((*iPFCand)->pdgId()) == 13 ||
0057           std::abs((*iPFCand)->pdgId()) == 11) {
0058         // *** Whether the charged hadron candidate will be selected or not depends on its rec. tk properties.
0059         const reco::Track* PFChargedHadrCand_rectk = (*iPFCand)->bestTrack();
0060         if (PFChargedHadrCand_rectk != nullptr) {
0061           if ((*PFChargedHadrCand_rectk).pt() >= ChargedHadrCand_tkminPt &&
0062               (*PFChargedHadrCand_rectk).normalizedChi2() <= ChargedHadrCand_tkmaxChi2 &&
0063               fabs((*PFChargedHadrCand_rectk).dxy(pv.position())) <= ChargedHadrCand_tkmaxipt &&
0064               (*PFChargedHadrCand_rectk).numberOfValidHits() >= ChargedHadrCand_tkminTrackerHitsn &&
0065               (*PFChargedHadrCand_rectk).hitPattern().numberOfValidPixelHits() >= ChargedHadrCand_tkminPixelHitsn)
0066             filteredPFChargedHadrCands.push_back(*iPFCand);
0067         }
0068       }
0069     }
0070     return filteredPFChargedHadrCands;
0071   }
0072   std::vector<reco::CandidatePtr> filteredPFChargedHadrCands(const std::vector<reco::CandidatePtr>& theInitialPFCands,
0073                                                              double ChargedHadrCand_tkminPt,
0074                                                              int ChargedHadrCand_tkminPixelHitsn,
0075                                                              int ChargedHadrCand_tkminTrackerHitsn,
0076                                                              double ChargedHadrCand_tkmaxipt,
0077                                                              double ChargedHadrCand_tkmaxChi2,
0078                                                              double ChargedHadrCand_tktorefpointmaxDZ,
0079                                                              const Vertex& pv,
0080                                                              double refpoint_Z) {
0081     if (pv.isFake())
0082       ChargedHadrCand_tktorefpointmaxDZ = 30.;
0083     std::vector<reco::CandidatePtr> filteredPFChargedHadrCands;
0084     for (std::vector<reco::CandidatePtr>::const_iterator iPFCand = theInitialPFCands.begin();
0085          iPFCand != theInitialPFCands.end();
0086          iPFCand++) {
0087       if (std::abs((*iPFCand)->pdgId()) == 211 || std::abs((*iPFCand)->pdgId()) == 13 ||
0088           std::abs((*iPFCand)->pdgId()) == 11) {
0089         // *** Whether the charged hadron candidate will be selected or not depends on its rec. tk properties.
0090         const reco::Track* PFChargedHadrCand_rectk = (*iPFCand)->bestTrack();
0091         if (PFChargedHadrCand_rectk != nullptr) {
0092           if ((*PFChargedHadrCand_rectk).pt() >= ChargedHadrCand_tkminPt &&
0093               (*PFChargedHadrCand_rectk).normalizedChi2() <= ChargedHadrCand_tkmaxChi2 &&
0094               fabs((*PFChargedHadrCand_rectk).dxy(pv.position())) <= ChargedHadrCand_tkmaxipt &&
0095               (*PFChargedHadrCand_rectk).numberOfValidHits() >= ChargedHadrCand_tkminTrackerHitsn &&
0096               (*PFChargedHadrCand_rectk).hitPattern().numberOfValidPixelHits() >= ChargedHadrCand_tkminPixelHitsn &&
0097               fabs((*PFChargedHadrCand_rectk).dz(pv.position())) <= ChargedHadrCand_tktorefpointmaxDZ)
0098             filteredPFChargedHadrCands.push_back(*iPFCand);
0099         }
0100       }
0101     }
0102     return filteredPFChargedHadrCands;
0103   }
0104 
0105   std::vector<reco::CandidatePtr> filteredPFNeutrHadrCands(const std::vector<reco::CandidatePtr>& theInitialPFCands,
0106                                                            double NeutrHadrCand_HcalclusMinEt) {
0107     std::vector<reco::CandidatePtr> filteredPFNeutrHadrCands;
0108     for (std::vector<reco::CandidatePtr>::const_iterator iPFCand = theInitialPFCands.begin();
0109          iPFCand != theInitialPFCands.end();
0110          iPFCand++) {
0111       if (std::abs((*iPFCand)->pdgId()) == 130) {
0112         // *** Whether the neutral hadron candidate will be selected or not depends on its rec. HCAL cluster properties.
0113         if ((**iPFCand).et() >= NeutrHadrCand_HcalclusMinEt) {
0114           filteredPFNeutrHadrCands.push_back(*iPFCand);
0115         }
0116       }
0117     }
0118     return filteredPFNeutrHadrCands;
0119   }
0120 
0121   std::vector<reco::CandidatePtr> filteredPFGammaCands(const std::vector<reco::CandidatePtr>& theInitialPFCands,
0122                                                        double GammaCand_EcalclusMinEt) {
0123     std::vector<reco::CandidatePtr> filteredPFGammaCands;
0124     for (std::vector<reco::CandidatePtr>::const_iterator iPFCand = theInitialPFCands.begin();
0125          iPFCand != theInitialPFCands.end();
0126          iPFCand++) {
0127       if (std::abs((*iPFCand)->pdgId()) == 22) {
0128         // *** Whether the gamma candidate will be selected or not depends on its rec. ECAL cluster properties.
0129         if ((**iPFCand).et() >= GammaCand_EcalclusMinEt) {
0130           filteredPFGammaCands.push_back(*iPFCand);
0131         }
0132       }
0133     }
0134     return filteredPFGammaCands;
0135   }
0136 
0137 }  // namespace tautagtools