CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
ElectronLikelihood Class Reference

#include <ElectronLikelihood.h>

Public Member Functions

 ElectronLikelihood ()
 ctor, not used for this algo (need initialization from ES) More...
 
 ElectronLikelihood (const ElectronLikelihoodCalibration *calibration, LikelihoodSwitches eleIDSwitches, std::string signalWeightSplitting, std::string backgroundWeightSplitting, bool splitSignalPdfs, bool splitBackgroundPdfs)
 ctor More...
 
float result (const reco::GsfElectron &electron, EcalClusterLazyTools) const
 get the result of the algorithm More...
 
void setup (const edm::ParameterSet &conf)
 not used for this algo More...
 
virtual ~ElectronLikelihood ()
 dtor More...
 

Private Member Functions

void getInputVar (const reco::GsfElectron &electron, std::vector< float > &measuremnts, EcalClusterLazyTools) const
 get the input variables from the electron and the e-Setup More...
 
void Setup (const ElectronLikelihoodCalibration *calibration, std::string signalWeightSplitting, std::string backgroundWeightSplitting, bool splitSignalPdfs, bool splitBackgroundPdfs)
 

Private Attributes

LikelihoodPdfProduct_EBgt15lh
 likelihood above 15GeV/c More...
 
LikelihoodPdfProduct_EBlt15lh
 likelihood below 15GeV/c More...
 
LikelihoodPdfProduct_EEgt15lh
 
LikelihoodPdfProduct_EElt15lh
 
std::string m_backgroundWeightSplitting
 
LikelihoodSwitches m_eleIDSwitches
 general parameters of all the ele id algorithms More...
 
std::string m_signalWeightSplitting
 splitting rule for PDF's More...
 
bool m_splitBackgroundPdfs
 
bool m_splitSignalPdfs
 

Detailed Description

Definition at line 17 of file ElectronLikelihood.h.

Constructor & Destructor Documentation

ElectronLikelihood::ElectronLikelihood ( )
inline

ctor, not used for this algo (need initialization from ES)

Definition at line 22 of file ElectronLikelihood.h.

22 {} ;
ElectronLikelihood::ElectronLikelihood ( const ElectronLikelihoodCalibration calibration,
LikelihoodSwitches  eleIDSwitches,
std::string  signalWeightSplitting,
std::string  backgroundWeightSplitting,
bool  splitSignalPdfs,
bool  splitBackgroundPdfs 
)

ctor

Definition at line 8 of file ElectronLikelihood.cc.

References Setup().

13  :
14  _EBlt15lh (new LikelihoodPdfProduct ("electronID_EB_ptLt15_likelihood",0,0)) ,
15  _EElt15lh (new LikelihoodPdfProduct ("electronID_EE_ptLt15_likelihood",1,0)) ,
16  _EBgt15lh (new LikelihoodPdfProduct ("electronID_EB_ptGt15_likelihood",0,1)) ,
17  _EEgt15lh (new LikelihoodPdfProduct ("electronID_EE_ptGt15_likelihood",1,1)) ,
18  m_eleIDSwitches (eleIDSwitches) ,
19  m_signalWeightSplitting (signalWeightSplitting),
20  m_backgroundWeightSplitting (backgroundWeightSplitting),
21  m_splitSignalPdfs (splitSignalPdfs),
22  m_splitBackgroundPdfs (splitBackgroundPdfs)
23 {
24  Setup (calibration,
25  signalWeightSplitting, backgroundWeightSplitting,
26  splitSignalPdfs, splitBackgroundPdfs) ;
27 }
std::string m_backgroundWeightSplitting
LikelihoodPdfProduct * _EBgt15lh
likelihood above 15GeV/c
LikelihoodPdfProduct * _EElt15lh
LikelihoodSwitches m_eleIDSwitches
general parameters of all the ele id algorithms
std::string m_signalWeightSplitting
splitting rule for PDF&#39;s
void Setup(const ElectronLikelihoodCalibration *calibration, std::string signalWeightSplitting, std::string backgroundWeightSplitting, bool splitSignalPdfs, bool splitBackgroundPdfs)
LikelihoodPdfProduct * _EEgt15lh
LikelihoodPdfProduct * _EBlt15lh
likelihood below 15GeV/c
ElectronLikelihood::~ElectronLikelihood ( )
virtual

dtor

Definition at line 35 of file ElectronLikelihood.cc.

