CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronLikelihood.cc
Go to the documentation of this file.
5 
7 #include <iostream>
8 
9 
11  LikelihoodSwitches eleIDSwitches,
12  std::string signalWeightSplitting,
13  std::string backgroundWeightSplitting,
14  bool splitSignalPdfs,
15  bool splitBackgroundPdfs) :
16  _EB0lt15lh (new LikelihoodPdfProduct ("electronID_EB0_ptLt15_likelihood",0,0)) ,
17  _EB1lt15lh (new LikelihoodPdfProduct ("electronID_EB1_ptLt15_likelihood",1,0)) ,
18  _EElt15lh (new LikelihoodPdfProduct ("electronID_EE_ptLt15_likelihood",2,0)) ,
19  _EB0gt15lh (new LikelihoodPdfProduct ("electronID_EB0_ptGt15_likelihood",0,1)) ,
20  _EB1gt15lh (new LikelihoodPdfProduct ("electronID_EB1_ptGt15_likelihood",1,1)) ,
21  _EEgt15lh (new LikelihoodPdfProduct ("electronID_EE_ptGt15_likelihood",2,1)) ,
22  m_eleIDSwitches (eleIDSwitches) ,
23  m_signalWeightSplitting (signalWeightSplitting),
24  m_backgroundWeightSplitting (backgroundWeightSplitting),
25  m_splitSignalPdfs (splitSignalPdfs),
26  m_splitBackgroundPdfs (splitBackgroundPdfs)
27 {
28  Setup (calibration,
29  signalWeightSplitting, backgroundWeightSplitting,
30  splitSignalPdfs, splitBackgroundPdfs) ;
31 }
32 
33 
34 
35 // --------------------------------------------------------
36 
37 
38 
40  delete _EB0lt15lh ;
41  delete _EB1lt15lh ;
42  delete _EElt15lh ;
43  delete _EB0gt15lh ;
44  delete _EB1gt15lh ;
45  delete _EEgt15lh ;
46 }
47 
48 
49 
50 // --------------------------------------------------------
51 
52 
53 void
55  std::string signalWeightSplitting,
56  std::string backgroundWeightSplitting,
57  bool splitSignalPdfs,
58  bool splitBackgroundPdfs)
59 {
60 
61  // ECAL BARREL0 (|eta|<1.0) LIKELIHOOD - Pt < 15 GeV region
62  _EB0lt15lh->initFromDB (calibration) ;
63 
64  _EB0lt15lh->addSpecies ("electrons") ;
65  _EB0lt15lh->addSpecies ("hadrons") ;
66 
67  if(signalWeightSplitting.compare("class")==0) {
68  _EB0lt15lh->setSplitFrac ("electrons", "class0") ;
69  _EB0lt15lh->setSplitFrac ("electrons", "class1") ;
70  }
71  else {
72  throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
73  << " and fullclass (golden / bigbrem / narrow / showering)"
74  << " splitting is implemented right now";
75  }
76 
77  if (m_eleIDSwitches.m_useDeltaPhi) _EB0lt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
78  if (m_eleIDSwitches.m_useDeltaEta) _EB0lt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
79  if (m_eleIDSwitches.m_useEoverP) _EB0lt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
80  if (m_eleIDSwitches.m_useHoverE) _EB0lt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
81  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB0lt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
82  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB0lt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
83  if (m_eleIDSwitches.m_useFBrem) _EB0lt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
84  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB0lt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
85 
86  if(backgroundWeightSplitting.compare("class")==0) {
87  _EB0lt15lh->setSplitFrac ("hadrons", "class0") ;
88  _EB0lt15lh->setSplitFrac ("hadrons", "class1") ;
89  }
90  else {
91  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
92  << " splitting is implemented right now";
93  }
94 
95  if (m_eleIDSwitches.m_useDeltaPhi) _EB0lt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
96  if (m_eleIDSwitches.m_useDeltaEta) _EB0lt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
97  if (m_eleIDSwitches.m_useEoverP) _EB0lt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
98  if (m_eleIDSwitches.m_useHoverE) _EB0lt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
99  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB0lt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
100  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB0lt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
101  if (m_eleIDSwitches.m_useFBrem) _EB0lt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
102  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB0lt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
103 
104  // ECAL BARREL0 (|eta|<1.0) LIKELIHOOD - Pt >= 15 GeV region
105  _EB0gt15lh->initFromDB (calibration) ;
106 
107  _EB0gt15lh->addSpecies ("electrons") ;
108  _EB0gt15lh->addSpecies ("hadrons") ;
109 
110  if(signalWeightSplitting.compare("class")==0) {
111  _EB0gt15lh->setSplitFrac ("electrons", "class0") ;
112  _EB0gt15lh->setSplitFrac ("electrons", "class1") ;
113  }
114  else {
115  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
116  << " splitting is implemented right now";
117  }
118 
119  if (m_eleIDSwitches.m_useDeltaPhi) _EB0gt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
120  if (m_eleIDSwitches.m_useDeltaEta) _EB0gt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
121  if (m_eleIDSwitches.m_useEoverP) _EB0gt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
122  if (m_eleIDSwitches.m_useHoverE) _EB0gt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
123  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB0gt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
124  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB0gt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
125  if (m_eleIDSwitches.m_useFBrem) _EB0gt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
126  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB0gt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
127 
128  if(backgroundWeightSplitting.compare("class")==0) {
129  _EB0gt15lh->setSplitFrac ("hadrons", "class0") ;
130  _EB0gt15lh->setSplitFrac ("hadrons", "class1") ;
131  }
132  else {
133  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
134  << " splitting is implemented right now";
135  }
136 
137  if (m_eleIDSwitches.m_useDeltaPhi) _EB0gt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
138  if (m_eleIDSwitches.m_useDeltaEta) _EB0gt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
139  if (m_eleIDSwitches.m_useEoverP) _EB0gt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
140  if (m_eleIDSwitches.m_useHoverE) _EB0gt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
141  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB0gt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
142  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB0gt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
143  if (m_eleIDSwitches.m_useFBrem) _EB0gt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
144  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB0gt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
145 
146 
147  // ECAL BARREL1 (|eta|>1.0) LIKELIHOOD - Pt < 15 GeV region
148  _EB1lt15lh->initFromDB (calibration) ;
149 
150  _EB1lt15lh->addSpecies ("electrons") ;
151  _EB1lt15lh->addSpecies ("hadrons") ;
152 
153  if(signalWeightSplitting.compare("class")==0) {
154  _EB1lt15lh->setSplitFrac ("electrons", "class0") ;
155  _EB1lt15lh->setSplitFrac ("electrons", "class1") ;
156  }
157  else {
158  throw cms::Exception("BadConfig") << "Only class (non-showering / showering)"
159  << " and fullclass (golden / bigbrem / narrow / showering)"
160  << " splitting is implemented right now";
161  }
162 
163  if (m_eleIDSwitches.m_useDeltaPhi) _EB1lt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
164  if (m_eleIDSwitches.m_useDeltaEta) _EB1lt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
165  if (m_eleIDSwitches.m_useEoverP) _EB1lt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
166  if (m_eleIDSwitches.m_useHoverE) _EB1lt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
167  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB1lt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
168  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB1lt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
169  if (m_eleIDSwitches.m_useFBrem) _EB1lt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
170  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB1lt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
171 
172  if(backgroundWeightSplitting.compare("class")==0) {
173  _EB1lt15lh->setSplitFrac ("hadrons", "class0") ;
174  _EB1lt15lh->setSplitFrac ("hadrons", "class1") ;
175  }
176  else {
177  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
178  << " splitting is implemented right now";
179  }
180 
181  if (m_eleIDSwitches.m_useDeltaPhi) _EB1lt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
182  if (m_eleIDSwitches.m_useDeltaEta) _EB1lt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
183  if (m_eleIDSwitches.m_useEoverP) _EB1lt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
184  if (m_eleIDSwitches.m_useHoverE) _EB1lt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
185  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB1lt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
186  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB1lt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
187  if (m_eleIDSwitches.m_useFBrem) _EB1lt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
188  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB1lt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
189 
190  // ECAL BARREL1 (|eta|>1.0) LIKELIHOOD - Pt >= 15 GeV region
191  _EB1gt15lh->initFromDB (calibration) ;
192 
193  _EB1gt15lh->addSpecies ("electrons") ;
194  _EB1gt15lh->addSpecies ("hadrons") ;
195 
196  if(signalWeightSplitting.compare("class")==0) {
197  _EB1gt15lh->setSplitFrac ("electrons", "class0") ;
198  _EB1gt15lh->setSplitFrac ("electrons", "class1") ;
199  }
200  else {
201  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
202  << " splitting is implemented right now";
203  }
204 
205  if (m_eleIDSwitches.m_useDeltaPhi) _EB1gt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
206  if (m_eleIDSwitches.m_useDeltaEta) _EB1gt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
207  if (m_eleIDSwitches.m_useEoverP) _EB1gt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
208  if (m_eleIDSwitches.m_useHoverE) _EB1gt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
209  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB1gt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
210  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB1gt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
211  if (m_eleIDSwitches.m_useFBrem) _EB1gt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
212  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB1gt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
213 
214  if(backgroundWeightSplitting.compare("class")==0) {
215  _EB1gt15lh->setSplitFrac ("hadrons", "class0") ;
216  _EB1gt15lh->setSplitFrac ("hadrons", "class1") ;
217  }
218  else {
219  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
220  << " splitting is implemented right now";
221  }
222 
223  if (m_eleIDSwitches.m_useDeltaPhi) _EB1gt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
224  if (m_eleIDSwitches.m_useDeltaEta) _EB1gt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
225  if (m_eleIDSwitches.m_useEoverP) _EB1gt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
226  if (m_eleIDSwitches.m_useHoverE) _EB1gt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
227  if (m_eleIDSwitches.m_useSigmaEtaEta) _EB1gt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
228  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EB1gt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
229  if (m_eleIDSwitches.m_useFBrem) _EB1gt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
230  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EB1gt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
231 
232  // ECAL ENDCAP LIKELIHOOD - Pt < 15 GeV
233  _EElt15lh->initFromDB (calibration) ;
234 
235  _EElt15lh->addSpecies ("electrons") ;
236  _EElt15lh->addSpecies ("hadrons") ;
237 
238  if(signalWeightSplitting.compare("class")==0) {
239  _EElt15lh->setSplitFrac ("electrons", "class0") ;
240  _EElt15lh->setSplitFrac ("electrons", "class1") ;
241  }
242  else {
243  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
244  << " splitting is implemented right now";
245  }
246 
247  if (m_eleIDSwitches.m_useDeltaPhi) _EElt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
248  if (m_eleIDSwitches.m_useDeltaEta) _EElt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
249  if (m_eleIDSwitches.m_useEoverP) _EElt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
250  if (m_eleIDSwitches.m_useHoverE) _EElt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
251  if (m_eleIDSwitches.m_useSigmaEtaEta) _EElt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
252  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EElt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
253  if (m_eleIDSwitches.m_useFBrem) _EElt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
254  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EElt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
255 
256  if(backgroundWeightSplitting.compare("class")==0) {
257  _EElt15lh->setSplitFrac ("hadrons", "class0") ;
258  _EElt15lh->setSplitFrac ("hadrons", "class1") ;
259  }
260  else {
261  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
262  << " splitting is implemented right now";
263  }
264 
265  if (m_eleIDSwitches.m_useDeltaPhi) _EElt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
266  if (m_eleIDSwitches.m_useDeltaEta) _EElt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
267  if (m_eleIDSwitches.m_useEoverP) _EElt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
268  if (m_eleIDSwitches.m_useHoverE) _EElt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
269  if (m_eleIDSwitches.m_useSigmaEtaEta) _EElt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
270  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EElt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
271  if (m_eleIDSwitches.m_useFBrem) _EElt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
272  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EElt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
273 
274  // ECAL ENDCAP LIKELIHOOD - Pt >= 15 GeV
275  _EEgt15lh->initFromDB (calibration) ;
276 
277  _EEgt15lh->addSpecies ("electrons") ;
278  _EEgt15lh->addSpecies ("hadrons") ;
279 
280  if(signalWeightSplitting.compare("class")==0) {
281  _EEgt15lh->setSplitFrac ("electrons", "class0") ;
282  _EEgt15lh->setSplitFrac ("electrons", "class1") ;
283  }
284  else {
285  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
286  << " splitting is implemented right now";
287  }
288 
289  if (m_eleIDSwitches.m_useDeltaPhi) _EEgt15lh->addPdf ("electrons", "dPhi", splitSignalPdfs) ;
290  if (m_eleIDSwitches.m_useDeltaEta) _EEgt15lh->addPdf ("electrons", "dEta", splitSignalPdfs) ;
291  if (m_eleIDSwitches.m_useEoverP) _EEgt15lh->addPdf ("electrons", "EoP", splitSignalPdfs) ;
292  if (m_eleIDSwitches.m_useHoverE) _EEgt15lh->addPdf ("electrons", "HoE", splitSignalPdfs) ;
293  if (m_eleIDSwitches.m_useSigmaEtaEta) _EEgt15lh->addPdf ("electrons", "sigmaIEtaIEta", splitSignalPdfs) ;
294  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EEgt15lh->addPdf ("electrons", "sigmaIPhiIPhi", splitSignalPdfs) ;
295  if (m_eleIDSwitches.m_useFBrem) _EEgt15lh->addPdf ("electrons", "fBrem", splitSignalPdfs) ;
296  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EEgt15lh->addPdf ("electrons", "OneOverEMinusOneOverP", splitSignalPdfs) ;
297 
298  if(backgroundWeightSplitting.compare("class")==0) {
299  _EEgt15lh->setSplitFrac ("hadrons", "class0") ;
300  _EEgt15lh->setSplitFrac ("hadrons", "class1") ;
301  }
302  else {
303  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
304  << " splitting is implemented right now";
305  }
306 
307  if (m_eleIDSwitches.m_useDeltaPhi) _EEgt15lh->addPdf ("hadrons", "dPhi", splitBackgroundPdfs) ;
308  if (m_eleIDSwitches.m_useDeltaEta) _EEgt15lh->addPdf ("hadrons", "dEta", splitBackgroundPdfs) ;
309  if (m_eleIDSwitches.m_useEoverP) _EEgt15lh->addPdf ("hadrons", "EoP", splitBackgroundPdfs) ;
310  if (m_eleIDSwitches.m_useHoverE) _EEgt15lh->addPdf ("hadrons", "HoE", splitBackgroundPdfs) ;
311  if (m_eleIDSwitches.m_useSigmaEtaEta) _EEgt15lh->addPdf ("hadrons", "sigmaIEtaIEta", splitBackgroundPdfs) ;
312  if (m_eleIDSwitches.m_useSigmaPhiPhi) _EEgt15lh->addPdf ("hadrons", "sigmaIPhiIPhi", splitBackgroundPdfs) ;
313  if (m_eleIDSwitches.m_useFBrem) _EEgt15lh->addPdf ("hadrons", "fBrem", splitBackgroundPdfs) ;
314  if (m_eleIDSwitches.m_useOneOverEMinusOneOverP) _EEgt15lh->addPdf ("hadrons", "OneOverEMinusOneOverP", splitBackgroundPdfs) ;
315 }
316 
317 
318 
319 // --------------------------------------------------------
320 
321 
322 
323 void
325  std::vector<float> &measurements,
326  const EcalClusterLazyTools& _myEcalCluster) const
327 {
328  EcalClusterLazyTools myEcalCluster = _myEcalCluster;
329  // the variables entering the likelihood
330  if (m_eleIDSwitches.m_useDeltaPhi) measurements.push_back ( electron.deltaPhiSuperClusterTrackAtVtx () ) ;
331  if (m_eleIDSwitches.m_useDeltaEta) measurements.push_back ( electron.deltaEtaSuperClusterTrackAtVtx () ) ;
332  if (m_eleIDSwitches.m_useEoverP) measurements.push_back ( electron.eSuperClusterOverP () ) ;
333  if (m_eleIDSwitches.m_useHoverE) measurements.push_back ( electron.hadronicOverEm () ) ;
334  std::vector<float> vCov = myEcalCluster.localCovariances(*(electron.superCluster()->seed())) ;
335  if (m_eleIDSwitches.m_useSigmaEtaEta) measurements.push_back ( sqrt (vCov[0]) );
336  if (m_eleIDSwitches.m_useSigmaPhiPhi) measurements.push_back ( sqrt (vCov[2]) );
337  if(m_eleIDSwitches.m_useFBrem) measurements.push_back( electron.fbrem() );
338  // 1/E - 1/P calculated consistently with the variables used to make the PDFs
339  reco::GsfTrackRef trkRef = electron.get<reco::GsfTrackRef>();
340  float OneOverEMinusOneOverP = 1.0/(electron.eSuperClusterOverP() * trkRef->p()) - 1.0/trkRef->p();
341  if(m_eleIDSwitches.m_useOneOverEMinusOneOverP) measurements.push_back( OneOverEMinusOneOverP );
342 
343 }
344 
345 
346 
347 // --------------------------------------------------------
348 
349 
350 
351 float
353  const EcalClusterLazyTools& myEcalCluster) const
354 {
355 
356  //=======================================================
357  // used classification:
358  // nbrem clusters = 0 => 0
359  // nbrem clusters >= 1 => 1
360  //=======================================================
361 
362  std::vector<float> measurements ;
363  getInputVar (electron, measurements, myEcalCluster) ;
364 
365  // Split using only the number of brem clusters
366  int bitVal=(electron.numberOfBrems()==0) ? 0 : 1 ;
367 
368  char className[20];
369  if(m_signalWeightSplitting.compare("class")==0) {
370  snprintf(className, 20, "class%d", bitVal);
371  } else {
372  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
373  << " splitting is implemented right now";
374  }
375 
376  reco::SuperClusterRef sclusRef = electron.superCluster() ;
377  EcalSubdetector subdet = EcalSubdetector (sclusRef->hitsAndFractions()[0].first.subdetId ()) ;
378  float thisPt = electron.pt();
379 
380  if (subdet==EcalBarrel && fabs(electron.eta())<=1.0 && thisPt<15.)
381  return _EB0lt15lh->getRatio ("electrons",measurements,std::string (className)) ;
382  else if (subdet==EcalBarrel && fabs(electron.eta())<=1.0 && thisPt>=15.)
383  return _EB0gt15lh->getRatio ("electrons",measurements,std::string (className)) ;
384  else if (subdet==EcalBarrel && fabs(electron.eta())>1.0 && thisPt<15.)
385  return _EB1lt15lh->getRatio ("electrons",measurements,std::string (className)) ;
386  else if (subdet==EcalBarrel && fabs(electron.eta())>1.0 && thisPt>=15.)
387  return _EB1gt15lh->getRatio ("electrons",measurements,std::string (className)) ;
388  else if (subdet==EcalEndcap && thisPt<15.)
389  return _EElt15lh->getRatio ("electrons",measurements,std::string (className)) ;
390  else if (subdet==EcalEndcap && thisPt>=15.)
391  return _EEgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
392  else return -999. ;
393 }
394 
395 float
397  const EcalClusterLazyTools& myEcalCluster) const
398 {
399 
400  //=======================================================
401  // used classification:
402  // nbrem clusters = 0 => 0
403  // nbrem clusters >= 1 => 1
404  //=======================================================
405 
406  std::vector<float> measurements ;
407  getInputVar (electron, measurements, myEcalCluster) ;
408 
409  // Split using only the number of brem clusters
410  int bitVal=(electron.numberOfBrems()==0) ? 0 : 1 ;
411 
412  char className[20];
413  if(m_signalWeightSplitting.compare("class")==0) {
414  snprintf(className, 20, "class%d", bitVal);
415  } else {
416  throw cms::Exception("BadConfig") << "Only class (0 brem clusters / >=1 brem clusters)"
417  << " splitting is implemented right now";
418  }
419 
420  reco::SuperClusterRef sclusRef = electron.superCluster() ;
421  EcalSubdetector subdet = EcalSubdetector (sclusRef->hitsAndFractions()[0].first.subdetId ()) ;
422  float thisPt = electron.pt();
423 
424  float lh=-999.;
425 
426  if (subdet==EcalBarrel && fabs(electron.eta())<=1.0 && thisPt<15.)
427  lh = _EB0lt15lh->getRatio ("electrons",measurements,std::string (className)) ;
428  else if (subdet==EcalBarrel && fabs(electron.eta())<=1.0 && thisPt>=15.)
429  lh = _EB0gt15lh->getRatio ("electrons",measurements,std::string (className)) ;
430  else if (subdet==EcalBarrel && fabs(electron.eta())>1.0 && thisPt<15.)
431  lh = _EB1lt15lh->getRatio ("electrons",measurements,std::string (className)) ;
432  else if (subdet==EcalBarrel && fabs(electron.eta())>1.0 && thisPt>=15.)
433  lh = _EB1gt15lh->getRatio ("electrons",measurements,std::string (className)) ;
434  else if (subdet==EcalEndcap && thisPt<15.)
435  lh = _EElt15lh->getRatio ("electrons",measurements,std::string (className)) ;
436  else if (subdet==EcalEndcap && thisPt>=15.)
437  lh = _EEgt15lh->getRatio ("electrons",measurements,std::string (className)) ;
438  else lh = -999. ;
439 
440  if(lh<=0) return -20.;
441  else if(lh==1) return 20.;
442  else return log(lh/(1.0-lh));
443 
444 }
445 
virtual ~ElectronLikelihood()
dtor
LikelihoodPdfProduct * _EB0gt15lh
likelihood above 15GeV/c
void initFromDB(const ElectronLikelihoodCalibration *calibration)
initialize the PDFs from CondDB
float result(const reco::GsfElectron &electron, const EcalClusterLazyTools &) const
get the result of the algorithm
float eSuperClusterOverP() const
Definition: GsfElectron.h:230
void setSplitFrac(const char *specname, const char *catName, float frac=1.0)
set the fraction of one category for a given species
ElectronLikelihood()
ctor, not used for this algo (need initialization from ES)
void addSpecies(const char *name, float priorWeight=1.)
add a species (hypothesis) to the likelihood, with a priori probability
LikelihoodPdfProduct * _EB0lt15lh
likelihood below 15GeV/c
float fbrem() const
Definition: GsfElectron.h:653
bool int lh
Definition: SIMDVec.h:19
LikelihoodPdfProduct * _EElt15lh
LikelihoodPdfProduct * _EB1gt15lh
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:234
float hadronicOverEm() const
Definition: GsfElectron.h:410
LikelihoodSwitches m_eleIDSwitches
general parameters of all the ele id algorithms
T sqrt(T t)
Definition: SSEVec.h:48
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:169
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:237
int numberOfBrems() const
Definition: GsfElectron.h:652
T get() const
get a component
float resultLog(const reco::GsfElectron &electron, const EcalClusterLazyTools &) const
get the log-expanded result of the algorithm
void getInputVar(const reco::GsfElectron &electron, std::vector< float > &measuremnts, const EcalClusterLazyTools &) const
get the input variables from the electron and the e-Setup
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
std::string m_signalWeightSplitting
splitting rule for PDF&#39;s
void addPdf(const char *specname, const char *name, bool splitPdf=false)
add a pdf for a species, splitted or not
LikelihoodPdfProduct * _EB1lt15lh
float getRatio(const char *specName, const std::vector< float > &measurements, std::string)
get the likelihood ratio p(a priori) * L(specName) / L_tot
void Setup(const ElectronLikelihoodCalibration *calibration, std::string signalWeightSplitting, std::string backgroundWeightSplitting, bool splitSignalPdfs, bool splitBackgroundPdfs)
LikelihoodPdfProduct * _EEgt15lh
virtual float pt() const GCC11_FINAL
transverse momentum
EcalSubdetector
std::string className(const T &t)
Definition: ClassName.h:30