CMS 3D CMS Logo

L2TauPixelTrackMatch.cc
Go to the documentation of this file.
1 
10 
11 
12 // all these debug printouts will need to be removed at some point
13 //#define DBG_PRINT(arg) (arg)
14 #define DBG_PRINT(arg)
15 
16 
18 {
19  m_jetSrc = consumes<trigger::TriggerFilterObjectWithRefs>(conf.getParameter<edm::InputTag>("JetSrc"));
20  m_jetMinPt = conf.getParameter<double>("JetMinPt");
21  m_jetMaxEta = conf.getParameter<double>("JetMaxEta");
22  //m_jetMinN = conf.getParameter<int>("JetMinN");
23  m_trackSrc = consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("TrackSrc"));
24  m_trackMinPt = conf.getParameter<double>("TrackMinPt");
25  m_deltaR = conf.getParameter<double>("deltaR");
26  m_beamSpotTag = consumes<reco::BeamSpot>(conf.getParameter<edm::InputTag>("BeamSpotSrc"));
27 
28  produces<reco::CaloJetCollection>();
29 }
30 
31 
33 
34 
36 {
37  using namespace std;
38  using namespace reco;
39 
40  // *** Pick up beam spot ***
41 
42  // use beam spot for vertex x,y
43  edm::Handle<BeamSpot> bsHandle;
44  ev.getByToken( m_beamSpotTag, bsHandle);
45  const reco::BeamSpot & bs = *bsHandle;
46  math::XYZPoint beam_spot(bs.x0(), bs.y0(), bs.z0());
47  DBG_PRINT(cout<<endl<<"beamspot "<<beam_spot<<endl);
48 
49  // *** Pick up pixel tracks ***
50 
51  edm::Handle<TrackCollection> tracksHandle;
52  ev.getByToken(m_trackSrc,tracksHandle);
53 
54  // *** Pick up L2 tau jets that were previously selected by some other filter ***
55 
56  // first, get L2 object refs by the label
58  ev.getByToken(m_jetSrc, jetsHandle);
59 
60  // now we can get pre-selected L2 tau jets
61  std::vector<CaloJetRef> tau_jets;
62  jetsHandle->getObjects(trigger::TriggerTau, tau_jets);
63  const size_t n_jets = tau_jets.size();
64 
65  // *** Selects interesting tracks ***
66 
67  vector<TinyTrack> good_tracks;
68  for(TrackCollection::const_iterator itrk = tracksHandle->begin(); itrk != tracksHandle->end(); ++itrk)
69  {
70  if (itrk->pt() < m_trackMinPt) continue;
71  if ( std::abs(itrk->eta()) > m_jetMaxEta + m_deltaR ) continue;
72 
73  TinyTrack trk;
74  trk.pt = itrk->pt();
75  trk.phi = itrk->phi();
76  trk.eta = itrk->eta();
77  double dz = itrk->dz(beam_spot);
78  trk.vtx = math::XYZPoint(bs.x(dz), bs.y(dz), dz);
79  good_tracks.push_back(trk);
80  }
81  DBG_PRINT(cout<<"got "<<good_tracks.size()<<" good tracks; "<<n_jets<<" tau jets"<<endl);
82 
83 
84  // *** Match tau jets to intertesting tracks and assign them new vertices ***
85 
86  // the new product
87  std::unique_ptr<CaloJetCollection> new_tau_jets(new CaloJetCollection);
88  int n_uniq = 0;
89  if (!good_tracks.empty()) for (size_t i=0; i < n_jets; ++i)
90  {
91  reco::CaloJetRef jet = tau_jets[i];
92  if ( jet->pt() < m_jetMinPt || std::abs(jet->eta()) > m_jetMaxEta ) continue;
93 
94  DBG_PRINT(cout<<i<<" :"<<endl);
95 
96  size_t n0 = new_tau_jets->size();
97 
98  for(vector<TinyTrack>::const_iterator itrk = good_tracks.begin(); itrk != good_tracks.end(); ++itrk)
99  {
100  DBG_PRINT(cout<<" trk pt,eta,phi,z: "<<itrk->pt<<" "<<itrk->eta<<" "<<itrk->phi<<" "<<itrk->vtx.z()<<" \t\t ");
101 
102  math::XYZTLorentzVector new_jet_dir = Jet::physicsP4(itrk->vtx, *jet, itrk->vtx);
103  float dphi = reco::deltaPhi(new_jet_dir.phi(), itrk->phi);
104  float deta = new_jet_dir.eta() - itrk->eta;
105 
106  DBG_PRINT(cout<<" jet pt,deta,dphi,dr: "<<jet->pt()<<" "<<deta<<" "<<dphi<<" "<<sqrt(dphi*dphi + deta*deta)<<endl);
107 
108  if ( dphi*dphi + deta*deta > m_deltaR*m_deltaR ) continue;
109 
110  DBG_PRINT(cout<<" jet-trk match!"<<endl);
111 
112  // create a jet copy and assign a new vertex to it
113  CaloJet new_jet = *jet;
114  new_jet.setVertex(itrk->vtx);
115 
116  new_tau_jets->push_back(new_jet);
117  }
118  DBG_PRINT(cout<<" nmatchedjets "<<new_tau_jets->size() - n0<<endl);
119  if (new_tau_jets->size() - n0 > 0 ) n_uniq++;
120 
122  }
123  DBG_PRINT(cout<<"n_uniq_matched_jets "<<n_uniq<<endl<<"storing njets "<<new_tau_jets->size()<<endl);
124 
125  // store the result
126  ev.put(std::move(new_tau_jets));
127 }
L2TauPixelTrackMatch(const edm::ParameterSet &)
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Jets made from CaloTowers.
Definition: CaloJet.h:29
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
bool ev
#define DBG_PRINT(arg)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_jetSrc
void setVertex(const Point &vertex) override
set vertex
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int n0
Definition: AMPTWrapper.h:34
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< reco::TrackCollection > m_trackSrc
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double x(const double z) const
x coordinate of the beeam spot position at a given z value (it takes into account the dxdz slope) ...
Definition: BeamSpot.h:71
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< reco::BeamSpot > m_beamSpotTag
fixed size matrix
double y(const double z) const
y coordinate of the beeam spot position at a given z value (it takes into account the dydz slope) ...
Definition: BeamSpot.h:73
double y0() const
y coordinate
Definition: BeamSpot.h:66
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:510
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
double x0() const
x coordinate
Definition: BeamSpot.h:64