Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:43

0001 
0002 #include "RecoTauTag/HLTProducers/interface/L2TauPixelTrackMatch.h"
0003 #include "FWCore/Utilities/interface/Exception.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "DataFormats/Common/interface/RefToBase.h"
0006 #include "DataFormats/JetReco/interface/CaloJet.h"
0007 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0008 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
0009 #include "DataFormats/Math/interface/deltaPhi.h"
0010 
0011 L2TauPixelTrackMatch::L2TauPixelTrackMatch(const edm::ParameterSet& conf) {
0012   m_jetSrc = consumes<trigger::TriggerFilterObjectWithRefs>(conf.getParameter<edm::InputTag>("JetSrc"));
0013   m_jetMinPt = conf.getParameter<double>("JetMinPt");
0014   m_jetMaxEta = conf.getParameter<double>("JetMaxEta");
0015   //m_jetMinN        = conf.getParameter<int>("JetMinN");
0016   m_trackSrc = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackSrc"));
0017   m_trackMinPt = conf.getParameter<double>("TrackMinPt");
0018   m_deltaR = conf.getParameter<double>("deltaR");
0019   m_beamSpotTag = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("BeamSpotSrc"));
0020 
0021   produces<reco::CaloJetCollection>();
0022 }
0023 
0024 L2TauPixelTrackMatch::~L2TauPixelTrackMatch() {}
0025 
0026 void L2TauPixelTrackMatch::produce(edm::StreamID, edm::Event& ev, const edm::EventSetup& es) const {
0027   using namespace std;
0028   using namespace reco;
0029 
0030   // *** Pick up beam spot ***
0031 
0032   // use beam spot for vertex x,y
0033   edm::Handle<BeamSpot> bsHandle;
0034   ev.getByToken(m_beamSpotTag, bsHandle);
0035   const reco::BeamSpot& bs = *bsHandle;
0036   math::XYZPoint beam_spot(bs.x0(), bs.y0(), bs.z0());
0037 
0038   // *** Pick up pixel tracks ***
0039 
0040   edm::Handle<TrackCollection> tracksHandle;
0041   ev.getByToken(m_trackSrc, tracksHandle);
0042 
0043   // *** Pick up L2 tau jets that were previously selected by some other filter ***
0044 
0045   // first, get L2 object refs by the label
0046   edm::Handle<trigger::TriggerFilterObjectWithRefs> jetsHandle;
0047   ev.getByToken(m_jetSrc, jetsHandle);
0048 
0049   // now we can get pre-selected L2 tau jets
0050   std::vector<CaloJetRef> tau_jets;
0051   jetsHandle->getObjects(trigger::TriggerTau, tau_jets);
0052   const size_t n_jets = tau_jets.size();
0053 
0054   // *** Selects interesting tracks ***
0055 
0056   vector<TinyTrack> good_tracks;
0057   for (TrackCollection::const_iterator itrk = tracksHandle->begin(); itrk != tracksHandle->end(); ++itrk) {
0058     if (itrk->pt() < m_trackMinPt)
0059       continue;
0060     if (std::abs(itrk->eta()) > m_jetMaxEta + m_deltaR)
0061       continue;
0062 
0063     TinyTrack trk;
0064     trk.pt = itrk->pt();
0065     trk.phi = itrk->phi();
0066     trk.eta = itrk->eta();
0067     double dz = itrk->dz(beam_spot);
0068     trk.vtx = math::XYZPoint(bs.x(dz), bs.y(dz), dz);
0069     good_tracks.push_back(trk);
0070   }
0071 
0072   // *** Match tau jets to intertesting tracks  and assign them new vertices ***
0073 
0074   // the new product
0075   std::unique_ptr<CaloJetCollection> new_tau_jets(new CaloJetCollection);
0076   if (!good_tracks.empty())
0077     for (size_t i = 0; i < n_jets; ++i) {
0078       reco::CaloJetRef jet = tau_jets[i];
0079       if (jet->pt() < m_jetMinPt || std::abs(jet->eta()) > m_jetMaxEta)
0080         continue;
0081 
0082       for (vector<TinyTrack>::const_iterator itrk = good_tracks.begin(); itrk != good_tracks.end(); ++itrk) {
0083         math::XYZTLorentzVector new_jet_dir = Jet::physicsP4(itrk->vtx, *jet, itrk->vtx);
0084         float dphi = reco::deltaPhi(new_jet_dir.phi(), itrk->phi);
0085         float deta = new_jet_dir.eta() - itrk->eta;
0086 
0087         if (dphi * dphi + deta * deta > m_deltaR * m_deltaR)
0088           continue;
0089 
0090         // create a jet copy and assign a new vertex to it
0091         CaloJet new_jet = *jet;
0092         new_jet.setVertex(itrk->vtx);
0093         new_tau_jets->push_back(new_jet);
0094       }
0095       ///if (jet_with_vertices.size()) new_tau_jets->push_back(jet_with_vertices);
0096     }
0097 
0098   // store the result
0099   ev.put(std::move(new_tau_jets));
0100 }