CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
L1TrackJetEmulationProducer Class Reference
Inheritance diagram for L1TrackJetEmulationProducer:
edm::stream::EDProducer<>

Public Types

typedef vector< L1TTTrackTypeL1TTTrackCollectionType
 
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

virtual TrackJetEmulationEtaPhiBinL1_cluster (TrackJetEmulationEtaPhiBin *phislice)
 
 L1TrackJetEmulationProducer (const ParameterSet &)
 
void L2_cluster (vector< Ptr< L1TTTrackType >> L1TrkPtrs_, TrackJetEmulationMaxZBin &mzb)
 
bool trackQualityCuts (int trk_nstub, float trk_chi2, float trk_bendchi2)
 
 ~L1TrackJetEmulationProducer () 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 (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void beginStream (StreamID) override
 
void endStream () override
 
void produce (Event &, const EventSetup &) override
 

Private Attributes

float d0CutNStubs4_
 
float d0CutNStubs5_
 
bool displaced_
 
int etaBins_
 
glbeta_intern etaStep_
 
float highpTJetMinpT_
 
int highpTJetMinTrackMultiplicity_
 
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
 
float lowpTJetMinpT_
 
int lowpTJetMinTrackMultiplicity_
 
const float MaxDzTrackPV
 
const double minTrkJetpT_
 
int nDisplacedTracks_
 
float nStubs4DisplacedBend_
 
float nStubs4DisplacedChi2_
 
float nStubs5DisplacedBend_
 
float nStubs5DisplacedChi2_
 
int phiBins_
 
glbphi_intern phiStep_
 
const EDGetTokenT< l1t::VertexWordCollectionPVtxToken_
 
float PVz
 
vector< int > tdtrk_
 
const EDGetTokenT< vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
 
const float trkBendChi2Max_
 
const float trkChi2dofMax_
 
const float trkEtaMax_
 
const int trkNPSStubMin_
 
const float trkPtMax_
 
const float trkPtMin_
 
const float trkZMax_
 
ESGetToken< TrackerTopology, TrackerTopologyRcdtTopoToken_
 
vector< int > zBinCount_
 
int zBins_
 
z0_intern zStep_
 

Detailed Description

Definition at line 54 of file L1TrackJetEmulationProducer.cc.

Member Typedef Documentation

◆ L1TTTrackCollectionType

Definition at line 59 of file L1TrackJetEmulationProducer.cc.

◆ L1TTTrackType

Definition at line 58 of file L1TrackJetEmulationProducer.cc.

Constructor & Destructor Documentation

◆ L1TrackJetEmulationProducer()

L1TrackJetEmulationProducer::L1TrackJetEmulationProducer ( const ParameterSet iConfig)
explicit

Definition at line 111 of file L1TrackJetEmulationProducer.cc.

References displaced_, etaBins_, etaStep_, M_PI, convert::makeGlbEta(), convert::makeGlbPhi(), convert::makeZ0(), phiBins_, phiStep_, trkEtaMax_, trkZMax_, zBins_, and zStep_.

112  : MaxDzTrackPV((float)iConfig.getParameter<double>("MaxDzTrackPV")),
113  trkZMax_((float)iConfig.getParameter<double>("trk_zMax")),
114  trkPtMax_((float)iConfig.getParameter<double>("trk_ptMax")),
115  trkPtMin_((float)iConfig.getParameter<double>("trk_ptMin")),
116  trkEtaMax_((float)iConfig.getParameter<double>("trk_etaMax")),
117  trkChi2dofMax_((float)iConfig.getParameter<double>("trk_chi2dofMax")),
118  trkBendChi2Max_((float)iConfig.getParameter<double>("trk_bendChi2Max")),
119  trkNPSStubMin_((int)iConfig.getParameter<int>("trk_nPSStubMin")),
120  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
121  etaBins_((int)iConfig.getParameter<int>("etaBins")),
122  phiBins_((int)iConfig.getParameter<int>("phiBins")),
123  zBins_((int)iConfig.getParameter<int>("zBins")),
124  d0CutNStubs4_((float)iConfig.getParameter<double>("d0_cutNStubs4")),
125  d0CutNStubs5_((float)iConfig.getParameter<double>("d0_cutNStubs5")),
126  lowpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
127  lowpTJetMinpT_((float)iConfig.getParameter<double>("lowpTJetMinpT")),
128  highpTJetMinTrackMultiplicity_((int)iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
129  highpTJetMinpT_((float)iConfig.getParameter<double>("highpTJetMinpT")),
130  displaced_(iConfig.getParameter<bool>("displaced")),
131  nStubs4DisplacedChi2_((float)iConfig.getParameter<double>("nStubs4DisplacedChi2")),
132  nStubs5DisplacedChi2_((float)iConfig.getParameter<double>("nStubs5DisplacedChi2")),
133  nStubs4DisplacedBend_((float)iConfig.getParameter<double>("nStubs4Displacedbend")),
134  nStubs5DisplacedBend_((float)iConfig.getParameter<double>("nStubs5Displacedbend")),
135  nDisplacedTracks_((int)iConfig.getParameter<int>("nDisplacedTracks")),
136  trackToken_(consumes<vector<TTTrack<Ref_Phase2TrackerDigi_>>>(iConfig.getParameter<InputTag>("L1TrackInputTag"))),
137  PVtxToken_(consumes<l1t::VertexWordCollection>(iConfig.getParameter<InputTag>("VertexInputTag"))),
138  tTopoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>(edm::ESInputTag("", ""))) {
139  zStep_ = convert::makeZ0(2.0 * trkZMax_ / (zBins_ + 1));
140  etaStep_ = convert::makeGlbEta(2.0 * trkEtaMax_ / etaBins_); //etaStep is the width of an etabin
142 
143  if (displaced_)
144  produces<l1t::TkJetWordCollection>("L1TrackJetsExtended");
145  else
146  produces<l1t::TkJetWordCollection>("L1TrackJets");
147 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
glbeta_intern makeGlbEta(float eta)
z0_intern makeZ0(float z0)
glbphi_intern makeGlbPhi(float phi)
const EDGetTokenT< l1t::VertexWordCollection > PVtxToken_
#define M_PI
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29
const EDGetTokenT< vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_

◆ ~L1TrackJetEmulationProducer()

L1TrackJetEmulationProducer::~L1TrackJetEmulationProducer ( )
override

Definition at line 149 of file L1TrackJetEmulationProducer.cc.

149 {}

Member Function Documentation

◆ beginStream()

void L1TrackJetEmulationProducer::beginStream ( StreamID  )
overrideprivate

Definition at line 656 of file L1TrackJetEmulationProducer.cc.

656 {}

◆ endStream()

void L1TrackJetEmulationProducer::endStream ( )
overrideprivate

Definition at line 658 of file L1TrackJetEmulationProducer.cc.

658 {}

◆ fillDescriptions()

void L1TrackJetEmulationProducer::fillDescriptions ( ConfigurationDescriptions descriptions)
static

Definition at line 671 of file L1TrackJetEmulationProducer.cc.

References edm::ConfigurationDescriptions::addDefault(), and submitPVResolutionJobs::desc.

671  {
672  //The following says we do not know what parameters are allowed so do no validation
673  // Please change this to state exactly what you do use, even if it is no parameters
675  desc.setUnknown();
676  descriptions.addDefault(desc);
677 }
void addDefault(ParameterSetDescription const &psetDescription)

◆ L1_cluster()

TrackJetEmulationEtaPhiBin * L1TrackJetEmulationProducer::L1_cluster ( TrackJetEmulationEtaPhiBin phislice)
virtual

Definition at line 568 of file L1TrackJetEmulationProducer.cc.

References bsc_activity_cfg::clusters, PVValHelper::eta, etaBins_, etaStep_, mps_fire::i, visualization-live-secondInstance_cfg::m, DeadROCCounter::nclust, TrackJetEmulationEtaPhiBin::ntracks, TrackJetEmulationEtaPhiBin::nxtracks, TrackJetEmulationEtaPhiBin::pTtot, and TrackJetEmulationEtaPhiBin::used.

Referenced by L2_cluster().

568  {
570  for (int etabin = 0; etabin < etaBins_ / 2; ++etabin) {
571  clusters[etabin].pTtot = 0;
572  clusters[etabin].ntracks = 0;
573  clusters[etabin].nxtracks = 0;
574  clusters[etabin].phi = 0;
575  clusters[etabin].eta = 0;
576  clusters[etabin].used = false;
577  }
578 
579  if (clusters == nullptr)
580  edm::LogWarning("L1TrackJetEmulationProducer") << "Clusters memory not assigned!\n";
581 
582  // Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used
583  pt_intern my_pt, left_pt, right_pt, right2pt;
584  int nclust = 0;
585  right2pt = 0;
586  for (int etabin = 0; etabin < etaBins_; ++etabin) {
587  // assign values for my pT and neighbors' pT
588  if (phislice[etabin].used)
589  continue;
590  my_pt = phislice[etabin].pTtot;
591  if (etabin > 0 && !phislice[etabin - 1].used) {
592  left_pt = phislice[etabin - 1].pTtot;
593  } else
594  left_pt = 0;
595  if (etabin < etaBins_ - 1 && !phislice[etabin + 1].used) {
596  right_pt = phislice[etabin + 1].pTtot;
597  if (etabin < etaBins_ - 2 && !phislice[etabin + 2].used) {
598  right2pt = phislice[etabin + 2].pTtot;
599  } else
600  right2pt = 0;
601  } else
602  right_pt = 0;
603 
604  // if I'm not a cluster, move on
605  if (my_pt < left_pt || my_pt <= right_pt) {
606  // if unused pT in the left neighbor, spit it out as a cluster
607  if (left_pt > 0) {
608  clusters[nclust] = phislice[etabin - 1];
609  phislice[etabin - 1].used = true;
610  nclust++;
611  }
612  continue;
613  }
614 
615  // I guess I'm a cluster-- should I use my right neighbor?
616  // Note: left neighbor will definitely be used because if it
617  // didn't belong to me it would have been used already
618  clusters[nclust] = phislice[etabin];
619  phislice[etabin].used = true;
620  if (left_pt > 0) {
621  clusters[nclust].pTtot += left_pt;
622  clusters[nclust].ntracks += phislice[etabin - 1].ntracks;
623  clusters[nclust].nxtracks += phislice[etabin - 1].nxtracks;
624  }
625  if (my_pt >= right2pt && right_pt > 0) {
626  clusters[nclust].pTtot += right_pt;
627  clusters[nclust].ntracks += phislice[etabin + 1].ntracks;
628  clusters[nclust].nxtracks += phislice[etabin + 1].nxtracks;
629  phislice[etabin + 1].used = true;
630  }
631  nclust++;
632  } // for each etabin
633 
634  // Now merge clusters, if necessary
635  for (int m = 0; m < nclust - 1; ++m) {
636  if (clusters[m + 1].eta - clusters[m].eta < 3 * etaStep_ / 2 &&
637  clusters[m].eta - clusters[m + 1].eta < 3 * etaStep_ / 2) {
638  if (clusters[m + 1].pTtot > clusters[m].pTtot) {
639  clusters[m].eta = clusters[m + 1].eta;
640  }
641  clusters[m].pTtot += clusters[m + 1].pTtot;
642  clusters[m].ntracks += clusters[m + 1].ntracks; // Previous version didn't add tracks when merging
643  clusters[m].nxtracks += clusters[m + 1].nxtracks;
644  for (int m1 = m + 1; m1 < nclust - 1; ++m1)
645  clusters[m1] = clusters[m1 + 1];
646  nclust--;
647  m = -1;
648  } // end if clusters neighbor in eta
649  } // end for (m) loop
650 
651  for (int i = nclust; i < etaBins_ / 2; ++i) // zero out remaining unused clusters
652  clusters[i].pTtot = 0;
653  return clusters;
654 }
ap_ufixed< 16+PT_EXTRABITS, 11, AP_TRN, AP_SAT > pt_intern
Log< level::Warning, false > LogWarning

◆ L2_cluster()

void L1TrackJetEmulationProducer::L2_cluster ( vector< Ptr< L1TTTrackType >>  L1TrkPtrs_,
TrackJetEmulationMaxZBin mzb 
)

Definition at line 254 of file L1TrackJetEmulationProducer.cc.

References TrackJetEmulationMaxZBin::clusters, TrackJetEmulationEtaPhiBin::eta, PVValHelper::eta, convert::ETA_LSB, etaBins_, muonTiming_cfi::etamax, muonTiming_cfi::etamin, etaStep_, highpTJetMinpT_, highpTJetMinTrackMultiplicity_, TrackJetEmulationMaxZBin::ht, mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, l1tmetemu::kPtMagSize, L1_cluster(), L1TrkPtrs_, lowpTJetMinpT_, lowpTJetMinTrackMultiplicity_, visualization-live-secondInstance_cfg::m, M_PI, convert::makeGlbEta(), convert::makeGlbPhi(), convert::makePtFromFloat(), convert::makeZ0(), minTrkJetpT_, dqmiodumpmetadata::n, DeadROCCounter::nclust, TrackJetEmulationMaxZBin::nclust, vertices_cff::ntracks, TrackJetEmulationEtaPhiBin::ntracks, TrackJetEmulationEtaPhiBin::nxtracks, phi, TrackJetEmulationEtaPhiBin::phi, convert::PHI_LSB, phiBins_, phimax, phimin, phiStep_, TrackJetEmulationEtaPhiBin::pTtot, TTTrack_TrackWord::stepPhi0, TTTrack_TrackWord::stepZ0, trkEtaMax_, trkPtMax_, trkZMax_, TrackJetEmulationEtaPhiBin::used, convert::Z0_LSB, TrackJetEmulationMaxZBin::zbincenter, zBinCount_, zBins_, SiStripMonitorCluster_cfi::zmax, SiStripMonitorCluster_cfi::zmin, TrackJetEmulationMaxZBin::znum, and zStep_.

Referenced by produce().

254  {
255  enum TrackBitWidths {
256  kEtaSize = 16, // Width of z-position (40cm / 0.1)
257  kEtaMagSize = 3, // Width of z-position magnitude (signed)
258  kPtSize = 14, // Width of pt
259  kPtMagSize = 9, // Width of pt magnitude (unsigned)
260  kPhiSize = 12,
261  kPhiMagSize = 1,
262  };
263 
264  const int nz = zBins_ + 1;
265  //values are initialized by 0 as struct assigns default values for members
266  TrackJetEmulationMaxZBin all_zBins[nz];
267 
269  z0_intern zmax = zmin + 2 * zStep_;
270 
271  TrackJetEmulationEtaPhiBin epbins[phiBins_][etaBins_]; // create grid of phiBins
275  for (int i = 0; i < phiBins_; ++i) {
277  for (int j = 0; j < etaBins_; ++j) {
278  phimin = phi;
279  phimax = phi + phiStep_;
280  etamin = eta;
281  eta = eta + etaStep_;
282  etamax = eta;
283  epbins[i][j].phi = (phimin + phimax) / 2;
284  epbins[i][j].eta = (etamin + etamax) / 2;
285  epbins[i][j].pTtot = 0;
286  epbins[i][j].ntracks = 0;
287  epbins[i][j].nxtracks = 0;
288  } // for each etabin
289  phi = phi + phiStep_;
290  } // for each phibin (finished creating epbins)
291 
292  mzb = all_zBins[0];
293  mzb.ht = 0;
294  int ntracks = L1TrkPtrs_.size();
295  // uninitalized arrays
298 
299  for (int zbin = 0; zbin < zBins_; ++zbin) {
300  for (int i = 0; i < phiBins_; ++i) { //First initialize pT, numtracks, used to 0 (or false)
301  for (int j = 0; j < etaBins_; ++j) {
302  epbins[i][j].pTtot = 0;
303  epbins[i][j].used = false;
304  epbins[i][j].ntracks = 0;
305  epbins[i][j].nxtracks = 0;
306  } //for each etabin
307  L1clusters[i] = epbins[i];
308  } //for each phibin
309 
310  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
311  ap_ufixed<TrackBitWidths::kPtSize, TrackBitWidths::kPtMagSize, AP_RND_CONV, AP_SAT> inputTrkPt = 0;
312  inputTrkPt.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kRinvMSB - 1,
313  TTTrack_TrackWord::TrackBitLocations::kRinvLSB);
314  pt_intern trkpt = inputTrkPt;
315 
316  ap_fixed<TrackBitWidths::kEtaSize, TrackBitWidths::kEtaMagSize, AP_RND_CONV, AP_SAT> trketainput = 0;
317  trketainput.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kTanlMSB,
318  TTTrack_TrackWord::TrackBitLocations::kTanlLSB);
319  ap_ufixed<64 + ETA_EXTRABITS, 32 + ETA_EXTRABITS> eta_conv =
320  1.0 / convert::ETA_LSB; //conversion factor from input eta format to internal format
321  glbeta_intern trketa = eta_conv * trketainput;
322 
323  int Sector = L1TrkPtrs_[k]->phiSector();
324  glbphi_intern trkphiSector = 0;
325  if (Sector < 5) {
326  trkphiSector = Sector * convert::makeGlbPhi(2.0 * M_PI / 9.0);
327  } else {
328  trkphiSector =
329  convert::makeGlbPhi(-1.0 * M_PI + M_PI / 9.0) + (Sector - 5) * convert::makeGlbPhi(2.0 * M_PI / 9.0);
330  }
331 
332  ap_fixed<TrackBitWidths::kPhiSize, TrackBitWidths::kPhiMagSize, AP_RND_CONV, AP_SAT> trkphiinput = 0;
333  trkphiinput.V = L1TrkPtrs_[k]->getTrackWord()(TTTrack_TrackWord::TrackBitLocations::kPhiMSB,
334  TTTrack_TrackWord::TrackBitLocations::kPhiLSB);
335  ap_ufixed<64 + PHI_EXTRABITS, 32 + PHI_EXTRABITS> phi_conv =
337  (1 << (TrackBitWidths::kPhiSize - TrackBitWidths::kPhiMagSize));
338  //phi_conv is a conversion factor from input phi format to the internal format
339  glbeta_intern localphi = phi_conv * trkphiinput;
340  glbeta_intern trkphi = localphi + trkphiSector;
341 
342  ap_int<TTTrack_TrackWord::TrackBitWidths::kZ0Size> inputTrkZ0 = L1TrkPtrs_[k]->getTrackWord()(
343  TTTrack_TrackWord::TrackBitLocations::kZ0MSB, TTTrack_TrackWord::TrackBitLocations::kZ0LSB);
344  ap_ufixed<32 + Z0_EXTRABITS, 1 + Z0_EXTRABITS> z0_conv =
345  TTTrack_TrackWord::stepZ0 / convert::Z0_LSB; //conversion factor from input z format to internal format
346  z0_intern trkZ = z0_conv * inputTrkZ0;
347 
348  for (int i = 0; i < phiBins_; ++i) {
349  for (int j = 0; j < etaBins_; ++j) {
350  L2cluster[k] = epbins[i][j];
351  if ((zmin <= trkZ && zmax >= trkZ) &&
352  ((epbins[i][j].eta - etaStep_ / 2 < trketa && epbins[i][j].eta + etaStep_ / 2 >= trketa) &&
353  epbins[i][j].phi - phiStep_ / 2 < trkphi && epbins[i][j].phi + phiStep_ / 2 >= trkphi &&
354  (zBinCount_[k] != 2))) {
355  zBinCount_.at(k) = zBinCount_.at(k) + 1;
356 
357  if (trkpt < convert::makePtFromFloat(trkPtMax_))
358  epbins[i][j].pTtot += trkpt;
359  else
361  ++epbins[i][j].ntracks;
362  //x-bit is currently not used in firmware, so we leave nxtracks = 0 for now
363  } // if right bin
364  } // for each phibin: j loop
365  } // for each phibin: i loop
366  } // end loop over tracks
367 
368  for (int phislice = 0; phislice < phiBins_; ++phislice) {
369  L1clusters[phislice] = L1_cluster(epbins[phislice]);
370  for (int ind = 0; L1clusters[phislice][ind].pTtot != 0; ++ind) {
371  L1clusters[phislice][ind].used = false;
372  }
373  }
374 
375  //Create clusters array to hold output cluster data for Layer2; can't have more clusters than tracks.
376  //Find eta-phibin with maxpT, make center of cluster, add neighbors if not already used.
377  pt_intern hipT = 0;
378  int nclust = 0;
379  int phibin = 0;
380  int imax = -1;
381  int index1; //index of clusters array for each phislice
382  pt_intern E1 = 0;
383  pt_intern E0 = 0;
384  pt_intern E2 = 0;
385  l1t::TkJetWord::nt_t ntrk1, ntrk2;
386  l1t::TkJetWord::nx_t nxtrk1, nxtrk2;
387  int used1, used2, used3, used4;
388 
389  for (phibin = 0; phibin < phiBins_; ++phibin) { //Find eta-phibin with highest pT
390  while (true) {
391  hipT = 0;
392  for (index1 = 0; L1clusters[phibin][index1].pTtot > 0; ++index1) {
393  if (!L1clusters[phibin][index1].used && L1clusters[phibin][index1].pTtot >= hipT) {
394  hipT = L1clusters[phibin][index1].pTtot;
395  imax = index1;
396  }
397  } // for each index within the phibin
398 
399  if (hipT == 0)
400  break; // If highest pT is 0, all bins are used
401  E0 = hipT; // E0 is pT of first phibin of the cluster
402  E1 = 0;
403  E2 = 0;
404  ntrk1 = 0;
405  ntrk2 = 0;
406  nxtrk1 = 0;
407  nxtrk2 = 0;
408  L2cluster[nclust] = L1clusters[phibin][imax];
409 
410  L1clusters[phibin][imax].used = true;
411  // Add pT of upper neighbor
412  // E1 is pT of the middle phibin (should be highest pT)
413  if (phibin != phiBins_ - 1) {
414  used1 = -1;
415  used2 = -1;
416  for (index1 = 0; L1clusters[phibin + 1][index1].pTtot != 0; ++index1) {
417  if (L1clusters[phibin + 1][index1].used)
418  continue;
419  if (L1clusters[phibin + 1][index1].eta - L1clusters[phibin][imax].eta <= 3 * etaStep_ / 2 &&
420  L1clusters[phibin][imax].eta - L1clusters[phibin + 1][index1].eta <= 3 * etaStep_ / 2) {
421  E1 += L1clusters[phibin + 1][index1].pTtot;
422  ntrk1 += L1clusters[phibin + 1][index1].ntracks;
423  nxtrk1 += L1clusters[phibin + 1][index1].nxtracks;
424  if (used1 < 0)
425  used1 = index1;
426  else
427  used2 = index1;
428  } // if cluster is within one phibin
429  } // for each cluster in above phibin
430 
431  if (E1 < E0) { // if E1 isn't higher, E0 and E1 are their own cluster
432  L2cluster[nclust].pTtot += E1;
433  L2cluster[nclust].ntracks += ntrk1;
434  L2cluster[nclust].nxtracks += nxtrk1;
435  if (used1 >= 0)
436  L1clusters[phibin + 1][used1].used = true;
437  if (used2 >= 0)
438  L1clusters[phibin + 1][used2].used = true;
439  nclust++;
440  continue;
441  }
442 
443  if (phibin != phiBins_ - 2) { // E2 will be the pT of the third phibin (should be lower than E1)
444  used3 = -1;
445  used4 = -1;
446  for (index1 = 0; L1clusters[phibin + 2][index1].pTtot != 0; ++index1) {
447  if (L1clusters[phibin + 2][index1].used)
448  continue;
449  if (L1clusters[phibin + 2][index1].eta - L1clusters[phibin][imax].eta <= 3 * etaStep_ / 2 &&
450  L1clusters[phibin][imax].eta - L1clusters[phibin + 2][index1].eta <= 3 * etaStep_ / 2) {
451  E2 += L1clusters[phibin + 2][index1].pTtot;
452  ntrk2 += L1clusters[phibin + 2][index1].ntracks;
453  nxtrk2 += L1clusters[phibin + 2][index1].nxtracks;
454  if (used3 < 0)
455  used3 = index1;
456  else
457  used4 = index1;
458  }
459  }
460  // if indeed E2 < E1, add E1 and E2 to E0, they're all a cluster together
461  // otherwise, E0 is its own cluster
462  if (E2 < E1) {
463  L2cluster[nclust].pTtot += E1 + E2;
464  L2cluster[nclust].ntracks += ntrk1 + ntrk2;
465  L2cluster[nclust].nxtracks += nxtrk1 + nxtrk2;
466  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
467  if (used1 >= 0)
468  L1clusters[phibin + 1][used1].used = true;
469  if (used2 >= 0)
470  L1clusters[phibin + 1][used2].used = true;
471  if (used3 >= 0)
472  L1clusters[phibin + 2][used3].used = true;
473  if (used4 >= 0)
474  L1clusters[phibin + 2][used4].used = true;
475  }
476  nclust++;
477  continue;
478  } // end Not phiBins-2
479  else {
480  L2cluster[nclust].pTtot += E1;
481  L2cluster[nclust].ntracks += ntrk1;
482  L2cluster[nclust].nxtracks += nxtrk1;
483  L2cluster[nclust].phi = L1clusters[phibin + 1][used1].phi;
484  if (used1 >= 0)
485  L1clusters[phibin + 1][used1].used = true;
486  if (used2 >= 0)
487  L1clusters[phibin + 1][used2].used = true;
488  nclust++;
489  continue;
490  }
491  } //End not last phibin(23)
492  else { //if it is phibin 23
493  L1clusters[phibin][imax].used = true;
494  nclust++;
495  }
496  } // while hipT not 0
497  } // for each phibin
498 
499  for (phibin = 0; phibin < phiBins_; ++phibin)
500  delete[] L1clusters[phibin];
501 
502  // Now merge clusters, if necessary
503  for (int m = 0; m < nclust - 1; ++m) {
504  for (int n = m + 1; n < nclust; ++n) {
505  if (L2cluster[n].eta == L2cluster[m].eta &&
506  ((L2cluster[n].phi - L2cluster[m].phi < 3 * phiStep_ / 2 &&
507  L2cluster[m].phi - L2cluster[n].phi < 3 * phiStep_ / 2) ||
508  (L2cluster[n].phi - L2cluster[m].phi > 2 * convert::makeGlbPhi(M_PI) - phiStep_ ||
509  L2cluster[m].phi - L2cluster[n].phi > 2 * convert::makeGlbPhi(M_PI) - phiStep_))) {
510  if (L2cluster[n].pTtot > L2cluster[m].pTtot) {
511  L2cluster[m].phi = L2cluster[n].phi;
512  }
513  L2cluster[m].pTtot += L2cluster[n].pTtot;
514  L2cluster[m].ntracks += L2cluster[n].ntracks;
515  L2cluster[m].nxtracks += L2cluster[n].nxtracks;
516  for (int m1 = n; m1 < nclust - 1; ++m1) {
517  L2cluster[m1] = L2cluster[m1 + 1];
518  }
519  nclust--;
520  m = -1;
521  break; //?????
522  } // end if clusters neighbor in eta
523  }
524  } // end for (m) loop
525  // sum up all pTs in this zbin to find ht
526 
527  pt_intern ht = 0;
528  for (int k = 0; k < nclust; ++k) {
529  if (L2cluster[k].pTtot > convert::makePtFromFloat(lowpTJetMinpT_) &&
531  continue;
532  if (L2cluster[k].pTtot > convert::makePtFromFloat(highpTJetMinpT_) &&
534  continue;
535  if (L2cluster[k].pTtot > convert::makePtFromFloat(minTrkJetpT_))
536  ht += L2cluster[k].pTtot;
537  }
538  // if ht is larger than previous max, this is the new vertex zbin
539  all_zBins[zbin].znum = zbin;
540  all_zBins[zbin].clusters = new TrackJetEmulationEtaPhiBin[nclust];
541  all_zBins[zbin].nclust = nclust;
542  all_zBins[zbin].zbincenter = (zmin + zmax) / 2;
543  for (int k = 0; k < nclust; ++k) {
544  all_zBins[zbin].clusters[k].phi = L2cluster[k].phi;
545  all_zBins[zbin].clusters[k].eta = L2cluster[k].eta;
546  all_zBins[zbin].clusters[k].pTtot = L2cluster[k].pTtot;
547  all_zBins[zbin].clusters[k].ntracks = L2cluster[k].ntracks;
548  all_zBins[zbin].clusters[k].nxtracks = L2cluster[k].nxtracks;
549  }
550  all_zBins[zbin].ht = ht;
551  if (ht >= mzb.ht) {
552  mzb = all_zBins[zbin];
553  mzb.zbincenter = (zmin + zmax) / 2;
554  }
555  // Prepare for next zbin!
556  zmin = zmin + zStep_;
557  zmax = zmax + zStep_;
558  } // for each zbin
559 
560  for (int zbin = 0; zbin < zBins_; ++zbin) {
561  if (zbin == mzb.znum) {
562  continue;
563  }
564  delete[] all_zBins[zbin].clusters;
565  }
566 }
pt_intern makePtFromFloat(float pt)
glbeta_intern makeGlbEta(float eta)
ap_int< 14+PHI_EXTRABITS > glbphi_intern
static const float ETA_LSB
z0_intern makeZ0(float z0)
ap_ufixed< 16+PT_EXTRABITS, 11, AP_TRN, AP_SAT > pt_intern
ap_int< 10+Z0_EXTRABITS > z0_intern
ap_int< 14+ETA_EXTRABITS > glbeta_intern
glbphi_intern makeGlbPhi(float phi)
static const float PHI_LSB
ap_uint< kXtSize > nx_t
Definition: TkJetWord.h:59
#define M_PI
static constexpr double stepPhi0
const unsigned int kPtMagSize
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
ap_uint< kNtSize > nt_t
Definition: TkJetWord.h:58
TrackJetEmulationEtaPhiBin * clusters
static constexpr double stepZ0
static const float Z0_LSB
virtual TrackJetEmulationEtaPhiBin * L1_cluster(TrackJetEmulationEtaPhiBin *phislice)

◆ produce()

void L1TrackJetEmulationProducer::produce ( Event iEvent,
const EventSetup iSetup 
)
overrideprivate

Definition at line 151 of file L1TrackJetEmulationProducer.cc.

References funct::abs(), TTTrack< T >::chi2Red(), TrackJetEmulationMaxZBin::clusters, TTTrack< T >::d0(), d0CutNStubs4_, d0CutNStubs5_, displaced_, TrackJetEmulationEtaPhiBin::eta, PV3DBase< T, PVType, FrameType >::eta(), convert::ETA_LSB_POW, edm::EventSetup::getData(), TTTrack< T >::getStubRefs(), iEvent, createfilelist::int, dqmiolumiharvest::j, reco::btau::jetEta, reco::btau::jetPhi, reco::btau::jetPt, L1TrkPtrs_, L2_cluster(), convert::makePtFromFloat(), MaxDzTrackPV, TTTrack< T >::momentum(), eostools::move(), TrackJetEmulationMaxZBin::nclust, TrackJetEmulationEtaPhiBin::ntracks, TrackJetEmulationEtaPhiBin::nxtracks, PV3DBase< T, PVType, FrameType >::perp(), TrackJetEmulationEtaPhiBin::phi, convert::PHI_LSB_POW, TrackJetEmulationEtaPhiBin::pTtot, PVtxToken_, PVz, TTTrack< T >::stubPtConsistency(), tdtrk_, StripSubdetector::TID, TrackerTopology::tidRing(), StripSubdetector::TOB, TrackerTopology::tobLayer(), align::Tracker, trackQualityCuts(), trackToken_, trkEtaMax_, trkNPSStubMin_, trkPtMin_, trkZMax_, tTopoToken_, HLTMuonOfflineAnalyzer_cfi::z0, TTTrack< T >::z0(), convert::Z0_LSB_POW, TrackJetEmulationMaxZBin::zbincenter, and zBinCount_.

151  {
152  unique_ptr<l1t::TkJetWordCollection> L1L1TrackJetProducer(new l1t::TkJetWordCollection);
153 
154  // For TTStubs
155  const TrackerTopology &tTopo = iSetup.getData(tTopoToken_);
156 
158  iEvent.getByToken(trackToken_, TTTrackHandle);
159  vector<TTTrack<Ref_Phase2TrackerDigi_>>::const_iterator iterL1Track;
160 
162  iEvent.getByToken(PVtxToken_, PVtx);
163  float PVz = (PVtx->at(0)).z0();
164 
165  L1TrkPtrs_.clear();
166  zBinCount_.clear();
167  tdtrk_.clear();
168  unsigned int this_l1track = 0;
169  for (iterL1Track = TTTrackHandle->begin(); iterL1Track != TTTrackHandle->end(); iterL1Track++) {
170  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
171  this_l1track++;
172  float trk_pt = trkPtr->momentum().perp();
173  int trk_nstubs = (int)trkPtr->getStubRefs().size();
174  float trk_chi2dof = trkPtr->chi2Red();
175  float trk_d0 = trkPtr->d0();
176  float trk_bendchi2 = trkPtr->stubPtConsistency();
177  float trk_z0 = trkPtr->z0();
178 
179  int trk_nPS = 0;
180  for (int istub = 0; istub < trk_nstubs; istub++) { // loop over the stubs
181  DetId detId(trkPtr->getStubRefs().at(istub)->getDetId());
182  if (detId.det() == DetId::Detector::Tracker) {
183  if ((detId.subdetId() == StripSubdetector::TOB && tTopo.tobLayer(detId) <= 3) ||
184  (detId.subdetId() == StripSubdetector::TID && tTopo.tidRing(detId) <= 9))
185  trk_nPS++;
186  }
187  }
188 
189  if (trk_nPS < trkNPSStubMin_)
190  continue;
191  if (!trackQualityCuts(trk_nstubs, trk_chi2dof, trk_bendchi2))
192  continue;
193  if (std::abs(trk_z0 - PVz) > MaxDzTrackPV)
194  continue;
195  if (std::abs(trk_z0) > trkZMax_)
196  continue;
197  if (std::abs(trkPtr->momentum().eta()) > trkEtaMax_)
198  continue;
199  if (trk_pt < trkPtMin_)
200  continue;
201  L1TrkPtrs_.push_back(trkPtr);
202  zBinCount_.push_back(0);
203 
204  if ((std::abs(trk_d0) > d0CutNStubs5_ && trk_nstubs >= 5 && d0CutNStubs5_ >= 0) ||
205  (trk_nstubs == 4 && std::abs(trk_d0) > d0CutNStubs4_ && d0CutNStubs4_ >= 0))
206  tdtrk_.push_back(1); //displaced track
207  else
208  tdtrk_.push_back(0); // not displaced track
209  }
210 
211  if (!L1TrkPtrs_.empty()) {
213 
214  L2_cluster(L1TrkPtrs_, mzb);
215  if (mzb.clusters != nullptr) {
216  for (int j = 0; j < mzb.nclust; ++j) {
217  //FILL Two Layer Jets for Jet Collection
219  continue; //protects against reading bad memory
220  if (mzb.clusters[j].ntracks < 1)
221  continue;
222  if (mzb.clusters[j].ntracks > 5000)
223  continue;
228  l1t::TkJetWord::nt_t totalntracks_ = mzb.clusters[j].ntracks;
229  l1t::TkJetWord::nx_t totalxtracks_ = mzb.clusters[j].nxtracks;
230  l1t::TkJetWord::tkjetunassigned_t unassigned_ = 0;
231 
232  l1t::TkJetWord trkJet(jetPt, jetEta, jetPhi, jetZ0, totalntracks_, totalxtracks_, unassigned_);
233  //trkJet.setDispCounters(DispCounters);
234  L1L1TrackJetProducer->push_back(trkJet);
235  }
236  } else if (mzb.clusters == nullptr) {
237  edm::LogWarning("L1TrackJetEmulationProducer") << "mzb.clusters Not Assigned!\n";
238  }
239 
240  if (displaced_)
241  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
242  else
243  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
244  delete[] mzb.clusters;
245  } else if (L1TrkPtrs_.empty()) {
246  edm::LogWarning("L1TrackJetEmulationProducer") << "L1TrkPtrs Not Assigned!\n";
247  if (displaced_)
248  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJetsExtended");
249  else
250  iEvent.put(std::move(L1L1TrackJetProducer), "L1TrackJets");
251  }
252 }
static const float ETA_LSB_POW
pt_intern makePtFromFloat(float pt)
unsigned int tobLayer(const DetId &id) const
std::vector< l1t::TkJetWord > TkJetWordCollection
Definition: TkJetWord.h:156
static const float Z0_LSB_POW
ap_int< kGlbEtaSize > glbeta_t
Definition: TkJetWord.h:55
ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
ap_int< kGlbPhiSize > glbphi_t
Definition: TkJetWord.h:56
void L2_cluster(vector< Ptr< L1TTTrackType >> L1TrkPtrs_, TrackJetEmulationMaxZBin &mzb)
bool trackQualityCuts(int trk_nstub, float trk_chi2, float trk_bendchi2)
ap_ufixed< kPtSize, kPtMagSize, AP_TRN, AP_SAT > pt_t
Definition: TkJetWord.h:54
static const float PHI_LSB_POW
ap_uint< TkJetBitWidths::kUnassignedSize > tkjetunassigned_t
Definition: TkJetWord.h:60
ap_int< kZ0Size > z0_t
Definition: TkJetWord.h:57
int iEvent
Definition: GenABIO.cc:224
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr auto TOB
const EDGetTokenT< l1t::VertexWordCollection > PVtxToken_
ap_uint< kXtSize > nx_t
Definition: TkJetWord.h:59
Definition: DetId.h:17
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
ap_uint< kNtSize > nt_t
Definition: TkJetWord.h:58
TrackJetEmulationEtaPhiBin * clusters
unsigned int tidRing(const DetId &id) const
const EDGetTokenT< vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Log< level::Warning, false > LogWarning
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511

◆ trackQualityCuts()

bool L1TrackJetEmulationProducer::trackQualityCuts ( int  trk_nstub,
float  trk_chi2,
float  trk_bendchi2 
)

Definition at line 660 of file L1TrackJetEmulationProducer.cc.

References displaced_, nStubs4DisplacedBend_, nStubs4DisplacedChi2_, nStubs5DisplacedBend_, and trkBendChi2Max_.

Referenced by produce().

660  {
661  bool PassQuality = false;
662  if (trk_bendchi2 < trkBendChi2Max_ && trk_chi2 < trkChi2dofMax_ && trk_nstub >= 4 && !displaced_)
663  PassQuality = true;
664  if (displaced_ && trk_bendchi2 < nStubs4DisplacedBend_ && trk_chi2 < nStubs4DisplacedChi2_ && trk_nstub == 4)
665  PassQuality = true;
666  if (displaced_ && trk_bendchi2 < nStubs5DisplacedBend_ && trk_chi2 < nStubs5DisplacedChi2_ && trk_nstub > 4)
667  PassQuality = true;
668  return PassQuality;
669 }

Member Data Documentation

◆ d0CutNStubs4_

float L1TrackJetEmulationProducer::d0CutNStubs4_
private

Definition at line 88 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ d0CutNStubs5_

float L1TrackJetEmulationProducer::d0CutNStubs5_
private

Definition at line 89 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ displaced_

bool L1TrackJetEmulationProducer::displaced_
private

◆ etaBins_

int L1TrackJetEmulationProducer::etaBins_
private

◆ etaStep_

glbeta_intern L1TrackJetEmulationProducer::etaStep_
private

◆ highpTJetMinpT_

float L1TrackJetEmulationProducer::highpTJetMinpT_
private

Definition at line 93 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ highpTJetMinTrackMultiplicity_

int L1TrackJetEmulationProducer::highpTJetMinTrackMultiplicity_
private

Definition at line 92 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ L1TrkPtrs_

vector<Ptr<L1TTTrackType> > L1TrackJetEmulationProducer::L1TrkPtrs_
private

Definition at line 73 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster(), and produce().

◆ lowpTJetMinpT_

float L1TrackJetEmulationProducer::lowpTJetMinpT_
private

Definition at line 91 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ lowpTJetMinTrackMultiplicity_

int L1TrackJetEmulationProducer::lowpTJetMinTrackMultiplicity_
private

Definition at line 90 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ MaxDzTrackPV

const float L1TrackJetEmulationProducer::MaxDzTrackPV
private

Definition at line 76 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ minTrkJetpT_

const double L1TrackJetEmulationProducer::minTrkJetpT_
private

Definition at line 84 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ nDisplacedTracks_

int L1TrackJetEmulationProducer::nDisplacedTracks_
private

Definition at line 99 of file L1TrackJetEmulationProducer.cc.

◆ nStubs4DisplacedBend_

float L1TrackJetEmulationProducer::nStubs4DisplacedBend_
private

Definition at line 97 of file L1TrackJetEmulationProducer.cc.

Referenced by trackQualityCuts().

◆ nStubs4DisplacedChi2_

float L1TrackJetEmulationProducer::nStubs4DisplacedChi2_
private

Definition at line 95 of file L1TrackJetEmulationProducer.cc.

Referenced by trackQualityCuts().

◆ nStubs5DisplacedBend_

float L1TrackJetEmulationProducer::nStubs5DisplacedBend_
private

Definition at line 98 of file L1TrackJetEmulationProducer.cc.

Referenced by trackQualityCuts().

◆ nStubs5DisplacedChi2_

float L1TrackJetEmulationProducer::nStubs5DisplacedChi2_
private

Definition at line 96 of file L1TrackJetEmulationProducer.cc.

◆ phiBins_

int L1TrackJetEmulationProducer::phiBins_
private

Definition at line 86 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), and L2_cluster().

