CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
HLTHPDFilter Class Reference

#include <HLTHPDFilter.h>

Inheritance diagram for HLTHPDFilter:
edm::stream::EDFilter<> edm::stream::EDFilterBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

bool filter (edm::Event &, const edm::EventSetup &) override
 
 HLTHPDFilter (const edm::ParameterSet &)
 
 ~HLTHPDFilter () override
 
- Public Member Functions inherited from edm::stream::EDFilter<>
 EDFilter ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
- Public Member Functions inherited from edm::stream::EDFilterBase
 EDFilterBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDFilterBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::stream::EDFilterBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Attributes

edm::EDGetTokenT< HBHERecHitCollectionm_theRecHitCollectionToken
 
double mEnergyThreshold
 
double mHPDSpikeEnergyThreshold
 
double mHPDSpikeIsolationEnergyThreshold
 
edm::InputTag mInputTag
 
double mRBXSpikeEnergyThreshold
 
double mRBXSpikeUnbalanceThreshold
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDFilter<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDFilterBase
typedef EDFilterAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Author
Fedor Ratnikov (UMd)

Definition at line 21 of file HLTHPDFilter.h.

Constructor & Destructor Documentation

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

Definition at line 76 of file HLTHPDFilter.cc.

References m_theRecHitCollectionToken, mInputTag, and ~HLTHPDFilter().

76  :
77  mInputTag (iConfig.getParameter <edm::InputTag> ("inputTag")),
78  mEnergyThreshold (iConfig.getParameter <double> ("energy")),
79  mHPDSpikeEnergyThreshold (iConfig.getParameter <double> ("hpdSpikeEnergy")),
80  mHPDSpikeIsolationEnergyThreshold (iConfig.getParameter <double> ("hpdSpikeIsolationEnergy")),
81  mRBXSpikeEnergyThreshold (iConfig.getParameter <double> ("rbxSpikeEnergy")),
82  mRBXSpikeUnbalanceThreshold (iConfig.getParameter <double> ("rbxSpikeUnbalance"))
83 {
84  m_theRecHitCollectionToken = consumes<HBHERecHitCollection>(mInputTag);
85 }
T getParameter(std::string const &) const
edm::EDGetTokenT< HBHERecHitCollection > m_theRecHitCollectionToken
Definition: HLTHPDFilter.h:30
double mEnergyThreshold
Definition: HLTHPDFilter.h:32
double mRBXSpikeEnergyThreshold
Definition: HLTHPDFilter.h:35
double mHPDSpikeIsolationEnergyThreshold
Definition: HLTHPDFilter.h:34
double mRBXSpikeUnbalanceThreshold
Definition: HLTHPDFilter.h:36
double mHPDSpikeEnergyThreshold
Definition: HLTHPDFilter.h:33
edm::InputTag mInputTag
Definition: HLTHPDFilter.h:31
HLTHPDFilter::~HLTHPDFilter ( )
overridedefault

Referenced by HLTHPDFilter().

Member Function Documentation

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

Definition at line 90 of file HLTHPDFilter.cc.

References edm::ConfigurationDescriptions::add(), and edm::ParameterSetDescription::add().

