CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SoftElectronMVAEstimator.cc
Go to the documentation of this file.
7 
9  std::vector<std::string> weightsfiles;
10  std::string path_mvaWeightFileEleID;
11  for(unsigned ifile=0 ; ifile < cfg_.vweightsfiles.size() ; ++ifile) {
12  path_mvaWeightFileEleID = edm::FileInPath ( cfg_.vweightsfiles[ifile].c_str() ).fullPath();
13  weightsfiles.push_back(path_mvaWeightFileEleID);
14  }
15 
16  for (unsigned int i=0;i<fmvaReader.size(); ++i) {
17  if (fmvaReader[i]) delete fmvaReader[i];
18  }
19  fmvaReader.clear();
20 
21  //initialize
22  //Define expected number of bins
23  UInt_t ExpectedNBins = 1;
24 
25  //Check number of weight files given
26  if (ExpectedNBins != cfg_.vweightsfiles.size() ) {
27  std::cout << "Error: Expected Number of bins = " << ExpectedNBins << " does not equal to weightsfiles.size() = "
28  << cfg_.vweightsfiles.size() << std::endl;
29 
30  assert(ExpectedNBins == cfg_.vweightsfiles.size());
31  }
32 
33 
34  for (unsigned int i=0;i<ExpectedNBins; ++i) {
35  tmvaReader_ = new TMVA::Reader("!Color:Silent");
36  tmvaReader_->AddVariable("fbrem", &fbrem);
37  tmvaReader_->AddVariable("EtotOvePin", &EtotOvePin);
38  tmvaReader_->AddVariable("EClusOverPout", &eleEoPout);
39  tmvaReader_->AddVariable("EBremOverDeltaP", &EBremOverDeltaP);
40  tmvaReader_->AddVariable("logSigmaEtaEta", &logSigmaEtaEta);
41  tmvaReader_->AddVariable("DeltaEtaTrackEcalSeed", &DeltaEtaTrackEcalSeed);
42  tmvaReader_->AddVariable("HoE", &HoE);
43  tmvaReader_->AddVariable("gsfchi2", &gsfchi2);
44  tmvaReader_->AddVariable("kfchi2", &kfchi2);
45  tmvaReader_->AddVariable("kfhits", &kfhits);
46  tmvaReader_->AddVariable("SigmaPtOverPt", &SigmaPtOverPt);
47  tmvaReader_->AddVariable("deta", &deta);
48  tmvaReader_->AddVariable("dphi", &dphi);
49  tmvaReader_->AddVariable("detacalo", &detacalo);
50  tmvaReader_->AddVariable("see", &see);
51  tmvaReader_->AddVariable("spp", &spp);
52  tmvaReader_->AddVariable("R9", &R9);
53  tmvaReader_->AddVariable("etawidth", &etawidth);
54  tmvaReader_->AddVariable("phiwidth", &phiwidth);
55  tmvaReader_->AddVariable("e1x5e5x5", &OneMinusE1x5E5x5);
56  tmvaReader_->AddVariable("IoEmIoP", &IoEmIoP);
57  tmvaReader_->AddVariable("PreShowerOverRaw", &PreShowerOverRaw);
58  tmvaReader_->AddVariable("nPV", &nPV);
59 
60  tmvaReader_->AddVariable( "pt", &pt);
61  tmvaReader_->AddVariable( "eta", &eta);
62 
63  tmvaReader_->AddSpectator( "pt", &pt);
64  tmvaReader_->AddSpectator( "eta", &eta);
65 
66 // tmvaReader_->AddSpectator( "nPV", &nPV);
67 
68  // Taken from Daniele (his mail from the 30/11)
69  // tmvaReader_->BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm");
70  // training of the 7/12 with Nvtx added
71  tmvaReader_->BookMVA("BDT",weightsfiles[i]);
72  fmvaReader.push_back(tmvaReader_);
73 // delete tmvaReader_;
74  }
75 
76 }
77 
78 
80 {
81  for (unsigned int i=0;i<fmvaReader.size(); ++i) {
82  if (fmvaReader[i]) delete fmvaReader[i];
83  }
84 }
85 
86 
87 UInt_t SoftElectronMVAEstimator::GetMVABin(int pu, double eta, double pt) const {
88 
89  //Default is to return the first bin
90  unsigned int bin = 0;
91 
92  bool ptrange[3],etarange[3],purange[2];
93  ptrange[0]=pt > 2 && pt < 5;
94  ptrange[1]=pt > 5 && pt < 10;
95  ptrange[2]=pt > 10;
96  etarange[0]=fabs(eta) < 0.8;
97  etarange[1]=fabs(eta) > 0.8 && fabs(eta) <1.4;
98  etarange[2]=fabs(eta) > 1.4;
99  purange[0]=nPV<=20;
100  purange[1]=nPV>20;
101 
102  int index=0;
103  for(int kPU=0;kPU<2;kPU++)
104  for(int kETA=0;kETA<3;kETA++)
105  for(int kPT=0;kPT<3;kPT++){
106  if (purange[kPU] && ptrange[kPT] && etarange[kETA]) bin=index;
107  index++;
108  }
109  return bin;
110 }
111 
112 
113 
114 double SoftElectronMVAEstimator::mva(const reco::GsfElectron& myElectron,const edm::Event & evt) {
115 
116  edm::Handle<reco::VertexCollection> FullprimaryVertexCollection;
117  evt.getByLabel("offlinePrimaryVertices", FullprimaryVertexCollection);
118  const reco::VertexCollection pvc = *(FullprimaryVertexCollection.product());
119 
120  fbrem =myElectron.fbrem();
121  EtotOvePin =myElectron.eSuperClusterOverP();
122  eleEoPout =myElectron.eEleClusterOverPout();
123  float etot =myElectron.eSuperClusterOverP()*myElectron.trackMomentumAtVtx().R();
124  float eEcal =myElectron.eEleClusterOverPout()*myElectron.trackMomentumAtEleClus().R();
125  float dP =myElectron.trackMomentumAtVtx().R()-myElectron.trackMomentumAtEleClus().R();
126  EBremOverDeltaP =(etot-eEcal)/dP;
127  logSigmaEtaEta =log(myElectron.sigmaEtaEta());
129  HoE =myElectron.hcalOverEcalBc();
130 
131  bool validKF= false;
132  reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef();
133  validKF = (myTrackRef.isAvailable() && myTrackRef.isNonnull());
134  kfchi2 =(validKF) ? myTrackRef->normalizedChi2() : 0 ;
135  kfhits =(validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1. ;
136  gsfchi2 =myElectron.gsfTrack()->normalizedChi2();
137  SigmaPtOverPt =myElectron.gsfTrack().get()->ptModeError()/myElectron.gsfTrack().get()->ptMode() ;
138  deta =myElectron.deltaEtaSuperClusterTrackAtVtx();
139  dphi =myElectron.deltaPhiSuperClusterTrackAtVtx();
141  see =myElectron.sigmaIetaIeta();
142  spp =myElectron.sigmaIphiIphi();
143  R9 =myElectron.r9();
144  IoEmIoP = (1.0/myElectron.ecalEnergy()) - (1.0 / myElectron.p());
145  etawidth =myElectron.superCluster()->etaWidth();
146  phiwidth =myElectron.superCluster()->phiWidth();
147  OneMinusE1x5E5x5 =(myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ;
148  pt =myElectron.pt();
149  eta =myElectron.eta();
150  nPV=pvc.size();
151  PreShowerOverRaw=myElectron.superCluster()->preshowerEnergy() / myElectron.superCluster()->rawEnergy();
152 
153 /*
154  std::cout<<"fbrem "<<fbrem<<std::endl;
155  std::cout<<"EtotOvePin "<<EtotOvePin<<std::endl;
156  std::cout<<"eleEoPout "<<eleEoPout<<std::endl;
157  std::cout<<"EBremOverDeltaP "<<EBremOverDeltaP<<std::endl;
158  std::cout<<"logSigmaEtaEta "<<logSigmaEtaEta<<std::endl;
159  std::cout<<"DeltaEtaTrackEcalSeed "<<DeltaEtaTrackEcalSeed<<std::endl;
160  std::cout<<"HoE "<<HoE<<std::endl;
161  std::cout<<"kfchi2 "<<kfchi2<<std::endl;
162  std::cout<<"kfhits "<<kfhits<<std::endl;
163  std::cout<<"gsfchi2 "<<gsfchi2<<std::endl;
164  std::cout<<"SigmaPtOverPt "<<SigmaPtOverPt<<std::endl;
165  std::cout<<"deta "<<deta<<std::endl;
166  std::cout<<"dphi "<<dphi<<std::endl;
167  std::cout<<"detacalo "<<detacalo<<std::endl;
168  std::cout<<"see "<<see<<std::endl;
169  std::cout<< "spp " << spp<< std::endl;
170  std::cout<< "R9 " << R9<< std::endl;
171  std::cout<< "IoEmIoP " << IoEmIoP<< std::endl;
172  std::cout<<"etawidth "<<etawidth<<std::endl;
173  std::cout<<"phiwidth "<<phiwidth<<std::endl;
174  std::cout<<"OneMinusE1x5E5x5 "<<OneMinusE1x5E5x5<<std::endl;
175  std::cout<<"PreShowerOverRaw "<<PreShowerOverRaw<<std::endl;
176 */
177  bindVariables();
178 // double result= fmvaReader[GetMVABin(nPV,eta,pt)]->EvaluateMVA("BDT");
179  double result= fmvaReader[0]->EvaluateMVA("BDT");
180 // double result = tmvaReader_->EvaluateMVA("BDT");
181  return result;
182 }
183 
184 
186  if(fbrem < -1.)
187  fbrem = -1.;
188 
189  deta = fabs(deta);
190  if(deta > 0.06)
191  deta = 0.06;
192 
193 
194  dphi = fabs(dphi);
195  if(dphi > 0.6)
196  dphi = 0.6;
197 
198 
199  if(EoP > 20.)
200  EoP = 20.;
201 
202  if(eleEoPout > 20.)
203  eleEoPout = 20.;
204 
205 
206  detacalo = fabs(detacalo);
207  if(detacalo > 0.2)
208  detacalo = 0.2;
209 
210  if(OneMinusE1x5E5x5 < -1.)
211  OneMinusE1x5E5x5 = -1;
212 
213  if(OneMinusE1x5E5x5 > 2.)
214  OneMinusE1x5E5x5 = 2.;
215 
216 
217 
218  if(gsfchi2 > 200.)
219  gsfchi2 = 200;
220 
221 
222  if(kfchi2 > 10.)
223  kfchi2 = 10.;
224 
225 }
float sigmaIphiIphi() const
Definition: GsfElectron.h:400
int i
Definition: DBlmapReader.cc:9
virtual double p() const
magnitude of momentum vector
static std::vector< std::string > checklist log
virtual float pt() const
transverse momentum
float eSuperClusterOverP() const
Definition: GsfElectron.h:243
SoftElectronMVAEstimator(const Configuration &)
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:286
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
float fbrem() const
Definition: GsfElectron.h:684
bool isAvailable() const
Definition: Ref.h:276
T eta() const
TrackRef closestCtfTrackRef() const
Definition: GsfElectron.h:199
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:247
float sigmaIetaIeta() const
Definition: GsfElectron.h:399
UInt_t GetMVABin(int pu, double eta, double pt) const
virtual float eta() const
momentum pseudorapidity
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:182
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:250
tuple result
Definition: query.py:137
float eEleClusterOverPout() const
Definition: GsfElectron.h:246
float hcalOverEcalBc() const
Definition: GsfElectron.h:411
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
std::vector< TMVA::Reader * > fmvaReader
float deltaEtaEleClusterTrackAtCalo() const
Definition: GsfElectron.h:249
float e1x5() const
Definition: GsfElectron.h:401
T const * product() const
Definition: Handle.h:81
float ecalEnergy() const
Definition: GsfElectron.h:772
math::XYZVectorF trackMomentumAtEleClus() const
Definition: GsfElectron.h:289
double mva(const reco::GsfElectron &myElectron, const edm::Event &evt)
float e5x5() const
Definition: GsfElectron.h:403
float r9() const
Definition: GsfElectron.h:404
float deltaEtaSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:248
tuple cout
Definition: gather_cfg.py:121
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:242
float sigmaEtaEta() const
Definition: GsfElectron.h:398
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:183