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
26 
27 // L1 objects
32 
33 // geometry
37 
38 //mc
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  float trkZMax_; // in [cm]
73  float trkChi2dofMax_; // maximum track chi2dof
74  double trkBendChi2Max_; // maximum track bendchi2
75  float trkPtMin_; // in [GeV]
76  float trkEtaMax_; // in [rad]
77  int trkNStubMin_; // minimum number of stubs
78  int trkNPSStubMin_; // minimum number of PS stubs
79  double deltaZ0Cut_; // save with |L1z-z0| < maxZ0
80  double coneSize_; // Use anti-kt with this cone size
84  bool displaced_; //use prompt/displaced tracks
88 
93 };
94 
95 // constructor
97  : trackToken_(consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>>(
98  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
99  pvToken_(consumes<std::vector<l1t::Vertex>>(iConfig.getParameter<edm::InputTag>("L1PrimaryVertexTag"))),
100  genToken_(
101  consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>>(iConfig.getParameter<edm::InputTag>("GenInfo"))),
102  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
103  trkZMax_ = (float)iConfig.getParameter<double>("trk_zMax");
104  trkChi2dofMax_ = (float)iConfig.getParameter<double>("trk_chi2dofMax");
105  trkBendChi2Max_ = iConfig.getParameter<double>("trk_bendChi2Max");
106  trkPtMin_ = (float)iConfig.getParameter<double>("trk_ptMin");
107  trkEtaMax_ = (float)iConfig.getParameter<double>("trk_etaMax");
108  trkNStubMin_ = (int)iConfig.getParameter<int>("trk_nStubMin");
109  trkNPSStubMin_ = (int)iConfig.getParameter<int>("trk_nPSStubMin");
110  deltaZ0Cut_ = (float)iConfig.getParameter<double>("deltaZ0Cut");
111  coneSize_ = (float)iConfig.getParameter<double>("coneSize");
112  doTightChi2_ = iConfig.getParameter<bool>("doTightChi2");
113  trkPtTightChi2_ = (float)iConfig.getParameter<double>("trk_ptTightChi2");
114  trkChi2dofTightChi2_ = (float)iConfig.getParameter<double>("trk_chi2dofTightChi2");
115  displaced_ = iConfig.getParameter<bool>("displaced");
116  selectTrkMatchGenTight_ = iConfig.getParameter<bool>("selectTrkMatchGenTight");
117  selectTrkMatchGenLoose_ = iConfig.getParameter<bool>("selectTrkMatchGenLoose");
118  selectTrkMatchGenOrPU_ = iConfig.getParameter<bool>("selectTrkMatchGenOrPU");
119  if (displaced_)
120  produces<TkJetCollection>("L1FastTrackingJetsExtended");
121  else
122  produces<TkJetCollection>("L1FastTrackingJets");
123 }
124 
125 // destructor
127 
128 // producer
130  std::unique_ptr<TkJetCollection> L1FastTrackingJets(new TkJetCollection);
131 
132  // L1 tracks
134  iEvent.getByToken(trackToken_, TTTrackHandle);
135  std::vector<TTTrack<Ref_Phase2TrackerDigi_>>::const_iterator iterL1Track;
136 
137  // Gen
139  if (!iEvent.isRealData())
140  iEvent.getByToken(genToken_, MCTrkAssociation);
141 
142  // Tracker Topology
143  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
144 
145  edm::Handle<std::vector<l1t::Vertex>> L1VertexHandle;
146  iEvent.getByToken(pvToken_, L1VertexHandle);
147 
148  fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
149  std::vector<fastjet::PseudoJet> JetInputs;
150 
151  float recoVtx = L1VertexHandle->begin()->z0();
152  unsigned int this_l1track = 0;
153  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
154  this_l1track++;
155  float trk_pt = iterL1Track->momentum().perp();
156  float trk_z0 = iterL1Track->z0();
157  float trk_chi2dof = iterL1Track->chi2Red();
158  float trk_bendchi2 = iterL1Track->stubPtConsistency();
159  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
160  theStubs = iterL1Track->getStubRefs();
161  int trk_nstub = (int)theStubs.size();
162 
163  if (std::abs(trk_z0) > trkZMax_)
164  continue;
165  if (std::abs(iterL1Track->momentum().eta()) > trkEtaMax_)
166  continue;
167  if (trk_pt < trkPtMin_)
168  continue;
169  if (trk_nstub < trkNStubMin_)
170  continue;
171  if (trk_chi2dof > trkChi2dofMax_)
172  continue;
173  if (trk_bendchi2 > trkBendChi2Max_)
174  continue;
175  if (doTightChi2_ && (trk_pt > trkPtTightChi2_ && trk_chi2dof > trkChi2dofTightChi2_))
176  continue;
177 
178  int trk_nPS = 0;
179  for (int istub = 0; istub < trk_nstub; istub++) {
180  DetId detId(theStubs.at(istub)->getDetId());
181  bool tmp_isPS = false;
182  if (detId.det() == DetId::Detector::Tracker) {
183  if (detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3)
184  tmp_isPS = true;
185  else if (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9)
186  tmp_isPS = true;
187  }
188  if (tmp_isPS)
189  trk_nPS++;
190  }
191  if (trk_nPS < trkNPSStubMin_)
192  continue;
193  if (std::abs(recoVtx - trk_z0) > deltaZ0Cut_)
194  continue;
195  if (!iEvent.isRealData()) {
196  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_>> trk_ptr(TTTrackHandle, this_l1track);
197 
198  if (!(MCTrkAssociation->isGenuine(trk_ptr)) && selectTrkMatchGenTight_)
199  continue;
200  if (!(MCTrkAssociation->isLooselyGenuine(trk_ptr) || MCTrkAssociation->isGenuine(trk_ptr)) &&
202  continue;
203  if (!(MCTrkAssociation->isLooselyGenuine(trk_ptr) || MCTrkAssociation->isGenuine(trk_ptr) ||
204  MCTrkAssociation->isCombinatoric(trk_ptr)) &&
206  continue;
207  }
208 
209  fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
210  iterL1Track->momentum().y(),
211  iterL1Track->momentum().z(),
212  iterL1Track->momentum().mag());
213  JetInputs.push_back(psuedoJet); // input tracks for clustering
214  JetInputs.back().set_user_index(this_l1track - 1); // save track index in the collection
215  } // end loop over tracks
216 
217  fastjet::ClusterSequence cs(JetInputs, jet_def); // define the output jet collection
218  std::vector<fastjet::PseudoJet> JetOutputs =
219  fastjet::sorted_by_pt(cs.inclusive_jets(0)); // output jet collection, pT-ordered
220 
221  for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
223  JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
224  float sumpt = 0;
225  float avgZ = 0;
226  std::vector<edm::Ptr<L1TTTrackType>> L1TrackPtrs;
227  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
228 
229  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
230  auto index = fjConstituents[i].user_index();
231  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
232  L1TrackPtrs.push_back(trkPtr); // L1Tracks in the jet
233  sumpt = sumpt + trkPtr->momentum().perp();
234  avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
235  }
236  avgZ = avgZ / sumpt;
238  TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
239  L1FastTrackingJets->push_back(trkJet);
240  } //end loop over Jet Outputs
241 
242  if (displaced_)
243  iEvent.put(std::move(L1FastTrackingJets), "L1FastTrackingJetsExtended");
244  else
245  iEvent.put(std::move(L1FastTrackingJets), "L1FastTrackingJets");
246 }
247 
249 
251 
253  {
254  // L1FastTrackingJets
256  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
257  desc.add<std::string>("L1PrimaryVertexTag", "l1vertices");
258  desc.add<edm::InputTag>("GenInfo", edm::InputTag("TTTrackAssociatorFromPixelDigis", "Level1TTTracks"));
259  desc.add<double>("trk_zMax", 15.0);
260  desc.add<double>("trk_chi2dofMax", 10.0);
261  desc.add<double>("trk_bendChi2Max", 2.2);
262  desc.add<double>("trk_ptMin", 2.0);
263  desc.add<double>("trk_etaMax", 2.5);
264  desc.add<int>("trk_nStubMin", 4);
265  desc.add<int>("trk_nPSStubMin", -1);
266  desc.add<double>("deltaZ0Cut", 0.5);
267  desc.add<bool>("doTightChi2", true);
268  desc.add<double>("trk_ptTightChi2", 20.0);
269  desc.add<double>("trk_chi2dofTightChi2", 5.0);
270  desc.add<double>("coneSize", 0.4);
271  desc.add<bool>("displaced", false);
272  desc.add<bool>("selectTrkMatchGenTight", true);
273  desc.add<bool>("selectTrkMatchGenLoose", false);
274  desc.add<bool>("selectTrkMatchGenOrPU", false);
275  descriptions.add("L1FastTrackingJets", desc);
276  // or use the following to generate the label from the module's C++ type
277  //descriptions.addWithDefaultLabel(desc);
278  }
279  /*{
280  // L1FastTrackingJetsExtended
281  desc.add<double>("trk_bendChi2Max", 2.4);
282  desc.add<double>("trk_ptMin", 3.0);
283  desc.add<double>("trk_etaMax", 2.5);
284  desc.add<int>("trk_nStubMin", 4);
285  desc.add<int>("trk_nPSStubMin", -1);
286  desc.add<double>("deltaZ0Cut", 3.0);
287  desc.add<bool>("doTightChi2", true);
288  desc.add<double>("trk_ptTightChi2", 20.0);
289  desc.add<double>("trk_chi2dofTightChi2", 5.0);
290  desc.add<double>("coneSize", 0.4);
291  desc.add<bool>("displaced", true);
292  desc.add<bool>("selectTrkMatchGenTight", true);
293  desc.add<bool>("selectTrkMatchGenLoose", false);
294  desc.add<bool>("selectTrkMatchGenOrPU", false);
295  descriptions.add("L1FastTrackingJetsExtended", desc);
296  // or use the following to generate the label from the module's C++ type
297  //descriptions.addWithDefaultLabel(desc);
298  }*/
299 }
300 
301 //define this as a plug-in
L1FastTrackingJetProducer::trkZMax_
float trkZMax_
Definition: L1FastTrackingJetProducer.cc:72
bk::beginJob
void beginJob()
Definition: Breakpoints.cc:14
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
edm::ESInputTag
Definition: ESInputTag.h:87
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
TrackerGeometry.h
TTTrack::momentum
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
L1FastTrackingJetProducer::doTightChi2_
bool doTightChi2_
Definition: L1FastTrackingJetProducer.cc:81
ESHandle.h
L1FastTrackingJetProducer::trkPtMin_
float trkPtMin_
Definition: L1FastTrackingJetProducer.cc:75
TTTypes.h
fwrapper::cs
unique_ptr< ClusterSequence > cs
Definition: fastjetfortran_madfks.cc:47
L1FastTrackingJetProducer::genToken_
const edm::EDGetTokenT< TTTrackAssociationMap< Ref_Phase2TrackerDigi_ > > genToken_
Definition: L1FastTrackingJetProducer.cc:91
l1t::Vertex
Definition: Vertex.h:15
edm::EDGetTokenT
Definition: EDGetToken.h:33
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackerTopology
Definition: TrackerTopology.h:16
L1FastTrackingJetProducer::L1TTTrackCollectionType
std::vector< L1TTTrackType > L1TTTrackCollectionType
Definition: L1FastTrackingJetProducer.cc:59
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
L1FastTrackingJetProducer_cfi.L1FastTrackingJets
L1FastTrackingJets
Definition: L1FastTrackingJetProducer_cfi.py:3
L1FastTrackingJetProducer::selectTrkMatchGenOrPU_
bool selectTrkMatchGenOrPU_
Definition: L1FastTrackingJetProducer.cc:87
TTTrack
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
EDProducer.h
L1FastTrackingJetProducer::displaced_
bool displaced_
Definition: L1FastTrackingJetProducer.cc:84
L1FastTrackingJetProducer::trackToken_
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Definition: L1FastTrackingJetProducer.cc:89
edm::Handle
Definition: AssociativeIterator.h:50
L1FastTrackingJetProducer::coneSize_
double coneSize_
Definition: L1FastTrackingJetProducer.cc:80
TrackerTopology::tidRing
unsigned int tidRing(const DetId &id) const
Definition: TrackerTopology.h:218
edm::Ref
Definition: AssociativeIterator.h:58
L1FastTrackingJetProducer::deltaZ0Cut_
double deltaZ0Cut_
Definition: L1FastTrackingJetProducer.cc:79
align::Tracker
Definition: StructureType.h:70
DetId
Definition: DetId.h:17
MakerMacros.h
TrackerTopology.h
L1FastTrackingJetProducer::beginJob
virtual void beginJob()
Definition: L1FastTrackingJetProducer.cc:248
TrackerTopologyRcd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
L1FastTrackingJetProducer::trkBendChi2Max_
double trkBendChi2Max_
Definition: L1FastTrackingJetProducer.cc:74
l1t::TkJetCollection
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TTStub
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
VirtualJetProducer.h
L1FastTrackingJetProducer::~L1FastTrackingJetProducer
~L1FastTrackingJetProducer() override
Definition: L1FastTrackingJetProducer.cc:126
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
TTTrackAssociationMap
Stores association of Truth Particles (TP) to L1 Track-Trigger Tracks.
Definition: TTTrackAssociationMap.h:46
L1FastTrackingJetProducer
Definition: L1FastTrackingJetProducer.cc:56
TrackerDigiGeometryRecord.h
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
L1FastTrackingJetProducer::selectTrkMatchGenTight_
bool selectTrkMatchGenTight_
Definition: L1FastTrackingJetProducer.cc:85
L1FastTrackingJetProducer::trkChi2dofTightChi2_
float trkChi2dofTightChi2_
Definition: L1FastTrackingJetProducer.cc:83
l1t
delete x;
Definition: CaloConfig.h:22
TkJet.h
TTTrack::z0
double z0() const
Track z0.
Definition: TTTrack.h:330
LorentzVector.h
fwrapper::jet_def
JetDefinition jet_def
Definition: fastjetfortran_madfks.cc:46
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
TkJetFwd.h
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
L1FastTrackingJetProducer::endJob
virtual void endJob()
Definition: L1FastTrackingJetProducer.cc:250
edm::stream::EDProducer
Definition: EDProducer.h:36
edm::EventSetup
Definition: EventSetup.h:58
L1FastTrackingJetProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1FastTrackingJetProducer.cc:252
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd >
L1FastTrackingJetProducer::trkNStubMin_
int trkNStubMin_
Definition: L1FastTrackingJetProducer.cc:77
InputTag.h
edm::Ptr
Definition: AssociationVector.h:31
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
L1FastTrackingJetProducer::L1TTTrackType
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
Definition: L1FastTrackingJetProducer.cc:58
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
L1FastTrackingJetProducer::tTopoToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
Definition: L1FastTrackingJetProducer.cc:92
TrackingParticle.h
TrackerTopology::tobLayer
unsigned int tobLayer(const DetId &id) const
Definition: TrackerTopology.h:147
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
l1t::TkJet
Definition: TkJet.h:18
TTTrackAssociationMap.h
Frameworkfwd.h
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
L1FastTrackingJetProducer::selectTrkMatchGenLoose_
bool selectTrkMatchGenLoose_
Definition: L1FastTrackingJetProducer.cc:86
L1FastTrackingJetProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: L1FastTrackingJetProducer.cc:129
PixelGeomDetUnit.h
L1FastTrackingJetProducer::L1FastTrackingJetProducer
L1FastTrackingJetProducer(const edm::ParameterSet &)
Definition: L1FastTrackingJetProducer.cc:96
StripSubdetector::TOB
static constexpr auto TOB
Definition: StripSubdetector.h:18
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
L1FastTrackingJetProducer::pvToken_
edm::EDGetTokenT< std::vector< l1t::Vertex > > pvToken_
Definition: L1FastTrackingJetProducer.cc:90
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
L1FastTrackingJetProducer::trkNPSStubMin_
int trkNPSStubMin_
Definition: L1FastTrackingJetProducer.cc:78
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrackerTopologyRcd
Definition: TrackerTopologyRcd.h:10
ParameterSet.h
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::Event
Definition: Event.h:73
Vertex.h
L1FastTrackingJetProducer::trkEtaMax_
float trkEtaMax_
Definition: L1FastTrackingJetProducer.cc:76
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
edm::InputTag
Definition: InputTag.h:15
L1FastTrackingJetProducer::trkChi2dofMax_
float trkChi2dofMax_
Definition: L1FastTrackingJetProducer.cc:73
StripSubdetector::TID
static constexpr auto TID
Definition: StripSubdetector.h:17
L1FastTrackingJetProducer::trkPtTightChi2_
float trkPtTightChi2_
Definition: L1FastTrackingJetProducer.cc:82