CMS 3D CMS Logo

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

Public Types

typedef std::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

 L1TkElectronTrackProducer (const edm::ParameterSet &)
 
 ~L1TkElectronTrackProducer () 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

double getPtScaledCut (double pt, std::vector< double > &parameters)
 
float isolation (const edm::Handle< L1TTTrackCollectionType > &trkHandle, int match_index)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
bool selectMatchedTrack (double &d_r, double &d_phi, double &d_eta, double &tk_pt, float &eg_eta)
 

Private Attributes

float deltaZ_
 
std::vector< double > dEtaCutoff_
 
std::vector< double > dPhiCutoff_
 
std::vector< double > dRCutoff_
 
const float dRMax_
 
const float dRMin_
 
const edm::EDGetTokenT< EGammaBxCollectionegToken_
 
const float etMin_
 
float isoCut_
 
std::string label
 
std::string matchType_
 
const float maxChi2IsoTracks_
 
const unsigned int minNStubsIsoTracks_
 
bool primaryVtxConstrain_
 
const float pTMinTra_
 
bool relativeIsolation_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtGeomToken_
 
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
 
float trkQualityChi2_
 
float trkQualityPtMin_
 
bool useClusterET_
 
bool useTwoStubsPT_
 

Detailed Description

Definition at line 53 of file L1TkElectronTrackProducer.cc.

Member Typedef Documentation

◆ L1TTTrackCollectionType

Definition at line 56 of file L1TkElectronTrackProducer.cc.

◆ L1TTTrackType

Definition at line 55 of file L1TkElectronTrackProducer.cc.

Constructor & Destructor Documentation

◆ L1TkElectronTrackProducer()

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

Definition at line 103 of file L1TkElectronTrackProducer.cc.

References label.

104  : label(iConfig.getParameter<std::string>("label")),
105  etMin_((float)iConfig.getParameter<double>("ETmin")),
106  pTMinTra_((float)iConfig.getParameter<double>("PTMINTRA")),
107  dRMin_((float)iConfig.getParameter<double>("DRmin")),
108  dRMax_((float)iConfig.getParameter<double>("DRmax")),
109  deltaZ_((float)iConfig.getParameter<double>("DeltaZ")),
110  maxChi2IsoTracks_(iConfig.getParameter<double>("maxChi2IsoTracks")),
111  minNStubsIsoTracks_(iConfig.getParameter<int>("minNStubsIsoTracks")),
112  isoCut_((float)iConfig.getParameter<double>("IsoCut")),
113  relativeIsolation_(iConfig.getParameter<bool>("RelativeIsolation")),
114  trkQualityChi2_((float)iConfig.getParameter<double>("TrackChi2")),
115  trkQualityPtMin_((float)iConfig.getParameter<double>("TrackMinPt")),
116  useTwoStubsPT_(iConfig.getParameter<bool>("useTwoStubsPT")),
117  useClusterET_(iConfig.getParameter<bool>("useClusterET")),
118  dPhiCutoff_(iConfig.getParameter<std::vector<double>>("TrackEGammaDeltaPhi")),
119  dRCutoff_(iConfig.getParameter<std::vector<double>>("TrackEGammaDeltaR")),
120  dEtaCutoff_(iConfig.getParameter<std::vector<double>>("TrackEGammaDeltaEta")),
121  matchType_(iConfig.getParameter<std::string>("TrackEGammaMatchType")),
122  egToken_(consumes<EGammaBxCollection>(iConfig.getParameter<edm::InputTag>("L1EGammaInputTag"))),
124  iConfig.getParameter<edm::InputTag>("L1TrackInputTag"))),
125  tGeomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()) {
126  produces<TkElectronCollection>(label);
127 }
const edm::EDGetTokenT< EGammaBxCollection > egToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomToken_
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
Class to store the L1 Track Trigger tracks.
Definition: TTTrack.h:29

◆ ~L1TkElectronTrackProducer()

L1TkElectronTrackProducer::~L1TkElectronTrackProducer ( )
override

Definition at line 129 of file L1TkElectronTrackProducer.cc.

129 {}

Member Function Documentation

◆ fillDescriptions()

void L1TkElectronTrackProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 274 of file L1TkElectronTrackProducer.cc.

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