References _EBgt15lh, _EBlt15lh, _EEgt15lh, and _EElt15lh.

35  {
36  delete _EBlt15lh ;
37  delete _EElt15lh ;
38  delete _EBgt15lh ;
39  delete _EEgt15lh ;
40 }
LikelihoodPdfProduct * _EBgt15lh
likelihood above 15GeV/c
LikelihoodPdfProduct * _EElt15lh
LikelihoodPdfProduct * _EEgt15lh
LikelihoodPdfProduct * _EBlt15lh
likelihood below 15GeV/c

Member Function Documentation

void ElectronLikelihood::getInputVar ( const reco::GsfElectron electron,
std::vector< float > &  measuremnts,
EcalClusterLazyTools  myEcalCluster 
) const
private

get the input variables from the electron and the e-Setup

Definition at line 224 of file ElectronLikelihood.cc.

References EcalClusterLazyTools::covariances(), reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::eSuperClusterOverP(), reco::GsfElectron::fbrem(), reco::GsfElectron::hadronicOverEm(), m_eleIDSwitches, LikelihoodSwitches::m_useDeltaEta, LikelihoodSwitches::m_useDeltaPhi, LikelihoodSwitches::m_useEoverP, LikelihoodSwitches::m_useFBrem, LikelihoodSwitches::m_useHoverE, LikelihoodSwitches::m_useSigmaEtaEta, LikelihoodSwitches::m_useSigmaPhiPhi, mathSSE::sqrt(), and reco::GsfElectron::superCluster().

Referenced by result().

227 {
228 
229  // the variables entering the likelihood
230  if (m_eleIDSwitches.m_useDeltaPhi) measurements.push_back ( electron.deltaPhiSuperClusterTrackAtVtx () ) ;
231  if (m_eleIDSwitches.m_useDeltaEta) measurements.push_back ( electron.deltaEtaSuperClusterTrackAtVtx () ) ;
232  if (m_eleIDSwitches.m_useEoverP) measurements.push_back ( electron.eSuperClusterOverP () ) ;
233  if (m_eleIDSwitches.m_useHoverE) measurements.push_back ( electron.hadronicOverEm () ) ;
234  std::vector<float> vCov = myEcalCluster.covariances(*(electron.superCluster()->seed())) ;
235  if (m_eleIDSwitches.m_useSigmaEtaEta) measurements.push_back ( sqrt (vCov[0]) );
236  if (m_eleIDSwitches.m_useSigmaPhiPhi) measurements.push_back ( sqrt (vCov[2]) );
237  if(m_eleIDSwitches.m_useFBrem) measurements.push_back( electron.fbrem() );
238 
239 }
float eSuperClusterOverP() const
Definition: GsfElectron.h:215
float fbrem() const
Definition: GsfElectron.h:600
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:219
float hadronicOverEm() const
Definition: GsfElectron.h:399
LikelihoodSwitches m_eleIDSwitches
general parameters of all the ele id algorithms
T sqrt(T t)
Definition: SSEVec.h:28
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:222
std::vector< float > covariances(const reco::BasicCluster &cluster, float w0=4.7)
float ElectronLikelihood::result ( const reco::GsfElectron electron,
EcalClusterLazyTools  myEcalCluster 
) const

get the result of the algorithm

Definition at line 248 of file ElectronLikelihood.cc.

References _EBgt15lh, _EBlt15lh, _EEgt15lh, _EElt15lh, className(), EcalBarrel, EcalEndcap, edm::hlt::Exception, getInputVar(), LikelihoodPdfProduct::getRatio(), m_signalWeightSplitting, reco::GsfElectron::numberOfBrems(), reco::LeafCandidate::pt(), and reco::GsfElectron::superCluster().

