CMS 3D CMS Logo

JetValidation Class Reference

#include <RecoJets/JetAnalyzers/interface/JetValidation.h>

Inheritance diagram for JetValidation:

edm::EDAnalyzer

List of all members.

Public Member Functions

 JetValidation (const edm::ParameterSet &)

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
void beginJob (const edm::EventSetup &)
void endJob ()

Private Attributes

TH2F AllGenEnergyVsEta2D
TProfile CaloEnergyVsEta
TProfile CaloErespVsEta
TProfile corRespVsPt
TH1F dEta
int diagPrintNum
TH1F dPhi
TH1F dR
TH1F dRcor
TProfile EBfractionVsEta
TProfile EEfractionVsEta
TH1F emEnergyFraction
TH1F emEnergyInEB
TH1F emEnergyInEE
TH1F emEnergyInHF
TProfile emEnergyVsEta
TProfile emErespVsEta
TH1F etaIC5cal
TH1F etaIC5gen
TH1F etaKT10cal
TH1F etaKT10gen
TH1F etaMC5cal
TH1F etaMC5gen
int evtCount
TProfile GenEnergyVsEta
TH2F GenEnergyVsEta2D
std::string GenType
TH1F hadEnergyInHB
TH1F hadEnergyInHE
TH1F hadEnergyInHF
TH1F hadEnergyInHO
TProfile hadEnergyVsEta
TProfile hadErespVsEta
TProfile HBfractionVsEta
TProfile HEfractionVsEta
TProfile HFfractionVsEta
TProfile HOfractionVsEta
TH1F m2jIC5cal
TH1F m2jIC5gen
TH1F m2jKT10cal
TH1F m2jKT10gen
TH1F m2jMC5cal
TH1F m2jMC5gen
TFile * m_file
int numJets
TH1F phiIC5cal
TH1F phiIC5gen
TH1F phiKT10cal
TH1F phiKT10gen
TH1F phiMC5cal
TH1F phiMC5gen
double PtHistMax
TH1F ptIC5cal
TH1F ptIC5gen
TH1F ptKT10cal
TH1F ptKT10gen
TH1F ptMC5cal
TH1F ptMC5gen
TProfile respVsPt
TProfile WindowCaloEnergyVsEta
TProfile WindowCaloErespVsEta
TProfile WindowEBfractionVsEta
TProfile WindowEEfractionVsEta
TProfile WindowEmEnergyVsEta
TProfile WindowEmErespVsEta
TProfile WindowGenEnergyVsEta
TProfile WindowHadEnergyVsEta
TProfile WindowHadErespVsEta
TProfile WindowHBfractionVsEta
TProfile WindowHEfractionVsEta
TProfile WindowHFfractionVsEta
TProfile WindowHOfractionVsEta
TProfile WindowMaxEmErespVsEta
TProfile WindowMaxHadErespVsEta
TProfile WindowMaxTowErespVsEta


Detailed Description

Definition at line 17 of file JetValidation.h.


Constructor & Destructor Documentation

JetValidation::JetValidation ( const edm::ParameterSet cfg  ) 

Definition at line 26 of file JetValidation.cc.

00026                                                        :
00027   PtHistMax( cfg.getParameter<double>( "PtHistMax" ) ),
00028   diagPrintNum( cfg.getParameter<int>( "diagPrintNum" ) ),
00029   GenType( cfg.getParameter<string>( "GenType" ) )
00030   {
00031 }


Member Function Documentation

