CMS 3D CMS Logo

L1TruthTrackFastJetProducer.cc
Go to the documentation of this file.
1 // //
3 // Producer of TruTrkFastJet, //
4 // Cluster L1 tracks with truth info using fastjet //
5 // //
6 // Created by: Claire Savard (Oct. 2023) //
7 // //
9 
10 // system include files
11 #include <memory>
12 
13 // user include files
26 
27 // L1 objects
32 
33 // MC
36 
37 // geometry
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  ~L1TruthTrackFastJetProducer() override = default;
64 
65  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
66 
67 private:
68  void produce(edm::Event&, const edm::EventSetup&) override;
69 
70  // track selection criteria
71  const float trkZMax_; // in [cm]
72  const float trkPtMin_; // in [GeV]
73  const float trkEtaMax_; // in [rad]
74  const int trkNStubMin_; // minimum number of stubs
75  const int trkNPSStubMin_; // minimum number of PS stubs
76  const double coneSize_; // Use anti-kt with this cone size
77  const bool displaced_; //use prompt/displaced tracks
78 
82 };
83 
84 // constructor
86  : trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
87  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
88  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
89  trkNStubMin_((int)iConfig.getParameter<int>("trk_nStubMin")),
90  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
91  coneSize_((float)iConfig.getParameter<double>("coneSize")),
92  displaced_(iConfig.getParameter<bool>("displaced")),
93  trackToken_(consumes<std::vector<TTTrack<Ref_Phase2TrackerDigi_> > >(
94  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
96  ttTrackMCTruthToken_(consumes<TTTrackAssociationMap<Ref_Phase2TrackerDigi_> >(
97  iConfig.getParameter<edm::InputTag>("MCTruthTrackInputTag"))) {
98  if (displaced_)
99  produces<TkJetCollection>("L1TruthTrackFastJetsExtended");
100  else
101  produces<TkJetCollection>("L1TruthTrackFastJets");
102 }
103 
104 // producer
106  std::unique_ptr<TkJetCollection> L1TrackFastJets(new TkJetCollection);
107 
108  // L1 tracks
110  iEvent.getByToken(trackToken_, TTTrackHandle);
111  std::vector<TTTrack<Ref_Phase2TrackerDigi_> >::const_iterator iterL1Track;
112 
113  // MC truth
115  iEvent.getByToken(ttTrackMCTruthToken_, MCTruthTTTrackHandle);
116 
117  // Tracker Topology
118  const TrackerTopology& tTopo = iSetup.getData(tTopoToken_);
119 
120  fastjet::JetDefinition jet_def(fastjet::antikt_algorithm, coneSize_);
121  std::vector<fastjet::PseudoJet> JetInputs;
122 
123  unsigned int this_l1track = 0;
124  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
125  edm::Ptr<TTTrack<Ref_Phase2TrackerDigi_> > l1track_ptr(TTTrackHandle, this_l1track);
126  this_l1track++;
127  std::vector<edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> >, TTStub<Ref_Phase2TrackerDigi_> > >
128  theStubs = iterL1Track->getStubRefs();
129 
130  // standard quality cuts
131  if (std::abs(iterL1Track->z0()) > trkZMax_)
132  continue;
133  if (std::abs(iterL1Track->momentum().eta()) > trkEtaMax_)
134  continue;
135  if (iterL1Track->momentum().perp() < trkPtMin_)
136  continue;
137  int trk_nstub = (int)theStubs.size();
138  if (trk_nstub < trkNStubMin_)
139  continue;
140 
141  int trk_nPS = 0;
142  for (int istub = 0; istub < trk_nstub; istub++) {
143  DetId detId(theStubs.at(istub)->getDetId());
144  bool tmp_isPS = false;
145  if (detId.det() == DetId::Detector::Tracker) {
146  if (detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3)
147  tmp_isPS = true;
148  else if (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9)
149  tmp_isPS = true;
150  }
151  if (tmp_isPS)
152  trk_nPS++;
153  }
154  if (trk_nPS < trkNPSStubMin_)
155  continue;
156 
157  // check that trk is real and from hard interaction
158  edm::Ptr<TrackingParticle> my_tp = MCTruthTTTrackHandle->findTrackingParticlePtr(l1track_ptr);
159  if (my_tp.isNull()) // there is no tp match so the track is fake
160  continue;
161  int tp_eventid = my_tp->eventId().event();
162  if (tp_eventid > 0) // matched tp is from pileup
163  continue;
164 
165  fastjet::PseudoJet psuedoJet(iterL1Track->momentum().x(),
166  iterL1Track->momentum().y(),
167  iterL1Track->momentum().z(),
168  iterL1Track->momentum().mag());
169  JetInputs.push_back(psuedoJet); // input tracks for clustering
170  JetInputs.back().set_user_index(this_l1track - 1); // save track index in the collection
171  } // end loop over tracks
172 
173  fastjet::ClusterSequence cs(JetInputs, jet_def); // define the output jet collection
174  std::vector<fastjet::PseudoJet> JetOutputs =
175  fastjet::sorted_by_pt(cs.inclusive_jets(0)); // output jet collection, pT-ordered
176 
177  for (unsigned int ijet = 0; ijet < JetOutputs.size(); ++ijet) {
179  JetOutputs[ijet].px(), JetOutputs[ijet].py(), JetOutputs[ijet].pz(), JetOutputs[ijet].modp());
180  float sumpt = 0;
181  float avgZ = 0;
182  std::vector<edm::Ptr<L1TTTrackType> > L1TrackPtrs;
183  std::vector<fastjet::PseudoJet> fjConstituents = fastjet::sorted_by_pt(cs.constituents(JetOutputs[ijet]));
184 
185  for (unsigned int i = 0; i < fjConstituents.size(); ++i) {
186  auto index = fjConstituents[i].user_index();
187  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, index);
188  L1TrackPtrs.push_back(trkPtr); // L1Tracks in the jet
189  sumpt = sumpt + trkPtr->momentum().perp();
190  avgZ = avgZ + trkPtr->momentum().perp() * trkPtr->z0();
191  }
192  avgZ = avgZ / sumpt;
194  TkJet trkJet(jetP4, jetRef, L1TrackPtrs, avgZ);
195  L1TrackFastJets->push_back(trkJet);
196  } //end loop over Jet Outputs
197 
198  if (displaced_)
199  iEvent.put(std::move(L1TrackFastJets), "L1TruthTrackFastJetsExtended");
200  else
201  iEvent.put(std::move(L1TrackFastJets), "L1TruthTrackFastJets");
202 }
203 
205  // The following says we do not know what parameters are allowed so do no validation
207  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("l1tTTTracksFromTrackletEmulation", "Level1TTTracks"));
208  desc.add<edm::InputTag>("MCTruthTrackInputTag", edm::InputTag("TTTrackAssociatorFromPixelDigis", "Level1TTTracks"));
209  desc.add<double>("trk_zMax", 15.);
210  desc.add<double>("trk_ptMin", 2.0);
211  desc.add<double>("trk_etaMax", 2.4);
212  desc.add<int>("trk_nStubMin", 4);
213  desc.add<int>("trk_nPSStubMin", -1);
214  desc.add<double>("coneSize", 0.4);
215  desc.add<bool>("displaced", false);
216  descriptions.add("l1tTruthTrackFastJets", desc);
217 }
218 
219 //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
EncodedEventId eventId() const
Signal source, crossing number.
int event() const
get the contents of the subdetector field (should be protected?)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::vector< L1TTTrackType > L1TTTrackCollectionType
L1TruthTrackFastJetProducer(const edm::ParameterSet &)
delete x;
Definition: CaloConfig.h:22
GlobalVector momentum() const
Track momentum.
Definition: TTTrack.h:295
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
bool isNull() const
Checks for null.
Definition: Ptr.h:142
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr auto TOB
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
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
TTTrack< Ref_Phase2TrackerDigi_ > L1TTTrackType
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< TTTrackAssociationMap< Ref_Phase2TrackerDigi_ > > ttTrackMCTruthToken_