◆ phiStep_

glbphi_intern L1TrackJetEmulationProducer::phiStep_
private

Definition at line 104 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), and L2_cluster().

◆ PVtxToken_

const EDGetTokenT<l1t::VertexWordCollection> L1TrackJetEmulationProducer::PVtxToken_
private

Definition at line 107 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ PVz

float L1TrackJetEmulationProducer::PVz
private

Definition at line 101 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ tdtrk_

vector<int> L1TrackJetEmulationProducer::tdtrk_
private

Definition at line 75 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ trackToken_

const EDGetTokenT<vector<TTTrack<Ref_Phase2TrackerDigi_> > > L1TrackJetEmulationProducer::trackToken_
private

Definition at line 106 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ trkBendChi2Max_

const float L1TrackJetEmulationProducer::trkBendChi2Max_
private

Definition at line 82 of file L1TrackJetEmulationProducer.cc.

Referenced by trackQualityCuts().

◆ trkChi2dofMax_

const float L1TrackJetEmulationProducer::trkChi2dofMax_
private

Definition at line 81 of file L1TrackJetEmulationProducer.cc.

◆ trkEtaMax_

const float L1TrackJetEmulationProducer::trkEtaMax_
private

Definition at line 80 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), L2_cluster(), and produce().

◆ trkNPSStubMin_

