File indexing completed on 2024-04-06 12:27:54
0001 #include "RecoTauTag/RecoTau/interface/PFRecoTauTagInfoAlgorithm.h"
0002
0003 #include "DataFormats/Math/interface/deltaR.h"
0004 #include "RecoTauTag/RecoTau/interface/TauTagTools.h"
0005
0006 using namespace reco;
0007
0008 PFRecoTauTagInfoAlgorithm::PFRecoTauTagInfoAlgorithm(const edm::ParameterSet& parameters) {
0009
0010 ChargedHadronsAssociationCone_ = parameters.getParameter<double>("ChargedHadrCand_AssociationCone");
0011 ChargedHadrCand_tkminPt_ = parameters.getParameter<double>("ChargedHadrCand_tkminPt");
0012 ChargedHadrCand_tkminPixelHitsn_ = parameters.getParameter<int>("ChargedHadrCand_tkminPixelHitsn");
0013 ChargedHadrCand_tkminTrackerHitsn_ = parameters.getParameter<int>("ChargedHadrCand_tkminTrackerHitsn");
0014 ChargedHadrCand_tkmaxipt_ = parameters.getParameter<double>("ChargedHadrCand_tkmaxipt");
0015 ChargedHadrCand_tkmaxChi2_ = parameters.getParameter<double>("ChargedHadrCand_tkmaxChi2");
0016
0017 NeutrHadrCand_HcalclusMinEt_ = parameters.getParameter<double>("NeutrHadrCand_HcalclusMinEt");
0018
0019 GammaCand_EcalclusMinEt_ = parameters.getParameter<double>("GammaCand_EcalclusMinEt");
0020
0021 tkminPt_ = parameters.getParameter<double>("tkminPt");
0022 tkminPixelHitsn_ = parameters.getParameter<int>("tkminPixelHitsn");
0023 tkminTrackerHitsn_ = parameters.getParameter<int>("tkminTrackerHitsn");
0024 tkmaxipt_ = parameters.getParameter<double>("tkmaxipt");
0025 tkmaxChi2_ = parameters.getParameter<double>("tkmaxChi2");
0026
0027 UsePVconstraint_ = parameters.getParameter<bool>("UsePVconstraint");
0028 ChargedHadrCand_tkPVmaxDZ_ = parameters.getParameter<double>("ChargedHadrCand_tkPVmaxDZ");
0029 tkPVmaxDZ_ = parameters.getParameter<double>("tkPVmaxDZ");
0030 }
0031
0032 PFTauTagInfo PFRecoTauTagInfoAlgorithm::buildPFTauTagInfo(const JetBaseRef& thePFJet,
0033 const std::vector<reco::CandidatePtr>& thePFCandsInEvent,
0034 const TrackRefVector& theTracks,
0035 const Vertex& thePV) const {
0036 PFTauTagInfo resultExtended;
0037 resultExtended.setpfjetRef(thePFJet);
0038
0039 std::vector<reco::CandidatePtr> thePFCands;
0040 const float jetPhi = (*thePFJet).phi();
0041 const float jetEta = (*thePFJet).eta();
0042 auto dr2 = [jetPhi, jetEta](float phi, float eta) { return reco::deltaR2(jetEta, jetPhi, eta, phi); };
0043 for (const auto& iPFCand : thePFCandsInEvent) {
0044 float delta = dr2((*iPFCand).phi(), (*iPFCand).eta());
0045 if (delta < ChargedHadronsAssociationCone_ * ChargedHadronsAssociationCone_)
0046 thePFCands.push_back(iPFCand);
0047 }
0048 bool pvIsFake = (thePV.z() < -500.);
0049
0050 std::vector<reco::CandidatePtr> theFilteredPFChargedHadrCands;
0051 if (UsePVconstraint_ && !pvIsFake) {
0052 theFilteredPFChargedHadrCands = tautagtools::filteredPFChargedHadrCands(thePFCands,
0053 ChargedHadrCand_tkminPt_,
0054 ChargedHadrCand_tkminPixelHitsn_,
0055 ChargedHadrCand_tkminTrackerHitsn_,
0056 ChargedHadrCand_tkmaxipt_,
0057 ChargedHadrCand_tkmaxChi2_,
0058 ChargedHadrCand_tkPVmaxDZ_,
0059 thePV,
0060 thePV.z());
0061 } else {
0062 theFilteredPFChargedHadrCands = tautagtools::filteredPFChargedHadrCands(thePFCands,
0063 ChargedHadrCand_tkminPt_,
0064 ChargedHadrCand_tkminPixelHitsn_,
0065 ChargedHadrCand_tkminTrackerHitsn_,
0066 ChargedHadrCand_tkmaxipt_,
0067 ChargedHadrCand_tkmaxChi2_,
0068 thePV);
0069 }
0070 resultExtended.setPFChargedHadrCands(theFilteredPFChargedHadrCands);
0071 resultExtended.setPFNeutrHadrCands(tautagtools::filteredPFNeutrHadrCands(thePFCands, NeutrHadrCand_HcalclusMinEt_));
0072 resultExtended.setPFGammaCands(tautagtools::filteredPFGammaCands(thePFCands, GammaCand_EcalclusMinEt_));
0073
0074 TrackRefVector theFilteredTracks;
0075 if (UsePVconstraint_ && !pvIsFake) {
0076 theFilteredTracks = tautagtools::filteredTracks(
0077 theTracks, tkminPt_, tkminPixelHitsn_, tkminTrackerHitsn_, tkmaxipt_, tkmaxChi2_, tkPVmaxDZ_, thePV, thePV.z());
0078 } else {
0079 theFilteredTracks = tautagtools::filteredTracks(
0080 theTracks, tkminPt_, tkminPixelHitsn_, tkminTrackerHitsn_, tkmaxipt_, tkmaxChi2_, thePV);
0081 }
0082 resultExtended.setTracks(theFilteredTracks);
0083
0084 return resultExtended;
0085 }