CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/TopQuarkAnalysis/Examples/bin/TopHypothesisFWLiteAnalyzer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <string>
00003 #include <cstdlib>
00004 #include <sstream>
00005 #include <fstream>
00006 #include <iostream>
00007 
00008 #include "FWCore/FWLite/interface/AutoLibraryLoader.h"
00009 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
00010 
00011 #include "TopQuarkAnalysis/Examples/bin/NiceStyle.cc"
00012 #include "TopQuarkAnalysis/Examples/interface/RootSystem.h"
00013 #include "TopQuarkAnalysis/Examples/interface/RootHistograms.h"
00014 #include "TopQuarkAnalysis/Examples/interface/RootPostScript.h"
00015 
00016 int main(int argc, char* argv[]) 
00017 {
00018   if( argc<4 ){
00019     // -------------------------------------------------  
00020     std::cerr << "ERROR: Wrong number of arguments!"     << std::endl 
00021               << "Please specify: * file name"    << std::endl
00022               << "                * process name" << std::endl
00023               << "                * HypoClassKey" << std::endl
00024               << "Example: TopHypothesisFWLiteAnalyzer ttSemiLepEvtBuilder.root TEST kGeom" << std::endl;
00025     // -------------------------------------------------  
00026     return -1;
00027   }
00028 
00029   // get HypoClassKey
00030   std::string hypoClassKey = argv[3];
00031 
00032   // load framework libraries
00033   gSystem->Load( "libFWCoreFWLite" );
00034   AutoLibraryLoader::enable();
00035   
00036   // set nice style for histograms
00037   setNiceStyle();
00038 
00039   // define some histograms 
00040   TH1F* hadWPt_    = new TH1F("hadWPt",     "p_{t} (W_{had}) [GeV]", 100,  0., 500.);
00041   TH1F* hadWMass_  = new TH1F("hadWMass",   "M (W_{had}) [GeV]",      50,  0., 150.);
00042   TH1F* hadTopPt_  = new TH1F("hadTopPt",   "p_{t} (t_{had}) [GeV]", 100,  0., 500.);
00043   TH1F* hadTopMass_= new TH1F("hadTopMass", "M (t_{had}) [GeV]",      50, 50., 250.);
00044 
00045   TH1F* lepWPt_    = new TH1F("lepWPt",     "p_{t} (W_{lep}) [GeV]", 100,  0., 500.);
00046   TH1F* lepWMass_  = new TH1F("lepWMass",   "M (W_{lep}) [GeV]",      50,  0., 150.);
00047   TH1F* lepTopPt_  = new TH1F("lepTopPt",   "p_{t} (t_{lep}) [GeV]", 100,  0., 500.);
00048   TH1F* lepTopMass_= new TH1F("lepTopMass", "M (t_{lep}) [GeV]",      50, 50., 250.);
00049 
00050   // -------------------------------------------------  
00051   std::cout << "open file: " << argv[1] << std::endl;
00052   // -------------------------------------------------
00053   TFile* inFile = TFile::Open(argv[1]);
00054   TTree* events_= 0;
00055   if( inFile ) inFile->GetObject("Events", events_); 
00056   if( events_==0 ){
00057     // -------------------------------------------------  
00058     std::cerr << "ERROR: Unable to retrieve TTree Events!"           << std::endl
00059               << "Either wrong file name or the tree doesn't exist." << std::endl;
00060     // -------------------------------------------------  
00061     return -1;
00062   }
00063 
00064   // acess branch of ttSemiLepEvent
00065   char decayName[50];
00066   sprintf(decayName, "recoGenParticles_decaySubset__%s.obj", argv[2]);
00067   TBranch* decay_   = events_->GetBranch( decayName ); // referred to from within TtGenEvent class
00068   assert( decay_ != 0 ); 
00069   char genEvtName[50];
00070   sprintf(genEvtName, "TtGenEvent_genEvt__%s.obj", argv[2]);
00071   TBranch* genEvt_  = events_->GetBranch( genEvtName ); // referred to from within TtSemiLeptonicEvent class
00072   assert( genEvt_ != 0 ); 
00073   char semiLepEvtName[50];
00074   sprintf(semiLepEvtName, "TtSemiLeptonicEvent_ttSemiLepEvent__%s.obj", argv[2]);
00075   TBranch* semiLepEvt_ = events_->GetBranch( semiLepEvtName ); 
00076   assert( semiLepEvt_ != 0 );
00077   
00078   // loop over events and fill histograms  
00079   int nevt = events_->GetEntries();
00080   TtSemiLeptonicEvent semiLepEvt;
00081   // -------------------------------------------------  
00082   std::cout << "start looping " << nevt << " events..." << std::endl;
00083   // -------------------------------------------------
00084   for(int evt=0; evt<nevt; ++evt){
00085     // set branch address
00086     semiLepEvt_->SetAddress( &semiLepEvt );
00087     // get event
00088     decay_     ->GetEntry( evt );
00089     genEvt_    ->GetEntry( evt );
00090     semiLepEvt_->GetEntry( evt );
00091     events_    ->GetEntry( evt, 0 );
00092 
00093     // -------------------------------------------------  
00094     if(evt>0 && !(evt%10)) std::cout << "  processing event: " << evt << std::endl;
00095     // -------------------------------------------------  
00096 
00097     // fill histograms
00098     if( !semiLepEvt.isHypoAvailable(hypoClassKey) ){
00099       std::cerr << "NonValidHyp:: " << "Hypothesis not available for this event" << std::endl;
00100       continue;
00101     }
00102     if( !semiLepEvt.isHypoValid(hypoClassKey) ){
00103       std::cerr << "NonValidHyp::" << "Hypothesis not valid for this event" << std::endl;
00104       continue;
00105     }
00106     
00107     const reco::Candidate* hadTop = semiLepEvt.hadronicDecayTop(hypoClassKey);
00108     const reco::Candidate* hadW   = semiLepEvt.hadronicDecayW  (hypoClassKey);
00109     const reco::Candidate* lepTop = semiLepEvt.leptonicDecayTop(hypoClassKey);
00110     const reco::Candidate* lepW   = semiLepEvt.leptonicDecayW  (hypoClassKey);
00111     
00112     if(hadTop && hadW && lepTop && lepW){
00113       hadWPt_    ->Fill( hadW->pt()    );
00114       hadWMass_  ->Fill( hadW->mass()  );
00115       hadTopPt_  ->Fill( hadTop->pt()  );
00116       hadTopMass_->Fill( hadTop->mass());
00117       
00118       lepWPt_    ->Fill( lepW->pt()    );
00119       lepWMass_  ->Fill( lepW->mass()  );
00120       lepTopPt_  ->Fill( lepTop->pt()  );
00121       lepTopMass_->Fill( lepTop->mass());
00122     }
00123   }
00124   // -------------------------------------------------  
00125   std::cout << "close file" << std::endl;
00126   // -------------------------------------------------
00127   inFile->Close();
00128   
00129   // save histograms to file
00130   TFile outFile( "analyzeHypothesis.root", "recreate" );
00131   // strip the leading "k" from the hypoClassKey to build directory name
00132   TString outDir = "analyze" + std::string(hypoClassKey, 1, hypoClassKey.length());
00133   outFile.mkdir(outDir);
00134   outFile.cd(outDir);
00135   hadWPt_    ->Write( );
00136   hadWMass_  ->Write( );
00137   hadTopPt_  ->Write( );
00138   hadTopMass_->Write( );
00139   lepWPt_    ->Write( );
00140   lepWMass_  ->Write( );
00141   lepTopPt_  ->Write( );
00142   lepTopMass_->Write( );
00143   outFile.Close();
00144   
00145   // free allocated space
00146   delete hadWPt_;
00147   delete hadWMass_;
00148   delete hadTopPt_;
00149   delete hadTopMass_;
00150   delete lepWPt_;
00151   delete lepWMass_;
00152   delete lepTopPt_;
00153   delete lepTopMass_;
00154 
00155   return 0;
00156 }