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