274  {
275  {
276  // L1TkElectrons
278  desc.add<std::string>("label", "EG");
279  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("simCaloStage2Digis"));
280  desc.add<double>("ETmin", -1.0);
281  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
282  desc.add<double>("TrackChi2", 10000000000.0);
283  desc.add<double>("TrackMinPt", 10.0);
284  desc.add<bool>("useTwoStubsPT", false);
285  desc.add<bool>("useClusterET", false);
286  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
287  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
288  {
289  0.07,
290  0.0,
291  0.0,
292  });
293  desc.add<std::vector<double>>("TrackEGammaDeltaR",
294  {
295  0.08,
296  0.0,
297  0.0,
298  });
299  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
300  {
301  10000000000.0,
302  0.0,
303  0.0,
304  });
305  desc.add<bool>("RelativeIsolation", true);
306  desc.add<double>("IsoCut", -0.1);
307  desc.add<double>("PTMINTRA", 2.0);
308  desc.add<double>("DRmin", 0.03);
309  desc.add<double>("DRmax", 0.2);
310  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
311  desc.add<int>("minNStubsIsoTracks", 0);
312  desc.add<double>("DeltaZ", 0.6);
313  descriptions.add("L1TkElectrons", desc);
314  // or use the following to generate the label from the module's C++ type
315  //descriptions.addWithDefaultLabel(desc);
316  }
317  {
318  // L1TkIsoElectrons
320  desc.add<std::string>("label", "EG");
321  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("simCaloStage2Digis"));
322  desc.add<double>("ETmin", -1.0);
323  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
324  desc.add<double>("TrackChi2", 10000000000.0);
325  desc.add<double>("TrackMinPt", 10.0);
326  desc.add<bool>("useTwoStubsPT", false);
327  desc.add<bool>("useClusterET", false);
328  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
329  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
330  {
331  0.07,
332  0.0,
333  0.0,
334  });
335  desc.add<std::vector<double>>("TrackEGammaDeltaR",
336  {
337  0.08,
338  0.0,
339  0.0,
340  });
341  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
342  {
343  10000000000.0,
344  0.0,
345  0.0,
346  });
347  desc.add<bool>("RelativeIsolation", true);
348  desc.add<double>("IsoCut", 0.1);
349  desc.add<double>("PTMINTRA", 2.0);
350  desc.add<double>("DRmin", 0.03);
351  desc.add<double>("DRmax", 0.2);
352  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
353  desc.add<int>("minNStubsIsoTracks", 0);
354  desc.add<double>("DeltaZ", 0.6);
355  descriptions.add("L1TkIsoElectrons", desc);
356  // or use the following to generate the label from the module's C++ type
357  //descriptions.addWithDefaultLabel(desc);
358  }
359  {
360  // L1TkElectronsLoose
362  desc.add<std::string>("label", "EG");
363  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("simCaloStage2Digis"));
364  desc.add<double>("ETmin", -1.0);
365  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
366  desc.add<double>("TrackChi2", 10000000000.0);
367  desc.add<double>("TrackMinPt", 3.0);
368  desc.add<bool>("useTwoStubsPT", false);
369  desc.add<bool>("useClusterET", false);
370  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
371  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
372  {
373  0.07,
374  0.0,
375  0.0,
376  });
377  desc.add<std::vector<double>>("TrackEGammaDeltaR",
378  {
379  0.12,
380  0.0,
381  0.0,
382  });
383  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
384  {
385  10000000000.0,
386  0.0,
387  0.0,
388  });
389  desc.add<bool>("RelativeIsolation", true);
390  desc.add<double>("IsoCut", -0.1);
391  desc.add<double>("PTMINTRA", 2.0);
392  desc.add<double>("DRmin", 0.03);
393  desc.add<double>("DRmax", 0.2);
394  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
395  desc.add<int>("minNStubsIsoTracks", 0);
396  desc.add<double>("DeltaZ", 0.6);
397  descriptions.add("L1TkElectronsLoose", desc);
398  // or use the following to generate the label from the module's C++ type
399  //descriptions.addWithDefaultLabel(desc);
400  }
401  {
402  // L1TkElectronsCrystal
404  desc.add<std::string>("label", "EG");
405  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("L1EGammaClusterEmuProducer"));
406  desc.add<double>("ETmin", -1.0);
407  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
408  desc.add<double>("TrackChi2", 10000000000.0);
409  desc.add<double>("TrackMinPt", 10.0);
410  desc.add<bool>("useTwoStubsPT", false);
411  desc.add<bool>("useClusterET", false);
412  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
413  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
414  {
415  0.07,
416  0.0,
417  0.0,
418  });
419  desc.add<std::vector<double>>("TrackEGammaDeltaR",
420  {
421  0.08,
422  0.0,
423  0.0,
424  });
425  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
426  {
427  10000000000.0,
428  0.0,
429  0.0,
430  });
431  desc.add<bool>("RelativeIsolation", true);
432  desc.add<double>("IsoCut", -0.1);
433  desc.add<double>("PTMINTRA", 2.0);
434  desc.add<double>("DRmin", 0.03);
435  desc.add<double>("DRmax", 0.2);
436  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
437  desc.add<int>("minNStubsIsoTracks", 0);
438  desc.add<double>("DeltaZ", 0.6);
439  descriptions.add("L1TkElectronsCrystal", desc);
440  // or use the following to generate the label from the module's C++ type
441  //descriptions.addWithDefaultLabel(desc);
442  }
443  {
444  // L1TkIsoElectronsCrystal
446  desc.add<std::string>("label", "EG");
447  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("L1EGammaClusterEmuProducer"));
448  desc.add<double>("ETmin", -1.0);
449  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
450  desc.add<double>("TrackChi2", 10000000000.0);
451  desc.add<double>("TrackMinPt", 10.0);
452  desc.add<bool>("useTwoStubsPT", false);
453  desc.add<bool>("useClusterET", false);
454  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
455  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
456  {
457  0.07,
458  0.0,
459  0.0,
460  });
461  desc.add<std::vector<double>>("TrackEGammaDeltaR",
462  {
463  0.08,
464  0.0,
465  0.0,
466  });
467  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
468  {
469  10000000000.0,
470  0.0,
471  0.0,
472  });
473  desc.add<bool>("RelativeIsolation", true);
474  desc.add<double>("IsoCut", 0.1);
475  desc.add<double>("PTMINTRA", 2.0);
476  desc.add<double>("DRmin", 0.03);
477  desc.add<double>("DRmax", 0.2);
478  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
479  desc.add<int>("minNStubsIsoTracks", 0);
480  desc.add<double>("DeltaZ", 0.6);
481  descriptions.add("L1TkIsoElectronsCrystal", desc);
482  // or use the following to generate the label from the module's C++ type
483  //descriptions.addWithDefaultLabel(desc);
484  }
485  {
486  // L1TkElectronsLooseCrystal
488  desc.add<std::string>("label", "EG");
489  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("L1EGammaClusterEmuProducer"));
490  desc.add<double>("ETmin", -1.0);
491  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
492  desc.add<double>("TrackChi2", 10000000000.0);
493  desc.add<double>("TrackMinPt", 3.0);
494  desc.add<bool>("useTwoStubsPT", false);
495  desc.add<bool>("useClusterET", false);
496  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
497  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
498  {
499  0.07,
500  0.0,
501  0.0,
502  });
503  desc.add<std::vector<double>>("TrackEGammaDeltaR",
504  {
505  0.12,
506  0.0,
507  0.0,
508  });
509  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
510  {
511  10000000000.0,
512  0.0,
513  0.0,
514  });
515  desc.add<bool>("RelativeIsolation", true);
516  desc.add<double>("IsoCut", -0.1);
517  desc.add<double>("PTMINTRA", 2.0);
518  desc.add<double>("DRmin", 0.03);
519  desc.add<double>("DRmax", 0.2);
520  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
521  desc.add<int>("minNStubsIsoTracks", 0);
522  desc.add<double>("DeltaZ", 0.6);
523  descriptions.add("L1TkElectronsLooseCrystal", desc);
524  // or use the following to generate the label from the module's C++ type
525  //descriptions.addWithDefaultLabel(desc);
526  }
527  {
528  // L1TkElectronsEllipticMatchCrystal
530  desc.add<std::string>("label", "EG");
531  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("L1EGammaClusterEmuProducer"));
532  desc.add<double>("ETmin", -1.0);
533  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
534  desc.add<double>("TrackChi2", 10000000000.0);
535  desc.add<double>("TrackMinPt", 10.0);
536  desc.add<bool>("useTwoStubsPT", false);
537  desc.add<bool>("useClusterET", false);
538  desc.add<std::string>("TrackEGammaMatchType", "EllipticalCut");
539  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
540  {
541  0.07,
542  0.0,
543  0.0,
544  });
545  desc.add<std::vector<double>>("TrackEGammaDeltaR",
546  {
547  0.08,
548  0.0,
549  0.0,
550  });
551  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
552  {
553  0.015,
554  0.025,
555  10000000000.0,
556  });
557  desc.add<bool>("RelativeIsolation", true);
558  desc.add<double>("IsoCut", -0.1);
559  desc.add<double>("PTMINTRA", 2.0);
560  desc.add<double>("DRmin", 0.03);
561  desc.add<double>("DRmax", 0.2);
562  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
563  desc.add<int>("minNStubsIsoTracks", 0);
564  desc.add<double>("DeltaZ", 0.6);
565  descriptions.add("L1TkElectronsEllipticMatchCrystal", desc);
566  // or use the following to generate the label from the module's C++ type
567  //descriptions.addWithDefaultLabel(desc);
568  }
569  {
570  // L1TkElectronsHGC
572  desc.add<std::string>("label", "EG");
573  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("l1EGammaEEProducer", "L1EGammaCollectionBXVWithCuts"));
574  desc.add<double>("ETmin", -1.0);
575  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
576  desc.add<double>("TrackChi2", 10000000000.0);
577  desc.add<double>("TrackMinPt", 10.0);
578  desc.add<bool>("useTwoStubsPT", false);
579  desc.add<bool>("useClusterET", false);
580  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
581  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
582  {
583  0.07,
584  0.0,
585  0.0,
586  });
587  desc.add<std::vector<double>>("TrackEGammaDeltaR",
588  {
589  0.08,
590  0.0,
591  0.0,
592  });
593  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
594  {
595  10000000000.0,
596  0.0,
597  0.0,
598  });
599  desc.add<bool>("RelativeIsolation", true);
600  desc.add<double>("IsoCut", -0.1);
601  desc.add<double>("PTMINTRA", 2.0);
602  desc.add<double>("DRmin", 0.03);
603  desc.add<double>("DRmax", 0.2);
604  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
605  desc.add<int>("minNStubsIsoTracks", 0);
606  desc.add<double>("DeltaZ", 0.6);
607  descriptions.add("L1TkElectronsHGC", desc);
608  // or use the following to generate the label from the module's C++ type
609  //descriptions.addWithDefaultLabel(desc);
610  }
611  {
612  // L1TkElectronsEllipticMatchHGC
614  desc.add<std::string>("label", "EG");
615  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("l1EGammaEEProducer", "L1EGammaCollectionBXVWithCuts"));
616  desc.add<double>("ETmin", -1.0);
617  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
618  desc.add<double>("TrackChi2", 10000000000.0);
619  desc.add<double>("TrackMinPt", 10.0);
620  desc.add<bool>("useTwoStubsPT", false);
621  desc.add<bool>("useClusterET", false);
622  desc.add<std::string>("TrackEGammaMatchType", "EllipticalCut");
623  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
624  {
625  0.07,
626  0.0,
627  0.0,
628  });
629  desc.add<std::vector<double>>("TrackEGammaDeltaR",
630  {
631  0.08,
632  0.0,
633  0.0,
634  });
635  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
636  {
637  0.01,
638  0.01,
639  10000000000.0,
640  });
641  desc.add<bool>("RelativeIsolation", true);
642  desc.add<double>("IsoCut", -0.1);
643  desc.add<double>("PTMINTRA", 2.0);
644  desc.add<double>("DRmin", 0.03);
645  desc.add<double>("DRmax", 0.2);
646  desc.add<double>("maxChi2IsoTracks", 100);
647  desc.add<int>("minNStubsIsoTracks", 4);
648  desc.add<double>("DeltaZ", 0.6);
649  descriptions.add("L1TkElectronsEllipticMatchHGC", desc);
650  // or use the following to generate the label from the module's C++ type
651  //descriptions.addWithDefaultLabel(desc);
652  }
653  {
654  // L1TkIsoElectronsHGC
656  desc.add<std::string>("label", "EG");
657  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("l1EGammaEEProducer", "L1EGammaCollectionBXVWithCuts"));
658  desc.add<double>("ETmin", -1.0);
659  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
660  desc.add<double>("TrackChi2", 10000000000.0);
661  desc.add<double>("TrackMinPt", 10.0);
662  desc.add<bool>("useTwoStubsPT", false);
663  desc.add<bool>("useClusterET", false);
664  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
665  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
666  {
667  0.07,
668  0.0,
669  0.0,
670  });
671  desc.add<std::vector<double>>("TrackEGammaDeltaR",
672  {
673  0.08,
674  0.0,
675  0.0,
676  });
677  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
678  {
679  10000000000.0,
680  0.0,
681  0.0,
682  });
683  desc.add<bool>("RelativeIsolation", true);
684  desc.add<double>("IsoCut", 0.1);
685  desc.add<double>("PTMINTRA", 2.0);
686  desc.add<double>("DRmin", 0.03);
687  desc.add<double>("DRmax", 0.4);
688  desc.add<double>("maxChi2IsoTracks", 100);
689  desc.add<int>("minNStubsIsoTracks", 4);
690  desc.add<double>("DeltaZ", 1.0);
691  descriptions.add("L1TkIsoElectronsHGC", desc);
692  // or use the following to generate the label from the module's C++ type
693  //descriptions.addWithDefaultLabel(desc);
694  }
695  {
696  // L1TkElectronsLooseHGC
698  desc.add<std::string>("label", "EG");
699  desc.add<edm::InputTag>("L1EGammaInputTag", edm::InputTag("l1EGammaEEProducer", "L1EGammaCollectionBXVWithCuts"));
700  desc.add<double>("ETmin", -1.0);
701  desc.add<edm::InputTag>("L1TrackInputTag", edm::InputTag("TTTracksFromTrackletEmulation", "Level1TTTracks"));
702  desc.add<double>("TrackChi2", 10000000000.0);
703  desc.add<double>("TrackMinPt", 3.0);
704  desc.add<bool>("useTwoStubsPT", false);
705  desc.add<bool>("useClusterET", false);
706  desc.add<std::string>("TrackEGammaMatchType", "PtDependentCut");
707  desc.add<std::vector<double>>("TrackEGammaDeltaPhi",
708  {
709  0.07,
710  0.0,
711  0.0,
712  });
713  desc.add<std::vector<double>>("TrackEGammaDeltaR",
714  {
715  0.12,
716  0.0,
717  0.0,
718  });
719  desc.add<std::vector<double>>("TrackEGammaDeltaEta",
720  {
721  10000000000.0,
722  0.0,
723  0.0,
724  });
725  desc.add<bool>("RelativeIsolation", true);
726  desc.add<double>("IsoCut", -0.1);
727  desc.add<double>("PTMINTRA", 2.0);
728  desc.add<double>("DRmin", 0.03);
729  desc.add<double>("DRmax", 0.2);
730  desc.add<double>("maxChi2IsoTracks", 10000000000.0);
731  desc.add<int>("minNStubsIsoTracks", 0);
732  desc.add<double>("DeltaZ", 0.6);
733  descriptions.add("L1TkElectronsLooseHGC", desc);
734  // or use the following to generate the label from the module's C++ type
735  //descriptions.addWithDefaultLabel(desc);
736  }
737 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ getPtScaledCut()

