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 Attributes
AntiElectronIDMVA Class Reference

#include <AntiElectronIDMVA.h>

Public Member Functions

 AntiElectronIDMVA ()
 
void Initialize (std::string methodName, std::string oneProng0Pi0_BL, std::string oneProng1pi0wGSF_BL, std::string oneProng1pi0woGSF_BL, std::string oneProng0Pi0_EC, std::string oneProng1pi0wGSF_EC, std::string oneProng1pi0woGSF_EC)
 
double MVAValue (Float_t TauEta, Float_t TauPt, Float_t TauSignalPFChargedCands, Float_t TauSignalPFGammaCands, Float_t TauLeadPFChargedHadrMva, Float_t TauLeadPFChargedHadrHoP, Float_t TauLeadPFChargedHadrEoP, Float_t TauHasGsf, Float_t TauVisMass, Float_t TauEmFraction, std::vector< Float_t > *GammasdEta, std::vector< Float_t > *GammasdPhi, std::vector< Float_t > *GammasPt)
 
double MVAValue (Float_t TauEta, Float_t TauPt, Float_t TauSignalPFChargedCands, Float_t TauSignalPFGammaCands, Float_t TauLeadPFChargedHadrMva, Float_t TauLeadPFChargedHadrHoP, Float_t TauLeadPFChargedHadrEoP, Float_t TauHasGsf, Float_t TauVisMass, Float_t TauEmFraction, Float_t GammaEtaMom, Float_t GammaPhiMom, Float_t GammaEnFrac)
 
double MVAValue (const reco::PFTauRef &thePFTauRef)
 
 ~AntiElectronIDMVA ()
 

Private Attributes

TMVA::Reader * fTMVAReader_ [6]
 
Float_t GammadEta_
 
Float_t GammadPhi_
 
Float_t GammadPt_
 
Bool_t isInitialized_
 
std::string methodName_
 
Float_t TauEmFraction_
 
Float_t TauLeadPFChargedHadrEoP_
 
Float_t TauLeadPFChargedHadrHoP_
 
Float_t TauLeadPFChargedHadrMva_
 
Float_t TauSignalPFGammaCands_
 
Float_t TauVisMass_
 
int verbosity_
 

Detailed Description

Definition at line 33 of file AntiElectronIDMVA.h.

Constructor & Destructor Documentation

AntiElectronIDMVA::AntiElectronIDMVA ( )

Definition at line 6 of file AntiElectronIDMVA.cc.

References fTMVAReader_, i, and verbosity_.

6  :
7  isInitialized_(kFALSE),
8  methodName_("BDT")
9 {
10  for(UInt_t i=0; i<6; ++i) {
11  fTMVAReader_[i] = 0;
12  }
13 
14  verbosity_ = 1;
15 }
int i
Definition: DBlmapReader.cc:9
TMVA::Reader * fTMVAReader_[6]
unsigned int UInt_t
Definition: FUTypes.h:12
AntiElectronIDMVA::~AntiElectronIDMVA ( )

Definition at line 18 of file AntiElectronIDMVA.cc.

References fTMVAReader_, and i.

19 {
20  for(UInt_t i=0; i<6; ++i) {
21  if (fTMVAReader_[i]) delete fTMVAReader_[i];
22  }
23 }
int i
Definition: DBlmapReader.cc:9
TMVA::Reader * fTMVAReader_[6]
unsigned int UInt_t
Definition: FUTypes.h:12

Member Function Documentation

void AntiElectronIDMVA::Initialize ( std::string  methodName,
std::string  oneProng0Pi0_BL,
std::string  oneProng1pi0wGSF_BL,
std::string  oneProng1pi0woGSF_BL,
std::string  oneProng0Pi0_EC,
std::string  oneProng1pi0wGSF_EC,
std::string  oneProng1pi0woGSF_EC 
)

Definition at line 25 of file AntiElectronIDMVA.cc.

References fTMVAReader_, GammadEta_, GammadPhi_, GammadPt_, i, isInitialized_, reco::details::loadTMVAWeights(), methodName_, TauEmFraction_, TauLeadPFChargedHadrEoP_, TauLeadPFChargedHadrHoP_, TauLeadPFChargedHadrMva_, TauSignalPFGammaCands_, TauVisMass_, and verbosity_.

