CMS 3D CMS Logo

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

Public Types

typedef vector< L1TTTrackTypeL1TTTrackCollectionType
 
typedef edm::RefVector< L1TTTrackCollectionTypeL1TTTrackRefCollectionType
 
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

 L1TrackJetProducer (const ParameterSet &)
 
 ~L1TrackJetProducer () override=default
 
- 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 produce (Event &, const EventSetup &) override
 

Private Attributes

const float d0CutNStubs4_
 
const float d0CutNStubs5_
 
const bool displaced_
 
const int etaBins_
 
float etaStep_
 
const int highpTJetMinTrackMultiplicity_
 
const float highpTJetThreshold_
 
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
 
const int lowpTJetMinTrackMultiplicity_
 
const float lowpTJetThreshold_
 
const double minTrkJetpT_
 
const int nDisplacedTracks_
 
const int phiBins_
 
float phiStep_
 
const EDGetTokenT< L1TTTrackRefCollectionTypetrackToken_
 
const float trkEtaMax_
 
const float trkPtMax_
 
const float trkZMax_
 
const int zBins_
 
float zStep_
 

Detailed Description

Definition at line 46 of file L1TrackJetProducer.cc.

Member Typedef Documentation

◆ L1TTTrackCollectionType

Definition at line 51 of file L1TrackJetProducer.cc.

◆ L1TTTrackRefCollectionType

Definition at line 52 of file L1TrackJetProducer.cc.

◆ L1TTTrackType

Definition at line 50 of file L1TrackJetProducer.cc.

Constructor & Destructor Documentation

◆ L1TrackJetProducer()

L1TrackJetProducer::L1TrackJetProducer ( const ParameterSet iConfig)
explicit

Definition at line 83 of file L1TrackJetProducer.cc.

References displaced_, etaBins_, etaStep_, M_PI, phiBins_, phiStep_, trkEtaMax_, trkZMax_, zBins_, and zStep_.

Referenced by produce().

84  : trkZMax_(iConfig.getParameter<double>("trk_zMax")),
85  trkPtMax_(iConfig.getParameter<double>("trk_ptMax")),
86  trkEtaMax_(iConfig.getParameter<double>("trk_etaMax")),
87  lowpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("lowpTJetMinTrackMultiplicity")),
88  lowpTJetThreshold_(iConfig.getParameter<double>("lowpTJetThreshold")),
89  highpTJetMinTrackMultiplicity_(iConfig.getParameter<int>("highpTJetMinTrackMultiplicity")),
90  highpTJetThreshold_(iConfig.getParameter<double>("highpTJetThreshold")),
91  zBins_(iConfig.getParameter<int>("zBins")),
92  etaBins_(iConfig.getParameter<int>("etaBins")),
93  phiBins_(iConfig.getParameter<int>("phiBins")),
94  minTrkJetpT_(iConfig.getParameter<double>("minTrkJetpT")),
95  displaced_(iConfig.getParameter<bool>("displaced")),
96  d0CutNStubs4_(iConfig.getParameter<double>("d0_cutNStubs4")),
97  d0CutNStubs5_(iConfig.getParameter<double>("d0_cutNStubs5")),
98  nDisplacedTracks_(iConfig.getParameter<int>("nDisplacedTracks")),
99  trackToken_(consumes<L1TTTrackRefCollectionType>(iConfig.getParameter<InputTag>("L1TrackInputTag"))) {
100  zStep_ = 2.0 * trkZMax_ / (zBins_ + 1); // added +1 in denom
101  etaStep_ = 2.0 * trkEtaMax_ / etaBins_; //etaStep is the width of an etabin
102  phiStep_ = 2 * M_PI / phiBins_;
103 
104  if (displaced_)
105  produces<TkJetCollection>("L1TrackJetsExtended");
106  else
107  produces<TkJetCollection>("L1TrackJets");
108 }
const int highpTJetMinTrackMultiplicity_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const float lowpTJetThreshold_
const float highpTJetThreshold_
const EDGetTokenT< L1TTTrackRefCollectionType > trackToken_
#define M_PI
const int lowpTJetMinTrackMultiplicity_

◆ ~L1TrackJetProducer()

L1TrackJetProducer::~L1TrackJetProducer ( )
overridedefault

Member Function Documentation

◆ fillDescriptions()

void L1TrackJetProducer::fillDescriptions ( ConfigurationDescriptions descriptions)
static

Definition at line 283 of file L1TrackJetProducer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

