CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/PhysicsTools/Utilities/src/Lumi3DReWeighting.cc

Go to the documentation of this file.
00001 #ifndef PhysicsTools_Utilities_interface_Lumi3DReWeighting_cc
00002 #define PhysicsTools_Utilities_interface_Lumi3DReWeighting_cc
00003 
00004 
00020 #include "TRandom1.h"
00021 #include "TRandom2.h"
00022 #include "TRandom3.h"
00023 #include "TStopwatch.h"
00024 #include "TH1.h"
00025 #include "TH3.h"
00026 #include "TFile.h"
00027 #include <string>
00028 #include <algorithm>
00029 #include <boost/shared_ptr.hpp>
00030 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h" 
00031 #include "PhysicsTools/Utilities/interface/Lumi3DReWeighting.h"
00032 #include "DataFormats/Common/interface/Handle.h"
00033 #include "DataFormats/Provenance/interface/ProcessHistory.h"
00034 #include "DataFormats/Provenance/interface/Provenance.h"
00035 #include "FWCore/Common/interface/EventBase.h"
00036 
00037 using namespace edm;
00038 
00039 Lumi3DReWeighting::Lumi3DReWeighting( std::string generatedFile,
00040                                       std::string dataFile,
00041                                       std::string GenHistName = "pileup",
00042                                       std::string DataHistName = "pileup", 
00043                                       std::string WeightOutputFile ="") :
00044       generatedFileName_( generatedFile), 
00045       dataFileName_     ( dataFile ), 
00046       GenHistName_      ( GenHistName ), 
00047       DataHistName_     ( DataHistName ),
00048       weightFileName_   (WeightOutputFile)
00049       {
00050         generatedFile_ = boost::shared_ptr<TFile>( new TFile(generatedFileName_.c_str()) ); //MC distribution
00051         dataFile_      = boost::shared_ptr<TFile>( new TFile(dataFileName_.c_str()) );      //Data distribution
00052 
00053         boost::shared_ptr<TH1> Data_temp = boost::shared_ptr<TH1>(  (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
00054 
00055         boost::shared_ptr<TH1> MC_temp = boost::shared_ptr<TH1>(  (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
00056 
00057 
00058         MC_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
00059         Data_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
00060 
00061         // MC * data/MC = data, so the weights are data/MC:
00062 
00063         // normalize both histograms first
00064 
00065         Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00066         MC_distr_->Scale( 1.0/ MC_distr_->Integral() );
00067 
00068 }
00069 
00070 Lumi3DReWeighting::Lumi3DReWeighting( std::vector< float > MC_distr, std::vector< float > Lumi_distr,
00071                                       std::string WeightOutputFile ="") {
00072 
00073   weightFileName_ = WeightOutputFile;
00074 
00075   // no histograms for input: use vectors
00076   
00077   // now, make histograms out of them:
00078 
00079   Int_t NMCBins = MC_distr.size();
00080 
00081   MC_distr_ = boost::shared_ptr<TH1> ( new TH1F("MC_distr","MC dist",NMCBins,0., float(NMCBins)) );
00082 
00083   Int_t NDBins = Lumi_distr.size();
00084 
00085   Data_distr_ = boost::shared_ptr<TH1> ( new TH1F("Data_distr","Data dist",NDBins,0., float(NDBins)) );
00086 
00087   for(int ibin = 1; ibin<NMCBins+1; ++ibin ) {
00088     MC_distr_->SetBinContent(ibin,MC_distr[ibin-1]);
00089   }
00090 
00091   for(int ibin = 1; ibin<NDBins+1; ++ibin ) {
00092     Data_distr_->SetBinContent(ibin, Lumi_distr[ibin-1]);
00093   }
00094 
00095   // check integrals, make sure things are normalized
00096 
00097   float deltaH = Data_distr_->Integral();
00098   if(fabs(1.0 - deltaH) > 0.001 ) { //*OOPS*...
00099     Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00100   }
00101   float deltaMC = MC_distr_->Integral();
00102   if(fabs(1.0 - deltaMC) > 0.001 ) {
00103     MC_distr_->Scale(1.0/ MC_distr_->Integral());
00104   }
00105 
00106 }
00107 
00108 
00109 double Lumi3DReWeighting::weight3D( int pv1, int pv2, int pv3 ) {
00110 
00111   using std::min;
00112 
00113   int npm1 = min(pv1,49);
00114   int np0 = min(pv2,49);
00115   int npp1 = min(pv3,49);
00116 
00117   return Weight3D_[npm1][np0][npp1];
00118 
00119 }
00120 
00121 // This version does all of the work for you, just taking an event pointer as input
00122 
00123 double Lumi3DReWeighting::weight3D( const edm::EventBase &e ) {
00124 
00125   using std::min;
00126 
00127   // get pileup summary information
00128 
00129   Handle<std::vector< PileupSummaryInfo > >  PupInfo;
00130   e.getByLabel(edm::InputTag("addPileupInfo"), PupInfo);
00131 
00132   std::vector<PileupSummaryInfo>::const_iterator PVI;
00133 
00134   int npm1=-1;
00135   int np0=-1;
00136   int npp1=-1;
00137  
00138   for(PVI = PupInfo->begin(); PVI != PupInfo->end(); ++PVI) {
00139 
00140     int BX = PVI->getBunchCrossing();
00141 
00142     if(BX == -1) { 
00143       npm1 = PVI->getPU_NumInteractions();
00144     }
00145     if(BX == 0) { 
00146       np0 = PVI->getPU_NumInteractions();
00147     }
00148     if(BX == 1) { 
00149       npp1 = PVI->getPU_NumInteractions();
00150     }
00151 
00152   }
00153 
00154   npm1 = min(npm1,49);
00155   np0 = min(np0,49);
00156   npp1 = min(npp1,49);
00157 
00158   //  std::cout << " vertices " << npm1 << " " << np0 << " " << npp1 << std::endl; 
00159 
00160   return Weight3D_[npm1][np0][npp1];
00161  
00162 }
00163 
00164 void Lumi3DReWeighting::weight3D_set( std::string generatedFile, std::string dataFile, std::string GenHistName, std::string DataHistName, std::string WeightOutputFile )
00165 {
00166  
00167   generatedFileName_ = generatedFile;
00168   dataFileName_ = dataFile ; 
00169   GenHistName_ = GenHistName ; 
00170   DataHistName_= DataHistName ;
00171   weightFileName_  = WeightOutputFile;  
00172     
00173   std::cout<< " seting values: " << generatedFileName_ << " " << dataFileName_ << " " << GenHistName_ << " " << DataHistName_ << std::endl;
00174 
00175   generatedFile_ = boost::shared_ptr<TFile>( new TFile(generatedFileName_.c_str()) ); //MC distribution
00176   dataFile_      = boost::shared_ptr<TFile>( new TFile(dataFileName_.c_str()) );      //Data distribution
00177   
00178   boost::shared_ptr<TH1> Data_temp = boost::shared_ptr<TH1>(  (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
00179   
00180   boost::shared_ptr<TH1> MC_temp = boost::shared_ptr<TH1>(  (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
00181   
00182   
00183   MC_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(generatedFile_->Get( GenHistName_.c_str() )->Clone() )) );
00184   Data_distr_ = boost::shared_ptr<TH1>(  (static_cast<TH1*>(dataFile_->Get( DataHistName_.c_str() )->Clone() )) );
00185   
00186   // MC * data/MC = data, so the weights are data/MC:
00187   
00188   // normalize both histograms first
00189   
00190   Data_distr_->Scale( 1.0/ Data_distr_->Integral() );
00191   MC_distr_->Scale( 1.0/ MC_distr_->Integral() );
00192 }
00193 
00194 void Lumi3DReWeighting::weight3D_init( float ScaleFactor ) { 
00195 
00196   // Scale factor is used to shift target distribution (i.e. luminosity scale)  1. = no shift
00197 
00198   //create histogram to write output weights, save pain of generating them again...
00199 
00200   TH3D* WHist = new TH3D("WHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
00201   TH3D* DHist = new TH3D("DHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
00202   TH3D* MHist = new TH3D("MHist","3D weights",50,0.,50.,50,0.,50.,50,0.,50. );
00203 
00204   using std::min;
00205 
00206   if( MC_distr_->GetEntries() == 0 ) {
00207     std::cout << " MC and Data distributions are not initialized! You must call the Lumi3DReWeighting constructor. " << std::endl;
00208   }
00209 
00210 
00211   // arrays for storing number of interactions
00212 
00213   double MC_ints[50][50][50];
00214   double Data_ints[50][50][50];
00215 
00216   for (int i=0; i<50; i++) {
00217     for(int j=0; j<50; j++) {
00218       for(int k=0; k<50; k++) {
00219         MC_ints[i][j][k] = 0.;
00220         Data_ints[i][j][k] = 0.;
00221       }
00222     }
00223   }
00224 
00225   double factorial[50];
00226   double PowerSer[50];
00227   double base = 1.;
00228 
00229   factorial[0] = 1.;
00230   PowerSer[0]=1.;
00231 
00232   for (int i = 1; i<50; ++i) {
00233     base = base*float(i);
00234     factorial[i] = base;
00235   }
00236 
00237 
00238   double x;
00239   double xweight;
00240   double probi, probj, probk;
00241   double Expval, mean;
00242   int xi;
00243 
00244   // Get entries for Data, MC, fill arrays:
00245 
00246   int NMCbin = MC_distr_->GetNbinsX();
00247 
00248   for (int jbin=1;jbin<NMCbin+1;jbin++) {       
00249     x =  MC_distr_->GetBinCenter(jbin);
00250     xweight = MC_distr_->GetBinContent(jbin); //use as weight for matrix
00251 
00252     //for Summer 11, we have this int feature:
00253     xi = int(x);
00254 
00255     // Generate Poisson distribution for each value of the mean
00256 
00257     mean = double(xi);
00258 
00259     if(mean<0.) {
00260       throw cms::Exception("BadInputValue") << " Your histogram generates MC luminosity values less than zero!"
00261                                                 << " Please Check.  Terminating." << std::endl;
00262     }
00263 
00264 
00265     if(mean==0.){
00266       Expval = 1.;
00267     }
00268     else {
00269       Expval = exp(-1.*mean);
00270     }
00271 
00272     base = 1.;
00273 
00274     for (int i = 1; i<50; ++i) {
00275       base = base*mean;
00276       PowerSer[i] = base; // PowerSer is mean^i
00277     }
00278 
00279     // compute poisson probability for each Nvtx in weight matrix
00280 
00281     for (int i=0; i<50; i++) {
00282       probi = PowerSer[i]/factorial[i]*Expval;
00283       for(int j=0; j<50; j++) {
00284         probj = PowerSer[j]/factorial[j]*Expval;
00285         for(int k=0; k<50; k++) {
00286           probk = PowerSer[k]/factorial[k]*Expval;
00287           // joint probability is product of event weights multiplied by weight of input distribution bin
00288           MC_ints[i][j][k] = MC_ints[i][j][k]+probi*probj*probk*xweight;
00289         }
00290       }
00291     }
00292 
00293   }
00294   
00295   int NDatabin = Data_distr_->GetNbinsX();
00296 
00297   for (int jbin=1;jbin<NDatabin+1;jbin++) {       
00298     mean =  (Data_distr_->GetBinCenter(jbin))*ScaleFactor;
00299     xweight = Data_distr_->GetBinContent(jbin);
00300 
00301     // Generate poisson distribution for each value of the mean
00302 
00303     if(mean<0.) {
00304       throw cms::Exception("BadInputValue") << " Your histogram generates Data luminosity values less than zero!"
00305                                                 << " Please Check.  Terminating." << std::endl;
00306     }
00307 
00308     if(mean==0.){
00309       Expval = 1.;
00310     }
00311     else {
00312       Expval = exp(-1.*mean);
00313     }
00314 
00315     base = 1.;
00316 
00317     for (int i = 1; i<50; ++i) {
00318       base = base*mean;
00319       PowerSer[i] = base;
00320     }
00321 
00322     // compute poisson probability for each Nvtx in weight matrix                                                                  
00323 
00324     for (int i=0; i<50; i++) {
00325       probi = PowerSer[i]/factorial[i]*Expval;
00326       for(int j=0; j<50; j++) {
00327         probj = PowerSer[j]/factorial[j]*Expval;
00328         for(int k=0; k<50; k++) {
00329           probk = PowerSer[k]/factorial[k]*Expval;
00330           // joint probability is product of event weights multiplied by weight of input distribution bin
00331           Data_ints[i][j][k] = Data_ints[i][j][k]+probi*probj*probk*xweight;
00332         }
00333       }
00334     }
00335 
00336   }
00337  
00338 
00339   for (int i=0; i<50; i++) {  
00340     //if(i<5) std::cout << "i = " << i << std::endl;
00341     for(int j=0; j<50; j++) {
00342       for(int k=0; k<50; k++) {
00343         if( (MC_ints[i][j][k])>0.) {
00344           Weight3D_[i][j][k]  =  Data_ints[i][j][k]/MC_ints[i][j][k];
00345         }
00346         else {
00347           Weight3D_[i][j][k]  = 0.;
00348         }
00349         WHist->SetBinContent( i+1,j+1,k+1,Weight3D_[i][j][k] );
00350         DHist->SetBinContent( i+1,j+1,k+1,Data_ints[i][j][k] );
00351         MHist->SetBinContent( i+1,j+1,k+1,MC_ints[i][j][k] );
00352         //      if(i<5 && j<5 && k<5) std::cout << Weight3D_[i][j][k] << " " ;
00353       }
00354       //      if(i<5 && j<5) std::cout << std::endl;
00355     }
00356   }
00357 
00358   if(! weightFileName_.empty() ) { 
00359     std::cout << " 3D Weight Matrix initialized! " << std::endl;
00360     std::cout << " Writing weights to file " << weightFileName_ << " for re-use...  " << std::endl;
00361 
00362 
00363     TFile * outfile = new TFile(weightFileName_.c_str(),"RECREATE");
00364     WHist->Write();
00365     MHist->Write();
00366     DHist->Write();
00367     outfile->Write();
00368     outfile->Close();
00369     outfile->Delete();              
00370   }
00371 
00372   return;
00373 
00374 
00375 }
00376 
00377 
00378 void Lumi3DReWeighting::weight3D_init( std::string WeightFileName ) { 
00379 
00380   TFile *infile = new TFile(WeightFileName.c_str());
00381   TH1F *WHist = (TH1F*)infile->Get("WHist");
00382 
00383   // Check if the histogram exists           
00384   if (!WHist) {
00385     throw cms::Exception("HistogramNotFound") << " Could not find the histogram WHist in the file "
00386                                               << "in the file " << WeightFileName << "." << std::endl;
00387   }
00388 
00389   for (int i=0; i<50; i++) {  
00390     //    if(i<5) std::cout << "i = " << i << std::endl;
00391     for(int j=0; j<50; j++) {
00392       for(int k=0; k<50; k++) {
00393         Weight3D_[i][j][k] = WHist->GetBinContent(i+1,j+1,k+1);
00394         //      if(i<5 && j<5 && k<5) std::cout << Weight3D_[i][j][k] << " ";
00395       }
00396       //      if(i<5 && j<5) std::cout << std::endl;
00397 
00398     }
00399   }
00400 
00401   std::cout << " 3D Weight Matrix initialized! " << std::endl;
00402 
00403   return;
00404 
00405 
00406 }
00407 
00408 void Lumi3DReWeighting::weight3D_init( std::string MCWeightFileName, std::string DataWeightFileName ) { 
00409 
00410   TFile *infileMC = new TFile(MCWeightFileName.c_str());
00411   TH3D *MHist = (TH3D*)infileMC->Get("MHist");
00412 
00413   // Check if the histogram exists           
00414   if (!MHist) {
00415     throw cms::Exception("HistogramNotFound") << " Could not find the histogram MHist in the file "
00416                                               << "in the file " << MCWeightFileName << "." << std::endl;
00417   }
00418 
00419   TFile *infileD = new TFile(DataWeightFileName.c_str());
00420   TH3D *DHist = (TH3D*)infileD->Get("DHist");
00421 
00422   // Check if the histogram exists           
00423   if (!DHist) {
00424     throw cms::Exception("HistogramNotFound") << " Could not find the histogram DHist in the file "
00425                                               << "in the file " << DataWeightFileName << "." << std::endl;
00426   }
00427 
00428   for (int i=0; i<50; i++) {  
00429     for(int j=0; j<50; j++) {
00430       for(int k=0; k<50; k++) {
00431         Weight3D_[i][j][k] = DHist->GetBinContent(i+1,j+1,k+1)/MHist->GetBinContent(i+1,j+1,k+1);
00432       }
00433     }
00434   }
00435 
00436   std::cout << " 3D Weight Matrix initialized! " << std::endl;
00437 
00438   return;
00439 
00440 
00441 }
00442 
00443 
00444 #endif