CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoBTag/ImpactParameterLearning/plugins/ImpactParameterCalibration.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    ImpactParameterCalibration
00004 // Class:      ImpactParameterCalibration
00005 // 
00013 //
00014 // Original Author:  Jeremy Andrea/Andrea Rizzi
00015 //         Created:  Mon Aug  6 16:10:38 CEST 2007
00016 // $Id: ImpactParameterCalibration.cc,v 1.14 2010/02/11 00:13:30 wmtan Exp $
00017 //
00018 //
00019 // system include files
00020 
00021 #include "FWCore/Framework/interface/Frameworkfwd.h"
00022 #include "FWCore/Framework/interface/EDAnalyzer.h"
00023 #include "FWCore/Framework/interface/Event.h"
00024 #include "FWCore/Framework/interface/MakerMacros.h"
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "FWCore/Utilities/interface/InputTag.h"
00027 #include "FWCore/Framework/interface/ESHandle.h"
00028 #include "FWCore/Framework/interface/EventSetup.h"
00029 #include "Utilities/General/interface/FileInPath.h"
00030 
00031 #include "FWCore/Framework/interface/IOVSyncValue.h"
00032 #include "FWCore/ServiceRegistry/interface/Service.h"
00033 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00034 
00035 
00036 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
00037 #include "CondFormats/BTauObjects/interface/TrackProbabilityCalibration.h"
00038 #include "CondFormats/BTauObjects/interface/CalibratedHistogram.h"
00039 
00040 #include "CondFormats/DataRecord/interface/BTagTrackProbability2DRcd.h"
00041 #include "CondFormats/DataRecord/interface/BTagTrackProbability3DRcd.h"
00042 
00043 #include "TClass.h"
00044 
00045 #include "TBufferFile.h"
00046 
00047 #include "TBufferXML.h"
00048 #include <iostream>
00049 #include <fstream>
00050 #include <sstream>
00051 #include <string>
00052 #include <vector>
00053 #include <memory>
00054 
00055 #include "TrackClassMatch.h"
00056 
00057 
00058 using namespace reco;
00059 using namespace std;
00060 //
00061 // class decleration
00062 //
00063 
00064 
00065 
00066 class ImpactParameterCalibration : public edm::EDAnalyzer {
00067    public:
00068       explicit ImpactParameterCalibration(const edm::ParameterSet&);
00069       ~ImpactParameterCalibration();
00070 
00071 
00072    private:
00073       virtual void beginJob() ;
00074       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00075       virtual void endJob() ;
00076 
00077       virtual void initFromFirstES(const edm::EventSetup&);
00078       edm::ParameterSet config;
00079       bool m_needInitFromES;
00080    TrackProbabilityCalibration * fromXml(edm::FileInPath xmlCalibration);
00081 
00082    static TrackProbabilityCategoryData createCategory(double  pmin,double  pmax,
00083                  double  etamin,  double  etamax,
00084                   int  nhitmin,  int  nhitmax,
00085                  int  npixelhitsmin, int  npixelhitsmax,
00086                   double cmin, double cmax, int withFirst)
00087   {
00088   TrackProbabilityCategoryData c;
00089   c.pMin=pmin;
00090   c.pMax=pmax;
00091   c.etaMin=etamin;
00092   c.etaMax=etamax;
00093   c.nHitsMin=nhitmin;
00094   c.nHitsMax=nhitmax;
00095   c.nPixelHitsMin=npixelhitsmin;
00096   c.nPixelHitsMax=npixelhitsmax;
00097   c.chiMin=cmin;
00098   c.chiMax=cmax;
00099   c.withFirstPixel=withFirst;
00100   return c;
00101  }
00102 
00103    TrackProbabilityCalibration * m_calibration[2];
00104    edm::InputTag m_iptaginfo;
00105    edm::InputTag m_pv;
00106   unsigned int minLoop, maxLoop;
00107 
00108 };
00109 
00110 ImpactParameterCalibration::ImpactParameterCalibration(const edm::ParameterSet& iConfig):config(iConfig)
00111 {
00112   m_needInitFromES = false;
00113   m_iptaginfo = iConfig.getParameter<edm::InputTag>("tagInfoSrc");
00114   m_pv = iConfig.getParameter<edm::InputTag>("primaryVertexSrc");
00115   bool createOnlyOne = iConfig.getUntrackedParameter<bool>("createOnlyOneCalibration", false);
00116   minLoop=0;
00117   maxLoop=1;
00118   if (createOnlyOne == true){
00119     int whichCalib = iConfig.getUntrackedParameter<int>("dimension", 2);
00120     if (whichCalib==2){
00121       std::cout <<" Writing only 2D calibrations"<<std::endl;
00122       minLoop=1;
00123       maxLoop=1;
00124     }else if (whichCalib==3){
00125       std::cout <<" Writing only 3D calibrations"<<std::endl;
00126       minLoop=0;
00127       maxLoop=0;
00128     }else {
00129       std::cout <<" Dimension not found: "<<whichCalib<<"; it must be either 2 or 3"<<std::endl;
00130     }
00131   }
00132 
00133 }
00134 
00135 
00136 ImpactParameterCalibration::~ImpactParameterCalibration()
00137 {
00138 }
00139 
00140 
00141 void
00142 ImpactParameterCalibration::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00143 {
00144   if(m_needInitFromES) initFromFirstES(iSetup);
00145   using namespace edm;
00146   using namespace reco;
00147 
00148   Handle<TrackIPTagInfoCollection> ipHandle;
00149   iEvent.getByLabel(m_iptaginfo, ipHandle);
00150   const TrackIPTagInfoCollection & ip = *(ipHandle.product());
00151 
00152 //  cout << "Found " << ip.size() << " TagInfo" << endl;
00153 
00154   Handle<reco::VertexCollection> primaryVertex;
00155   iEvent.getByLabel(m_pv,primaryVertex);
00156 
00157   vector<TrackProbabilityCalibration::Entry>::iterator found;
00158   vector<TrackProbabilityCalibration::Entry>::iterator it_begin;
00159   vector<TrackProbabilityCalibration::Entry>::iterator it_end;
00160 
00161       
00162    TrackIPTagInfoCollection::const_iterator it = ip.begin();
00163    for(; it != ip.end(); it++)
00164      {
00165       TrackRefVector selTracks=it->selectedTracks();
00166 //      if(it->primaryVertex().isNull()) continue;
00167       if(primaryVertex.product()->size() == 0) 
00168        {
00169        std::cout << "No PV in the event!!" << std::endl;
00170          continue;
00171         }
00172        const Vertex & pv = *(primaryVertex.product()->begin());
00173            
00174       for(unsigned int i=minLoop; i <= maxLoop;i++)
00175       { 
00176         it_begin=m_calibration[i]->data.begin();
00177         it_end=m_calibration[i]->data.end();
00178   
00179       for(unsigned int j=0;j<selTracks.size(); j++)
00180         {
00181           double ipsig;
00182           if (i==0) ipsig  = it->impactParameterData()[j].ip3d.significance();
00183           else  ipsig  = it->impactParameterData()[j].ip2d.significance();
00184           TrackClassMatch::Input input(*selTracks[j],*it->jet(),pv);
00185           if(ipsig < 0) 
00186            {
00187             found = std::find_if(it_begin,it_end,bind1st(TrackClassMatch(),input));
00188 //            std::cout << ip[j].significance() << std::endl; 
00189             if(found!=it_end) 
00190               found->histogram.fill(-ipsig);
00191             else
00192               {std::cout << "No category for this track!!" << std::endl;
00193               std::cout << "p       "  <<(*selTracks[j]).p ()  << std::endl;
00194               std::cout << "eta     " << (*selTracks[j]).eta() << std::endl;
00195               std::cout << "NHit    " << (*selTracks[j]).numberOfValidHits() << std::endl;
00196               std::cout << "NPixHit " << (*selTracks[j]).hitPattern().numberOfValidPixelHits() << std::endl;
00197               std::cout << "FPIXHIT " << (*selTracks[j]).hitPattern().hasValidHitInFirstPixelBarrel() << std::endl;}
00198               
00199            }
00200          }
00201       } 
00202      }  
00203       
00204          
00205   
00206 }
00207 
00208 
00209 
00210 
00211 
00212 
00213 
00214 void ImpactParameterCalibration::initFromFirstES(const edm::EventSetup& iSetup)
00215 {
00216   using namespace edm;
00217 
00218     CalibratedHistogram hist(config.getParameter<int>("nBins"),0,config.getParameter<double>("maxSignificance"));
00219     bool resetHistogram = config.getParameter<bool>("resetHistograms");
00220     ESHandle<TrackProbabilityCalibration> calib2DHandle;
00221     iSetup.get<BTagTrackProbability2DRcd>().get(calib2DHandle);
00222     ESHandle<TrackProbabilityCalibration> calib3DHandle;
00223     iSetup.get<BTagTrackProbability3DRcd>().get(calib3DHandle);
00224     const TrackProbabilityCalibration * ca[2];
00225     ca[0]  = calib3DHandle.product();
00226     ca[1]  = calib2DHandle.product();
00227     for(unsigned int i=minLoop;i <=maxLoop ;i++)
00228     for(unsigned int j=0;j<ca[i]->data.size() ; j++)
00229     {
00230      TrackProbabilityCalibration::Entry e;
00231      e.category=ca[i]->data[j].category;
00232 
00233      if(resetHistogram)
00234       e.histogram=hist;
00235      else
00236       e.histogram=ca[i]->data[j].histogram;
00237 
00238      m_calibration[i]->data.push_back(e);
00239     }
00240 
00241 
00242 }
00243 
00244 
00245 // ------------ method called once each job just before starting event loop  ------------
00246 void 
00247 ImpactParameterCalibration::beginJob()
00248 {
00249   using namespace edm;
00250   m_calibration[0] =   new TrackProbabilityCalibration();
00251   m_calibration[1] =   new TrackProbabilityCalibration();
00252 
00253   CalibratedHistogram hist(config.getParameter<int>("nBins"),0,config.getParameter<double>("maxSignificance"));
00254 
00255   std::string categories = config.getParameter<std::string>("inputCategories");
00256 
00257   if(categories == "HardCoded")
00258   {
00259    vector<TrackProbabilityCategoryData> v;
00260     //TrackProbabilityCategoryData {pMin, pMax, etaMin, etaMax,
00261     //nHitsMin, nHitsMax, nPixelHitsMin, nPixelHitsMax, chiMin,chiMax, withFirstPixel;
00262     //trackQuality;
00263   v.push_back(createCategory(0, 5000, 0  , 2.5, 8 , 50, 1, 1, 0  , 5  , 0));
00264   v.push_back(createCategory(0, 5000, 0  , 2.5, 8 , 50, 2, 8, 2.5, 5  , 0));
00265   v.push_back(createCategory(0, 8   , 0  , 0.8, 8 , 50, 3, 8, 0  , 2.5, 0));
00266   v.push_back(createCategory(0, 8   , 0.8, 1.6, 8 , 50, 3, 8, 0  , 2.5, 0));
00267   v.push_back(createCategory(0, 8   , 1.6, 2.5, 8 , 50, 3, 8, 0  , 2.5, 0));
00268   v.push_back(createCategory(0, 8   , 0  , 2.5, 8 , 50, 2, 8, 0  , 2.5, 0));
00269   v.push_back(createCategory(8, 5000, 0  , 0.8, 8 , 50, 3, 8, 0  , 2.5, 0));
00270   v.push_back(createCategory(8, 5000, 0.8, 1.6, 8 , 50, 3, 8, 0  , 2.5, 0));
00271   v.push_back(createCategory(8, 5000, 1.6, 2.5, 8 , 50, 3, 8, 0  , 2.5, 0));
00272   v.push_back(createCategory(8, 5000, 0  , 2.5, 8 , 50, 2 ,2, 0  , 2.5, 0));
00273   for(unsigned int i=minLoop;i <=maxLoop ;i++)
00274    for(unsigned int j=0;j<v.size() ; j++)
00275     {
00276      TrackProbabilityCalibration::Entry e;
00277      e.category=v[j];
00278      e.histogram=hist;
00279      m_calibration[i]->data.push_back(e);
00280     }
00281 
00282   }
00283   if(categories == "RootXML")
00284   {
00285     bool resetHistogram = config.getParameter<bool>("resetHistograms");
00286     const TrackProbabilityCalibration * ca[2];
00287     ca[0]  = fromXml(config.getParameter<edm::FileInPath>("calibFile3d"));
00288     ca[1]  = fromXml(config.getParameter<edm::FileInPath>("calibFile2d"));
00289   
00290     for(unsigned int i=minLoop;i <=maxLoop ;i++)
00291      for(unsigned int j=0;j<ca[i]->data.size() ; j++)
00292      {
00293       TrackProbabilityCalibration::Entry e;
00294       e.category=ca[i]->data[j].category;
00295 
00296       if(resetHistogram)
00297        e.histogram=hist;
00298       else
00299        e.histogram=ca[i]->data[j].histogram;
00300 
00301       m_calibration[i]->data.push_back(e);
00302      }
00303 
00304    delete ca[0];
00305    delete ca[1];
00306 
00307    }
00308   if(categories == "EventSetup")
00309    {
00310     m_needInitFromES=true;
00311    }
00312 
00313 
00314 
00315 /*  edm::FileInPath f2d(m_xmlfilename2D);
00316     edm::FileInPath f3d(m_xmlfilename3D);
00317     calibrationNew   =  new AlgorithmCalibration<TrackClassFilterCategory,CalibratedHistogramXML>((f3d.fullPath()).c_str());
00318     calibration2dNew =  new AlgorithmCalibration<TrackClassFilterCategory,CalibratedHistogramXML>((f2d.fullPath()).c_str());
00319     vector<float> * bins =0;
00320     if(m_resetData)
00321       {
00322         if(m_newBinning)  bins = new  vector<float>(CalibratedHistogram::constantBinning(m_nBin,0,m_range));
00323         vector<pair<TrackClassFilterCategory, CalibratedHistogramXML> > data = calibrationNew->categoriesWithData();
00324         vector<pair<TrackClassFilterCategory, CalibratedHistogramXML> > data2d = calibration2dNew->categoriesWithData();
00325         std::cout <<  data.size() <<  std::endl;
00326         for(unsigned int i = 0 ; i < data.size();i++)
00327           {
00328             data[i].second.reset();
00329             if(bins)  data[i].second.setUpperLimits(*bins);
00330           }
00331         for(unsigned int i = 0 ; i < data2d.size();i++)
00332           {
00333             data2d[i].second.reset();
00334             if(bins)  data2d[i].second.setUpperLimits(*bins);
00335           }
00336         
00337       }
00338     if(bins) delete bins;
00339     
00340 */
00341 
00342 
00343 }
00344 
00345 TrackProbabilityCalibration * ImpactParameterCalibration::fromXml(edm::FileInPath xmlCalibration)   
00346 {
00347      std::ifstream xmlFile(xmlCalibration.fullPath().c_str());
00348         if (!xmlFile.good())
00349                 throw cms::Exception("BTauFakeMVAJetTagConditions")
00350                         << "File \"" << xmlCalibration.fullPath()
00351                         << "\" could not be opened for reading."
00352                         << std::endl;
00353         std::ostringstream ss;
00354         ss << xmlFile.rdbuf();
00355         xmlFile.close();
00356         TClass *classType = 0;
00357         void *ptr = TBufferXML(TBuffer::kRead).ConvertFromXMLAny(
00358                                 ss.str().c_str(), &classType, kTRUE, kFALSE);
00359         if (!ptr)
00360                 throw cms::Exception("ImpactParameterCalibration")
00361                         << "Unknown error parsing XML serialization"
00362                         << std::endl;
00363 
00364         if (std::strcmp(classType->GetName(),
00365                 "TrackProbabilityCalibration")) {
00366                 classType->Destructor(ptr);
00367                 throw cms::Exception("ImpactParameterCalibration")
00368                         << "Serialized object has wrong C++ type."
00369                         << std::endl;
00370         }
00371 
00372         return static_cast<TrackProbabilityCalibration*>(ptr);
00373 }
00374 
00375 
00376 
00377 
00378 // ------------ method called once each job just after ending the event loop  ------------
00379 void 
00380 ImpactParameterCalibration::endJob() {
00381 
00382   if(config.getParameter<bool>("writeToDB"))
00383   {
00384     edm::Service<cond::service::PoolDBOutputService> mydbservice;
00385     if( !mydbservice.isAvailable() ) return;
00386     if(minLoop == 0 )  mydbservice->createNewIOV<TrackProbabilityCalibration>(m_calibration[0], mydbservice->beginOfTime(), mydbservice->endOfTime(),"BTagTrackProbability3DRcd");
00387     if(maxLoop == 1)   mydbservice->createNewIOV<TrackProbabilityCalibration>(m_calibration[1],  mydbservice->beginOfTime(), mydbservice->endOfTime(),"BTagTrackProbability2DRcd");
00388   } 
00389     
00390 
00391   if(config.getParameter<bool>("writeToRootXML"))
00392   {
00393     if(maxLoop == 1 ){
00394       std::ofstream of2("2d.xml");
00395       TBufferXML b2(TBuffer::kWrite);
00396       of2 << b2.ConvertToXML(const_cast<void*>(static_cast<const void*>(m_calibration[1])),
00397                              TClass::GetClass("TrackProbabilityCalibration"),
00398                              kTRUE, kFALSE);
00399       of2.close();
00400     }
00401     if(minLoop == 0 ){
00402       std::ofstream of3("3d.xml");
00403       TBufferXML b3(TBuffer::kWrite);
00404       of3 << b3.ConvertToXML(const_cast<void*>(static_cast<const void*>(m_calibration[0])),
00405                              TClass::GetClass("TrackProbabilityCalibration"),
00406                              kTRUE, kFALSE);
00407       of3.close();
00408     }
00409   }
00410 
00411  
00412   if(config.getParameter<bool>("writeToBinary"))
00413   {
00414     if(maxLoop == 1 ){
00415       
00416       std::ofstream ofile("2d.dat");
00417       TBufferFile buffer(TBuffer::kWrite);
00418       buffer.StreamObject(const_cast<void*>(static_cast<const void*>(m_calibration[1])),
00419                           TClass::GetClass("TrackProbabilityCalibration"));
00420       Int_t size = buffer.Length();
00421       ofile.write(buffer.Buffer(),size);
00422       ofile.close();
00423     }
00424     if(minLoop == 0 ){
00425       std::ofstream ofile3("3d.dat");
00426       TBufferFile buffer3(TBuffer::kWrite);
00427       buffer3.StreamObject(const_cast<void*>(static_cast<const void*>(m_calibration[0])),
00428                            TClass::GetClass("TrackProbabilityCalibration"));
00429       Int_t size3 = buffer3.Length();
00430       ofile3.write(buffer3.Buffer(),size3);
00431       ofile3.close();
00432     }
00433   }
00434   
00435   
00436 
00437 
00438 
00439 }
00440 
00441 //define this as a plug-in
00442 DEFINE_FWK_MODULE(ImpactParameterCalibration);