void JetValidation::analyze ( const edm::Event evt,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 149 of file JetValidation.cc.

References AllGenEnergyVsEta2D, CaloEnergyVsEta, CaloErespVsEta, corRespVsPt, GenMuonPlsPt100GeV_cfg::cout, deltaPhi(), deltaR(), dEta, diagPrintNum, dPhi, dR, dRcor, EBfractionVsEta, EEfractionVsEta, emEnergyFraction, reco::CaloJet::emEnergyInEB(), emEnergyInEB, reco::CaloJet::emEnergyInEE(), emEnergyInEE, emEnergyInHF, reco::CaloJet::emEnergyInHF(), emEnergyVsEta, emErespVsEta, lat::endl(), reco::Particle::energy(), reco::Particle::eta(), eta, etaIC5cal, etaIC5gen, etaKT10cal, etaKT10gen, etaMC5cal, etaMC5gen, evtCount, GenEnergyVsEta, GenEnergyVsEta2D, GenType, edm::Event::getByLabel(), reco::CaloJet::hadEnergyInHB(), hadEnergyInHB, hadEnergyInHE, reco::CaloJet::hadEnergyInHE(), hadEnergyInHF, reco::CaloJet::hadEnergyInHF(), reco::CaloJet::hadEnergyInHO(), hadEnergyInHO, hadEnergyVsEta, hadErespVsEta, HBfractionVsEta, HEfractionVsEta, HFfractionVsEta, HOfractionVsEta, m2jIC5cal, m2jIC5gen, m2jKT10cal, m2jKT10gen, m2jMC5cal, m2jMC5gen, reco::CaloJet::maxEInEmTowers(), reco::CaloJet::maxEInHadTowers(), numJets, phi, phiIC5cal, phiIC5gen, phiKT10cal, phiKT10gen, phiMC5cal, phiMC5gen, PtHistMax, ptIC5cal, ptIC5gen, ptKT10cal, ptKT10gen, ptMC5cal, ptMC5gen, respVsPt, WindowCaloEnergyVsEta, WindowCaloErespVsEta, WindowEBfractionVsEta, WindowEEfractionVsEta, WindowEmEnergyVsEta, WindowEmErespVsEta, WindowGenEnergyVsEta, WindowHadEnergyVsEta, WindowHadErespVsEta, WindowHBfractionVsEta, WindowHEfractionVsEta, WindowHFfractionVsEta, WindowHOfractionVsEta, WindowMaxEmErespVsEta, WindowMaxHadErespVsEta, and WindowMaxTowErespVsEta.

00149                                                                     {
00150 
00151   evtCount++;
00152   math::XYZTLorentzVector p4jet[2], p4gen[2], p4cal[2], p4cor[2];
00153   int jetInd;
00154   Handle<CaloJetCollection> caloJets;
00155   Handle<GenJetCollection> genJets;
00156 
00157   //Fill Simple Histos
00158   
00159   //MC5 cal
00160   evt.getByLabel( "midPointCone5CaloJets", caloJets );
00161   jetInd = 0;
00162   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<numJets; ++ cal ) {
00163     ptMC5cal.Fill( cal->pt() );   
00164     etaMC5cal.Fill( cal->eta() );
00165     phiMC5cal.Fill( cal->phi() );
00166     p4jet[jetInd] = cal->p4();
00167     jetInd++;
00168   }
00169   if(jetInd==2)m2jMC5cal.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00170 
00171   //MC5 gen
00172   evt.getByLabel( "midPointCone5GenJets", genJets );
00173   jetInd = 0;
00174   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<numJets; ++ gen ) {
00175     ptMC5gen.Fill( gen->pt() );   
00176     etaMC5gen.Fill( gen->eta() );
00177     phiMC5gen.Fill( gen->phi() );
00178     p4jet[jetInd] = gen->p4();
00179     jetInd++;
00180   }
00181   if(jetInd==2)m2jMC5gen.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00182 
00183   //IC5 cal
00184   evt.getByLabel( "iterativeCone5CaloJets", caloJets );
00185   jetInd = 0;
00186   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<numJets; ++ cal ) {
00187     ptIC5cal.Fill( cal->pt() );   
00188     etaIC5cal.Fill( cal->eta() );
00189     phiIC5cal.Fill( cal->phi() );
00190     p4jet[jetInd] = cal->p4();
00191     jetInd++;
00192   }
00193   if(jetInd==2)m2jIC5cal.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00194 
00195   //IC5 gen
00196   evt.getByLabel( "iterativeCone5GenJets", genJets );
00197   jetInd = 0;
00198   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<numJets; ++ gen ) {
00199     ptIC5gen.Fill( gen->pt() );   
00200     etaIC5gen.Fill( gen->eta() );
00201     phiIC5gen.Fill( gen->phi() );
00202     p4jet[jetInd] = gen->p4();
00203 
00204     GenEnergyVsEta.Fill(fabs(gen->eta()),gen->energy());
00205     GenEnergyVsEta2D.Fill(fabs(gen->eta()),gen->energy());
00206 
00207     jetInd++;
00208   }
00209   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end(); ++ gen ) {
00210     AllGenEnergyVsEta2D.Fill(fabs(gen->eta()),gen->energy());
00211     if(gen->energy()<0.5*PtHistMax && gen->energy()>0.3*PtHistMax){
00212       WindowGenEnergyVsEta.Fill(fabs(gen->eta()),gen->energy());
00213     }
00214 
00215   } 
00216  
00217   if(jetInd==2)m2jIC5gen.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00218 
00219   //KT10 cal
00220   evt.getByLabel( "ktCaloJets", caloJets );
00221   jetInd = 0;
00222   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<numJets; ++ cal ) {
00223     ptKT10cal.Fill( cal->pt() );   
00224     etaKT10cal.Fill( cal->eta() );
00225     phiKT10cal.Fill( cal->phi() );
00226     p4jet[jetInd] = cal->p4();
00227     jetInd++;
00228   }
00229   if(jetInd==2)m2jKT10cal.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00230 
00231   //KT10 gen
00232   evt.getByLabel( "ktGenJets", genJets );
00233   jetInd = 0;
00234   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<numJets; ++ gen ) {
00235     ptKT10gen.Fill( gen->pt() );   
00236     etaKT10gen.Fill( gen->eta() );
00237     phiKT10gen.Fill( gen->phi() );
00238     p4jet[jetInd] = gen->p4();
00239     jetInd++;
00240   }
00241   if(jetInd==2)m2jKT10gen.Fill( (p4jet[0]+p4jet[1]).mass() ); 
00242 
00243  //Calorimeter Sub-System Analysis Histograms for IC5 CaloJets only
00244   evt.getByLabel( "iterativeCone5CaloJets", caloJets );
00245   jetInd = 0;
00246   for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end() && jetInd<numJets; ++ cal ) {
00247     emEnergyFraction.Fill(cal->emEnergyFraction()); 
00248     emEnergyInEB.Fill(cal->emEnergyInEB()); 
00249     emEnergyInEE.Fill(cal->emEnergyInEE()); 
00250     emEnergyInHF.Fill(cal->emEnergyInHF()); 
00251     hadEnergyInHB.Fill(cal->hadEnergyInHB()); 
00252     hadEnergyInHE.Fill(cal->hadEnergyInHE()); 
00253     hadEnergyInHF.Fill(cal->hadEnergyInHF()); 
00254     hadEnergyInHO.Fill(cal->hadEnergyInHO()); 
00255 
00256     EBfractionVsEta.Fill(fabs(cal->eta()),cal->emEnergyInEB()/cal->energy());
00257     EEfractionVsEta.Fill(fabs(cal->eta()),cal->emEnergyInEE()/cal->energy());
00258     HBfractionVsEta.Fill(fabs(cal->eta()),cal->hadEnergyInHB()/cal->energy());
00259     HOfractionVsEta.Fill(fabs(cal->eta()),cal->hadEnergyInHO()/cal->energy());
00260     HEfractionVsEta.Fill(fabs(cal->eta()),cal->hadEnergyInHE()/cal->energy());
00261     HFfractionVsEta.Fill(fabs(cal->eta()), (cal->hadEnergyInHF()+cal->emEnergyInHF())/cal->energy());
00262     
00263     CaloEnergyVsEta.Fill(fabs(cal->eta()),cal->energy());
00264     emEnergyVsEta.Fill(fabs(cal->eta()),cal->emEnergyInEB()+cal->emEnergyInEE()+cal->emEnergyInHF());
00265     hadEnergyVsEta.Fill(fabs(cal->eta()),cal->hadEnergyInHB()+cal->hadEnergyInHO()+cal->hadEnergyInHE()+cal->hadEnergyInHF());
00266    
00267    jetInd++;
00268   }
00269 
00270   //Matching for IC5 Jets for window energy fraction study
00271   evt.getByLabel( "iterativeCone5GenJets", genJets );
00272   evt.getByLabel( "iterativeCone5CaloJets", caloJets );
00273   double dRminSave;
00274   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end(); ++ gen ) { //all genJets
00275     CaloJet MatchedJet;
00276     dRminSave=1000.0;
00277     if( gen->energy()<0.5*PtHistMax && gen->energy()>0.3*PtHistMax){
00278       for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) { //all CaloJets
00279          double delR = deltaR( cal->eta(), cal->phi(), gen->eta(), gen->phi() ); 
00280          if(delR<dRminSave){
00281            dRminSave=delR;           //delta R of match
00282            MatchedJet = *cal;       //Matched Jet
00283          }
00284       }
00285       if(dRminSave>0.5)cout << "Window Match Warning: dR=" << dRminSave << ", Gen eta=" << gen->eta() << ", Cal eta=" << MatchedJet.eta() << endl;
00286       else{
00287         WindowEBfractionVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.emEnergyInEB()/MatchedJet.energy());
00288         WindowEEfractionVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.emEnergyInEE()/MatchedJet.energy());
00289         WindowHBfractionVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.hadEnergyInHB()/MatchedJet.energy());
00290         WindowHOfractionVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.hadEnergyInHO()/MatchedJet.energy());
00291         WindowHEfractionVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.hadEnergyInHE()/MatchedJet.energy());
00292         WindowHFfractionVsEta.Fill(fabs(MatchedJet.eta()), (MatchedJet.hadEnergyInHF()+MatchedJet.emEnergyInHF())/MatchedJet.energy());      
00293         WindowCaloEnergyVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.energy());
00294         WindowEmEnergyVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.emEnergyInEB()+MatchedJet.emEnergyInEE()+MatchedJet.emEnergyInHF());
00295         WindowHadEnergyVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.hadEnergyInHB()+MatchedJet.hadEnergyInHO()+MatchedJet.hadEnergyInHE()+MatchedJet.hadEnergyInHF());
00296         WindowCaloErespVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.energy()/gen->energy());
00297         WindowEmErespVsEta.Fill(fabs(MatchedJet.eta()),(MatchedJet.emEnergyInEB()+MatchedJet.emEnergyInEE()+MatchedJet.emEnergyInHF())/gen->energy());
00298         WindowHadErespVsEta.Fill(fabs(MatchedJet.eta()),(MatchedJet.hadEnergyInHB()+MatchedJet.hadEnergyInHO()+MatchedJet.hadEnergyInHE()+MatchedJet.hadEnergyInHF())/gen->energy());
00299         WindowMaxTowErespVsEta.Fill(fabs(MatchedJet.eta()),(MatchedJet.maxEInEmTowers()+MatchedJet.maxEInHadTowers())/gen->energy());
00300         WindowMaxEmErespVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.maxEInEmTowers()/gen->energy());
00301         WindowMaxHadErespVsEta.Fill(fabs(MatchedJet.eta()),MatchedJet.maxEInHadTowers()/gen->energy());
00302       }
00303     }
00304     
00305   }
00306 
00307   //Matching for IC5 Jets for leading jets energy fraction
00308   evt.getByLabel( "iterativeCone5GenJets", genJets );
00309   evt.getByLabel( "iterativeCone5CaloJets", caloJets );
00310   jetInd = 0;
00311   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<numJets; ++ gen ) { //leading genJets
00312     CaloJet MatchedJet;
00313     dRminSave=1000.0;
00314     for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) { //all CaloJets
00315        double delR = deltaR( cal->eta(), cal->phi(), gen->eta(), gen->phi() ); 
00316        if(delR<dRminSave){
00317          dRminSave=delR;           //delta R of match
00318          MatchedJet = *cal;       //Matched Jet
00319        }
00320     }
00321     if(dRminSave>0.5)cout << "Leading Match Warning: dR=" << dRminSave << ", Gen eta=" << gen->eta() << ", Cal eta=" << MatchedJet.eta() << endl;       
00322     else{
00323       float abs_eta;
00324       if(GenType=="chargedPion"){abs_eta=fabs(gen->eta());} else {abs_eta=fabs(MatchedJet.eta());}
00325       CaloErespVsEta.Fill(abs_eta, MatchedJet.energy()/gen->energy());
00326       emErespVsEta.Fill(abs_eta,(MatchedJet.emEnergyInEB()+MatchedJet.emEnergyInEE()+MatchedJet.emEnergyInHF())/gen->energy());
00327       hadErespVsEta.Fill(abs_eta,(MatchedJet.hadEnergyInHB()+MatchedJet.hadEnergyInHO()+MatchedJet.hadEnergyInHE()+MatchedJet.hadEnergyInHF())/gen->energy());
00328     }
00329     jetInd++;
00330   }
00331 
00332 
00333   //Matching for MC5 Jets: leading genjets matched to any CaloJet
00334   evt.getByLabel( "midPointCone5GenJets", genJets );
00335   evt.getByLabel( "midPointCone5CaloJets", caloJets );
00336   jetInd = 0;
00337   double dRmin[2];
00338   for( GenJetCollection::const_iterator gen = genJets->begin(); gen != genJets->end() && jetInd<numJets; ++ gen ) { //leading genJets
00339     p4gen[jetInd] = gen->p4();    //Gen 4-vector
00340     dRmin[jetInd]=1000.0;
00341     for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) { //all CaloJets
00342        double delR = deltaR( cal->eta(), cal->phi(), gen->eta(), gen->phi() ); 
00343        if(delR<dRmin[jetInd]){
00344          dRmin[jetInd]=delR;           //delta R of match
00345          p4cal[jetInd] = cal->p4();  //Matched Cal 4-vector
00346        }
00347     }
00348     dR.Fill(dRmin[jetInd]);
00349     double dphi = deltaPhi(p4cal[jetInd].phi(), p4gen[jetInd].phi());
00350     dPhi.Fill(dphi);
00351     double deta = p4cal[jetInd].eta() - p4gen[jetInd].eta();
00352     dEta.Fill(deta);
00353     if(dRmin[jetInd]>0.5)cout << "MC5 Match Warning: dR=" <<dRmin[jetInd]<<", GenPt="<<p4gen[jetInd].Pt()<<", CalPt="<<p4cal[jetInd].Pt()<<endl;
00354     jetInd++;    
00355   }
00356   //Fill Resp vs Pt profile histogram with response of two leading gen jets
00357   for( jetInd=0; jetInd<numJets; ++jetInd ){
00358     if(fabs(p4gen[jetInd].eta())<1.){
00359       respVsPt.Fill(p4gen[jetInd].Pt(), p4cal[jetInd].Pt()/p4gen[jetInd].Pt() );
00360     }
00361   }
00362 
00363   //Find the Corrected CaloJets that match the two uncorrected CaloJets
00364   evt.getByLabel( "corJetMcone5", caloJets );
00365   for( jetInd=0; jetInd<numJets; ++jetInd ){
00366     bool found=kFALSE;
00367     for( CaloJetCollection::const_iterator cor = caloJets->begin(); cor != caloJets->end() && !found; ++ cor ) { //all corrected CaloJets
00368        double delR = deltaR( cor->eta(), cor->phi(), p4cal[jetInd].eta(),  p4cal[jetInd].phi()); 
00369        if(delR<0.01){
00370          dRmin[jetInd]=delR;           //delta R of match
00371          p4cor[jetInd] = cor->p4();  //Matched Cal 4-vector
00372          found=kTRUE;
00373          dRcor.Fill(dRmin[jetInd]);
00374        }
00375     }
00376     if(!found)cout << "Warning: corrected jet not found. jetInd=" << jetInd << endl;
00377   }
00378   //Fill Resp vs Pt profile histogram with corrected response of two leading gen jets
00379   for( jetInd=0; jetInd<numJets; ++jetInd ){
00380     if(fabs(p4gen[jetInd].eta())<1.){
00381      corRespVsPt.Fill(p4gen[jetInd].Pt(), p4cor[jetInd].Pt()/p4gen[jetInd].Pt() ); 
00382     }
00383   }
00384    
00385 
00386   //Diagnostic Printout
00387   if(evtCount<=diagPrintNum){
00388 
00389     // Printout the uncorrected calojets for midpoint cone size 0.5
00390     evt.getByLabel( "midPointCone5CaloJets", caloJets );
00391     jetInd=0;
00392     for( CaloJetCollection::const_iterator cal = caloJets->begin(); cal != caloJets->end(); ++ cal ) {
00393      cout <<"CaloJet "<< jetInd << ", pt=" <<cal->pt()<<", eta="<<cal->eta()<<", phi="<<cal->phi()<<endl;
00394      jetInd++;
00395     }
00396 
00397     // Printout the Corrected jets for midpoint cone size 0.5
00398     evt.getByLabel( "corJetMcone5", caloJets );
00399     jetInd=0;
00400     for( CaloJetCollection::const_iterator cor = caloJets->begin(); cor != caloJets->end(); ++ cor ) { 
00401      cout <<"CorJet "<< jetInd << ", pt=" <<cor->pt()<<", eta="<<cor->eta()<<", phi="<<cor->phi()<<endl;
00402      jetInd++;
00403     }
00404 
00405   }
00406 
00407 
00408 }