250 {
251 
252  //=======================================================
253  // used classification:
254  // nbrem clusters = 0 => 0
255  // nbrem clusters >= 1 => 1
256  //=======================================================
257 
258  std::vector<float> measurements ;
259  getInputVar (electron, measurements, myEcalCluster) ;
260 
261  // Split using only the number of brem clusters
262  int bitVal=(electron.numberOfBrems()==0) ? 0 : 1 ;
263 
264  char className[20] ;
265  if(m_signalWeightSplitting.compare("class")==0) {
266  sprintf (className,"class%d",bitVal);
267  } else {
268  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
269  << " splitting is implemented right now";
270  }
271 
272  reco::SuperClusterRef sclusRef = electron.superCluster() ;
273  EcalSubdetector subdet = EcalSubdetector (sclusRef->hitsAndFractions()[0].first.subdetId ()) ;
274  float thisPt = electron.pt();
275 
276  if (subdet==EcalBarrel && thisPt<15.)
277  return _EBlt15lh->getRatio ("electrons",measurements,std::string (className)) ;
278  else if (subdet==EcalBarrel && thisPt>=15.)
279  return _EBgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
280  else if (subdet==EcalEndcap && thisPt<15.)
281  return _EElt15lh->getRatio ("electrons",measurements,std::string (className)) ;
282  else if (subdet==EcalEndcap && thisPt>=15.)
283  return _EEgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
284  else return -999. ;
285 }
void getInputVar(const reco::GsfElectron &electron, std::vector< float > &measuremnts, EcalClusterLazyTools) const
get the input variables from the electron and the e-Setup
LikelihoodPdfProduct * _EBgt15lh
likelihood above 15GeV/c
float getRatio(const char *specName, std::vector< float > measurements, std::string)
get the likelihood ratio p(a priori) * L(specName) / L_tot
SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:168
LikelihoodPdfProduct * _EElt15lh
int numberOfBrems() const
Definition: GsfElectron.h:601
std::string m_signalWeightSplitting
splitting rule for PDF&#39;s
virtual double pt() const
transverse momentum
LikelihoodPdfProduct * _EEgt15lh
LikelihoodPdfProduct * _EBlt15lh
likelihood below 15GeV/c
EcalSubdetector
std::string className(const T &t)
Definition: ClassName.h:30
void ElectronLikelihood::setup ( const edm::ParameterSet conf)
inline

not used for this algo

Definition at line 36 of file ElectronLikelihood.h.

36 {} ;
void ElectronLikelihood::Setup ( const ElectronLikelihoodCalibration calibration,
std::string  signalWeightSplitting,
std::string  backgroundWeightSplitting,
bool  splitSignalPdfs,
bool  splitBackgroundPdfs 
)
private

build the likelihood model from histograms in Barrel file and Endcap file

Definition at line 48 of file ElectronLikelihood.cc.

References _EBgt15lh, _EBlt15lh, _EEgt15lh, _EElt15lh, LikelihoodPdfProduct::addPdf(), LikelihoodPdfProduct::addSpecies(), edm::hlt::Exception, LikelihoodPdfProduct::initFromDB(), m_eleIDSwitches, LikelihoodSwitches::m_useDeltaEta, LikelihoodSwitches::m_useDeltaPhi, LikelihoodSwitches::m_useEoverP, LikelihoodSwitches::m_useFBrem, LikelihoodSwitches::m_useHoverE, LikelihoodSwitches::m_useSigmaEtaEta, LikelihoodSwitches::m_useSigmaPhiPhi, and LikelihoodPdfProduct::setSplitFrac().

Referenced by ElectronLikelihood().