double L1TkElectronTrackProducer::getPtScaledCut ( double  pt,
std::vector< double > &  parameters 
)
private

◆ isolation()

float L1TkElectronTrackProducer::isolation ( const edm::Handle< L1TTTrackCollectionType > &  trkHandle,
int  match_index 
)
private

Definition at line 740 of file L1TkElectronTrackProducer.cc.

References funct::abs(), reco::deltaPhi(), HLT_2022v12_cff::dEta, HLT_2022v12_cff::dPhi, dRMax_, dRMin_, PV3DBase< T, PVType, FrameType >::eta(), maxChi2IsoTracks_, minNStubsIsoTracks_, TTTrack< T >::momentum(), PV3DBase< T, PVType, FrameType >::phi(), TTTrack< T >::POCA(), pTMinTra_, TtFullHadEvtBuilder_cfi::sumPt, and PV3DBase< T, PVType, FrameType >::z().

Referenced by produce().

740  {
741  edm::Ptr<L1TTTrackType> matchedTrkPtr(trkHandle, match_index);
742  L1TTTrackCollectionType::const_iterator trackIter;
743 
744  float sumPt = 0.0;
745  int itr = 0;
746 
747  float dRMin2 = dRMin_ * dRMin_;
748  float dRMax2 = dRMax_ * dRMax_;
749 
750  for (trackIter = trkHandle->begin(); trackIter != trkHandle->end(); ++trackIter) {
751  if (itr++ != match_index) {
752  if (trackIter->chi2() > maxChi2IsoTracks_ || trackIter->momentum().perp() <= pTMinTra_ ||
753  trackIter->getStubRefs().size() < minNStubsIsoTracks_) {
754  continue;
755  }
756 
757  float dZ = std::abs(trackIter->POCA().z() - matchedTrkPtr->POCA().z());
758 
759  float phi1 = trackIter->momentum().phi();
760  float phi2 = matchedTrkPtr->momentum().phi();
761  float dPhi = reco::deltaPhi(phi1, phi2);
762  float dEta = (trackIter->momentum().eta() - matchedTrkPtr->momentum().eta());
763  float dR2 = (dPhi * dPhi + dEta * dEta);
764 
765  if (dR2 > dRMin2 && dR2 < dRMax2 && dZ < deltaZ_ && trackIter->momentum().perp() > pTMinTra_) {
766  sumPt += trackIter->momentum().perp();
767  }
768  }
769  }
770 
771  return sumPt;
772 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ produce()

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

Definition at line 132 of file L1TkElectronTrackProducer.cc.

References BXVector< T >::begin(), HLT_2022v12_cff::dEta, L1TkElectronTrackMatchAlgo::doMatch(), L1TkElectronTrackMatchAlgo::doMatchClusterET(), HLT_2022v12_cff::dPhi, HGC3DClusterGenMatchSelector_cfi::dR, egToken_, BXVector< T >::end(), etMin_, Exception, edm::EventSetup::getData(), iEvent, isoCut_, isolation(), edm::HandleBase::isValid(), label, eostools::move(), edm::Handle< T >::product(), pTFrom2Stubs::pTFrom2(), relativeIsolation_, mps_fire::result, TTTrack< T >::rInv(), selectMatchedTrack(), l1t::TkElectron::setTrackCurvature(), tGeomToken_, trackToken_, trkQualityChi2_, trkQualityPtMin_, useClusterET_, and useTwoStubsPT_.

132  {
133  std::unique_ptr<TkElectronCollection> result(new TkElectronCollection);
134 
135  // geometry needed to call pTFrom2Stubs
136  const TrackerGeometry& tGeom = iSetup.getData(tGeomToken_);
137  const TrackerGeometry* tGeometry = &tGeom;
138 
139  // the L1EGamma objects
140  edm::Handle<EGammaBxCollection> eGammaHandle;
141  iEvent.getByToken(egToken_, eGammaHandle);
142  EGammaBxCollection eGammaCollection = (*eGammaHandle.product());
144 
145  // the L1Tracks
146  edm::Handle<L1TTTrackCollectionType> L1TTTrackHandle;
147  iEvent.getByToken(trackToken_, L1TTTrackHandle);
148  L1TTTrackCollectionType::const_iterator trackIter;
149 
150  if (!eGammaHandle.isValid()) {
151  throw cms::Exception("L1TkElectronTrackProducer")
152  << "\nWarning: L1EmCollection not found in the event. Exit" << std::endl;
153  return;
154  }
155  if (!L1TTTrackHandle.isValid()) {
156  throw cms::Exception("L1TkElectronTrackProducer")
157  << "\nWarning: L1TTTrackCollectionType not found in the event. Exit." << std::endl;
158  return;
159  }
160 
161  int ieg = 0;
162  for (egIter = eGammaCollection.begin(0); egIter != eGammaCollection.end(0); ++egIter) { // considering BX = only
163  edm::Ref<EGammaBxCollection> EGammaRef(eGammaHandle, ieg);
164  ieg++;
165 
166  float e_ele = egIter->energy();
167  float eta_ele = egIter->eta();
168  float et_ele = 0;
169  float cosh_eta_ele = cosh(eta_ele);
170  if (cosh_eta_ele > 0.0)
171  et_ele = e_ele / cosh_eta_ele;
172  else
173  et_ele = -1.0;
174  if (etMin_ > 0.0 && et_ele <= etMin_)
175  continue;
176  // match the L1EG object with a L1Track
177  float drmin = 999;
178  int itr = 0;
179  int itrack = -1;
180  for (trackIter = L1TTTrackHandle->begin(); trackIter != L1TTTrackHandle->end(); ++trackIter) {
181  edm::Ptr<L1TTTrackType> L1TrackPtr(L1TTTrackHandle, itr);
182  double trkPt_fit = trackIter->momentum().perp();
183  double trkPt_stubs = pTFrom2Stubs::pTFrom2(trackIter, tGeometry);
184  double trkPt = trkPt_fit;
185  if (useTwoStubsPT_)
186  trkPt = trkPt_stubs;
187 
188  if (trkPt > trkQualityPtMin_ && trackIter->chi2() < trkQualityChi2_) {
189  double dPhi = 99.;
190  double dR = 99.;
191  double dEta = 99.;
192  if (useClusterET_)
194  else
195  L1TkElectronTrackMatchAlgo::doMatch(egIter, L1TrackPtr, dPhi, dR, dEta);
196  if (dR < drmin && selectMatchedTrack(dR, dPhi, dEta, trkPt, eta_ele)) {
197  drmin = dR;
198  itrack = itr;
199  }
200  }
201  itr++;
202  }
203  if (itrack >= 0) {
204  edm::Ptr<L1TTTrackType> matchedL1TrackPtr(L1TTTrackHandle, itrack);
205 
206  const math::XYZTLorentzVector P4 = egIter->p4();
207  float trkisol = isolation(L1TTTrackHandle, itrack);
208  if (relativeIsolation_ && et_ele > 0.0) { // relative isolation
209  trkisol = trkisol / et_ele;
210  }
211 
212  TkElectron trkEm(P4, EGammaRef, matchedL1TrackPtr, trkisol);
213 
214  trkEm.setTrackCurvature(matchedL1TrackPtr->rInv()); // should this have npars? 4? 5?
215 
216  //std::cout<<matchedL1TrackPtr->rInv()<<" "<<matchedL1TrackPtr->rInv(4)<<" "<<matchedL1TrackPtr->rInv()<<std::endl;
217 
218  if (isoCut_ <= 0) {
219  // write the L1TkEm particle to the collection,
220  // irrespective of its relative isolation
221  result->push_back(trkEm);
222  } else {
223  // the object is written to the collection only
224  // if it passes the isolation cut
225  if (trkisol <= isoCut_)
226  result->push_back(trkEm);
227  }
228  }
229 
230  } // end loop over EGamma objects
231 
232  iEvent.put(std::move(result), label);
233 }
const edm::EDGetTokenT< EGammaBxCollection > egToken_
std::vector< TkElectron > TkElectronCollection
Definition: TkElectronFwd.h:15
edm::Ref< EGammaBxCollection > EGammaRef
Definition: TkEGTau.h:35
T const * product() const
Definition: Handle.h:70
void doMatchClusterET(BXVector< l1t::EGamma >::const_iterator egIter, const edm::Ptr< L1TTTrackType > &pTrk, double &dph, double &dr, double &deta)
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tGeomToken_
const_iterator begin(int bx) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool selectMatchedTrack(double &d_r, double &d_phi, double &d_eta, double &tk_pt, float &eg_eta)
std::vector< T >::const_iterator const_iterator
Definition: BXVector.h:18
const edm::EDGetTokenT< std::vector< TTTrack< Ref_Phase2TrackerDigi_ > > > trackToken_
int iEvent
Definition: GenABIO.cc:224
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void doMatch(BXVector< l1t::EGamma >::const_iterator egIter, const edm::Ptr< L1TTTrackType > &pTrk, double &dph, double &dr, double &deta)
float isolation(const edm::Handle< L1TTTrackCollectionType > &trkHandle, int match_index)
bool isValid() const
Definition: HandleBase.h:70
const_iterator end(int bx) const
float pTFrom2(std::vector< TTTrack< Ref_Phase2TrackerDigi_ > >::const_iterator trk, const TrackerGeometry *tkGeometry)
Definition: pTFrom2Stubs.cc:72
def move(src, dest)
Definition: eostools.py:511

◆ selectMatchedTrack()

bool L1TkElectronTrackProducer::selectMatchedTrack ( double &  d_r,
double &  d_phi,
double &  d_eta,
double &  tk_pt,
float &  eg_eta 
)
private

Definition at line 776 of file L1TkElectronTrackProducer.cc.

References funct::abs(), dEtaCutoff_, dPhiCutoff_, dRCutoff_, EB_MaxEta, getPtScaledCut(), and matchType_.

Referenced by produce().

777  {
778  if (matchType_ == "PtDependentCut") {
779  if (std::abs(d_phi) < getPtScaledCut(tk_pt, dPhiCutoff_) && d_r < getPtScaledCut(tk_pt, dRCutoff_))
780  return true;
781  } else {
782  double deta_max = dEtaCutoff_[0];
783  if (std::abs(eg_eta) < EB_MaxEta)
784  deta_max = dEtaCutoff_[1];
785  double dphi_max = dPhiCutoff_[0];
786  if ((d_eta / deta_max) * (d_eta / deta_max) + (d_phi / dphi_max) * (d_phi / dphi_max) < 1)
787  return true;
788  }
789  return false;
790 }
static constexpr float EB_MaxEta
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double getPtScaledCut(double pt, std::vector< double > &parameters)

Member Data Documentation

◆ deltaZ_

float L1TkElectronTrackProducer::deltaZ_
private

Definition at line 79 of file L1TkElectronTrackProducer.cc.

◆ dEtaCutoff_

std::vector<double> L1TkElectronTrackProducer::dEtaCutoff_
private

Definition at line 92 of file L1TkElectronTrackProducer.cc.

Referenced by selectMatchedTrack().

◆ dPhiCutoff_

std::vector<double> L1TkElectronTrackProducer::dPhiCutoff_
private

Definition at line 90 of file L1TkElectronTrackProducer.cc.

Referenced by selectMatchedTrack().

◆ dRCutoff_

std::vector<double> L1TkElectronTrackProducer::dRCutoff_
private

Definition at line 91 of file L1TkElectronTrackProducer.cc.

Referenced by selectMatchedTrack().

◆ dRMax_

const float L1TkElectronTrackProducer::dRMax_
private

Definition at line 78 of file L1TkElectronTrackProducer.cc.

Referenced by isolation().

◆ dRMin_

const float L1TkElectronTrackProducer::dRMin_
private

Definition at line 77 of file L1TkElectronTrackProducer.cc.

Referenced by isolation().

◆ egToken_

const edm::EDGetTokenT<EGammaBxCollection> L1TkElectronTrackProducer::egToken_
private

Definition at line 95 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ etMin_

const float L1TkElectronTrackProducer::etMin_
private

Definition at line 75 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ isoCut_

float L1TkElectronTrackProducer::isoCut_
private

Definition at line 83 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ label

std::string L1TkElectronTrackProducer::label
private

◆ matchType_

std::string L1TkElectronTrackProducer::matchType_
private

Definition at line 93 of file L1TkElectronTrackProducer.cc.

Referenced by selectMatchedTrack().

◆ maxChi2IsoTracks_

const float L1TkElectronTrackProducer::maxChi2IsoTracks_
private

Definition at line 81 of file L1TkElectronTrackProducer.cc.

Referenced by isolation().

◆ minNStubsIsoTracks_

const unsigned int L1TkElectronTrackProducer::minNStubsIsoTracks_
private

Definition at line 82 of file L1TkElectronTrackProducer.cc.

Referenced by isolation().

◆ primaryVtxConstrain_

bool L1TkElectronTrackProducer::primaryVtxConstrain_
private

Definition at line 85 of file L1TkElectronTrackProducer.cc.

◆ pTMinTra_

const float L1TkElectronTrackProducer::pTMinTra_
private

Definition at line 76 of file L1TkElectronTrackProducer.cc.

Referenced by isolation().

◆ relativeIsolation_

bool L1TkElectronTrackProducer::relativeIsolation_
private

Definition at line 84 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ tGeomToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> L1TkElectronTrackProducer::tGeomToken_
private

Definition at line 97 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ trackToken_

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

Definition at line 96 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ trkQualityChi2_

float L1TkElectronTrackProducer::trkQualityChi2_
private

Definition at line 86 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ trkQualityPtMin_

float L1TkElectronTrackProducer::trkQualityPtMin_
private

Definition at line 87 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ useClusterET_

bool L1TkElectronTrackProducer::useClusterET_
private

Definition at line 89 of file L1TkElectronTrackProducer.cc.

Referenced by produce().

◆ useTwoStubsPT_

bool L1TkElectronTrackProducer::useTwoStubsPT_
private

Definition at line 88 of file L1TkElectronTrackProducer.cc.

Referenced by produce().