void JetValidation::beginJob ( const edm::EventSetup  )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 33 of file JetValidation.cc.

References AllGenEnergyVsEta2D, CaloEnergyVsEta, CaloErespVsEta, corRespVsPt, GenMuonPlsPt100GeV_cfg::cout, dEta, dPhi, dR, dRcor, EBfractionVsEta, EEfractionVsEta, emEnergyFraction, emEnergyInEB, emEnergyInEE, emEnergyInHF, emEnergyVsEta, emErespVsEta, lat::endl(), etaIC5cal, etaIC5gen, etaKT10cal, etaKT10gen, etaMC5cal, etaMC5gen, evtCount, GenEnergyVsEta, GenEnergyVsEta2D, GenType, hadEnergyInHB, hadEnergyInHE, hadEnergyInHF, hadEnergyInHO, hadEnergyVsEta, hadErespVsEta, HBfractionVsEta, HEfractionVsEta, HFfractionVsEta, HOfractionVsEta, m2jIC5cal, m2jIC5gen, m2jKT10cal, m2jKT10gen, m2jMC5cal, m2jMC5gen, m_file, numJets, phiIC5cal, phiIC5gen, phiKT10cal, phiKT10gen, phiMC5cal, phiMC5gen, PtHistMax, ptIC5cal, ptIC5gen, ptKT10cal, ptKT10gen, ptMC5cal, ptMC5gen, respVsPt, WindowCaloEnergyVsEta, WindowCaloErespVsEta, WindowEBfractionVsEta, WindowEEfractionVsEta, WindowEmEnergyVsEta, WindowEmErespVsEta, WindowGenEnergyVsEta, WindowHadEnergyVsEta, WindowHadErespVsEta, WindowHBfractionVsEta, WindowHEfractionVsEta, WindowHFfractionVsEta, WindowHOfractionVsEta, WindowMaxEmErespVsEta, WindowMaxHadErespVsEta, and WindowMaxTowErespVsEta.

