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