32  {
34  TauVisMass_ = 0.;
35  GammadEta_ = 0.;
36  GammadPhi_ = 0.;
37  GammadPt_ = 0.;
41  TauEmFraction_ = 0.;
42 
43  for(UInt_t i=0; i<6; ++i) {
44  if (fTMVAReader_[i]) delete fTMVAReader_[i];
45  }
46 
47  isInitialized_ = kTRUE;
48  methodName_ = methodName;
49 
50  //TMVA::Tools::Instance();
51 
52  TMVA::Reader *readerX0BL = new TMVA::Reader( "!Color:Silent:Error" );
53  readerX0BL->AddVariable("HoP", &TauLeadPFChargedHadrHoP_);
54  readerX0BL->AddVariable("EoP", &TauLeadPFChargedHadrEoP_);
55  //readerX0BL->AddVariable("emFraction",&TauEmFraction_);
56  readerX0BL->SetVerbose(verbosity_);
57  reco::details::loadTMVAWeights(readerX0BL, methodName_, oneProng0Pi0_BL );
58 
59  TMVA::Reader *reader11BL = new TMVA::Reader( "!Color:Silent:Error" );
60  reader11BL->AddVariable("mva", &TauLeadPFChargedHadrMva_);
61  reader11BL->AddVariable("visMass", &TauVisMass_);
62  reader11BL->AddVariable("etaMom2*TMath::Sqrt(gammaFrac)*pt", &GammadEta_);
63  reader11BL->AddVariable("phiMom2*TMath::Sqrt(gammaFrac)*pt", &GammadPhi_);
64  reader11BL->AddVariable("gammaFrac", &GammadPt_);
65  reader11BL->SetVerbose(verbosity_);
66  reco::details::loadTMVAWeights(reader11BL, methodName_, oneProng1pi0wGSF_BL );
67 
68  TMVA::Reader *reader01BL = new TMVA::Reader( "!Color:Silent:Error" );
69  reader01BL->AddVariable("visMass", &TauVisMass_);
70  reader01BL->AddVariable("etaMom2*TMath::Sqrt(gammaFrac)*pt", &GammadEta_);
71  reader01BL->AddVariable("phiMom2*TMath::Sqrt(gammaFrac)*pt", &GammadPhi_);
72  reader01BL->AddVariable("gammaFrac", &GammadPt_);
73  reader01BL->SetVerbose(verbosity_);
74  reco::details::loadTMVAWeights(reader01BL, methodName_, oneProng1pi0woGSF_BL );
75 
77 
78  TMVA::Reader *readerX0EC = new TMVA::Reader( "!Color:Silent:Error" );
79  readerX0EC->AddVariable("HoP", &TauLeadPFChargedHadrHoP_);
80  readerX0EC->AddVariable("EoP", &TauLeadPFChargedHadrEoP_);
81  //readerX0EC->AddVariable("emFraction",&TauEmFraction_);
82  readerX0EC->SetVerbose(verbosity_);
83  reco::details::loadTMVAWeights(readerX0EC, methodName_, oneProng0Pi0_EC );
84 
85  TMVA::Reader *reader11EC = new TMVA::Reader( "!Color:Silent:Error" );
86  reader11EC->AddVariable("mva", &TauLeadPFChargedHadrMva_);
87  reader11EC->AddVariable("visMass", &TauVisMass_);
88  reader11EC->AddVariable("etaMom2*TMath::Sqrt(gammaFrac)*pt", &GammadEta_);
89  reader11EC->AddVariable("phiMom2*TMath::Sqrt(gammaFrac)*pt", &GammadPhi_);
90  reader11EC->AddVariable("gammaFrac", &GammadPt_);
91  reader11EC->SetVerbose(verbosity_);
92  reco::details::loadTMVAWeights(reader11EC, methodName_, oneProng1pi0wGSF_EC );
93 
94  TMVA::Reader *reader01EC = new TMVA::Reader( "!Color:Silent:Error" );
95  reader01EC->AddVariable("visMass", &TauVisMass_);
96  reader01EC->AddVariable("etaMom2*TMath::Sqrt(gammaFrac)*pt", &GammadEta_);
97  reader01EC->AddVariable("phiMom2*TMath::Sqrt(gammaFrac)*pt", &GammadPhi_);
98  reader01EC->AddVariable("gammaFrac", &GammadPt_);
99  reader01EC->SetVerbose(verbosity_);
100  reco::details::loadTMVAWeights(reader01EC, methodName_, oneProng1pi0woGSF_EC );
101 
102  fTMVAReader_[0] = readerX0BL;
103  fTMVAReader_[1] = reader11BL;
104  fTMVAReader_[2] = reader01BL;
105  fTMVAReader_[3] = readerX0EC;
106  fTMVAReader_[4] = reader11EC;
107  fTMVAReader_[5] = reader01EC;
108 }
int i
Definition: DBlmapReader.cc:9
Float_t TauLeadPFChargedHadrMva_
Float_t TauLeadPFChargedHadrHoP_
TMVA::Reader * fTMVAReader_[6]
Float_t TauLeadPFChargedHadrEoP_
unsigned int UInt_t
Definition: FUTypes.h:12
void loadTMVAWeights(TMVA::Reader *reader, const std::string &method, const std::string &weightFile, bool verbose=false)
double AntiElectronIDMVA::MVAValue ( Float_t  TauEta,
Float_t  TauPt,
Float_t  TauSignalPFChargedCands,
Float_t  TauSignalPFGammaCands,
Float_t  TauLeadPFChargedHadrMva,
Float_t  TauLeadPFChargedHadrHoP,
Float_t  TauLeadPFChargedHadrEoP,
Float_t  TauHasGsf,
Float_t  TauVisMass,
Float_t  TauEmFraction,
std::vector< Float_t > *  GammasdEta,
std::vector< Float_t > *  GammasdPhi,
std::vector< Float_t > *  GammasPt 
)