00033                                                  {
00034   cout << "JetValidation: Maximum bin edge for Pt Hists = " << PtHistMax << endl;
00035   if(GenType=="chargedPion"){
00036     numJets=1;
00037     cout << "charged pion analysis mode" << endl;
00038   }
00039   else
00040   {
00041     numJets=2;
00042   }
00043 
00044   //Initialize some stuff
00045   evtCount = 0;
00046 
00047   // Open the histogram file and book some associated histograms
00048   m_file=new TFile("JetValidationHistos.root","RECREATE"); 
00049   
00050   //Simple histos
00051   
00052   //MC5 cal
00053   ptMC5cal  = TH1F("ptMC5cal","p_{T} of leading CaloJets (MC5)",50,0.0,PtHistMax);
00054   etaMC5cal = TH1F("etaMC5cal","#eta of leading CaloJets (MC5)",100,-5.0,5.0);
00055   phiMC5cal = TH1F("phiMC5cal","#phi of leading CaloJets (MC5)",72,-M_PI, M_PI);
00056   m2jMC5cal = TH1F("m2jMC5cal","Dijet Mass of leading CaloJets (MC5)",100,0.0,2*PtHistMax);
00057 
00058   //MC5 gen
00059   ptMC5gen = TH1F("ptMC5gen","p_{T} of leading GenJets (MC5)",50,0.0,PtHistMax);
00060   etaMC5gen = TH1F("etaMC5gen","#eta of leading GenJets (MC5)",100,-5.0,5.0);
00061   phiMC5gen = TH1F("phiMC5gen","#phi of leading GenJets (MC5)",72,-M_PI, M_PI);
00062   m2jMC5gen = TH1F("m2jMC5gen","Dijet Mass of leading GenJets (MC5)",100,0.0,2*PtHistMax);
00063 
00064   //IC5 cal
00065   ptIC5cal  = TH1F("ptIC5cal","p_{T} of leading CaloJets (IC5)",50,0.0,PtHistMax);
00066   etaIC5cal = TH1F("etaIC5cal","#eta of leading CaloJets (IC5)",100,-5.0,5.0);
00067   phiIC5cal = TH1F("phiIC5cal","#phi of leading CaloJets (IC5)",72,-M_PI, M_PI);
00068   m2jIC5cal = TH1F("m2jIC5cal","Dijet Mass of leading CaloJets (IC5)",100,0.0,2*PtHistMax);
00069 
00070   //IC5 gen
00071   ptIC5gen = TH1F("ptIC5gen","p_{T} of leading GenJets (IC5)",50,0.0,PtHistMax);
00072   etaIC5gen = TH1F("etaIC5gen","#eta of leading GenJets (IC5)",100,-5.0,5.0);
00073   phiIC5gen = TH1F("phiIC5gen","#phi of leading GenJets (IC5)",72,-M_PI, M_PI);
00074   m2jIC5gen = TH1F("m2jIC5gen","Dijet Mass of leading GenJets (IC5)",100,0.0,2*PtHistMax);
00075 
00076   //KT10 cal
00077   ptKT10cal  = TH1F("ptKT10cal","p_{T} of leading CaloJets (KT10)",50,0.0,PtHistMax);
00078   etaKT10cal = TH1F("etaKT10cal","#eta of leading CaloJets (KT10)",100,-5.0,5.0);
00079   phiKT10cal = TH1F("phiKT10cal","#phi of leading CaloJets (KT10)",72,-M_PI, M_PI);
00080   m2jKT10cal = TH1F("m2jKT10cal","Dijet Mass of leading CaloJets (KT10)",100,0.0,2*PtHistMax);
00081 
00082   //KT10 gen
00083   ptKT10gen = TH1F("ptKT10gen","p_{T} of leading GenJets (KT10)",50,0.0,PtHistMax);
00084   etaKT10gen = TH1F("etaKT10gen","#eta of leading GenJets (KT10)",100,-5.0,5.0);
00085   phiKT10gen = TH1F("phiKT10gen","#phi of leading GenJets (KT10)",72,-M_PI, M_PI);
00086   m2jKT10gen = TH1F("m2jKT10gen","Dijet Mass of leading GenJets (KT10)",100,0.0,2*PtHistMax);
00087 
00088   //Calorimeter Sub-System Analysis Histograms for IC5 CaloJets only
00089   emEnergyFraction =  TH1F("emEnergyFraction","Leading Jets EM Fraction",100,0.0,1.00001);
00090   emEnergyInEB = TH1F("emEnergyInEB","Leading Jets emEnergyInEB",100,0.0,2*PtHistMax);
00091   emEnergyInEE = TH1F("emEnergyInEE","Leading Jets emEnergyInEE",100,0.0,2*PtHistMax);
00092   emEnergyInHF = TH1F("emEnergyInHF","Leading Jets emEnergyInHF",100,0.0,2*PtHistMax);
00093   hadEnergyInHB = TH1F("hadEnergyInHB","Leading Jets hadEnergyInHB",100,0.0,2*PtHistMax);
00094   hadEnergyInHE = TH1F("hadEnergyInHE","Leading Jets hadEnergyInHE",100,0.0,2*PtHistMax);
00095   hadEnergyInHF = TH1F("hadEnergyInHF","Leading Jets hadEnergyInHF",100,0.0,2*PtHistMax);
00096   hadEnergyInHO = TH1F("hadEnergyInHO","Leading Jets hadEnergyInHO",100,0.0,0.5*PtHistMax);
00097 
00098   EBfractionVsEta = TProfile("EBfractionVsEta","Leading Jets EBfraction vs #eta",50,0.0,5.0);
00099   EEfractionVsEta = TProfile("EEfractionVsEta","Leading Jets EEfraction vs #eta",50,0.0,5.0);
00100   HBfractionVsEta = TProfile("HBfractionVsEta","Leading Jets HBfraction vs #eta",50,0.0,5.0);
00101   HOfractionVsEta = TProfile("HOfractionVsEta","Leading Jets HOfraction vs #eta",50,0.0,5.0);
00102   HEfractionVsEta = TProfile("HEfractionVsEta","Leading Jets HEfraction vs #eta",50,0.0,5.0);
00103   HFfractionVsEta = TProfile("HFfractionVsEta","Leading Jets HFfraction vs #eta",50,0.0,5.0);
00104 
00105   CaloEnergyVsEta  = TProfile("CaloEnergyVsEta","Leading CaloJets Energy Vs. Eta",50,0.0,5.0);
00106   GenEnergyVsEta  = TProfile("GenEnergyVsEta","Leading GenJets Energy Vs. Eta",50,0.0,5.0);
00107   emEnergyVsEta  = TProfile("emEnergyVsEta","Leading Jets EM Energy Vs. Eta",50,0.0,5.0);
00108   hadEnergyVsEta  = TProfile("hadEnergyVsEta","Leading Jets HAD Energy Vs. Eta",50,0.0,5.0);
00109 
00110   CaloErespVsEta  = TProfile("CaloErespVsEta","Leading Jets Energy Response Vs. Eta",50,0.0,5.0);
00111   emErespVsEta  = TProfile("emErespVsEta","Leading Jets EM Energy Response Vs. Eta",50,0.0,5.0);
00112   hadErespVsEta  = TProfile("hadErespVsEta","Leading Jets HAD Enegy Response Vs. Eta",50,0.0,5.0);
00113 
00114 
00115   WindowEBfractionVsEta = TProfile("WindowEBfractionVsEta","Jets EBfraction vs #eta (E Window)",50,0.0,5.0);
00116   WindowEEfractionVsEta = TProfile("WindowEEfractionVsEta","Jets EEfraction vs #eta (E Window)",50,0.0,5.0);
00117   WindowHBfractionVsEta = TProfile("WindowHBfractionVsEta","Jets HBfraction vs #eta (E Window)",50,0.0,5.0);
00118   WindowHOfractionVsEta = TProfile("WindowHOfractionVsEta","Jets HOfraction vs #eta (E Window)",50,0.0,5.0);
00119   WindowHEfractionVsEta = TProfile("WindowHEfractionVsEta","Jets HEfraction vs #eta (E Window)",50,0.0,5.0);
00120   WindowHFfractionVsEta = TProfile("WindowHFfractionVsEta","Jets HFfraction vs #eta (E Window)",50,0.0,5.0);
00121 
00122 
00123   WindowCaloEnergyVsEta  = TProfile("WindowCaloEnergyVsEta","CaloJets Energy Vs. Eta (E Window)",50,0.0,5.0);
00124   WindowGenEnergyVsEta  = TProfile("WindowGenEnergyVsEta","GenJets Energy Vs. Eta (E Window)",50,0.0,5.0);
00125   WindowEmEnergyVsEta  = TProfile("WindowEmEnergyVsEta","Jets EM Energy Vs. Eta (E Window)",50,0.0,5.0);
00126   WindowHadEnergyVsEta  = TProfile("WindowHadEnergyVsEta","Jets HAD Energy Vs. Eta (E Window)",50,0.0,5.0);
00127 
00128   WindowCaloErespVsEta  = TProfile("WindowCaloErespVsEta","Jets Energy Response Vs. Eta (E Window)",50,0.0,5.0);
00129   WindowEmErespVsEta  = TProfile("WindowEmErespVsEta","Jets EM Energy Response Vs. Eta (E Window)",50,0.0,5.0);
00130   WindowHadErespVsEta  = TProfile("WindowHadErespVsEta","Jets HAD Enegy Response Vs. Eta (E Window)",50,0.0,5.0);
00131 
00132   WindowMaxTowErespVsEta  = TProfile("WindowMaxTowErespVsEta","Max EM + HAD Tower Response Vs. Jet Eta (E Window)",50,0.0,5.0);
00133   WindowMaxEmErespVsEta  = TProfile("WindowMaxEmErespVsEta","Max EM Tower Response Vs. Jet Eta (E Window)",50,0.0,5.0);
00134   WindowMaxHadErespVsEta  = TProfile("WindowMaxHadErespVsEta","Max Had Tower Response Vs. Jet Eta (E Window)",50,0.0,5.0);
00135 
00136   GenEnergyVsEta2D  = TH2F("GenEnergyVsEta2D","Leading GenJets Energy Vs. Eta",50,0.0,5.0,100,0.0,5*PtHistMax);
00137   AllGenEnergyVsEta2D  = TH2F("AllGenEnergyVsEta2D","All GenJets Energy Vs. Eta",50,0.0,5.0,100,0.0,5*PtHistMax);
00138 
00139 
00140   //Matched jets Analysis Histograms for MC5 CaloJets only
00141   dR = TH1F("dR","Leading GenJets dR with matched CaloJet",100,0,0.5);
00142   dPhi = TH1F("dPhi","Leading GenJets dPhi with matched CaloJet",200,-0.5,0.5);
00143   dEta = TH1F("dEta","Leading GenJets dEta with matched CaloJet",200,-0.5,0.5);
00144   respVsPt = TProfile("respVsPt","CaloJet Response of Leading GenJets in Barrel",100,0.0,PtHistMax/2); 
00145   dRcor = TH1F("dRcor","CorJets dR with matched CaloJet",100,0.0,0.01);
00146   corRespVsPt = TProfile("corRespVsPt","Corrected CaloJet Response of Leading GenJets in Barrel",100,0.0,PtHistMax/2); 
00147 }

