CMS 3D CMS Logo

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