53 {
54 
55  // ECAL BARREL LIKELIHOOD - Pt < 15 GeV region
56  _EBlt15lh->initFromDB (calibration) ;
57 
58  _EBlt15lh->addSpecies ("electrons") ;
59  _EBlt15lh->addSpecies ("hadrons") ;
60 
61  if(signalWeightSplitting.compare("class")==0) {
62  _EBlt15lh->setSplitFrac ("electrons", "class0") ;
63  _EBlt15lh->setSplitFrac ("electrons", "class1") ;
64  }
65  else {
66  throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
67  << " and fullclass (golden / bigbrem / narrow / showering)"
68  << " splitting is implemented right now";
69  }
70 
71  if (m_eleIDSwitches.m_useDeltaPhi) _EBlt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
72  if (m_eleIDSwitches.m_useDeltaEta) _EBlt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
73  if (m_eleIDSwitches.m_useEoverP) _EBlt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
74  if (m_eleIDSwitches.m_useHoverE) _EBlt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
75  if (m_eleIDSwitches.m_useSigmaEtaEta) _EBlt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
76  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EBlt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
77  if (m_eleIDSwitches.m_useFBrem) _EBlt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
78 
79  if(backgroundWeightSplitting.compare("class")==0) {
80  _EBlt15lh->setSplitFrac ("hadrons", "class0") ;
81  _EBlt15lh->setSplitFrac ("hadrons", "class1") ;
82  }
83  else {
84  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
85  << " splitting is implemented right now";
86  }
87 
88  if (m_eleIDSwitches.m_useDeltaPhi) _EBlt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
89  if (m_eleIDSwitches.m_useDeltaEta) _EBlt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
90  if (m_eleIDSwitches.m_useEoverP) _EBlt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
91  if (m_eleIDSwitches.m_useHoverE) _EBlt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
92  if (m_eleIDSwitches.m_useSigmaEtaEta) _EBlt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
93  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EBlt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
94  if (m_eleIDSwitches.m_useFBrem) _EBlt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
95 
96  // ECAL BARREL LIKELIHOOD - Pt >= 15 GeV region
97  _EBgt15lh->initFromDB (calibration) ;
98 
99  _EBgt15lh->addSpecies ("electrons") ;
100  _EBgt15lh->addSpecies ("hadrons") ;
101 
102  if(signalWeightSplitting.compare("class")==0) {
103  _EBgt15lh->setSplitFrac ("electrons", "class0") ;
104  _EBgt15lh->setSplitFrac ("electrons", "class1") ;
105  }
106  else {
107  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
108  << " splitting is implemented right now";
109  }
110 
111  if (m_eleIDSwitches.m_useDeltaPhi) _EBgt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
112  if (m_eleIDSwitches.m_useDeltaEta) _EBgt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
113  if (m_eleIDSwitches.m_useEoverP) _EBgt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
114  if (m_eleIDSwitches.m_useHoverE) _EBgt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
115  if (m_eleIDSwitches.m_useSigmaEtaEta) _EBgt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
116  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EBgt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
117  if (m_eleIDSwitches.m_useFBrem) _EBgt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
118 
119  if(backgroundWeightSplitting.compare("class")==0) {
120  _EBgt15lh->setSplitFrac ("hadrons", "class0") ;
121  _EBgt15lh->setSplitFrac ("hadrons", "class1") ;
122  }
123  else {
124  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
125  << " splitting is implemented right now";
126  }
127 
128  if (m_eleIDSwitches.m_useDeltaPhi) _EBgt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
129  if (m_eleIDSwitches.m_useDeltaEta) _EBgt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
130  if (m_eleIDSwitches.m_useEoverP) _EBgt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
131  if (m_eleIDSwitches.m_useHoverE) _EBgt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
132  if (m_eleIDSwitches.m_useSigmaEtaEta) _EBgt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
133  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EBgt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
134  if (m_eleIDSwitches.m_useFBrem) _EBgt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
135 
136  // ECAL ENDCAP LIKELIHOOD - Pt < 15 GeV
137  _EElt15lh->initFromDB (calibration) ;
138 
139  _EElt15lh->addSpecies ("electrons") ;
140  _EElt15lh->addSpecies ("hadrons") ;
141 
142  if(signalWeightSplitting.compare("class")==0) {
143  _EElt15lh->setSplitFrac ("electrons", "class0") ;
144  _EElt15lh->setSplitFrac ("electrons", "class1") ;
145  }
146  else {
147  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
148  << " splitting is implemented right now";
149  }
150 
151  if (m_eleIDSwitches.m_useDeltaPhi) _EElt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
152  if (m_eleIDSwitches.m_useDeltaEta) _EElt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
153  if (m_eleIDSwitches.m_useEoverP) _EElt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
154  if (m_eleIDSwitches.m_useHoverE) _EElt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
155  if (m_eleIDSwitches.m_useSigmaEtaEta) _EElt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
156  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EElt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
157  if (m_eleIDSwitches.m_useFBrem) _EElt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
158 
159  if(backgroundWeightSplitting.compare("class")==0) {
160  _EElt15lh->setSplitFrac ("hadrons", "class0") ;
161  _EElt15lh->setSplitFrac ("hadrons", "class1") ;
162  }
163  else {
164  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
165  << " splitting is implemented right now";
166  }
167 
168  if (m_eleIDSwitches.m_useDeltaPhi) _EElt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
169  if (m_eleIDSwitches.m_useDeltaEta) _EElt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
170  if (m_eleIDSwitches.m_useEoverP) _EElt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
171  if (m_eleIDSwitches.m_useHoverE) _EElt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
172  if (m_eleIDSwitches.m_useSigmaEtaEta) _EElt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
173  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EElt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
174  if (m_eleIDSwitches.m_useFBrem) _EElt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
175 
176  // ECAL ENDCAP LIKELIHOOD - Pt >= 15 GeV
177  _EEgt15lh->initFromDB (calibration) ;
178 
179  _EEgt15lh->addSpecies ("electrons") ;
180  _EEgt15lh->addSpecies ("hadrons") ;
181 
182  if(signalWeightSplitting.compare("class")==0) {
183  _EEgt15lh->setSplitFrac ("electrons", "class0") ;
184  _EEgt15lh->setSplitFrac ("electrons", "class1") ;
185  }
186  else {
187  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
188  << " splitting is implemented right now";
189  }
190 
191  if (m_eleIDSwitches.m_useDeltaPhi) _EEgt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
192  if (m_eleIDSwitches.m_useDeltaEta) _EEgt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
193  if (m_eleIDSwitches.m_useEoverP) _EEgt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
194  if (m_eleIDSwitches.m_useHoverE) _EEgt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
195  if (m_eleIDSwitches.m_useSigmaEtaEta) _EEgt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
196  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EEgt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
197  if (m_eleIDSwitches.m_useFBrem) _EEgt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
198 
199  if(backgroundWeightSplitting.compare("class")==0) {
200  _EEgt15lh->setSplitFrac ("hadrons", "class0") ;
201  _EEgt15lh->setSplitFrac ("hadrons", "class1") ;
202  }
203  else {
204  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
205  << " splitting is implemented right now";
206  }
207 
208  if (m_eleIDSwitches.m_useDeltaPhi) _EEgt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
209  if (m_eleIDSwitches.m_useDeltaEta) _EEgt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
210  if (m_eleIDSwitches.m_useEoverP) _EEgt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
211  if (m_eleIDSwitches.m_useHoverE) _EEgt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
212  if (m_eleIDSwitches.m_useSigmaEtaEta) _EEgt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
213  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EEgt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
214  if (m_eleIDSwitches.m_useFBrem) _EEgt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
215 }
void initFromDB(const ElectronLikelihoodCalibration *calibration)
initialize the PDFs from CondDB
void setSplitFrac(const char *specname, const char *catName, float frac=1.0)
set the fraction of one category for a given species
LikelihoodPdfProduct * _EBgt15lh
likelihood above 15GeV/c
void addSpecies(const char *name, float priorWeight=1.)
add a species (hypothesis) to the likelihood, with a priori probability
LikelihoodPdfProduct * _EElt15lh
LikelihoodSwitches m_eleIDSwitches
general parameters of all the ele id algorithms
void addPdf(const char *specname, const char *name, bool splitPdf=false)
add a pdf for a species, splitted or not
LikelihoodPdfProduct * _EEgt15lh
LikelihoodPdfProduct * _EBlt15lh
likelihood below 15GeV/c

