CMS 3D CMS Logo

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