CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoEgamma/ElectronIdentification/src/ElectronMVAEstimator.cc

Go to the documentation of this file.
00001 #include "RecoEgamma/ElectronIdentification/interface/ElectronMVAEstimator.h"
00002 #include "DataFormats/TrackReco/interface/Track.h"
00003 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00004 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00005 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00006 #include "FWCore/ParameterSet/interface/FileInPath.h"
00007 
00008 ElectronMVAEstimator::ElectronMVAEstimator(){}
00009 
00010 ElectronMVAEstimator::ElectronMVAEstimator(std::string fileName){
00011   init(fileName);
00012 }
00013 
00014 
00015 void ElectronMVAEstimator::init(std::string fileName) {
00016   tmvaReader_ = new TMVA::Reader("!Color:Silent");
00017   tmvaReader_->AddVariable("fbrem",&fbrem);
00018   tmvaReader_->AddVariable("detain", &detain);
00019   tmvaReader_->AddVariable("dphiin", &dphiin);
00020   tmvaReader_->AddVariable("sieie", &sieie);
00021   tmvaReader_->AddVariable("hoe", &hoe);
00022   tmvaReader_->AddVariable("eop", &eop);
00023   tmvaReader_->AddVariable("e1x5e5x5", &e1x5e5x5);
00024   tmvaReader_->AddVariable("eleopout", &eleopout);
00025   tmvaReader_->AddVariable("detaeleout", &detaeleout);
00026   tmvaReader_->AddVariable("kfchi2", &kfchi2);
00027   tmvaReader_->AddVariable("kfhits", &mykfhits);
00028   tmvaReader_->AddVariable("mishits",&mymishits);
00029   tmvaReader_->AddVariable("dist", &absdist);
00030   tmvaReader_->AddVariable("dcot", &absdcot);
00031   tmvaReader_->AddVariable("nvtx", &myNvtx);
00032 
00033   tmvaReader_->AddSpectator("eta",&eta);
00034   tmvaReader_->AddSpectator("pt",&pt);
00035   tmvaReader_->AddSpectator("ecalseed",&ecalseed);
00036   
00037   // Taken from Daniele (his mail from the 30/11)
00038   //  tmvaReader_->BookMVA("BDTSimpleCat","../Training/weights_Root527b_3Depth_DanVarConvRej_2PtBins_10Pt_800TPrune5_Min100Events_NoBjets_half/TMVA_BDTSimpleCat.weights.xm");
00039   // training of the 7/12 with Nvtx added
00040   tmvaReader_->BookMVA("BDTSimpleCat",fileName.c_str());
00041 }
00042 
00043 
00044 
00045 double ElectronMVAEstimator::mva(const reco::GsfElectron& myElectron, int nvertices )  {
00046   fbrem = myElectron.fbrem();
00047   detain = myElectron.deltaEtaSuperClusterTrackAtVtx();
00048   dphiin = myElectron.deltaPhiSuperClusterTrackAtVtx();
00049   sieie = myElectron.sigmaIetaIeta();
00050   hoe = myElectron.hcalOverEcal();
00051   eop = myElectron.eSuperClusterOverP();
00052   e1x5e5x5 = (myElectron.e5x5()) !=0. ? 1.-(myElectron.e1x5()/myElectron.e5x5()) : -1. ;
00053   eleopout = myElectron.eEleClusterOverPout();
00054   detaeleout = myElectron.deltaEtaEleClusterTrackAtCalo();
00055   
00056   bool validKF= false;
00057 
00058   reco::TrackRef myTrackRef = myElectron.closestCtfTrackRef();
00059   validKF = (myTrackRef.isAvailable());
00060   validKF = (myTrackRef.isNonnull());  
00061 
00062   kfchi2 = (validKF) ? myTrackRef->normalizedChi2() : 0 ;
00063   kfhits = (validKF) ? myTrackRef->hitPattern().trackerLayersWithMeasurement() : -1. ; 
00064   dist = myElectron.convDist();
00065   dcot = myElectron.convDcot();
00066   eta = myElectron.eta();
00067   pt = myElectron.pt();
00068   
00069   mishits = myElectron.gsfTrack()->trackerExpectedHitsInner().numberOfLostHits();
00070   ecalseed = myElectron.ecalDrivenSeed();
00071   
00072   Nvtx = nvertices;
00073 
00074   bindVariables();
00075   double result =  tmvaReader_->EvaluateMVA("BDTSimpleCat");
00076 //  
00077 //  std::cout << "fbrem" <<fbrem << std::endl;
00078 //  std::cout << "detain"<< detain << std::endl;
00079 //  std::cout << "dphiin"<< dphiin << std::endl;
00080 //  std::cout << "sieie"<< sieie << std::endl;
00081 //  std::cout << "hoe"<< hoe << std::endl;
00082 //  std::cout << "eop"<< eop << std::endl;
00083 //  std::cout << "e1x5e5x5"<< e1x5e5x5 << std::endl;
00084 //  std::cout << "eleopout"<< eleopout << std::endl;
00085 //  std::cout << "detaeleout"<< detaeleout << std::endl;
00086 //  std::cout << "kfchi2"<< kfchi2 << std::endl;
00087 //  std::cout << "kfhits"<< mykfhits << std::endl;
00088 //  std::cout << "mishits"<<mymishits << std::endl;
00089 //  std::cout << "dist"<< absdist << std::endl;
00090 //  std::cout << "dcot"<< absdcot << std::endl;
00091 //
00092 //  std::cout << "eta"<<eta << std::endl;
00093 //  std::cout << "pt"<<pt << std::endl;
00094 //  std::cout << "ecalseed"<<ecalseed << std::endl;
00095 //
00096 //  std::cout << " MVA " << result << std::endl;
00097   return result;
00098 }
00099 
00100 
00101 void ElectronMVAEstimator::bindVariables() {
00102   if(fbrem < -1.)
00103     fbrem = -1.;  
00104   
00105   detain = fabs(detain);
00106   if(detain > 0.06)
00107     detain = 0.06;
00108   
00109   
00110   dphiin = fabs(dphiin);
00111   if(dphiin > 0.6)
00112     dphiin = 0.6;
00113 
00114   
00115   if(eop > 20.)
00116     eop = 20.;
00117   
00118   
00119   if(eleopout > 20.)
00120     eleopout = 20;
00121   
00122   detaeleout = fabs(detaeleout);
00123   if(detaeleout > 0.2)
00124     detaeleout = 0.2;
00125   
00126   mykfhits = float(kfhits);
00127   mymishits = float(mishits);
00128   
00129   if(kfchi2 < 0.)
00130     kfchi2 = 0.;
00131   
00132   if(kfchi2 > 15.)
00133     kfchi2 = 15.;
00134   
00135   
00136   if(e1x5e5x5 < -1.)
00137     e1x5e5x5 = -1;
00138 
00139   if(e1x5e5x5 > 2.)
00140     e1x5e5x5 = 2.; 
00141   
00142   
00143   if(dist > 15.)
00144     dist = 15.;
00145   if(dist < -15.)
00146     dist = -15.;
00147   
00148   if(dcot > 3.)
00149     dcot = 3.;
00150   if(dcot < -3.)
00151     dcot = -3.;
00152   
00153   absdist = fabs(dist);
00154   absdcot = fabs(dcot);
00155   myNvtx = float(Nvtx);
00156 
00157 }