void JetValidation::endJob ( void   )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 410 of file JetValidation.cc.

References m_file.

00410                            {
00411 
00412   //Write out the histogram file.
00413   m_file->Write(); 
00414 
00415 }


Member Data Documentation

TH2F JetValidation::AllGenEnergyVsEta2D [private]

Definition at line 54 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::CaloEnergyVsEta [private]

Definition at line 47 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::CaloErespVsEta [private]

Definition at line 48 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::corRespVsPt [private]

Definition at line 58 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::dEta [private]

Definition at line 57 of file JetValidation.h.

Referenced by analyze(), and beginJob().

int JetValidation::diagPrintNum [private]

Definition at line 29 of file JetValidation.h.

Referenced by analyze().

TH1F JetValidation::dPhi [private]

Definition at line 57 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::dR [private]

Definition at line 57 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::dRcor [private]

Definition at line 57 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::EBfractionVsEta [private]

Definition at line 45 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::EEfractionVsEta [private]

Definition at line 45 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::emEnergyFraction [private]

Definition at line 43 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::emEnergyInEB [private]

Definition at line 43 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::emEnergyInEE [private]

Definition at line 43 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::emEnergyInHF [private]

Definition at line 43 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::emEnergyVsEta [private]

Definition at line 47 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::emErespVsEta [private]

