CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
L1RegionData< T1 > Class Template Reference
Inheritance diagram for L1RegionData< T1 >:
L1RegionDataBase

Public Member Functions

void getEtaPhiRegions (const edm::Event &, const edm::EventSetup &, std::vector< RectangularEtaPhiRegion > &) const override
 
template<>
void getEtaPhiRegions (const edm::Event &event, const edm::EventSetup &setup, std::vector< RectangularEtaPhiRegion > &regions) const
 
template<>
void getEtaPhiRegions (const edm::Event &event, const edm::EventSetup &setup, std::vector< RectangularEtaPhiRegion > &regions) const
 
template<typename T2 >
bool isEmpty (const T2 &coll) const
 
template<typename T2 >
bool isEmpty (const BXVector< T2 > &coll) const
 
 L1RegionData (const edm::ParameterSet &para, edm::ConsumesCollector &consumesColl)
 
- Public Member Functions inherited from L1RegionDataBase
virtual ~L1RegionDataBase ()
 

Static Public Member Functions

template<typename T2 >
static T2::const_iterator beginIt (const T2 &coll)
 
template<typename T2 >
static BXVector< T2 >::const_iterator beginIt (const BXVector< T2 > &coll)
 
template<typename T2 >
static T2::const_iterator endIt (const T2 &coll)
 
template<typename T2 >
static BXVector< T2 >::const_iterator endIt (const BXVector< T2 > &coll)
 

Private Member Functions

void eventSetupConsumes (edm::ConsumesCollector &consumesColl)
 
template<>
void eventSetupConsumes (edm::ConsumesCollector &consumesColl)
 
template<>
void eventSetupConsumes (edm::ConsumesCollector &consumesColl)
 

Private Attributes

edm::ESGetToken< L1CaloGeometry, L1CaloGeometryRecordl1CaloGeometryToken_
 
double const maxEt_
 
double const minEt_
 
double const regionEtaMargin_
 
double const regionPhiMargin_
 
edm::EDGetTokenT< T1 > const token_
 

Detailed Description

template<typename T1>
class L1RegionData< T1 >

Definition at line 59 of file HLTRecHitInAllL1RegionsProducer.cc.

Constructor & Destructor Documentation

◆ L1RegionData()

template<typename T1>
L1RegionData< T1 >::L1RegionData ( const edm::ParameterSet para,
edm::ConsumesCollector consumesColl 
)
inline

Definition at line 71 of file HLTRecHitInAllL1RegionsProducer.cc.

References L1RegionData< T1 >::eventSetupConsumes().

