CMS 3D CMS Logo

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