283  {
285  desc.add<edm::InputTag>(
286  "L1TrackInputTag", edm::InputTag("l1tTrackVertexAssociationProducerForJets", "Level1TTTracksSelectedAssociated"));
287  desc.add<double>("trk_zMax", 15.0);
288  desc.add<double>("trk_ptMax", 200.0);
289  desc.add<double>("trk_etaMax", 2.4);
290  desc.add<double>("minTrkJetpT", -1.0);
291  desc.add<int>("etaBins", 24);
292  desc.add<int>("phiBins", 27);
293  desc.add<int>("zBins", 1);
294  desc.add<double>("d0_cutNStubs4", -1);
295  desc.add<double>("d0_cutNStubs5", -1);
296  desc.add<int>("lowpTJetMinTrackMultiplicity", 2);
297  desc.add<double>("lowpTJetThreshold", 50.0);
298  desc.add<int>("highpTJetMinTrackMultiplicity", 3);
299  desc.add<double>("highpTJetThreshold", 100.0);
300  desc.add<bool>("displaced", false);
301  desc.add<int>("nDisplacedTracks", 2);
302  descriptions.add("l1tTrackJets", desc);
303 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

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

Definition at line 110 of file L1TrackJetProducer.cc.

References a, funct::abs(), b, l1ttrackjet::MaxZBin::clusters, filterCSVwithJSON::copy, funct::cos(), d0CutNStubs4_, d0CutNStubs5_, displaced_, l1ttrackjet::EtaPhiBin::eta, PVValHelper::eta, egammaIdentification::eta_max, egammaIdentification::eta_min, etaBins_, etaStep_, highpTJetMinTrackMultiplicity_, highpTJetThreshold_, l1ttrackjet::MaxZBin::ht, mps_fire::i, iEvent, createfilelist::int, dqmiolumiharvest::j, reco::btau::jetEta, reco::btau::jetPhi, reco::btau::jetPt, isotrackApplyRegressor::k, L1TrackJetProducer(), L1TrkPtrs_, lowpTJetMinTrackMultiplicity_, lowpTJetThreshold_, M_PI, minTrkJetpT_, eostools::move(), l1ttrackjet::MaxZBin::nclust, nDisplacedTracks_, vertices_cff::ntracks, l1ttrackjet::EtaPhiBin::ntracks, l1ttrackjet::EtaPhiBin::nxtracks, phi, l1ttrackjet::EtaPhiBin::phi, phiBins_, phiStep_, l1ttrackjet::EtaPhiBin::pTtot, funct::sin(), edm::RefVector< C, T, F >::size(), jetUpdater_cfi::sort, l1ttrackjet::EtaPhiBin::trackidx, trackToken_, trkEtaMax_, trkPtMax_, trkZMax_, l1ttrackjet::EtaPhiBin::used, l1ttrackjet::MaxZBin::zbincenter, zBins_, SiStripMonitorCluster_cfi::zmax, SiStripMonitorCluster_cfi::zmin, l1ttrackjet::MaxZBin::znum, and zStep_.

110  {
111  unique_ptr<TkJetCollection> L1TrackJetProducer(new TkJetCollection);
112 
113  // L1 tracks
115  iEvent.getByToken(trackToken_, TTTrackHandle);
116 
117  L1TrkPtrs_.clear();
118 
119  // track selection
120  for (unsigned int this_l1track = 0; this_l1track < TTTrackHandle->size(); this_l1track++) {
121  edm::Ptr<L1TTTrackType> trkPtr(TTTrackHandle, this_l1track);
122  L1TrkPtrs_.push_back(trkPtr);
123  }
124 
125  // if no tracks pass selection return empty containers
126  if (L1TrkPtrs_.empty()) {
127  if (displaced_)
128  iEvent.put(std::move(L1TrackJetProducer), "L1TrackJetsExtended");
129  else
130  iEvent.put(std::move(L1TrackJetProducer), "L1TrackJets");
131  return;
132  }
133 
134  MaxZBin mzb;
135  mzb.znum = -1;
136  mzb.nclust = -1;
137  mzb.ht = -1;
138 
139  // create 2D grid of eta phi bins
140  EtaPhiBin epbins_default[phiBins_][etaBins_];
141  float phi = -1.0 * M_PI;
142  for (int i = 0; i < phiBins_; ++i) {
143  float eta = -1.0 * trkEtaMax_;
144  for (int j = 0; j < etaBins_; ++j) {
145  epbins_default[i][j].phi = (phi + (phi + phiStep_)) / 2.0;
146  epbins_default[i][j].eta = (eta + (eta + etaStep_)) / 2.0;
147  eta += etaStep_;
148  } // for each etabin
149  phi += phiStep_;
150  } // for each phibin (finished creating bins)
151 
152  // create z-bins (might be useful for displaced if we run w/o dz cut)
153  std::vector<float> zmins, zmaxs;
154  for (int zbin = 0; zbin < zBins_; zbin++) {
155  zmins.push_back(-1.0 * trkZMax_ + zStep_ * zbin);
156  zmaxs.push_back(-1.0 * trkZMax_ + zStep_ * zbin + zStep_ * 2);
157  }
158 
159  // create vectors of clusters
160  std::vector<std::vector<EtaPhiBin>> L1clusters;
161  L1clusters.reserve(phiBins_);
162  std::vector<EtaPhiBin> L2clusters;
163 
164  // default (empty) grid
165  for (int i = 0; i < phiBins_; ++i) {
166  for (int j = 0; j < etaBins_; ++j) {
167  epbins_default[i][j].pTtot = 0;
168  epbins_default[i][j].used = false;
169  epbins_default[i][j].ntracks = 0;
170  epbins_default[i][j].nxtracks = 0;
171  epbins_default[i][j].trackidx.clear();
172  }
173  }
174 
175  for (unsigned int zbin = 0; zbin < zmins.size(); ++zbin) {
176  // initialize grid
177  float zmin = zmins[zbin];
178  float zmax = zmaxs[zbin];
179  EtaPhiBin epbins[phiBins_][etaBins_];
180  std::copy(&epbins_default[0][0], &epbins_default[0][0] + phiBins_ * etaBins_, &epbins[0][0]);
181 
182  //clear cluster containers
183  L1clusters.clear();
184  L2clusters.clear();
185 
186  // fill grid with tracks
187  for (unsigned int k = 0; k < L1TrkPtrs_.size(); ++k) {
188  float trkZ = L1TrkPtrs_[k]->z0();
189  if (zmax < trkZ)
190  continue;
191  if (zmin > trkZ)
192  continue;
193  if (zbin == 0 && zmin == trkZ)
194  continue;
195  float trkpt = L1TrkPtrs_[k]->momentum().perp();
196  float trketa = L1TrkPtrs_[k]->momentum().eta();
197  float trkphi = L1TrkPtrs_[k]->momentum().phi();
198  float trkd0 = L1TrkPtrs_[k]->d0();
199  int trknstubs = (int)L1TrkPtrs_[k]->getStubRefs().size();
200  for (int i = 0; i < phiBins_; ++i) {
201  for (int j = 0; j < etaBins_; ++j) {
202  float eta_min = epbins[i][j].eta - etaStep_ / 2.0; //eta min
203  float eta_max = epbins[i][j].eta + etaStep_ / 2.0; //eta max
204  float phi_min = epbins[i][j].phi - phiStep_ / 2.0; //phi min
205  float phi_max = epbins[i][j].phi + phiStep_ / 2.0; //phi max
206 
207  if ((trketa < eta_min) || (trketa > eta_max) || (trkphi < phi_min) || (trkphi > phi_max))
208  continue;
209 
210  if (trkpt < trkPtMax_)
211  epbins[i][j].pTtot += trkpt;
212  else
213  epbins[i][j].pTtot += trkPtMax_;
214  if ((std::abs(trkd0) > d0CutNStubs5_ && trknstubs >= 5 && d0CutNStubs5_ >= 0) ||
215  (trknstubs == 4 && std::abs(trkd0) > d0CutNStubs4_ && d0CutNStubs4_ >= 0))
216  epbins[i][j].nxtracks += 1;
217  epbins[i][j].trackidx.push_back(k);
218  ++epbins[i][j].ntracks;
219  } // for each etabin
220  } // for each phibin
221  } //end loop over tracks
222 
223  // cluster tracks in eta (first layer) using grid
224  for (int phibin = 0; phibin < phiBins_; ++phibin) {
225  L1clusters.push_back(L1_clustering<EtaPhiBin, float, float, float>(epbins[phibin], etaBins_, etaStep_));
226  }
227 
228  // second layer clustering in phi for using eta clusters
229  L2clusters = L2_clustering<EtaPhiBin, float, float, float>(L1clusters, phiBins_, phiStep_, etaStep_);
230 
231  // sum all cluster pTs in this zbin to find max
232  float sum_pt = 0;
233  for (unsigned int k = 0; k < L2clusters.size(); ++k) {
234  if (L2clusters[k].pTtot > lowpTJetThreshold_ && L2clusters[k].ntracks < lowpTJetMinTrackMultiplicity_)
235  continue;
236  if (L2clusters[k].pTtot > highpTJetThreshold_ && L2clusters[k].ntracks < highpTJetMinTrackMultiplicity_)
237  continue;
238  if (L2clusters[k].pTtot > minTrkJetpT_)
239  sum_pt += L2clusters[k].pTtot;
240  }
241 
242  if (sum_pt < mzb.ht)
243  continue;
244 
245  // if ht is larger than previous max, this is the new vertex zbin
246  mzb.ht = sum_pt;
247  mzb.znum = zbin;
248  mzb.clusters = L2clusters;
249  mzb.nclust = L2clusters.size();
250  mzb.zbincenter = (zmin + zmax) / 2.0;
251  } //zbin loop
252 
253  // output
254  vector<Ptr<L1TTTrackType>> L1TrackAssocJet;
255  for (unsigned int j = 0; j < mzb.clusters.size(); ++j) {
256  float jetEta = mzb.clusters[j].eta;
257  float jetPhi = mzb.clusters[j].phi;
258  float jetPt = mzb.clusters[j].pTtot;
259  float jetPx = jetPt * cos(jetPhi);
260  float jetPy = jetPt * sin(jetPhi);
261  float jetPz = jetPt * sinh(jetEta);
262  float jetP = jetPt * cosh(jetEta);
263  int totalDisptrk = mzb.clusters[j].nxtracks;
264  bool isDispJet = (totalDisptrk >= nDisplacedTracks_);
265 
266  math::XYZTLorentzVector jetP4(jetPx, jetPy, jetPz, jetP);
267  L1TrackAssocJet.clear();
268  for (unsigned int itrk = 0; itrk < mzb.clusters[j].trackidx.size(); itrk++)
269  L1TrackAssocJet.push_back(L1TrkPtrs_[mzb.clusters[j].trackidx[itrk]]);
270 
271  TkJet trkJet(jetP4, L1TrackAssocJet, mzb.zbincenter, mzb.clusters[j].ntracks, 0, totalDisptrk, 0, isDispJet);
272 
273  L1TrackJetProducer->push_back(trkJet);
274  }
275 
276  std::sort(L1TrackJetProducer->begin(), L1TrackJetProducer->end(), [](auto &a, auto &b) { return a.pt() > b.pt(); });
277  if (displaced_)
278  iEvent.put(std::move(L1TrackJetProducer), "L1TrackJetsExtended");
279  else
280  iEvent.put(std::move(L1TrackJetProducer), "L1TrackJets");
281 }
const int highpTJetMinTrackMultiplicity_
const float lowpTJetThreshold_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< unsigned int > trackidx
L1TrackJetProducer(const ParameterSet &)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
const float highpTJetThreshold_
vector< Ptr< L1TTTrackType > > L1TrkPtrs_
std::vector< EtaPhiBin > clusters
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const EDGetTokenT< L1TTTrackRefCollectionType > trackToken_
#define M_PI
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
double b
Definition: hdecay.h:120
const int lowpTJetMinTrackMultiplicity_
double a
Definition: hdecay.h:121
std::vector< TkJet > TkJetCollection
Definition: TkJetFwd.h:20
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ d0CutNStubs4_

const float L1TrackJetProducer::d0CutNStubs4_
private

Definition at line 76 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ d0CutNStubs5_

const float L1TrackJetProducer::d0CutNStubs5_
private

Definition at line 77 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ displaced_

const bool L1TrackJetProducer::displaced_
private

Definition at line 75 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ etaBins_

const int L1TrackJetProducer::etaBins_
private

Definition at line 69 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ etaStep_

float L1TrackJetProducer::etaStep_
private

Definition at line 73 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ highpTJetMinTrackMultiplicity_

const int L1TrackJetProducer::highpTJetMinTrackMultiplicity_
private

Definition at line 66 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ highpTJetThreshold_

const float L1TrackJetProducer::highpTJetThreshold_
private

Definition at line 67 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ L1TrkPtrs_

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

Definition at line 60 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ lowpTJetMinTrackMultiplicity_

const int L1TrackJetProducer::lowpTJetMinTrackMultiplicity_
private

Definition at line 64 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ lowpTJetThreshold_

const float L1TrackJetProducer::lowpTJetThreshold_
private

Definition at line 65 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ minTrkJetpT_

const double L1TrackJetProducer::minTrkJetpT_
private

Definition at line 71 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ nDisplacedTracks_

const int L1TrackJetProducer::nDisplacedTracks_
private

Definition at line 78 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ phiBins_

const int L1TrackJetProducer::phiBins_
private

Definition at line 70 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ phiStep_

float L1TrackJetProducer::phiStep_
private

Definition at line 74 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ trackToken_

const EDGetTokenT<L1TTTrackRefCollectionType> L1TrackJetProducer::trackToken_
private

Definition at line 80 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ trkEtaMax_

const float L1TrackJetProducer::trkEtaMax_
private

Definition at line 63 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ trkPtMax_

const float L1TrackJetProducer::trkPtMax_
private

Definition at line 62 of file L1TrackJetProducer.cc.

Referenced by produce().

◆ trkZMax_

const float L1TrackJetProducer::trkZMax_
private

Definition at line 61 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ zBins_

const int L1TrackJetProducer::zBins_
private

Definition at line 68 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().

◆ zStep_

float L1TrackJetProducer::zStep_
private

Definition at line 72 of file L1TrackJetProducer.cc.

Referenced by L1TrackJetProducer(), and produce().