Definition at line 48 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaIC5cal [private]

Definition at line 37 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaIC5gen [private]

Definition at line 38 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaKT10cal [private]

Definition at line 39 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaKT10gen [private]

Definition at line 40 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaMC5cal [private]

Definition at line 35 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::etaMC5gen [private]

Definition at line 36 of file JetValidation.h.

Referenced by analyze(), and beginJob().

int JetValidation::evtCount [private]

Definition at line 64 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::GenEnergyVsEta [private]

Definition at line 47 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH2F JetValidation::GenEnergyVsEta2D [private]

Definition at line 54 of file JetValidation.h.

Referenced by analyze(), and beginJob().

std::string JetValidation::GenType [private]

Definition at line 30 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::hadEnergyInHB [private]

Definition at line 44 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::hadEnergyInHE [private]

Definition at line 44 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::hadEnergyInHF [private]

Definition at line 44 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::hadEnergyInHO [private]

Definition at line 44 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::hadEnergyVsEta [private]

Definition at line 47 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::hadErespVsEta [private]

Definition at line 48 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::HBfractionVsEta [private]

Definition at line 45 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::HEfractionVsEta [private]

Definition at line 46 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::HFfractionVsEta [private]

Definition at line 46 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::HOfractionVsEta [private]

Definition at line 46 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jIC5cal [private]