90  {
92  desc.add<edm::InputTag>("inputTag",edm::InputTag("hltHbhereco"));
93  desc.add<double>("energy",-99.0);
94  desc.add<double>("hpdSpikeEnergy",10.0);
95  desc.add<double>("hpdSpikeIsolationEnergy",1.0);
96  desc.add<double>("rbxSpikeEnergy",50.0);
97  desc.add<double>("rbxSpikeUnbalance",0.2);
98  descriptions.add("hltHPDFilter",desc);
99 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool HLTHPDFilter::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 101 of file HLTHPDFilter.cc.

References edm::Event::getByToken(), photonIsolationHIProducer_cfi::hbhe, mps_fire::i, createfilelist::int, m_theRecHitCollectionToken, particleFlowClusterECALTimeSelected_cfi::maxEnergy, mEnergyThreshold, mHPDSpikeEnergyThreshold, mHPDSpikeIsolationEnergyThreshold, particleFlowClusterECALTimeSelected_cfi::minEnergy, mRBXSpikeEnergyThreshold, mRBXSpikeUnbalanceThreshold, and edm::SortedCollection< T, SORT >::size().

102 {
103  if (mHPDSpikeEnergyThreshold <= 0 && mRBXSpikeEnergyThreshold <= 0) return true; // nothing to filter
104  // get hits
107 
108  // collect energies
109  float hpdEnergy[4][73];
110  for (auto & i : hpdEnergy) for (size_t j = 0; j < 73; ++j) i[j] = 0;
111 
112  // select hist above threshold
113  for (unsigned i = 0; i < hbhe->size(); ++i) {
114  if ((*hbhe)[i].energy() > mEnergyThreshold) {
115  std::pair<Partition,int> hpd = hpdId ((*hbhe)[i].id());
116  hpdEnergy[int (hpd.first)][hpd.second] += (*hbhe)[i].energy ();
117  }
118  }
119 
120  // not single HPD spike
121  if (mHPDSpikeEnergyThreshold > 0) {
122  for (auto & partition : hpdEnergy) {
123  for (size_t i = 1; i < 73; ++i) {
124  if (partition[i] > mHPDSpikeEnergyThreshold) {
125  int hpdPlus = i + 1;
126  if (hpdPlus == 73) hpdPlus = 1;
127  int hpdMinus = i - 1;
128  if (hpdMinus == 0) hpdMinus = 72;
129  double maxNeighborEnergy = fmax (partition[hpdPlus], partition[hpdMinus]);
130  if (maxNeighborEnergy < mHPDSpikeIsolationEnergyThreshold) return false; // HPD spike found
131  }
132  }
133  }
134  }
135 
136  // not RBX flash
137  if (mRBXSpikeEnergyThreshold > 0) {
138  for (auto & partition : hpdEnergy) {
139  for (size_t rbx = 1; rbx < 19; ++rbx) {
140  int ifirst = (rbx-1)*4-1;
141  int iend = (rbx-1)*4+3;
142  double minEnergy = 0;
143  double maxEnergy = -1;
144  double totalEnergy = 0;
145  for (int irm = ifirst; irm < iend; ++irm) {
146  int hpd = irm;
147  if (hpd <= 0) hpd = 72 + hpd;
148  totalEnergy += partition[hpd];
149  if (minEnergy > maxEnergy) {
150  minEnergy = maxEnergy = partition[hpd];
151  }
152  else {
153  if (partition[hpd] < minEnergy) minEnergy = partition[hpd];
154  if (partition[hpd] > maxEnergy) maxEnergy = partition[hpd];
155  }
156  }
157  if (totalEnergy > mRBXSpikeEnergyThreshold) {
158  if (minEnergy / maxEnergy > mRBXSpikeUnbalanceThreshold) return false; // likely HPF flash
159  }
160  }
161  }
162  }
163  return true;
164 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< HBHERecHitCollection > m_theRecHitCollectionToken
Definition: HLTHPDFilter.h:30
double mEnergyThreshold
Definition: HLTHPDFilter.h:32
double mRBXSpikeEnergyThreshold
Definition: HLTHPDFilter.h:35
double mHPDSpikeIsolationEnergyThreshold
Definition: HLTHPDFilter.h:34
double mRBXSpikeUnbalanceThreshold
Definition: HLTHPDFilter.h:36
double mHPDSpikeEnergyThreshold
Definition: HLTHPDFilter.h:33
size_type size() const

Member Data Documentation

edm::EDGetTokenT<HBHERecHitCollection> HLTHPDFilter::m_theRecHitCollectionToken
private

Definition at line 30 of file HLTHPDFilter.h.

Referenced by filter(), and HLTHPDFilter().

double HLTHPDFilter::mEnergyThreshold
private

Definition at line 32 of file HLTHPDFilter.h.

Referenced by filter().

double HLTHPDFilter::mHPDSpikeEnergyThreshold
private

Definition at line 33 of file HLTHPDFilter.h.

Referenced by filter().

double HLTHPDFilter::mHPDSpikeIsolationEnergyThreshold
private

Definition at line 34 of file HLTHPDFilter.h.

Referenced by filter().

edm::InputTag HLTHPDFilter::mInputTag
private

Definition at line 31 of file HLTHPDFilter.h.

Referenced by HLTHPDFilter().

double HLTHPDFilter::mRBXSpikeEnergyThreshold
private

Definition at line 35 of file HLTHPDFilter.h.

Referenced by filter().

double HLTHPDFilter::mRBXSpikeUnbalanceThreshold
private

Definition at line 36 of file HLTHPDFilter.h.

Referenced by filter().