Definition at line 111 of file AntiElectronIDMVA.cc.

References gather_cfg::cout, dPhi(), fTMVAReader_, GammadEta_, GammadPhi_, GammadPt_, isInitialized_, gen::k, siStripFEDMonitor_P5_cff::Max, methodName_, Pi, TauEmFraction_, TauLeadPFChargedHadrEoP_, TauLeadPFChargedHadrHoP_, TauLeadPFChargedHadrMva_, and TauVisMass_.

117  {
118 
119  if (!isInitialized_) {
120  std::cout << "Error: AntiElectronMVA with method 1 not properly initialized.\n";
121  return -999;
122  }
123 
124  double mva;
125 
126  TauVisMass_ = TauVisMass;
127  TauLeadPFChargedHadrMva_ = TMath::Max(TauLeadPFChargedHadrMva,float(-1.0));
128  TauLeadPFChargedHadrHoP_ = TauLeadPFChargedHadrHoP;
129  TauLeadPFChargedHadrEoP_ = TauLeadPFChargedHadrEoP;
130  TauEmFraction_ = TMath::Max(TauEmFraction,float(0.0));
131 
132  float sumPt = 0;
133  float dEta = 0;
134  float dEta2 = 0;
135  float dPhi = 0;
136  float dPhi2 = 0;
137  float sumPt2 = 0;
138 
139  for(unsigned int k = 0 ; k < GammasPt->size() ; k++){
140  float pt_k = (*GammasPt)[k];
141  float phi_k = (*GammasdPhi)[k];
142  if ((*GammasdPhi)[k] > TMath::Pi()) phi_k = (*GammasdPhi)[k] - 2*TMath::Pi();
143  else if((*GammasdPhi)[k] < -TMath::Pi()) phi_k = (*GammasdPhi)[k] + 2*TMath::Pi();
144  float eta_k = (*GammasdEta)[k];
145  sumPt += pt_k;
146  sumPt2 += (pt_k*pt_k);
147  dEta += (pt_k*eta_k);
148  dEta2 += (pt_k*eta_k*eta_k);
149  dPhi += (pt_k*phi_k);
150  dPhi2 += (pt_k*phi_k*phi_k);
151  }
152 
153  GammadPt_ = sumPt/TauPt;
154 
155  if(sumPt>0){
156  dEta /= sumPt;
157  dPhi /= sumPt;
158  dEta2 /= sumPt;
159  dPhi2 /= sumPt;
160 
161  }
162 
163  //GammadEta_ = dEta;
164  //GammadPhi_ = dPhi;
165 
166  GammadEta_ = TMath::Sqrt(dEta2)*TMath::Sqrt(GammadPt_)*TauPt;
167  GammadPhi_ = TMath::Sqrt(dPhi2)*TMath::Sqrt(GammadPt_)*TauPt;
168 
169 
170  if( TauSignalPFChargedCands==3 )
171  mva = 1.0;
172  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands==0){
173  if(TMath::Abs(TauEta)<1.5)
174  mva = fTMVAReader_[0]->EvaluateMVA( methodName_ );
175  else
176  mva = fTMVAReader_[3]->EvaluateMVA( methodName_ );
177  }
178  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands>0 && TauHasGsf>0.5){
179  if(TMath::Abs(TauEta)<1.5)
180  mva = fTMVAReader_[1]->EvaluateMVA( methodName_ );
181  else
182  mva = fTMVAReader_[4]->EvaluateMVA( methodName_ );
183  }
184  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands>0 && TauHasGsf<0.5){
185  if(TMath::Abs(TauEta)<1.5)
186  mva = fTMVAReader_[2]->EvaluateMVA( methodName_ );
187  else
188  mva = fTMVAReader_[5]->EvaluateMVA( methodName_ );
189  }
190  else{
191  mva = -99;
192  }
193 
194  return mva;
195 
196 }
const double Pi
Float_t TauLeadPFChargedHadrMva_
Float_t TauLeadPFChargedHadrHoP_
TMVA::Reader * fTMVAReader_[6]
Float_t TauLeadPFChargedHadrEoP_
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
int k[5][pyjets_maxn]
tuple cout
Definition: gather_cfg.py:121
double AntiElectronIDMVA::MVAValue ( Float_t  TauEta,
Float_t  TauPt,
Float_t  TauSignalPFChargedCands,
Float_t  TauSignalPFGammaCands,
Float_t  TauLeadPFChargedHadrMva,
Float_t  TauLeadPFChargedHadrHoP,
Float_t  TauLeadPFChargedHadrEoP,
Float_t  TauHasGsf,
Float_t  TauVisMass,
Float_t  TauEmFraction,
Float_t  GammaEtaMom,
Float_t  GammaPhiMom,
Float_t  GammaEnFrac 
)