Definition at line 37 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jIC5gen [private]

Definition at line 38 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jKT10cal [private]

Definition at line 39 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jKT10gen [private]

Definition at line 40 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jMC5cal [private]

Definition at line 35 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::m2jMC5gen [private]

Definition at line 36 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TFile* JetValidation::m_file [private]

Definition at line 61 of file JetValidation.h.

Referenced by beginJob(), and endJob().

int JetValidation::numJets [private]

Definition at line 65 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiIC5cal [private]

Definition at line 37 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiIC5gen [private]

Definition at line 38 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiKT10cal [private]

Definition at line 39 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiKT10gen [private]

Definition at line 40 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiMC5cal [private]

Definition at line 35 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::phiMC5gen [private]

Definition at line 36 of file JetValidation.h.

Referenced by analyze(), and beginJob().

double JetValidation::PtHistMax [private]

Definition at line 28 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptIC5cal [private]

Definition at line 37 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptIC5gen [private]

Definition at line 38 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptKT10cal [private]

Definition at line 39 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptKT10gen [private]

Definition at line 40 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptMC5cal [private]

Definition at line 35 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TH1F JetValidation::ptMC5gen [private]

Definition at line 36 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::respVsPt [private]

Definition at line 58 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowCaloEnergyVsEta [private]

Definition at line 52 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowCaloErespVsEta [private]

Definition at line 51 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowEBfractionVsEta [private]

Definition at line 49 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowEEfractionVsEta [private]

Definition at line 49 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowEmEnergyVsEta [private]

Definition at line 52 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowEmErespVsEta [private]

Definition at line 51 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowGenEnergyVsEta [private]

Definition at line 52 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHadEnergyVsEta [private]

Definition at line 52 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHadErespVsEta [private]

Definition at line 51 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHBfractionVsEta [private]

Definition at line 49 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHEfractionVsEta [private]

Definition at line 50 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHFfractionVsEta [private]

Definition at line 50 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowHOfractionVsEta [private]

Definition at line 50 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowMaxEmErespVsEta [private]

Definition at line 53 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowMaxHadErespVsEta [private]

Definition at line 53 of file JetValidation.h.

Referenced by analyze(), and beginJob().

TProfile JetValidation::WindowMaxTowErespVsEta [private]

Definition at line 53 of file JetValidation.h.

Referenced by analyze(), and beginJob().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:26:14 2009 for CMSSW by  doxygen 1.5.4