CMS 3D CMS Logo

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

#include <HLTEgammaGenericQuadraticFilter.h>

Inheritance diagram for HLTEgammaGenericQuadraticFilter:
HLTFilter edm::global::EDFilter<> edm::global::EDFilterBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 HLTEgammaGenericQuadraticFilter (const edm::ParameterSet &)
 
bool hltFilter (edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
 
 ~HLTEgammaGenericQuadraticFilter () override
 
- Public Member Functions inherited from HLTFilter
 HLTFilter (const edm::ParameterSet &config)
 
int module (edm::Event const &) const
 
const std::string * moduleLabel () const
 
int path (edm::Event const &) const
 
const std::string * pathName (edm::Event const &) const
 
std::pair< int, int > pmid (edm::Event const &) const
 
bool saveTags () const
 
 ~HLTFilter () override
 
- Public Member Functions inherited from edm::global::EDFilter<>
 EDFilter ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::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 HLTFilter
static void makeHLTFilterDescription (edm::ParameterSetDescription &desc)
 
- Static Public Member Functions inherited from edm::global::EDFilterBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Attributes

std::vector< double > absEtaLowEdges_
 
edm::InputTag candTag_
 
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefscandToken_
 
bool doRhoCorrection_
 
std::vector< double > effectiveAreas_
 
std::vector< double > energyLowEdges_
 
edm::InputTag l1EGTag_
 
bool lessThan_
 
int ncandcut_
 
double rhoMax_
 
double rhoScale_
 
edm::InputTag rhoTag_
 
edm::EDGetTokenT< double > rhoToken_
 
std::vector< double > thrOverE2EB_
 
std::vector< double > thrOverE2EE_
 
std::vector< double > thrOverEEB_
 
std::vector< double > thrOverEEE_
 
std::vector< double > thrRegularEB_
 
std::vector< double > thrRegularEE_
 
bool useEt_
 
edm::InputTag varTag_
 
edm::EDGetTokenT< reco::RecoEcalCandidateIsolationMapvarToken_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDFilterBase
typedef EDFilterBase 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
Roberto Covarelli (CERN) modified by Chris Tully (Princeton)

Definition at line 23 of file HLTEgammaGenericQuadraticFilter.h.

Constructor & Destructor Documentation

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

Definition at line 25 of file HLTEgammaGenericQuadraticFilter.cc.

References absEtaLowEdges_, candTag_, candToken_, doRhoCorrection_, effectiveAreas_, energyLowEdges_, Exception, edm::ParameterSet::getParameter(), l1EGTag_, lessThan_, ncandcut_, or, rhoMax_, rhoScale_, rhoTag_, rhoToken_, thrOverE2EB_, thrOverE2EE_, thrOverEEB_, thrOverEEE_, thrRegularEB_, thrRegularEE_, useEt_, varTag_, and varToken_.

25  : HLTFilter(iConfig) {
26  candTag_ = iConfig.getParameter< edm::InputTag > ("candTag");
27  varTag_ = iConfig.getParameter< edm::InputTag > ("varTag");
28  rhoTag_ = iConfig.getParameter< edm::InputTag > ("rhoTag");
29 
30  energyLowEdges_ = iConfig.getParameter<std::vector<double> > ("energyLowEdges");
31  lessThan_ = iConfig.getParameter<bool> ("lessThan");
32  useEt_ = iConfig.getParameter<bool> ("useEt");
33 
34  thrRegularEB_ = iConfig.getParameter<std::vector<double> > ("thrRegularEB");
35  thrRegularEE_ = iConfig.getParameter<std::vector<double> > ("thrRegularEE");
36  thrOverEEB_ = iConfig.getParameter<std::vector<double> > ("thrOverEEB");
37  thrOverEEE_ = iConfig.getParameter<std::vector<double> > ("thrOverEEE");
38  thrOverE2EB_ = iConfig.getParameter<std::vector<double> > ("thrOverE2EB");
39  thrOverE2EE_ = iConfig.getParameter<std::vector<double> > ("thrOverE2EE");
40 
41  ncandcut_ = iConfig.getParameter<int> ("ncandcut");
42 
43  doRhoCorrection_ = iConfig.getParameter<bool> ("doRhoCorrection");
44  rhoMax_ = iConfig.getParameter<double> ("rhoMax");
45  rhoScale_ = iConfig.getParameter<double> ("rhoScale");
46  effectiveAreas_ = iConfig.getParameter<std::vector<double> > ("effectiveAreas");
47  absEtaLowEdges_ = iConfig.getParameter<std::vector<double> > ("absEtaLowEdges");
48 
49  l1EGTag_= iConfig.getParameter< edm::InputTag > ("l1EGCand");
50 
51  candToken_ = consumes<trigger::TriggerFilterObjectWithRefs>(candTag_);
52  varToken_ = consumes<reco::RecoEcalCandidateIsolationMap>(varTag_);
53 
54  if (energyLowEdges_.size() != thrRegularEB_.size() or energyLowEdges_.size() != thrRegularEE_.size() or
55  energyLowEdges_.size() != thrOverEEB_.size() or energyLowEdges_.size() != thrOverEEE_.size() or
56  energyLowEdges_.size() != thrOverE2EB_.size() or energyLowEdges_.size() != thrOverE2EE_.size())
57  throw cms::Exception("IncompatibleVects") << "energyLowEdges and threshold vectors should be of the same size. \n";
58 
59  if (energyLowEdges_.at(0) != 0.0)
60  throw cms::Exception("IncompleteCoverage") << "energyLowEdges should start from 0. \n";
61 
62  for (unsigned int aIt = 0; aIt < energyLowEdges_.size() - 1; aIt++) {
63  if ( !(energyLowEdges_.at( aIt ) < energyLowEdges_.at( aIt + 1 )) )
64  throw cms::Exception("ImproperBinning") << "energyLowEdges entries should be in increasing order. \n";
65  }
66 
67  if (doRhoCorrection_) {
68  rhoToken_ = consumes<double> (rhoTag_);
69  if (absEtaLowEdges_.size() != effectiveAreas_.size())
70  throw cms::Exception("IncompatibleVects") << "absEtaLowEdges and effectiveAreas should be of the same size. \n";
71 
72  if (absEtaLowEdges_.at(0) != 0.0)
73  throw cms::Exception("IncompleteCoverage") << "absEtaLowEdges should start from 0. \n";
74 
75  for (unsigned int bIt = 0; bIt < absEtaLowEdges_.size() - 1; bIt++) {
76  if ( !(absEtaLowEdges_.at( bIt ) < absEtaLowEdges_.at( bIt + 1 )) )
77  throw cms::Exception("ImproperBinning") << "absEtaLowEdges entries should be in increasing order. \n";
78  }
79  }
80 }
T getParameter(std::string const &) const
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
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
edm::EDGetTokenT< reco::RecoEcalCandidateIsolationMap > varToken_
HLTFilter(const edm::ParameterSet &config)
Definition: HLTFilter.cc:20
HLTEgammaGenericQuadraticFilter::~HLTEgammaGenericQuadraticFilter ( )
override

Definition at line 108 of file HLTEgammaGenericQuadraticFilter.cc.

108 {}

Member Function Documentation

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

Definition at line 83 of file HLTEgammaGenericQuadraticFilter.cc.

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

83  {
86  desc.add<edm::InputTag>("candTag", edm::InputTag("hltSingleEgammaEtFilter"));
87  desc.add<edm::InputTag>("varTag", edm::InputTag("hltSingleEgammaHcalIsol"));
88  desc.add<edm::InputTag>("rhoTag", edm::InputTag("")); // No rho correction by default
89  desc.add<std::vector<double> >("energyLowEdges", {0.0}); // No energy-dependent cuts by default
90  desc.add<bool>("lessThan", true);
91  desc.add<bool>("useEt", false);
92  desc.add<std::vector<double> >("thrRegularEB", {0.0});
93  desc.add<std::vector<double> >("thrRegularEE", {0.0});
94  desc.add<std::vector<double> >("thrOverEEB", {-1.0});
95  desc.add<std::vector<double> >("thrOverEEE", {-1.0});
96  desc.add<std::vector<double> >("thrOverE2EB", {-1.0});
97  desc.add<std::vector<double> >("thrOverE2EE", {-1.0});
98  desc.add<int>("ncandcut",1);
99  desc.add<bool>("doRhoCorrection", false);
100  desc.add<double>("rhoMax", 9.9999999E7);
101  desc.add<double>("rhoScale", 1.0);
102  desc.add<std::vector<double> >("effectiveAreas", {0.0, 0.0});
103  desc.add<std::vector<double> >("absEtaLowEdges", {0.0, 1.479}); // EB, EE
104  desc.add<edm::InputTag>("l1EGCand", edm::InputTag("hltL1IsoRecoEcalCandidate"));
105  descriptions.add("hltEgammaGenericQuadraticFilter", desc);
106 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool HLTEgammaGenericQuadraticFilter::hltFilter ( edm::Event iEvent,
const edm::EventSetup iSetup,
trigger::TriggerFilterObjectWithRefs filterproduct 
) const
overridevirtual

Implements HLTFilter.

Definition at line 113 of file HLTEgammaGenericQuadraticFilter.cc.

References funct::abs(), absEtaLowEdges_, accept(), trigger::TriggerFilterObjectWithRefs::addCollectionTag(), trigger::TriggerRefsCollections::addObject(), candToken_, DEFINE_FWK_MODULE, SoftLeptonByDistance_cfi::distance, doRhoCorrection_, effectiveAreas_, energyLowEdges_, JetChargeProducer_cfi::exp, edm::Event::getByToken(), trigger::TriggerRefsCollections::getObjects(), mps_fire::i, l1EGTag_, lessThan_, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, gen::n, ncandcut_, edm::Handle< T >::product(), rho, rhoMax_, rhoScale_, rhoToken_, HLTFilter::saveTags(), funct::sin(), thrOverE2EB_, thrOverE2EE_, thrOverEEB_, thrOverEEE_, thrRegularEB_, thrRegularEE_, trigger::TriggerCluster, trigger::TriggerPhoton, useEt_, edm::helpers::KeyVal< K, V >::val, and varToken_.

114 {
115  using namespace trigger;
116  if (saveTags()) {
117  filterproduct.addCollectionTag(l1EGTag_);
118  }
119 
120  // Ref to Candidate object to be recorded in filter object
122 
123  // Set output format
124  int trigger_type = trigger::TriggerCluster;
125  if (saveTags()) trigger_type = trigger::TriggerPhoton;
126 
128  iEvent.getByToken (candToken_, PrevFilterOutput);
129 
130  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands;
131  PrevFilterOutput->getObjects(TriggerCluster, recoecalcands);
132  if(recoecalcands.empty()) PrevFilterOutput->getObjects(TriggerPhoton,recoecalcands); //we dont know if its type trigger cluster or trigger photon
133 
134  //get hold of isolated association map
136  iEvent.getByToken (varToken_,depMap);
137 
138  // Get rho if needed
139  edm::Handle<double> rhoHandle;
140  double rho = 0.0;
141  if (doRhoCorrection_) {
142  iEvent.getByToken(rhoToken_, rhoHandle);
143  rho = *(rhoHandle.product());
144  }
145 
146  if (rho > rhoMax_)
147  rho = rhoMax_;
148  rho = rho * rhoScale_;
149 
150  // look at all photons, check cuts and add to filter object
151  int n = 0;
152  for (unsigned int i=0; i<recoecalcands.size(); i++) {
153 
154  ref = recoecalcands[i];
155  reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*depMap).find( ref );
156 
157  float vali = mapi->val;
158  float EtaSC = ref->eta();
159 
160  // Pick the right EA and do rhoCorr
161  if (doRhoCorrection_) {
162  auto cIt = std::lower_bound( absEtaLowEdges_.begin(), absEtaLowEdges_.end(), std::abs(EtaSC) ) - 1;
163  vali = vali - (rho * effectiveAreas_.at( std::distance( absEtaLowEdges_.begin(), cIt ) ));
164  }
165 
166  float energy = ref->superCluster()->energy();
167  if (useEt_) energy = energy * sin (2*atan(exp(-EtaSC)));
168  if (energy < 0.) energy = 0.; /* first and second order terms assume non-negative energies */
169 
170  // Pick the right cut threshold
171  double cutRegularEB_ = 9999., cutRegularEE_ = 9999.;
172  double cutOverEEB_ = 9999., cutOverEEE_ = 9999.;
173  double cutOverE2EB_ = 9999., cutOverE2EE_ = 9999.;
174 
175  auto dIt = std::lower_bound( energyLowEdges_.begin(), energyLowEdges_.end(), energy ) - 1;
176  unsigned iEn = std::distance( energyLowEdges_.begin(), dIt );
177 
178  cutRegularEB_ = thrRegularEB_.at(iEn);
179  cutRegularEE_ = thrRegularEE_.at(iEn);
180  cutOverEEB_ = thrOverEEB_.at(iEn);
181  cutOverEEE_ = thrOverEEE_.at(iEn);
182  cutOverE2EB_ = thrOverE2EB_.at(iEn);
183  cutOverE2EE_ = thrOverE2EE_.at(iEn);
184 
185  if ( lessThan_ ) {
186  if ((std::abs(EtaSC) < 1.479 && vali <= cutRegularEB_ + energy*cutOverEEB_ + energy*energy*cutOverE2EB_) || (std::abs(EtaSC) >= 1.479 && vali <= cutRegularEE_ + energy*cutOverEEE_ + energy*energy*cutOverE2EE_) ) {
187  n++;
188  filterproduct.addObject(trigger_type, ref);
189  continue;
190  }
191  } else {
192  if ((std::abs(EtaSC) < 1.479 && vali >= cutRegularEB_ + energy*cutOverEEB_ + energy*energy*cutOverE2EB_) || (std::abs(EtaSC) >= 1.479 && vali >= cutRegularEE_ + energy*cutOverEEE_ + energy*energy*cutOverE2EE_) ) {
193  n++;
194  filterproduct.addObject(trigger_type, ref);
195  continue;
196  }
197  }
198  }
199 
200  // filter decision
201  bool accept(n>=ncandcut_);
202 
203  return accept;
204 }
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:30
void addObject(int id, const reco::RecoEcalCandidateRef &ref)
setters for L3 collections: (id=physics type, and Ref<C>)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > candToken_
edm::EDGetTokenT< reco::RecoEcalCandidateIsolationMap > varToken_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T const * product() const
Definition: Handle.h:74
void addCollectionTag(const edm::InputTag &collectionTag)
collectionTags
bool saveTags() const
Definition: HLTFilter.h:45

Member Data Documentation

std::vector<double> HLTEgammaGenericQuadraticFilter::absEtaLowEdges_
private

Definition at line 61 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

edm::InputTag HLTEgammaGenericQuadraticFilter::candTag_
private

Definition at line 32 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter().

edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> HLTEgammaGenericQuadraticFilter::candToken_
private

Definition at line 34 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

bool HLTEgammaGenericQuadraticFilter::doRhoCorrection_
private

Definition at line 57 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::effectiveAreas_
private

Definition at line 60 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::energyLowEdges_
private

Definition at line 44 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

edm::InputTag HLTEgammaGenericQuadraticFilter::l1EGTag_
private

Definition at line 53 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

bool HLTEgammaGenericQuadraticFilter::lessThan_
private

Definition at line 37 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

int HLTEgammaGenericQuadraticFilter::ncandcut_
private

Definition at line 51 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

double HLTEgammaGenericQuadraticFilter::rhoMax_
private

Definition at line 58 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

double HLTEgammaGenericQuadraticFilter::rhoScale_
private

Definition at line 59 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

edm::InputTag HLTEgammaGenericQuadraticFilter::rhoTag_
private

Definition at line 55 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter().

edm::EDGetTokenT<double> HLTEgammaGenericQuadraticFilter::rhoToken_
private

Definition at line 56 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrOverE2EB_
private

Definition at line 49 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrOverE2EE_
private

Definition at line 50 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrOverEEB_
private

Definition at line 47 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrOverEEE_
private

Definition at line 48 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrRegularEB_
private

Definition at line 45 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

std::vector<double> HLTEgammaGenericQuadraticFilter::thrRegularEE_
private

Definition at line 46 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

bool HLTEgammaGenericQuadraticFilter::useEt_
private

Definition at line 38 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().

edm::InputTag HLTEgammaGenericQuadraticFilter::varTag_
private

Definition at line 33 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter().

edm::EDGetTokenT<reco::RecoEcalCandidateIsolationMap> HLTEgammaGenericQuadraticFilter::varToken_
private

Definition at line 35 of file HLTEgammaGenericQuadraticFilter.h.

Referenced by HLTEgammaGenericQuadraticFilter(), and hltFilter().