Member Data Documentation

LikelihoodPdfProduct* ElectronLikelihood::_EBgt15lh
private

likelihood above 15GeV/c

Definition at line 61 of file ElectronLikelihood.h.

Referenced by result(), Setup(), and ~ElectronLikelihood().

LikelihoodPdfProduct* ElectronLikelihood::_EBlt15lh
private

likelihood below 15GeV/c

Definition at line 59 of file ElectronLikelihood.h.

Referenced by result(), Setup(), and ~ElectronLikelihood().

LikelihoodPdfProduct * ElectronLikelihood::_EEgt15lh
private

Definition at line 61 of file ElectronLikelihood.h.

Referenced by result(), Setup(), and ~ElectronLikelihood().

LikelihoodPdfProduct * ElectronLikelihood::_EElt15lh
private

Definition at line 59 of file ElectronLikelihood.h.

Referenced by result(), Setup(), and ~ElectronLikelihood().

std::string ElectronLikelihood::m_backgroundWeightSplitting
private

Definition at line 68 of file ElectronLikelihood.h.

LikelihoodSwitches ElectronLikelihood::m_eleIDSwitches
private

general parameters of all the ele id algorithms

Definition at line 64 of file ElectronLikelihood.h.

Referenced by getInputVar(), and Setup().

std::string ElectronLikelihood::m_signalWeightSplitting
private

splitting rule for PDF's

Definition at line 67 of file ElectronLikelihood.h.

Referenced by result().

bool ElectronLikelihood::m_splitBackgroundPdfs
private

Definition at line 70 of file ElectronLikelihood.h.

bool ElectronLikelihood::m_splitSignalPdfs
private

Definition at line 69 of file ElectronLikelihood.h.