CMS 3D CMS Logo

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