CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
L1FastTrackingJetProducer Class Reference
Inheritance diagram for L1FastTrackingJetProducer:
edm::stream::EDProducer<>

Public Types

typedef std::vector
< L1TTTrackType
L1TTTrackCollectionType
 
typedef TTTrack
< Ref_Phase2TrackerDigi_
L1TTTrackType
 
- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T...>
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T...>
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Public Member Functions

 L1FastTrackingJetProducer (const edm::ParameterSet &)
 
 ~L1FastTrackingJetProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

virtual void beginJob ()
 
virtual void endJob ()
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

double coneSize_
 
double deltaZ0Cut_
 
bool displaced_
 
bool doTightChi2_
 
const edm::EDGetTokenT
< TTTrackAssociationMap
< Ref_Phase2TrackerDigi_ > > 
genToken_
 
edm::EDGetTokenT< std::vector
< l1t::Vertex > > 
pvToken_
 
bool selectTrkMatchGenLoose_
 
bool selectTrkMatchGenOrPU_
 
bool selectTrkMatchGenTight_
 
const edm::EDGetTokenT
< std::vector< TTTrack
< Ref_Phase2TrackerDigi_ > > > 
trackToken_
 
double trkBendChi2Max_
 
float trkChi2dofMax_
 
float trkChi2dofTightChi2_
 
float trkEtaMax_
 
int trkNPSStubMin_
 
int trkNStubMin_
 
float trkPtMin_
 
float trkPtTightChi2_
 
float trkZMax_
 
edm::ESGetToken
< TrackerTopology,
TrackerTopologyRcd
tTopoToken_
 

Detailed Description

Definition at line 57 of file L1FastTrackingJetProducer.cc.

Member Typedef Documentation

Definition at line 60 of file L1FastTrackingJetProducer.cc.

Definition at line 59 of file L1FastTrackingJetProducer.cc.

Constructor & Destructor Documentation

L1FastTrackingJetProducer::L1FastTrackingJetProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 97 of file L1FastTrackingJetProducer.cc.

References coneSize_, deltaZ0Cut_, displaced_, doTightChi2_, edm::ParameterSet::getParameter(), selectTrkMatchGenLoose_, selectTrkMatchGenOrPU_, selectTrkMatchGenTight_, trkBendChi2Max_, trkChi2dofMax_, trkChi2dofTightChi2_, trkEtaMax_, trkNPSStubMin_, trkNStubMin_, trkPtMin_, trkPtTightChi2_, and trkZMax_.

99  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
100  pvToken_(consumes<std::vector<l1t::Vertex>>(iConfig.getParameter<edm::InputTag>("L1PrimaryVertexTag"))),
101  genToken_(
103  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
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 }
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const edm::EDGetTokenT< TTTrackAssociationMap< Ref_Phase2TrackerDigi_ > > genToken_
edm::EDGetTokenT< std::vector< l1t::Vertex > > pvToken_
Stores association of Truth Particles (TP) to L1 Track-Trigger Tracks.
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
L1FastTrackingJetProducer::~L1FastTrackingJetProducer ( )
override

Definition at line 127 of file L1FastTrackingJetProducer.cc.

127 {}

Member Function Documentation

void L1FastTrackingJetProducer::beginJob ( void  )
privatevirtual

Definition at line 249 of file L1FastTrackingJetProducer.cc.

249 {}
void L1FastTrackingJetProducer::endJob ( void  )
privatevirtual

Definition at line 251 of file L1FastTrackingJetProducer.cc.

251 {}
void L1FastTrackingJetProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 253 of file L1FastTrackingJetProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), submitPVResolutionJobs::desc, HLT_FULL_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

253  {
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 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void L1FastTrackingJetProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 130 of file L1FastTrackingJetProducer.cc.

References funct::abs(), coneSize_, fwrapper::cs, deltaZ0Cut_, displaced_, doTightChi2_, genToken_, edm::Event::getByToken(), edm::EventSetup::getData(), mps_fire::i, edm::EventBase::isRealData(), fwrapper::jet_def, L1FastTrackingJetProducer_cfi::L1FastTrackingJets, eostools::move(), edm::Event::put(), pvToken_, selectTrkMatchGenLoose_, selectTrkMatchGenOrPU_, selectTrkMatchGenTight_, StripSubdetector::TID, TrackerTopology::tidRing(), StripSubdetector::TOB, TrackerTopology::tobLayer(), align::Tracker, trackToken_, trkBendChi2Max_, trkChi2dofMax_, trkChi2dofTightChi2_, trkEtaMax_, trkNPSStubMin_, trkNStubMin_, trkPtMin_, trkPtTightChi2_, trkZMax_, and tTopoToken_.

130  {
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 }
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_
unsigned int tidRing(const DetId &id) const
JetDefinition jet_def
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
unique_ptr< ClusterSequence > cs
edm::EDGetTokenT< std::vector< l1t::Vertex > > pvToken_
bool isRealData() const
Definition: EventBase.h:62
bool getData(T &iHolder) const
Definition: EventSetup.h:128
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr auto TOB
Class to store the L1 Track Trigger stubs.
Definition: TTStub.h:22
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Definition: DetId.h:17
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
static constexpr auto TID
unsigned int tobLayer(const DetId &id) const

Member Data Documentation

double L1FastTrackingJetProducer::coneSize_
private

Definition at line 81 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

double L1FastTrackingJetProducer::deltaZ0Cut_
private

Definition at line 80 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

bool L1FastTrackingJetProducer::displaced_
private

Definition at line 85 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

bool L1FastTrackingJetProducer::doTightChi2_
private

Definition at line 82 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

const edm::EDGetTokenT<TTTrackAssociationMap<Ref_Phase2TrackerDigi_> > L1FastTrackingJetProducer::genToken_
private

Definition at line 92 of file L1FastTrackingJetProducer.cc.

Referenced by produce().

edm::EDGetTokenT<std::vector<l1t::Vertex> > L1FastTrackingJetProducer::pvToken_
private

Definition at line 91 of file L1FastTrackingJetProducer.cc.

Referenced by produce().

bool L1FastTrackingJetProducer::selectTrkMatchGenLoose_
private

Definition at line 87 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

bool L1FastTrackingJetProducer::selectTrkMatchGenOrPU_
private

Definition at line 88 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

bool L1FastTrackingJetProducer::selectTrkMatchGenTight_
private

Definition at line 86 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

const edm::EDGetTokenT<std::vector<TTTrack<Ref_Phase2TrackerDigi_> > > L1FastTrackingJetProducer::trackToken_
private

Definition at line 90 of file L1FastTrackingJetProducer.cc.

Referenced by produce().

double L1FastTrackingJetProducer::trkBendChi2Max_
private

Definition at line 75 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkChi2dofMax_
private

Definition at line 74 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkChi2dofTightChi2_
private

Definition at line 84 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkEtaMax_
private

Definition at line 77 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

int L1FastTrackingJetProducer::trkNPSStubMin_
private

Definition at line 79 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

int L1FastTrackingJetProducer::trkNStubMin_
private

Definition at line 78 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkPtMin_
private

Definition at line 76 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkPtTightChi2_
private

Definition at line 83 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

float L1FastTrackingJetProducer::trkZMax_
private

Definition at line 73 of file L1FastTrackingJetProducer.cc.

Referenced by L1FastTrackingJetProducer(), and produce().

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> L1FastTrackingJetProducer::tTopoToken_
private

Definition at line 93 of file L1FastTrackingJetProducer.cc.

Referenced by produce().