Definition at line 198 of file AntiElectronIDMVA.cc.

References gather_cfg::cout, fTMVAReader_, GammadEta_, GammadPhi_, GammadPt_, isInitialized_, siStripFEDMonitor_P5_cff::Max, methodName_, TauEmFraction_, TauLeadPFChargedHadrEoP_, TauLeadPFChargedHadrHoP_, TauLeadPFChargedHadrMva_, and TauVisMass_.

204  {
205 
206  if (!isInitialized_) {
207  std::cout << "Error: AntiElectronMVA with method 2 not properly initialized.\n";
208  return -999;
209  }
210 
211 
212  double mva;
213 
214  TauVisMass_ = TauVisMass;
215  TauLeadPFChargedHadrMva_ = TMath::Max(TauLeadPFChargedHadrMva,float(-1.0));
216  TauLeadPFChargedHadrHoP_ = TauLeadPFChargedHadrHoP;
217  TauLeadPFChargedHadrEoP_ = TauLeadPFChargedHadrEoP;
218  TauEmFraction_ = TMath::Max(TauEmFraction,float(0.0));
219  GammadPt_ = GammaEnFrac;
220  GammadEta_ = GammaEtaMom;
221  GammadPhi_ = GammaPhiMom;
222 
223  if( TauSignalPFChargedCands==3 )
224  mva = 1.0;
225  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands==0){
226  if(TMath::Abs(TauEta)<1.5)
227  mva = fTMVAReader_[0]->EvaluateMVA( methodName_ );
228  else
229  mva = fTMVAReader_[3]->EvaluateMVA( methodName_ );
230  }
231  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands>0 && TauHasGsf>0.5){
232  if(TMath::Abs(TauEta)<1.5)
233  mva = fTMVAReader_[1]->EvaluateMVA( methodName_ );
234  else
235  mva = fTMVAReader_[4]->EvaluateMVA( methodName_ );
236  }
237  else if( TauSignalPFChargedCands==1 && TauSignalPFGammaCands>0 && TauHasGsf<0.5){
238  if(TMath::Abs(TauEta)<1.5)
239  mva = fTMVAReader_[2]->EvaluateMVA( methodName_ );
240  else
241  mva = fTMVAReader_[5]->EvaluateMVA( methodName_ );
242  }
243  else{
244  mva = -99.;
245  }
246 
247  return mva;
248 
249 }
Float_t TauLeadPFChargedHadrMva_
Float_t TauLeadPFChargedHadrHoP_
TMVA::Reader * fTMVAReader_[6]
Float_t TauLeadPFChargedHadrEoP_
tuple cout
Definition: gather_cfg.py:121
double AntiElectronIDMVA::MVAValue ( const reco::PFTauRef thePFTauRef)