const int L1TrackJetEmulationProducer::trkNPSStubMin_
private

Definition at line 83 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ trkPtMax_

const float L1TrackJetEmulationProducer::trkPtMax_
private

Definition at line 78 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster().

◆ trkPtMin_

const float L1TrackJetEmulationProducer::trkPtMin_
private

Definition at line 79 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ trkZMax_

const float L1TrackJetEmulationProducer::trkZMax_
private

Definition at line 77 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), L2_cluster(), and produce().

◆ tTopoToken_

ESGetToken<TrackerTopology, TrackerTopologyRcd> L1TrackJetEmulationProducer::tTopoToken_
private

Definition at line 108 of file L1TrackJetEmulationProducer.cc.

Referenced by produce().

◆ zBinCount_

vector<int> L1TrackJetEmulationProducer::zBinCount_
private

Definition at line 74 of file L1TrackJetEmulationProducer.cc.

Referenced by L2_cluster(), and produce().

◆ zBins_

int L1TrackJetEmulationProducer::zBins_
private

Definition at line 87 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), and L2_cluster().

◆ zStep_

z0_intern L1TrackJetEmulationProducer::zStep_
private

Definition at line 102 of file L1TrackJetEmulationProducer.cc.

Referenced by L1TrackJetEmulationProducer(), and L2_cluster().