CMS 3D CMS Logo

L1FastTrackingJetProducer.cc
Go to the documentation of this file.
1 //
3 // Producer of L1FastTrackingJets
4 //
5 // FastTracking Jets: Jets created by running the FastJet clustering algorithm on L1 tracks that have been matched to a Tracking Particle
6 // Author: G.Karathanasis , CU Boulder
8 
9 // system include files
10 #include <memory>
11 
12 // user include files
25 
26 // L1 objects
31 
32 // geometry
35 
36 //mc
39 
40 #include <fastjet/JetDefinition.hh>
41 
42 #include <string>
43 #include "TMath.h"
44 #include "TH1.h"
45 
46 using namespace l1t;
47 using namespace edm;
48 using namespace std;
49 
51 // //
52 // CLASS DEFINITION //
53 // //
55 
57 public:
59  typedef std::vector<L1TTTrackType> L1TTTrackCollectionType;
60 
62  ~L1FastTrackingJetProducer() override;
63 
64  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
65 
66 private:
67  virtual void beginJob();
68  void produce(edm::Event&, const edm::EventSetup&) override;
69  virtual void endJob();
70 
71  // track selection criteria
72  const float trkZMax_; // in [cm]
73  const float trkChi2dofMax_; // maximum track chi2dof
74  const double trkBendChi2Max_; // maximum track bendchi2
75  const float trkPtMin_; // in [GeV]
76  const float trkEtaMax_; // in [rad]
77  const int trkNStubMin_; // minimum number of stubs
78  const int trkNPSStubMin_; // minimum number of PS stubs
79  const double deltaZ0Cut_; // save with |L1z-z0| < maxZ0
80  const double coneSize_; // Use anti-kt with this cone size
81  const bool doTightChi2_;
82  const float trkPtTightChi2_;
83  const float trkChi2dofTightChi2_;
84  const bool displaced_; //use prompt/displaced tracks
88 
94 };
95 
96 // constructor
98  : trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
99  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
100  trkBendChi2Max_(iConfig.getParameter<double>("trk_bendChi2Max")),
101  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
102  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
103  trkNStubMin_((int)iConfig.getParameter<int>("trk_nStubMin")),
104  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
105  deltaZ0Cut_((float)iConfig.getParameter<double>("deltaZ0Cut")),
106  coneSize_((float)iConfig.getParameter<double>("coneSize")),
107  doTightChi2_(iConfig.getParameter<bool>("doTightChi2")),
108  trkPtTightChi2_((float)iConfig.getParameter<double>("trk_ptTightChi2")),
109  trkChi2dofTightChi2_((float)iConfig.getParameter<double>("trk_chi2dofTightChi2")),
110  displaced_(iConfig.getParameter<bool>("displaced")),
111  selectTrkMatchGenTight_(iConfig.getParameter<bool>("selectTrkMatchGenTight")),
112  selectTrkMatchGenLoose_(iConfig.getParameter<bool>("selectTrkMatchGenLoose")),
113  selectTrkMatchGenOrPU_(iConfig.getParameter<bool>("selectTrkMatchGenOrPU")),
114  trackToken_(consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>>(
115  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
116  pvToken_(consumes<std::vector<l1t::Vertex>>(iConfig.getParameter<edm::InputTag>("L1PrimaryVertexTag"))),
117  genToken_(
118  consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>(iConfig.getParameter<edm::InputTag>("GenInfo"))),
121  if (displaced_)
122  produces<TkJetCollection>("L1FastTrackingJetsExtended");
123  else
124  produces<TkJetCollection>("L1FastTrackingJets");
125 }
126 
127 // destructor
129 
130 // producer
132  std::unique_ptr<TkJetCollection> L1FastTrackingJets(new TkJetCollection);
133 
134  // L1 tracks
136  iEvent.getByToken(trackToken_, TTTrackHandle);
137  std::vector<TTTrack<Ref_Phase2TrackerDigi_>>::const_iterator iterL1Track;
138 
139  // Gen
141  if (!iEvent.isRealData())
142  iEvent.getByToken(genToken_, MCTrkAssociation);
143 
144  // Tracker Topology
145  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
146  //const TrackerGeometry& tGeom = iSetup.getData(tGeomToken_);
147 
148  edm::Handle<std::vector<l1t::Vertex>> L1VertexHandle;
149  iEvent.getByToken(pvToken_, L1VertexHandle);
150 
151  fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
152  std::vector<fastjet::PseudoJet> JetInputs;
153 
154  float recoVtx = L1VertexHandle->begin()->z0();
155  unsigned int this_l1track = 0;
156  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
157  this_l1track++;
158  float trk_pt = iterL1Track->momentum().perp();
159  float trk_z0 = iterL1Track->z0();
160  float trk_chi2dof = iterL1Track->chi2Red();
161  float trk_bendchi2 = iterL1Track->stubPtConsistency();
162  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
163  theStubs = iterL1Track->getStubRefs();
164  int trk_nstub = (int)theStubs.size();
165 
166  if (std::abs(trk_z0) > trkZMax_)
167  continue;
168  if (std::abs(iterL1Track->momentum().eta()) > trkEtaMax_)
169  continue;
170  if (trk_pt < trkPtMin_)
171  continue;
172  if (trk_nstub < trkNStubMin_)
173  continue;
174  if (trk_chi2dof > trkChi2dofMax_)
175  continue;
176  if (trk_bendchi2 > trkBendChi2Max_)
177  continue;
178  if (doTightChi2_ && (trk_pt > trkPtTightChi2_ && trk_chi2dof > trkChi2dofTightChi2_))
179  continue;
180 
181  int trk_nPS = 0;
182  for (int istub = 0; istub < trk_nstub; istub++) {
183  DetId detId(theStubs.at(istub)->getDetId());
184  bool tmp_isPS = false;
185  if (detId.det() == DetId::Detector::Tracker) {
186  if (detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3)
187  tmp_isPS = true;
188  else if (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9)
189  tmp_isPS = true;
190  }
191  if (tmp_isPS)
192  trk_nPS++;
193  }
194  if (trk_nPS < trkNPSStubMin_)
195  continue;
196  if (std::abs(recoVtx - trk_z0) > deltaZ0Cut_)
197  continue;
198  if (!iEvent.isRealData()) {
199  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>> trk_ptr(TTTrackHandle, this_l1track);
200 
201  if (!(MCTrkAssociation->isGenuine(trk_ptr)) && selectTrkMatchGenTight_)
202  continue;
203  if (!(MCTrkAssociation->isLooselyGenuine(trk_ptr) || MCTrkAssociation->isGenuine(trk_ptr)) &&
205  continue;
206  if (!(MCTrkAssociation->isLooselyGenuine(trk_ptr) || MCTrkAssociation->isGenuine(trk_ptr) ||
207  MCTrkAssociation->isCombinatoric(trk_ptr)) &&
209  continue;
210  }
211 
212  fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
213  iterL1Track->momentum().y(),
214  iterL1Track->momentum().z(),
215  iterL1Track->momentum().mag());
216  JetInputs.push_back(psuedoJet); // input tracks for clustering
217  JetInputs.back().set_user_index(this_l1track - 1); // save track index in the collection
218  } // end loop over tracks
219 
220  fastjet::ClusterSequence cs(JetInputs, jet_def); // define the output jet collection
221  std::vector<fastjet::PseudoJet> JetOutputs =
222  fastjet::sorted_by_pt(cs.inclusive_jets(0)); // output jet collection, pT-ordered
223 
224  for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
226  JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
227  float sumpt = 0;
228  float avgZ = 0;
229  std::vector<edm::Ptr<L1TTTrackType>> L1TrackPtrs;
230  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
231 
232  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
233  auto index = fjConstituents[i].user_index();
234  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
235  L1TrackPtrs.push_back(trkPtr); // L1Tracks in the jet
236  sumpt = sumpt + trkPtr->momentum().perp();
237  avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
238  }
239  avgZ = avgZ / sumpt;
241  TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
242  L1FastTrackingJets->push_back(trkJet);
243  } //end loop over Jet Outputs
244 
245  if (displaced_)
246  iEvent.put(std::move(L1FastTrackingJets), "L1FastTrackingJetsExtended");
247  else
248  iEvent.put(std::move(L1FastTrackingJets), "L1FastTrackingJets");
249 }
250 
252 
254 
256  {
257  // L1FastTrackingJets
259  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
260  desc.add<std::string>("L1PrimaryVertexTag", "l1vertices");
261  desc.add<edm::InputTag>("GenInfo", edm::InputTag("TTTrackAssociatorFromPixelDigis", "Level1TTTracks"));
262  desc.add<double>("trk_zMax", 15.0);
263  desc.add<double>("trk_chi2dofMax", 10.0);
264  desc.add<double>("trk_bendChi2Max", 2.2);
265  desc.add<double>("trk_ptMin", 2.0);
266  desc.add<double>("trk_etaMax", 2.5);
267  desc.add<int>("trk_nStubMin", 4);
268  desc.add<int>("trk_nPSStubMin", -1);
269  desc.add<double>("deltaZ0Cut", 0.5);
270  desc.add<bool>("doTightChi2", true);
271  desc.add<double>("trk_ptTightChi2", 20.0);
272  desc.add<double>("trk_chi2dofTightChi2", 5.0);
273  desc.add<double>("coneSize", 0.4);
274  desc.add<bool>("displaced", false);
275  desc.add<bool>("selectTrkMatchGenTight", true);
276  desc.add<bool>("selectTrkMatchGenLoose", false);
277  desc.add<bool>("selectTrkMatchGenOrPU", false);
278  descriptions.add("L1FastTrackingJets", desc);
279  // or use the following to generate the label from the module's C++ type
280  //descriptions.addWithDefaultLabel(desc);
281  }
282  /*{
283  // L1FastTrackingJetsExtended
284  desc.add<double>("trk_bendChi2Max", 2.4);
285  desc.add<double>("trk_ptMin", 3.0);
286  desc.add<double>("trk_etaMax", 2.5);
287  desc.add<int>("trk_nStubMin", 4);
288  desc.add<int>("trk_nPSStubMin", -1);
289  desc.add<double>("deltaZ0Cut", 3.0);
290  desc.add<bool>("doTightChi2", true);
291  desc.add<double>("trk_ptTightChi2", 20.0);
292  desc.add<double>("trk_chi2dofTightChi2", 5.0);
293  desc.add<double>("coneSize", 0.4);
294  desc.add<bool>("displaced", true);
295  desc.add<bool>("selectTrkMatchGenTight", true);
296  desc.add<bool>("selectTrkMatchGenLoose", false);
297  desc.add<bool>("selectTrkMatchGenOrPU", false);
298  descriptions.add("L1FastTrackingJetsExtended", desc);
299  // or use the following to generate the label from the module's C++ type
300  //descriptions.addWithDefaultLabel(desc);
301  }*/
302 }
303 
304 //define this as a plug-in
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
unsigned int tobLayer(const DetId &id) const
T perp() const
Definition: PV3DBase.h:69
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const edm::EDGetTokenT< TTTrackAssociationMap< Ref_Phase2TrackerDigi_ > > genToken_
L1FastTrackingJetProducer(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< std::vector< l1t::Vertex > > pvToken_
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
delete x;
Definition: CaloConfig.h:22
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
void beginJob()
Definition: Breakpoints.cc:14
Stores association of Truth Particles (TP) to L1 Track-Trigger Tracks.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
std::vector< L1TTTrackType > L1TTTrackCollectionType
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
void produce(edm::Event &, const edm::EventSetup &) override
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Definition: DetId.h:17
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
double z0() const
Track z0.
Definition: TTTrack.h:330
unsigned int tidRing(const DetId &id) const
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511