72  : minEt_(para.getParameter<double>("minEt")),
73  maxEt_(para.getParameter<double>("maxEt")),
74  regionEtaMargin_(para.getParameter<double>("regionEtaMargin")),
75  regionPhiMargin_(para.getParameter<double>("regionPhiMargin")),
76  token_(consumesColl.consumes<T1>(para.getParameter<edm::InputTag>("inputColl"))) {
77  eventSetupConsumes(consumesColl);
78  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void eventSetupConsumes(edm::ConsumesCollector &consumesColl)
edm::EDGetTokenT< T1 > const token_

Member Function Documentation

◆ beginIt() [1/2]

template<typename T1>
template<typename T2 >
static T2::const_iterator L1RegionData< T1 >::beginIt ( const T2 &  coll)
inlinestatic

Definition at line 88 of file HLTRecHitInAllL1RegionsProducer.cc.

88  {
89  return coll.begin();
90  }

◆ beginIt() [2/2]

template<typename T1>
template<typename T2 >
static BXVector<T2>::const_iterator L1RegionData< T1 >::beginIt ( const BXVector< T2 > &  coll)
inlinestatic

Definition at line 100 of file HLTRecHitInAllL1RegionsProducer.cc.

References BXVector< T >::begin().

100  {
101  return coll.begin(0);
102  }
const_iterator begin(int bx) const

◆ endIt() [1/2]

template<typename T1>
template<typename T2 >
static T2::const_iterator L1RegionData< T1 >::endIt ( const T2 &  coll)
inlinestatic

Definition at line 92 of file HLTRecHitInAllL1RegionsProducer.cc.

92  {
93  return coll.end();
94  }

◆ endIt() [2/2]

template<typename T1>
template<typename T2 >
static BXVector<T2>::const_iterator L1RegionData< T1 >::endIt ( const BXVector< T2 > &  coll)
inlinestatic

Definition at line 104 of file HLTRecHitInAllL1RegionsProducer.cc.

References BXVector< T >::end().

104  {
105  return coll.end(0);
106  }
const_iterator end(int bx) const

◆ eventSetupConsumes() [1/3]

template<typename L1CollType >
void L1RegionData< L1CollType >::eventSetupConsumes ( edm::ConsumesCollector consumesColl)
private

Definition at line 289 of file HLTRecHitInAllL1RegionsProducer.cc.

Referenced by L1RegionData< T1 >::L1RegionData().

289 {}

◆ eventSetupConsumes() [2/3]

template<>
void L1RegionData< l1extra::L1JetParticleCollection >::eventSetupConsumes ( edm::ConsumesCollector consumesColl)
private

Definition at line 318 of file HLTRecHitInAllL1RegionsProducer.cc.

References edm::ConsumesCollector::esConsumes().

318  {
319  l1CaloGeometryToken_ = consumesColl.esConsumes();
320 }
edm::ESGetToken< L1CaloGeometry, L1CaloGeometryRecord > l1CaloGeometryToken_

◆ eventSetupConsumes() [3/3]

template<>
void L1RegionData< l1extra::L1EmParticleCollection >::eventSetupConsumes ( edm::ConsumesCollector consumesColl)
private

Definition at line 353 of file HLTRecHitInAllL1RegionsProducer.cc.

References edm::ConsumesCollector::esConsumes().

353  {
354  l1CaloGeometryToken_ = consumesColl.esConsumes();
355 }
edm::ESGetToken< L1CaloGeometry, L1CaloGeometryRecord > l1CaloGeometryToken_

◆ getEtaPhiRegions() [1/3]

template<typename L1CollType >
void L1RegionData< L1CollType >::getEtaPhiRegions ( const edm::Event event,
const edm::EventSetup ,
std::vector< RectangularEtaPhiRegion > &  regions 
) const
overridevirtual

Implements L1RegionDataBase.

Definition at line 292 of file HLTRecHitInAllL1RegionsProducer.cc.

References LogDebug, Skims_PA_cff::name, l1tPhase1JetProducer_cfi::phiLow, and edm::typeDemangle().

294  {
295  edm::Handle<L1CollType> l1Cands;
296  event.getByToken(token_, l1Cands);
297 
298  if (isEmpty(*l1Cands)) {
299  LogDebug("HLTRecHitInAllL1RegionsProducerL1RegionData")
300  << "The input collection of L1T candidates is empty (L1CollType = \""
301  << edm::typeDemangle(typeid(L1CollType).name()) << "\"). No regions selected.";
302  return;
303  }
304 
305  for (auto l1CandIt = beginIt(*l1Cands); l1CandIt != endIt(*l1Cands); ++l1CandIt) {
306  if (l1CandIt->et() >= minEt_ && l1CandIt->et() < maxEt_) {
307  double etaLow = l1CandIt->eta() - regionEtaMargin_;
308  double etaHigh = l1CandIt->eta() + regionEtaMargin_;
309  double phiLow = l1CandIt->phi() - regionPhiMargin_;
310  double phiHigh = l1CandIt->phi() + regionPhiMargin_;
311 
312  regions.push_back(RectangularEtaPhiRegion(etaLow, etaHigh, phiLow, phiHigh));
313  }
314  }
315 }
bool isEmpty(const T2 &coll) const
static T2::const_iterator endIt(const T2 &coll)
edm::EDGetTokenT< T1 > const token_
std::string typeDemangle(char const *mangledName)
static T2::const_iterator beginIt(const T2 &coll)
#define LogDebug(id)

◆ getEtaPhiRegions() [2/3]

template<>
void L1RegionData< l1extra::L1JetParticleCollection >::getEtaPhiRegions ( const edm::Event event,
const edm::EventSetup setup,
std::vector< RectangularEtaPhiRegion > &  regions 
) const
virtual

Implements L1RegionDataBase.

Definition at line 323 of file HLTRecHitInAllL1RegionsProducer.cc.

References l1tPhase1JetProducer_cfi::phiLow, and singleTopDQM_cfi::setup.

324  {
326  event.getByToken(token_, l1Cands);
327 
328  auto const& l1CaloGeom = setup.getData(l1CaloGeometryToken_);
329 
330  for (const auto& l1Cand : *l1Cands) {
331  if (l1Cand.et() >= minEt_ && l1Cand.et() < maxEt_) {
332  // Access the GCT hardware object corresponding to the L1Extra EM object.
333  int etaIndex = l1Cand.gctJetCand()->etaIndex();
334  int phiIndex = l1Cand.gctJetCand()->phiIndex();
335 
336  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
337  double etaLow = l1CaloGeom.etaBinLowEdge(etaIndex);
338  double etaHigh = l1CaloGeom.etaBinHighEdge(etaIndex);
339  double phiLow = l1CaloGeom.emJetPhiBinLowEdge(phiIndex);
340  double phiHigh = l1CaloGeom.emJetPhiBinHighEdge(phiIndex);
341 
342  etaLow -= regionEtaMargin_;
343  etaHigh += regionEtaMargin_;
345  phiHigh += regionPhiMargin_;
346 
347  regions.push_back(RectangularEtaPhiRegion(etaLow, etaHigh, phiLow, phiHigh));
348  }
349  }
350 }
edm::EDGetTokenT< T1 > const token_
edm::ESGetToken< L1CaloGeometry, L1CaloGeometryRecord > l1CaloGeometryToken_

◆ getEtaPhiRegions() [3/3]

template<>
void L1RegionData< l1extra::L1EmParticleCollection >::getEtaPhiRegions ( const edm::Event event,
const edm::EventSetup setup,
std::vector< RectangularEtaPhiRegion > &  regions 
) const
virtual

Implements L1RegionDataBase.

Definition at line 358 of file HLTRecHitInAllL1RegionsProducer.cc.

References l1tPhase1JetProducer_cfi::phiLow, and singleTopDQM_cfi::setup.

359  {
361  event.getByToken(token_, l1Cands);
362 
363  auto const& l1CaloGeom = setup.getData(l1CaloGeometryToken_);
364 
365  for (const auto& l1Cand : *l1Cands) {
366  if (l1Cand.et() >= minEt_ && l1Cand.et() < maxEt_) {
367  // Access the GCT hardware object corresponding to the L1Extra EM object.
368  int etaIndex = l1Cand.gctEmCand()->etaIndex();
369  int phiIndex = l1Cand.gctEmCand()->phiIndex();
370 
371  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
372  double etaLow = l1CaloGeom.etaBinLowEdge(etaIndex);
373  double etaHigh = l1CaloGeom.etaBinHighEdge(etaIndex);
374  double phiLow = l1CaloGeom.emJetPhiBinLowEdge(phiIndex);
375  double phiHigh = l1CaloGeom.emJetPhiBinHighEdge(phiIndex);
376 
377  etaLow -= regionEtaMargin_;
378  etaHigh += regionEtaMargin_;
380  phiHigh += regionPhiMargin_;
381 
382  regions.push_back(RectangularEtaPhiRegion(etaLow, etaHigh, phiLow, phiHigh));
383  }
384  }
385 }
edm::EDGetTokenT< T1 > const token_
edm::ESGetToken< L1CaloGeometry, L1CaloGeometryRecord > l1CaloGeometryToken_

◆ isEmpty() [1/2]

template<typename T1>
template<typename T2 >
bool L1RegionData< T1 >::isEmpty ( const T2 &  coll) const
inline

Definition at line 84 of file HLTRecHitInAllL1RegionsProducer.cc.

Referenced by plotting.Plot::clone().

84  {
85  return coll.empty();
86  }

◆ isEmpty() [2/2]

template<typename T1>
template<typename T2 >
bool L1RegionData< T1 >::isEmpty ( const BXVector< T2 > &  coll) const
inline

Definition at line 96 of file HLTRecHitInAllL1RegionsProducer.cc.

References BXVector< T >::isEmpty(), or, and BXVector< T >::size().

Referenced by plotting.Plot::clone().

96  {
97  return (coll.size() == 0 or coll.isEmpty(0));
98  }
bool isEmpty(int bx) const
unsigned size(int bx) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12

Member Data Documentation

◆ l1CaloGeometryToken_

template<typename T1>
edm::ESGetToken<L1CaloGeometry, L1CaloGeometryRecord> L1RegionData< T1 >::l1CaloGeometryToken_
private

Definition at line 66 of file HLTRecHitInAllL1RegionsProducer.cc.

◆ maxEt_

template<typename T1>
double const L1RegionData< T1 >::maxEt_
private

Definition at line 62 of file HLTRecHitInAllL1RegionsProducer.cc.

◆ minEt_

template<typename T1>
double const L1RegionData< T1 >::minEt_
private

Definition at line 61 of file HLTRecHitInAllL1RegionsProducer.cc.

◆ regionEtaMargin_

template<typename T1>
double const L1RegionData< T1 >::regionEtaMargin_
private

Definition at line 63 of file HLTRecHitInAllL1RegionsProducer.cc.

◆ regionPhiMargin_

template<typename T1>
double const L1RegionData< T1 >::regionPhiMargin_
private

Definition at line 64 of file HLTRecHitInAllL1RegionsProducer.cc.

◆ token_

template<typename T1>
edm::EDGetTokenT<T1> const L1RegionData< T1 >::token_
private

Definition at line 65 of file HLTRecHitInAllL1RegionsProducer.cc.