File indexing completed on 2024-04-06 12:27:43
0001 #include "RecoTauTag/HLTProducers/interface/DQMTauProducer.h"
0002 #include "DataFormats/TrackReco/interface/Track.h"
0003
0004
0005
0006
0007 DQMTauProducer::DQMTauProducer(const edm::ParameterSet& iConfig) {
0008 trackIsolatedJets_ =
0009 consumes<reco::IsolatedTauTagInfoCollection>(iConfig.getParameter<edm::InputTag>("TrackIsoJets"));
0010 matchingCone_ = iConfig.getParameter<double>("MatchingCone");
0011 signalCone_ = iConfig.getParameter<double>("SignalCone");
0012 ptMin_ = iConfig.getParameter<double>("MinPtTracks");
0013
0014 isolationCone_ = iConfig.getParameter<double>("IsolationCone");
0015 produces<reco::HLTTauCollection>();
0016 }
0017
0018 DQMTauProducer::~DQMTauProducer() {}
0019
0020 void DQMTauProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iES) const {
0021 using namespace reco;
0022 using namespace edm;
0023 using namespace std;
0024
0025 auto jetCollection = std::make_unique<HLTTauCollection>();
0026
0027 edm::Handle<IsolatedTauTagInfoCollection> tauL25Jets;
0028 iEvent.getByToken(trackIsolatedJets_, tauL25Jets);
0029
0030 IsolatedTauTagInfoCollection tau = *(tauL25Jets.product());
0031
0032 float eta_, phi_, pt_;
0033 float ptLeadTk = 0.;
0034 int trackIsolation = 1000.;
0035 int nTracks = 1000.;
0036
0037 for (unsigned int i = 0; i < tau.size(); i++) {
0038 JetTracksAssociationRef jetTracks = tau[i].jtaRef();
0039 math::XYZVector jetDir(jetTracks->first->px(), jetTracks->first->py(), jetTracks->first->pz());
0040 eta_ = jetDir.eta();
0041 phi_ = jetDir.phi();
0042 pt_ = jetTracks->first->pt();
0043
0044 const TrackRef leadTk = tau[i].leadingSignalTrack(jetDir, matchingCone_, 1.);
0045 if (!leadTk) {
0046 } else {
0047 trackIsolation = (int)tau[i].discriminator(jetDir, matchingCone_, signalCone_, isolationCone_, 1., 1., 0);
0048 ptLeadTk = (*leadTk).pt();
0049 nTracks = (tau[i].tracksInCone((*leadTk).momentum(), isolationCone_, ptMin_)).size() -
0050 (tau[i].tracksInCone((*leadTk).momentum(), signalCone_, ptMin_)).size();
0051 }
0052 HLTTau pippo(eta_, phi_, pt_, -1, trackIsolation, ptLeadTk, trackIsolation, ptLeadTk);
0053 pippo.setNL25TrackIsolation(nTracks);
0054 pippo.setNL3TrackIsolation(nTracks);
0055 jetCollection->push_back(pippo);
0056 }
0057
0058 iEvent.put(std::move(jetCollection));
0059 }