Definition at line 253 of file AntiElectronIDMVA.cc.

References asciidump::at, gather_cfg::cout, dPhi(), fTMVAReader_, GammadEta_, GammadPhi_, GammadPt_, isInitialized_, gen::k, siStripFEDMonitor_P5_cff::Max, methodName_, Pi, findQualityFiles::size, TauEmFraction_, TauLeadPFChargedHadrEoP_, TauLeadPFChargedHadrHoP_, TauLeadPFChargedHadrMva_, and TauVisMass_.

253  {
254 
255  if (!isInitialized_) {
256  std::cout << "Error: AntiElectronMVA with method 3 not properly initialized.\n";
257  return -999;
258  }
259 
260  double mva;
261 
262  TauVisMass_ = (*thePFTauRef).mass();
263  TauLeadPFChargedHadrMva_ = TMath::Max((*thePFTauRef).electronPreIDOutput(),float(-1.0));
264  TauLeadPFChargedHadrHoP_ = ((*thePFTauRef).leadPFChargedHadrCand())->hcalEnergy()/(*thePFTauRef).leadPFChargedHadrCand()->p();
265  TauLeadPFChargedHadrEoP_ = ((*thePFTauRef).leadPFChargedHadrCand())->ecalEnergy()/(*thePFTauRef).leadPFChargedHadrCand()->p();
266  TauEmFraction_ = TMath::Max((*thePFTauRef).emFraction(),float(0.0));
267 
268  std::vector<float> GammasdEta;
269  std::vector<float> GammasdPhi;
270  std::vector<float> GammasPt;
271 
272  for(unsigned int k = 0 ; k < ((*thePFTauRef).signalPFGammaCands()).size() ; k++){
273  reco::PFCandidateRef gamma = ((*thePFTauRef).signalPFGammaCands()).at(k);
274  if( ((*thePFTauRef).leadPFChargedHadrCand()).isNonnull() ){
275  GammasdEta.push_back( gamma->eta() - (*thePFTauRef).leadPFChargedHadrCand()->eta() );
276  GammasdPhi.push_back( gamma->phi() - (*thePFTauRef).leadPFChargedHadrCand()->phi() );
277  }
278  else{
279  GammasdEta.push_back( gamma->eta() - (*thePFTauRef).eta() );
280  GammasdPhi.push_back( gamma->phi() - (*thePFTauRef).phi() );
281  }
282  GammasPt.push_back( gamma->pt() );
283  }
284 
285  float sumPt = 0;
286  float dEta = 0;
287  float dEta2 = 0;
288  float dPhi = 0;
289  float dPhi2 = 0;
290  float sumPt2 = 0;
291 
292  for(unsigned int k = 0 ; k < GammasPt.size() ; k++){
293  float pt_k = GammasPt[k];
294  float phi_k = GammasdPhi[k];
295  if (GammasdPhi[k] > TMath::Pi()) phi_k = GammasdPhi[k] - 2*TMath::Pi();
296  else if(GammasdPhi[k] < -TMath::Pi()) phi_k = GammasdPhi[k] + 2*TMath::Pi();
297  float eta_k = GammasdEta[k];
298  sumPt += pt_k;
299  sumPt2 += (pt_k*pt_k);
300  dEta += (pt_k*eta_k);
301  dEta2 += (pt_k*eta_k*eta_k);
302  dPhi += (pt_k*phi_k);
303  dPhi2 += (pt_k*phi_k*phi_k);
304  }
305 
306  GammadPt_ = sumPt/(*thePFTauRef).pt();
307 
308  if(sumPt>0){
309  dEta /= sumPt;
310  dPhi /= sumPt;
311  dEta2 /= sumPt;
312  dPhi2 /= sumPt;
313  }
314 
315  //GammadEta_ = dEta;
316  //GammadPhi_ = dPhi;
317 
318  GammadEta_ = TMath::Sqrt(dEta2)*TMath::Sqrt(GammadPt_)*(*thePFTauRef).pt();
319  GammadPhi_ = TMath::Sqrt(dPhi2)*TMath::Sqrt(GammadPt_)*(*thePFTauRef).pt();
320 
321  if( ((*thePFTauRef).signalPFChargedHadrCands()).size() == 3)
322  mva = 1.0;
323  else if( ((*thePFTauRef).signalPFChargedHadrCands()).size()==1 && ((*thePFTauRef).signalPFGammaCands()).size()==0){
324  if(TMath::Abs((*thePFTauRef).eta())<1.5)
325  mva = fTMVAReader_[0]->EvaluateMVA( methodName_ );
326  else
327  mva = fTMVAReader_[3]->EvaluateMVA( methodName_ );
328  }
329  else if( ((*thePFTauRef).signalPFChargedHadrCands()).size()==1 && ((*thePFTauRef).signalPFGammaCands()).size()>0 && (((*thePFTauRef).leadPFChargedHadrCand())->gsfTrackRef()).isNonnull()){
330  if(TMath::Abs((*thePFTauRef).eta())<1.5)
331  mva = fTMVAReader_[1]->EvaluateMVA( methodName_ );
332  else
333  mva = fTMVAReader_[4]->EvaluateMVA( methodName_ );
334  }
335  else if( ((*thePFTauRef).signalPFChargedHadrCands()).size()==1 && ((*thePFTauRef).signalPFGammaCands()).size()>0 && !(((*thePFTauRef).leadPFChargedHadrCand())->gsfTrackRef()).isNonnull()){
336  if(TMath::Abs((*thePFTauRef).eta())<1.5)
337  mva = fTMVAReader_[2]->EvaluateMVA( methodName_ );
338  else
339  mva = fTMVAReader_[5]->EvaluateMVA( methodName_ );
340  }
341  else{
342  mva = -99;
343  }
344 
345  return mva;
346 
347 }
const double Pi
Float_t TauLeadPFChargedHadrMva_
Float_t TauLeadPFChargedHadrHoP_
TMVA::Reader * fTMVAReader_[6]
Float_t TauLeadPFChargedHadrEoP_
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
int k[5][pyjets_maxn]
tuple cout
Definition: gather_cfg.py:121
tuple size
Write out results.
list at
Definition: asciidump.py:428

