CMS 3D CMS Logo

L1TrackFastJetProducer.cc
Go to the documentation of this file.
1 // //
3 // Producer of TkJet, //
4 // Cluster L1 tracks using fastjet //
5 // //
7 
8 // system include files
9 #include <memory>
10 
11 // user include files
24 
25 // L1 objects
30 
31 // geometry
34 
35 #include <fastjet/JetDefinition.hh>
36 
37 #include <string>
38 #include "TMath.h"
39 #include "TH1.h"
40 
41 using namespace l1t;
42 using namespace edm;
43 using namespace std;
44 
46 // //
47 // CLASS DEFINITION //
48 // //
50 
52 public:
54  typedef std::vector<L1TTTrackType> L1TTTrackCollectionType;
55 
57  ~L1TrackFastJetProducer() override;
58 
59  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
60 
61 private:
62  //virtual void beginJob();
63  void produce(edm::Event&, const edm::EventSetup&) override;
64  //virtual void endJob();
65 
66  // track selection criteria
67  const float trkZMax_; // in [cm]
68  const float trkChi2dofMax_; // maximum track chi2dof
69  const double trkBendChi2Max_; // maximum track bendchi2
70  const float trkPtMin_; // in [GeV]
71  const float trkEtaMax_; // in [rad]
72  const int trkNStubMin_; // minimum number of stubs
73  const int trkNPSStubMin_; // minimum number of PS stubs
74  const double deltaZ0Cut_; // save with |L1z-z0| < maxZ0
75  const double coneSize_; // Use anti-kt with this cone size
76  const bool doTightChi2_;
77  const float trkPtTightChi2_;
78  const float trkChi2dofTightChi2_;
79  const bool displaced_; //use prompt/displaced tracks
80 
84 };
85 
86 // constructor
88  : trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
89  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
90  trkBendChi2Max_(iConfig.getParameter<double>("trk_bendChi2Max")),
91  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
92  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
93  trkNStubMin_((int)iConfig.getParameter<int>("trk_nStubMin")),
94  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
95  deltaZ0Cut_((float)iConfig.getParameter<double>("deltaZ0Cut")),
96  coneSize_((float)iConfig.getParameter<double>("coneSize")),
97  doTightChi2_(iConfig.getParameter<bool>("doTightChi2")),
98  trkPtTightChi2_((float)iConfig.getParameter<double>("trk_ptTightChi2")),
99  trkChi2dofTightChi2_((float)iConfig.getParameter<double>("trk_chi2dofTightChi2")),
100  displaced_(iConfig.getParameter<bool>("displaced")),
101  trackToken_(consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_> > >(
102  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
103  pvToken_(consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("L1PrimaryVertexTag"))),
104  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
105  if (displaced_)
106  produces<TkJetCollection>("L1TrackFastJetsExtended");
107  else
108  produces<TkJetCollection>("L1TrackFastJets");
109 }
110 
111 // destructor
113 
114 // producer
116  std::unique_ptr<TkJetCollection> L1TrackFastJets(new TkJetCollection);
117 
118  // L1 tracks
120  iEvent.getByToken(trackToken_, TTTrackHandle);
121  std::vector<TTTrack<Ref_Phase2TrackerDigi_> >::const_iterator iterL1Track;
122 
123  // Tracker Topology
124  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
125 
126  edm::Handle<l1t::VertexCollection> TkPrimaryVertexHandle;
127  iEvent.getByToken(pvToken_, TkPrimaryVertexHandle);
128 
129  fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
130  std::vector<fastjet::PseudoJet> JetInputs;
131 
132  float recoVtx = TkPrimaryVertexHandle->begin()->z0();
133  unsigned int this_l1track = 0;
134  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
135  this_l1track++;
136  float trk_pt = iterL1Track->momentum().perp();
137  float trk_z0 = iterL1Track->z0();
138  float trk_chi2dof = iterL1Track->chi2Red();
139  float trk_bendchi2 = iterL1Track->stubPtConsistency();
140  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> >, TTStub<Ref_Phase2TrackerDigi_> > >
141  theStubs = iterL1Track->getStubRefs();
142  int trk_nstub = (int)theStubs.size();
143 
144  if (std::abs(trk_z0) > trkZMax_)
145  continue;
146  if (std::abs(iterL1Track->momentum().eta()) > trkEtaMax_)
147  continue;
148  if (trk_pt < trkPtMin_)
149  continue;
150  if (trk_nstub < trkNStubMin_)
151  continue;
152  if (trk_chi2dof > trkChi2dofMax_)
153  continue;
154  if (trk_bendchi2 > trkBendChi2Max_)
155  continue;
156  if (doTightChi2_ && (trk_pt > trkPtTightChi2_ && trk_chi2dof > trkChi2dofTightChi2_))
157  continue;
158 
159  int trk_nPS = 0;
160  for (int istub = 0; istub < trk_nstub; istub++) {
161  DetId detId(theStubs.at(istub)->getDetId());
162  bool tmp_isPS = false;
163  if (detId.det() == DetId::Detector::Tracker) {
164  if (detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3)
165  tmp_isPS = true;
166  else if (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9)
167  tmp_isPS = true;
168  }
169  if (tmp_isPS)
170  trk_nPS++;
171  }
172  if (trk_nPS < trkNPSStubMin_)
173  continue;
174  if (std::abs(recoVtx - trk_z0) > deltaZ0Cut_)
175  continue;
176 
177  fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
178  iterL1Track->momentum().y(),
179  iterL1Track->momentum().z(),
180  iterL1Track->momentum().mag());
181  JetInputs.push_back(psuedoJet); // input tracks for clustering
182  JetInputs.back().set_user_index(this_l1track - 1); // save track index in the collection
183  } // end loop over tracks
184 
185  fastjet::ClusterSequence cs(JetInputs, jet_def); // define the output jet collection
186  std::vector<fastjet::PseudoJet> JetOutputs =
187  fastjet::sorted_by_pt(cs.inclusive_jets(0)); // output jet collection, pT-ordered
188 
189  for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
191  JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
192  float sumpt = 0;
193  float avgZ = 0;
194  std::vector<edm::Ptr<L1TTTrackType> > L1TrackPtrs;
195  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
196 
197  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
198  auto index = fjConstituents[i].user_index();
199  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
200  L1TrackPtrs.push_back(trkPtr); // L1Tracks in the jet
201  sumpt = sumpt + trkPtr->momentum().perp();
202  avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
203  }
204  avgZ = avgZ / sumpt;
206  TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
207  L1TrackFastJets->push_back(trkJet);
208  } //end loop over Jet Outputs
209 
210  if (displaced_)
211  iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJetsExtended");
212  else
213  iEvent.put(std::move(L1TrackFastJets), "L1TrackFastJets");
214 }
215 
216 //void L1TrackFastJetProducer::beginJob() {}
217 
218 //void L1TrackFastJetProducer::endJob() {}
219 
221  //The following says we do not know what parameters are allowed so do no validation
222  // Please change this to state exactly what you do use, even if it is no parameters
224  desc.setUnknown();
225  descriptions.addDefault(desc);
226 }
227 
228 //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
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
delete x;
Definition: CaloConfig.h:22
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
edm::EDGetTokenT< VertexCollection > pvToken_
void produce(edm::Event &, const edm::EventSetup &) override
L1TrackFastJetProducer(const edm::ParameterSet &)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
std::vector< L1TTTrackType > L1TTTrackCollectionType
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getData(T &iHolder) const
Definition: EventSetup.h:122
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
Definition: DetId.h:17
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
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