CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/QCDAnalysis/UEAnalysis/src/AnalysisRootpleProducer.cc

Go to the documentation of this file.
00001 // Authors: F. Ambroglini, L. Fano'
00002 #include <QCDAnalysis/UEAnalysis/interface/AnalysisRootpleProducer.h>
00003 #include "FWCore/Common/interface/TriggerNames.h"
00004  
00005 using namespace edm;
00006 using namespace std;
00007 using namespace reco;
00008 
00009 class GreaterPt{
00010 public:
00011   bool operator()( const math::XYZTLorentzVector& a, const math::XYZTLorentzVector& b) {
00012     return a.pt() > b.pt();
00013   }
00014 };
00015 
00016 class GenJetSort{
00017 public:
00018   bool operator()(const GenJet& a, const GenJet& b) {
00019     return a.pt() > b.pt();
00020   }
00021 };
00022 
00023 class BasicJetSort{
00024 public:
00025   bool operator()(const BasicJet& a, const BasicJet& b) {
00026     return a.pt() > b.pt();
00027   }
00028 };
00029 
00030 class CaloJetSort{
00031 public:
00032   bool operator()(const CaloJet& a, const CaloJet& b) {
00033     return a.pt() > b.pt();
00034   }
00035 };
00036  
00037 
00038 void AnalysisRootpleProducer::store(){
00039 
00040   AnalysisTree->Fill();
00041 
00042   NumberMCParticles=0;
00043   NumberTracks=0;
00044   NumberInclusiveJet=0;
00045   NumberChargedJet=0;
00046   NumberTracksJet=0;
00047   NumberCaloJet=0;
00048 }
00049 
00050 void AnalysisRootpleProducer::fillEventInfo(int e){
00051   EventKind = e;
00052 }
00053 
00054 void AnalysisRootpleProducer::fillMCParticles(float p, float pt, float eta, float phi){
00055   MomentumMC[NumberMCParticles]=p;
00056   TransverseMomentumMC[NumberMCParticles]=pt;
00057   EtaMC[NumberMCParticles]=eta;
00058   PhiMC[NumberMCParticles]=phi;
00059   NumberMCParticles++;
00060 }
00061 
00062 void AnalysisRootpleProducer::fillTracks(float p, float pt, float eta, float phi){
00063   MomentumTK[NumberTracks]=p;
00064   TransverseMomentumTK[NumberTracks]=pt;
00065   EtaTK[NumberTracks]=eta;
00066   PhiTK[NumberTracks]=phi;
00067   NumberTracks++;
00068 }
00069 
00070 void AnalysisRootpleProducer::fillInclusiveJet(float p, float pt, float eta, float phi){
00071   MomentumIJ[NumberInclusiveJet]=p;
00072   TransverseMomentumIJ[NumberInclusiveJet]=pt;
00073   EtaIJ[NumberInclusiveJet]=eta;
00074   PhiIJ[NumberInclusiveJet]=phi;
00075   NumberInclusiveJet++;
00076 }
00077 
00078 void AnalysisRootpleProducer::fillChargedJet(float p, float pt, float eta, float phi){
00079   MomentumCJ[NumberChargedJet]=p;
00080   TransverseMomentumCJ[NumberChargedJet]=pt;
00081   EtaCJ[NumberChargedJet]=eta;
00082   PhiCJ[NumberChargedJet]=phi;
00083   NumberChargedJet++;
00084 }
00085 
00086 void AnalysisRootpleProducer::fillTracksJet(float p, float pt, float eta, float phi){
00087   MomentumTJ[NumberTracksJet]=p;
00088   TransverseMomentumTJ[NumberTracksJet]=pt;
00089   EtaTJ[NumberTracksJet]=eta;
00090   PhiTJ[NumberTracksJet]=phi;
00091   NumberTracksJet++;
00092 }
00093 
00094 void AnalysisRootpleProducer::fillCaloJet(float p, float pt, float eta, float phi){
00095   MomentumEHJ[NumberCaloJet]=p;
00096   TransverseMomentumEHJ[NumberCaloJet]=pt;
00097   EtaEHJ[NumberCaloJet]=eta;
00098   PhiEHJ[NumberCaloJet]=phi;
00099   NumberCaloJet++;
00100 }
00101 
00102 AnalysisRootpleProducer::AnalysisRootpleProducer( const ParameterSet& pset )
00103 {
00104   // flag to ignore gen-level analysis
00105   onlyRECO = pset.getUntrackedParameter<bool>("OnlyRECO",false);
00106 
00107   // particle, track and jet collections
00108   mcEvent = pset.getUntrackedParameter<InputTag>("MCEvent",std::string(""));
00109   genJetCollName = pset.getUntrackedParameter<InputTag>("GenJetCollectionName",std::string(""));
00110   chgJetCollName = pset.getUntrackedParameter<InputTag>("ChgGenJetCollectionName",std::string(""));
00111   tracksJetCollName = pset.getUntrackedParameter<InputTag>("TracksJetCollectionName",std::string(""));
00112   recoCaloJetCollName = pset.getUntrackedParameter<InputTag>("RecoCaloJetCollectionName",std::string(""));
00113   chgGenPartCollName = pset.getUntrackedParameter<InputTag>("ChgGenPartCollectionName",std::string(""));
00114   tracksCollName = pset.getUntrackedParameter<InputTag>("TracksCollectionName",std::string(""));
00115 
00116   // trigger results
00117   triggerResultsTag = pset.getParameter<InputTag>("triggerResults");
00118   //   hltFilterTag      = pset.getParameter<InputTag>("hltFilter");
00119   //   triggerName       = pset.getParameter<InputTag>("triggerName");
00120 
00121   piG = acos(-1.);
00122   NumberMCParticles=0;
00123   NumberTracks=0;
00124   NumberInclusiveJet=0;
00125   NumberChargedJet=0;
00126   NumberTracksJet=0;
00127   NumberCaloJet=0;
00128 }
00129 
00130 void AnalysisRootpleProducer::beginJob()
00131 {
00132  
00133   // use TFileService for output to root file
00134   AnalysisTree = fs->make<TTree>("AnalysisTree","MBUE Analysis Tree ");
00135 
00136   AnalysisTree->Branch("EventKind",&EventKind,"EventKind/I");
00137 
00138   // store p, pt, eta, phi for particles and jets
00139 
00140   // GenParticles at hadron level  
00141   AnalysisTree->Branch("NumberMCParticles",&NumberMCParticles,"NumberMCParticles/I");
00142   AnalysisTree->Branch("MomentumMC",MomentumMC,"MomentumMC[NumberMCParticles]/F");
00143   AnalysisTree->Branch("TransverseMomentumMC",TransverseMomentumMC,"TransverseMomentumMC[NumberMCParticles]/F");
00144   AnalysisTree->Branch("EtaMC",EtaMC,"EtaMC[NumberMCParticles]/F");
00145   AnalysisTree->Branch("PhiMC",PhiMC,"PhiMC[NumberMCParticles]/F");
00146   
00147   // tracks
00148   AnalysisTree->Branch("NumberTracks",&NumberTracks,"NumberTracks/I");
00149   AnalysisTree->Branch("MomentumTK",MomentumTK,"MomentumTK[NumberTracks]/F");
00150   AnalysisTree->Branch("TrasverseMomentumTK",TransverseMomentumTK,"TransverseMomentumTK[NumberTracks]/F");
00151   AnalysisTree->Branch("EtaTK",EtaTK,"EtaTK[NumberTracks]/F");
00152   AnalysisTree->Branch("PhiTK",PhiTK,"PhiTK[NumberTracks]/F");
00153   
00154   // GenJets
00155   AnalysisTree->Branch("NumberInclusiveJet",&NumberInclusiveJet,"NumberInclusiveJet/I");
00156   AnalysisTree->Branch("MomentumIJ",MomentumIJ,"MomentumIJ[NumberInclusiveJet]/F");
00157   AnalysisTree->Branch("TrasverseMomentumIJ",TransverseMomentumIJ,"TransverseMomentumIJ[NumberInclusiveJet]/F");
00158   AnalysisTree->Branch("EtaIJ",EtaIJ,"EtaIJ[NumberInclusiveJet]/F");
00159   AnalysisTree->Branch("PhiIJ",PhiIJ,"PhiIJ[NumberInclusiveJet]/F");
00160   
00161   // jets from charged GenParticles
00162   AnalysisTree->Branch("NumberChargedJet",&NumberChargedJet,"NumberChargedJet/I");
00163   AnalysisTree->Branch("MomentumCJ",MomentumCJ,"MomentumCJ[NumberChargedJet]/F");
00164   AnalysisTree->Branch("TrasverseMomentumCJ",TransverseMomentumCJ,"TransverseMomentumCJ[NumberChargedJet]/F");
00165   AnalysisTree->Branch("EtaCJ",EtaCJ,"EtaCJ[NumberChargedJet]/F");
00166   AnalysisTree->Branch("PhiCJ",PhiCJ,"PhiCJ[NumberChargedJet]/F");
00167   
00168   // jets from tracks
00169   AnalysisTree->Branch("NumberTracksJet",&NumberTracksJet,"NumberTracksJet/I");
00170   AnalysisTree->Branch("MomentumTJ",MomentumTJ,"MomentumTJ[NumberTracksJet]/F");
00171   AnalysisTree->Branch("TrasverseMomentumTJ",TransverseMomentumTJ,"TransverseMomentumTJ[NumberTracksJet]/F");
00172   AnalysisTree->Branch("EtaTJ",EtaTJ,"EtaTJ[NumberTracksJet]/F");
00173   AnalysisTree->Branch("PhiTJ",PhiTJ,"PhiTJ[NumberTracksJet]/F");
00174   
00175   // jets from calorimeter towers
00176   AnalysisTree->Branch("NumberCaloJet",&NumberCaloJet,"NumberCaloJet/I");
00177   AnalysisTree->Branch("MomentumEHJ",MomentumEHJ,"MomentumEHJ[NumberCaloJet]/F");
00178   AnalysisTree->Branch("TrasverseMomentumEHJ",TransverseMomentumEHJ,"TransverseMomentumEHJ[NumberCaloJet]/F");
00179   AnalysisTree->Branch("EtaEHJ",EtaEHJ,"EtaEHJ[NumberCaloJet]/F");
00180   AnalysisTree->Branch("PhiEHJ",PhiEHJ,"PhiEHJ[NumberCaloJet]/F");
00181   
00182 
00183   // alternative storage method:
00184   // save TClonesArrays of TLorentzVectors
00185   // i.e. store 4-vectors of particles and jets
00186 
00187   MonteCarlo = new TClonesArray("TLorentzVector", 10000);
00188   AnalysisTree->Branch("MonteCarlo", "TClonesArray", &MonteCarlo, 128000, 0);
00189 
00190   Track = new TClonesArray("TLorentzVector", 10000);
00191   AnalysisTree->Branch("Track", "TClonesArray", &Track, 128000, 0);
00192 
00193   InclusiveJet = new TClonesArray("TLorentzVector", 10000);
00194   AnalysisTree->Branch("InclusiveJet", "TClonesArray", &InclusiveJet, 128000, 0);
00195 
00196   ChargedJet = new TClonesArray("TLorentzVector", 10000);
00197   AnalysisTree->Branch("ChargedJet", "TClonesArray", &ChargedJet, 128000, 0);
00198 
00199   TracksJet = new TClonesArray("TLorentzVector", 10000);
00200   AnalysisTree->Branch("TracksJet", "TClonesArray", &TracksJet, 128000, 0);
00201 
00202   CalorimeterJet = new TClonesArray("TLorentzVector", 10000);
00203   AnalysisTree->Branch("CalorimeterJet", "TClonesArray", &CalorimeterJet, 128000, 0);
00204 
00205   acceptedTriggers = new TClonesArray("TObjString", 10000);
00206   AnalysisTree->Branch("acceptedTriggers", "TClonesArray", &acceptedTriggers, 128000, 0);
00207 
00208 }
00209 
00210   
00211 void AnalysisRootpleProducer::analyze( const Event& e, const EventSetup& )
00212 {
00213   
00214   e.getByLabel( triggerResultsTag, triggerResults );
00215   const edm::TriggerNames & triggerNames = e.triggerNames(*triggerResults);
00216 
00217   acceptedTriggers->Clear();
00218   unsigned int iAcceptedTriggers( 0 ); 
00219   if ( triggerResults.product()->wasrun() )
00220     {
00221       //cout << "at least one path out of " << triggerResults.product()->size() << " ran? " << triggerResults.product()->wasrun() << endl;
00222   
00223       if ( triggerResults.product()->accept() ) 
00224         {
00225           //cout << endl << "at least one path accepted? " << triggerResults.product()->accept() << endl;
00226 
00227           const unsigned int n_TriggerResults( triggerResults.product()->size() );
00228           for ( unsigned int itrig( 0 ); itrig < n_TriggerResults; ++itrig )
00229             {
00230               if ( triggerResults.product()->accept( itrig ) )
00231                 {
00232                   //cout << "path " << triggerNames.triggerName( itrig );
00233                   //cout << ", module index " << triggerResults.product()->index( itrig );
00234                   //cout << ", state (Ready = 0, Pass = 1, Fail = 2, Exception = 3) " << triggerResults.product()->state( itrig );
00235                   //cout << ", accept " << triggerResults.product()->accept( itrig );
00236                   //cout << endl;
00237 
00238                   // save name of accepted trigger path
00239                   new((*acceptedTriggers)[iAcceptedTriggers]) TObjString( (triggerNames.triggerName( itrig )).c_str() );
00240                   ++iAcceptedTriggers;
00241                 }
00242             }
00243         }
00244     }
00245 
00246   // gen level analysis
00247   // skipped, if onlyRECO flag set to true
00248 
00249   if(!onlyRECO){
00250     
00251     e.getByLabel( mcEvent           , EvtHandle        );
00252     e.getByLabel( chgGenPartCollName, CandHandleMC     );
00253     e.getByLabel( chgJetCollName    , ChgGenJetsHandle );
00254     e.getByLabel( genJetCollName    , GenJetsHandle    );
00255     
00256     const HepMC::GenEvent* Evt = EvtHandle->GetEvent() ;
00257     
00258     EventKind = Evt->signal_process_id();
00259 
00260     std::vector<math::XYZTLorentzVector> GenPart;
00261     std::vector<GenJet> ChgGenJetContainer;
00262     std::vector<GenJet> GenJetContainer;
00263     
00264     GenPart.clear();
00265     ChgGenJetContainer.clear();
00266     GenJetContainer.clear();
00267     MonteCarlo->Clear();
00268     InclusiveJet->Clear();
00269     ChargedJet->Clear();
00270 
00271     // jets from charged particles at hadron level
00272     if (ChgGenJetsHandle->size()){
00273 
00274       for ( GenJetCollection::const_iterator it(ChgGenJetsHandle->begin()), itEnd(ChgGenJetsHandle->end());
00275             it!=itEnd; ++it)
00276         {
00277           ChgGenJetContainer.push_back(*it);
00278         }
00279 
00280       std::stable_sort(ChgGenJetContainer.begin(),ChgGenJetContainer.end(),GenJetSort());
00281 
00282       std::vector<GenJet>::const_iterator it(ChgGenJetContainer.begin()), itEnd(ChgGenJetContainer.end());
00283       for ( int iChargedJet(0); it != itEnd; ++it, ++iChargedJet)
00284         {
00285           fillChargedJet(it->p(),it->pt(),it->eta(),it->phi());
00286           new((*ChargedJet)[iChargedJet]) TLorentzVector(it->px(), it->py(), it->pz(), it->energy());
00287         }
00288     }
00289 
00290 
00291     // GenJets
00292     if (GenJetsHandle->size()){
00293 
00294       for ( GenJetCollection::const_iterator it(GenJetsHandle->begin()), itEnd(GenJetsHandle->end());
00295             it!=itEnd; ++it )
00296         {
00297           GenJetContainer.push_back(*it);
00298         }
00299 
00300       std::stable_sort(GenJetContainer.begin(),GenJetContainer.end(),GenJetSort());
00301 
00302       std::vector<GenJet>::const_iterator it(GenJetContainer.begin()), itEnd(GenJetContainer.end());
00303       for ( int iInclusiveJet(0); it != itEnd; ++it, ++iInclusiveJet)
00304         {
00305           fillInclusiveJet(it->p(),it->pt(),it->eta(),it->phi());
00306           new((*InclusiveJet)[iInclusiveJet]) TLorentzVector(it->px(), it->py(), it->pz(), it->energy());
00307         }
00308     }
00309 
00310 
00311     // hadron level particles
00312     if (CandHandleMC->size()){
00313       
00314       for (vector<GenParticle>::const_iterator it(CandHandleMC->begin()), itEnd(CandHandleMC->end());
00315            it != itEnd;it++)
00316         {
00317           GenPart.push_back(it->p4());
00318         }
00319 
00320       std::stable_sort(GenPart.begin(),GenPart.end(),GreaterPt());
00321 
00322       std::vector<math::XYZTLorentzVector>::const_iterator it(GenPart.begin()), itEnd(GenPart.end());
00323       for( int iMonteCarlo(0); it != itEnd; ++it, ++iMonteCarlo )
00324         {
00325           fillMCParticles(it->P(),it->Pt(),it->Eta(),it->Phi());
00326           new((*MonteCarlo)[iMonteCarlo]) TLorentzVector(it->Px(), it->Py(), it->Pz(), it->E());
00327         }
00328     }
00329 
00330   } 
00331 
00332   
00333   // reco level analysis
00334   
00335   e.getByLabel( tracksCollName     , CandHandleRECO     );
00336   e.getByLabel( recoCaloJetCollName, RecoCaloJetsHandle );
00337   e.getByLabel( tracksJetCollName  , TracksJetsHandle   );
00338   
00339   std::vector<math::XYZTLorentzVector> Tracks;
00340   std::vector<BasicJet> TracksJetContainer;
00341   std::vector<CaloJet> RecoCaloJetContainer;
00342   
00343   Tracks.clear();
00344   TracksJetContainer.clear();
00345   RecoCaloJetContainer.clear();
00346   
00347   Track->Clear();
00348   TracksJet->Clear();
00349   CalorimeterJet->Clear();
00350 
00351   if(RecoCaloJetsHandle->size())
00352     {
00353     for(CaloJetCollection::const_iterator it(RecoCaloJetsHandle->begin()), itEnd(RecoCaloJetsHandle->end());
00354         it!=itEnd;++it)
00355       {
00356         RecoCaloJetContainer.push_back(*it);
00357       }
00358     std::stable_sort(RecoCaloJetContainer.begin(),RecoCaloJetContainer.end(),CaloJetSort());
00359 
00360     std::vector<CaloJet>::const_iterator it(RecoCaloJetContainer.begin()), itEnd(RecoCaloJetContainer.end());
00361     for( int iCalorimeterJet(0); it != itEnd; ++it, ++iCalorimeterJet)
00362       {
00363         fillCaloJet(it->p(),it->pt(),it->eta(),it->phi());
00364         new((*CalorimeterJet)[iCalorimeterJet]) TLorentzVector(it->px(), it->py(), it->pz(), it->energy());
00365       }
00366     }
00367     
00368   if(TracksJetsHandle->size())
00369     {
00370       for(BasicJetCollection::const_iterator it(TracksJetsHandle->begin()), itEnd(TracksJetsHandle->end());
00371           it!=itEnd;++it)
00372         {
00373           TracksJetContainer.push_back(*it);
00374         }
00375       std::stable_sort(TracksJetContainer.begin(),TracksJetContainer.end(),BasicJetSort());
00376     
00377       std::vector<BasicJet>::const_iterator it(TracksJetContainer.begin()), itEnd(TracksJetContainer.end());
00378       for(int iTracksJet(0); it != itEnd; ++it, ++iTracksJet)
00379         {
00380           fillTracksJet(it->p(),it->pt(),it->eta(),it->phi());
00381           new((*TracksJet)[iTracksJet]) TLorentzVector(it->px(), it->py(), it->pz(), it->energy());
00382         }
00383     }
00384   
00385   if(CandHandleRECO->size())
00386     {
00387       for(CandidateCollection::const_iterator it(CandHandleRECO->begin()), itEnd(CandHandleRECO->end());
00388           it!=itEnd;++it)
00389         {
00390           Tracks.push_back(it->p4());
00391         }
00392       std::stable_sort(Tracks.begin(),Tracks.end(),GreaterPt());
00393     
00394       std::vector<math::XYZTLorentzVector>::const_iterator it( Tracks.begin()), itEnd(Tracks.end());
00395       for(int iTracks(0); it != itEnd; ++it, ++iTracks)
00396         {
00397           fillTracks(it->P(),it->Pt(),it->Eta(),it->Phi());
00398           new ((*Track)[iTracks]) TLorentzVector(it->Px(), it->Py(), it->Pz(), it->E());
00399         }
00400     }
00401   
00402   store();
00403 }
00404 
00405 void AnalysisRootpleProducer::endJob()
00406 {
00407 }
00408