Member Data Documentation

TMVA::Reader* AntiElectronIDMVA::fTMVAReader_[6]
private

Definition at line 107 of file AntiElectronIDMVA.h.

Referenced by AntiElectronIDMVA(), Initialize(), MVAValue(), and ~AntiElectronIDMVA().

Float_t AntiElectronIDMVA::GammadEta_
private

Definition at line 110 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::GammadPhi_
private

Definition at line 111 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::GammadPt_
private

Definition at line 112 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Bool_t AntiElectronIDMVA::isInitialized_
private

Definition at line 105 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

std::string AntiElectronIDMVA::methodName_
private

Definition at line 106 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::TauEmFraction_
private

Definition at line 116 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::TauLeadPFChargedHadrEoP_
private

Definition at line 115 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::TauLeadPFChargedHadrHoP_
private

Definition at line 114 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::TauLeadPFChargedHadrMva_
private

Definition at line 113 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

Float_t AntiElectronIDMVA::TauSignalPFGammaCands_
private

Definition at line 108 of file AntiElectronIDMVA.h.

Referenced by Initialize().

Float_t AntiElectronIDMVA::TauVisMass_
private

Definition at line 109 of file AntiElectronIDMVA.h.

Referenced by Initialize(), and MVAValue().

int AntiElectronIDMVA::verbosity_
private

Definition at line 118 of file AntiElectronIDMVA.h.

Referenced by AntiElectronIDMVA(), and Initialize().