CMS 3D CMS Logo

JPTAnalyzer Class Reference

Inheritance diagram for JPTAnalyzer:

edm::EDAnalyzer

List of all members.

Public Member Functions

 JPTAnalyzer (const edm::ParameterSet &)
 ~JPTAnalyzer ()

Private Member Functions

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

Private Attributes

string calojetsSrc
double DRMAXgjet1
double DRMAXgjet2
double EtaGen1
double EtaGen2
double EtaRaw1
double EtaRaw2
double EtGen1
double EtGen2
double EtJPT1
double EtJPT2
double EtMCJ1
double EtMCJ2
double EtRaw1
double EtRaw2
double EtZSP1
double EtZSP2
string fOutputFileName
string genjetsSrc
TFile * hOutputFile
string JetCorrectionJPT
double PhiGen1
double PhiGen2
double PhiRaw1
double PhiRaw2
TTree * t1
string zspjetsSrc


Detailed Description

Definition at line 78 of file JPTAnalyzer.cc.


Constructor & Destructor Documentation

JPTAnalyzer::JPTAnalyzer ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 171 of file JPTAnalyzer.cc.

References calojetsSrc, fOutputFileName, genjetsSrc, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), JetCorrectionJPT, and zspjetsSrc.

00173 {
00174    //now do what ever initialization is needed
00175   using namespace edm;
00176   // 
00177   // get name of output file with histogramms
00178   fOutputFileName = iConfig.getUntrackedParameter<string>("HistOutFile");
00179   //
00180   // get names of input object collections
00181   // raw calo jets
00182   calojetsSrc   = iConfig.getParameter< std::string > ("calojets");
00183   // calo jets after zsp corrections
00184   zspjetsSrc    = iConfig.getParameter< std::string > ("zspjets");
00185   genjetsSrc    = iConfig.getParameter< std::string > ("genjets");
00186   //
00187   // MC jet energy corrections
00188   //  JetCorrectionMCJ = iConfig.getParameter< std::string > ("JetCorrectionMCJ");
00189   // ZSP jet energy corrections
00190   //  JetCorrectionZSP = iConfig.getParameter< std::string > ("JetCorrectionZSP");
00191   // Jet+tracks energy corrections
00192   JetCorrectionJPT = iConfig.getParameter< std::string > ("JetCorrectionJPT");
00193 }

JPTAnalyzer::~JPTAnalyzer (  ) 

Definition at line 196 of file JPTAnalyzer.cc.

00197 {
00198  
00199    // do anything here that needs to be done at desctruction time
00200    // (e.g. close files, deallocate resources etc.)
00201 
00202 }


Member Function Documentation

void JPTAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 211 of file JPTAnalyzer.cc.

References calojetsSrc, JetCorrector::correction(), DRMAXgjet1, DRMAXgjet2, EtaGen1, EtaGen2, EtaRaw1, EtaRaw2, EtGen1, EtGen2, EtJPT1, EtJPT2, EtMCJ1, EtMCJ2, EtRaw1, EtRaw2, EtZSP1, EtZSP2, genjetsSrc, edm::Event::getByLabel(), JetCorrector::getJetCorrector(), metsig::jet, JetCorrectionJPT, l1, l2, PhiGen1, PhiGen2, PhiRaw1, PhiRaw2, t1, and zspjetsSrc.

00212 {
00213    using namespace edm;
00214 
00215    // initialize vector containing two highest Et gen jets > 20 GeV
00216    // in this example they are checked not to be leptons from Z->ll decay (DR match)
00217    vector<HepLorentzVector> gjets;
00218    gjets.clear();
00219 
00220    // initialize tree variables
00221    EtaGen1 = 0.;
00222    PhiGen1 = 0.;
00223    EtaRaw1 = 0.;
00224    PhiRaw1 = 0.;
00225    EtGen1  = 0.;
00226    EtRaw1  = 0.;
00227    EtMCJ1  = 0.;
00228    EtZSP1  = 0.;
00229    EtJPT1  = 0.;
00230    DRMAXgjet1 = 1000.;
00231 
00232    EtaGen2 = 0.;
00233    PhiGen2 = 0.;
00234    EtaRaw2 = 0.;
00235    PhiRaw2 = 0.;
00236    EtGen2  = 0.;
00237    EtRaw2  = 0.;
00238    EtMCJ2  = 0.;
00239    EtZSP2  = 0.;
00240    EtJPT2  = 0.;
00241    DRMAXgjet2 = 1000.;
00242 
00243    //   edm::ESHandle<CaloGeometry> geometry;
00244    //   iSetup.get<IdealGeometryRecord>().get(geometry);
00245 
00246 
00247    // get MC info
00248    edm::Handle<HepMCProduct> EvtHandle ;
00249    iEvent.getByLabel( "source", EvtHandle ) ;
00250    //  iEvent.getByLabel( "VtxSmeared", EvtHandle ) ;
00251 
00252    // l1 and l2 are leptons from Z->ll to be checked they are not gen jets (DR match)
00253    HepLorentzVector l1(0.,0.,1.,1.);
00254    HepLorentzVector l2(0.,0.,1.,1.);
00255  
00256   //
00257      
00258    // get collection of towers
00259    /*
00260    edm::Handle<CandidateCollection> calotowers;
00261    iEvent.getByLabel(calotowersSrc, calotowers);   
00262    const CandidateCollection* inputCol = calotowers.product();
00263    CandidateCollection::const_iterator candidate;
00264    for( candidate = inputCol->begin(); candidate != inputCol->end(); ++candidate )
00265      {
00266        double phi   = candidate->phi();
00267        double theta = candidate->theta();
00268        double eta   = candidate->eta();
00269        double e     = candidate->energy();
00270        double et    = e*sin(theta);
00271        cout <<" towers: phi = " << phi
00272             <<" eta = " << eta
00273             <<" et = " << et << endl;
00274     }
00275    */
00276 
00277    /*
00278    vector<HepLorentzVector> cjetsRaw;
00279    cjetsRaw.clear();
00280 
00281    vector<HepLorentzVector> cjetsMCJ;
00282    cjetsMCJ.clear();
00283 
00284    vector<HepLorentzVector> cjetsZSP;
00285    cjetsZSP.clear();
00286 
00287    vector<HepLorentzVector> cjetsJPT;
00288    cjetsJPT.clear();
00289    */
00290 
00291    // get gen jets collection
00292    Handle<GenJetCollection> genjets;
00293    iEvent.getByLabel(genjetsSrc, genjets);
00294    int jg = 0;
00295    for(GenJetCollection::const_iterator gjet = genjets->begin(); 
00296        gjet != genjets->end(); ++gjet ) {
00297      if(gjet->pt() >= 20.) {
00298        HepLorentzVector jet(gjet->px(), gjet->py(), gjet->pz(), gjet->energy());
00299        double drjl1 = l1.deltaR(jet);
00300        double drjl2 = l2.deltaR(jet);
00301        /*
00302        cout <<" Gen Jet " << jg
00303             <<" pt = " << gjet->pt()
00304             <<" px = " << gjet->px()
00305             <<" py = " << gjet->py()
00306             <<" pz = " << gjet->pz()
00307             <<" energy = " << gjet->energy()
00308             <<" j eta = " << gjet->eta()
00309             <<" j phi = " << gjet->phi() 
00310             <<" l1 eta = " << l1.eta() 
00311             <<" l1 phi = " << l1.phi() 
00312             <<" l2 eta = " << l2.eta() 
00313             <<" l2 phi = " << l2.phi() 
00314             <<" dr1 = " << drjl1 
00315             <<" dr2 = " << drjl2 << endl;
00316        */
00317        if(drjl1 > 1.0 && drjl2 > 1.0) 
00318          {
00319            jg++;
00320            if(jg <= 2) {
00321              gjets.push_back(jet);
00322            }
00323          }
00324      }
00325    }
00326 
00327    //   cout <<" ==> NUMBER OF GOOD GEN JETS " << gjets.size() << endl;
00328 
00329    if(gjets.size() > 0) {
00330    // get calo jet collection
00331      edm::Handle<CaloJetCollection> calojets;
00332      iEvent.getByLabel(calojetsSrc, calojets);
00333    // get calo jet after zsp collection
00334      edm::Handle<CaloJetCollection> zspjets;
00335      iEvent.getByLabel(zspjetsSrc, zspjets);
00336      /*
00337      cout << "====> number of calo jets "<< calojets->size() 
00338           << " number if zsp jets = " << zspjets->size() << endl;
00339      */
00340      if(calojets->size() > 0) {
00341 
00342        // MC jet energy corrections
00343        //       const JetCorrector* correctorMCJ = JetCorrector::getJetCorrector (JetCorrectionMCJ, iSetup);
00344        // ZSP jet energy corrections
00345        //       const JetCorrector* correctorZSP = JetCorrector::getJetCorrector (JetCorrectionZSP, iSetup);
00346        // Jet+tracks energy corrections
00347 
00348        const JetCorrector* correctorJPT = JetCorrector::getJetCorrector (JetCorrectionJPT, iSetup);
00349        
00350        // loop over jets and do matching with gen jets
00351        int jc = 0;
00352        
00353        for( CaloJetCollection::const_iterator cjet = calojets->begin(); 
00354             cjet != calojets->end(); ++cjet ){ 
00355          //
00356          HepLorentzVector cjetc(cjet->px(), cjet->py(), cjet->pz(), cjet->energy());
00357          /*
00358          cout <<" ==> calo jet Et = " << cjet->pt()
00359               <<" eta = " << cjet->eta()
00360               <<" phi = " << cjet->phi() << endl;
00361          */
00362          CaloJetCollection::const_iterator zspjet;
00363          for( zspjet = zspjets->begin(); 
00364               zspjet != zspjets->end(); ++zspjet ){ 
00365            HepLorentzVector zspjetc(zspjet->px(), zspjet->py(), zspjet->pz(), zspjet->energy());
00366            double dr = zspjetc.deltaR(cjetc);
00367            /*
00368            cout <<"      zspjet Et = " << zspjet->pt()
00369                 <<" eta = " << zspjet->eta()
00370                 <<" phi = " << zspjet->phi()
00371                 <<" dr = " << dr << endl;
00372            */
00373            if(dr < 0.001) break;
00374          }
00375          /*
00376          cout <<" --> matched zsp jet found Et = " << zspjet->pt()
00377               <<" eta = " << zspjet->eta()
00378               <<" phi = " << zspjet->phi() << endl;
00379          */
00380          // JPT corrections
00381          double scaleJPT = correctorJPT->correction ((*zspjet),iEvent,iSetup);
00382          Jet::LorentzVector jetscaleJPT(zspjet->px()*scaleJPT, zspjet->py()*scaleJPT,
00383                                         zspjet->pz()*scaleJPT, zspjet->energy()*scaleJPT);
00384          /*
00385          cout <<" ..> corrected jpt jet Et = " << jetscaleJPT.pt()
00386               <<" eta = " << jetscaleJPT.eta()
00387               <<" phi = " << jetscaleJPT.phi() << endl;
00388          */
00389          CaloJet cjetJPT(jetscaleJPT, cjet->getSpecific(), cjet->getJetConstituents());
00390 
00391          double DRgjet1 = gjets[0].deltaR(cjetc);
00392 
00393          if(DRgjet1 < DRMAXgjet1) {
00394            DRMAXgjet1 = DRgjet1;
00395  
00396            EtaGen1 = gjets[0].eta();
00397            PhiGen1 = gjets[0].phi();
00398            EtGen1  = gjets[0].perp();
00399 
00400            EtaRaw1 = cjet->eta(); 
00401            PhiRaw1 = cjet->phi();
00402            EtRaw1  = cjet->pt();
00403            //      EtMCJ1  = cjetMCJ.pt(); 
00404            EtZSP1  = zspjet->pt(); 
00405            EtJPT1  = cjetJPT.pt(); 
00406          }
00407          if(gjets.size() == 2) {
00408            double DRgjet2 = gjets[1].deltaR(cjetc);
00409            if(DRgjet2 < DRMAXgjet2) { 
00410              DRMAXgjet2 = DRgjet2;
00411 
00412              EtaGen2 = gjets[1].eta();
00413              PhiGen2 = gjets[1].phi();
00414              EtGen2  = gjets[1].perp();
00415 
00416              EtaRaw2 = cjet->eta(); 
00417              PhiRaw2 = cjet->phi();
00418              EtRaw2  = cjet->pt();
00419              EtZSP2  = zspjet->pt(); 
00420              EtJPT2  = cjetJPT.pt(); 
00421            }
00422          }
00423          jc++;
00424        }
00425        /*
00426         cout <<" best match to 1st gen get = " << DRMAXgjet1
00427             <<" raw jet pt = " << EtRaw1 <<" eta = " << EtaRaw1 <<" phi " << PhiRaw1 
00428             <<" mcj pt = " << EtMCJ1 << " zsp pt = " << EtZSP1 <<" jpt = " << EtJPT1 << endl; 
00429        if(gjets.size() == 2) {
00430          cout <<" best match to 2st gen get = " << DRMAXgjet2
00431               <<" raw jet pt = " << EtRaw2 <<" eta = " << EtaRaw2 <<" phi " << PhiRaw2 
00432               <<" mcj pt = " << EtMCJ2 << " zsp pt = " << EtZSP2 <<" jpt = " << EtJPT2 << endl; 
00433        }
00434        */
00435      }
00436    }
00437    // fill tree
00438    t1->Fill();
00439 }

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

Reimplemented from edm::EDAnalyzer.

Definition at line 124 of file JPTAnalyzer.cc.

References DRMAXgjet1, DRMAXgjet2, EtaGen1, EtaGen2, EtaRaw1, EtaRaw2, EtGen1, EtGen2, EtJPT1, EtJPT2, EtMCJ1, EtMCJ2, EtRaw1, EtRaw2, EtZSP1, EtZSP2, fOutputFileName, hOutputFile, PhiGen1, PhiGen2, PhiRaw1, PhiRaw2, and t1.

00125 {
00126   using namespace edm;
00127   // creating a simple tree
00128 
00129   hOutputFile   = new TFile( fOutputFileName.c_str(), "RECREATE" ) ;
00130 
00131   t1 = new TTree("t1","analysis tree");
00132 
00133   t1->Branch("EtaGen1",&EtaGen1,"EtaGen1/D");
00134   t1->Branch("PhiGen1",&PhiGen1,"PhiGen1/D");
00135   t1->Branch("EtaRaw1",&EtaRaw1,"EtaRaw1/D");
00136   t1->Branch("PhiRaw1",&PhiRaw1,"PhiRaw1/D");
00137   t1->Branch("EtGen1",&EtGen1,"EtGen1/D");
00138   t1->Branch("EtRaw1",&EtRaw1,"EtRaw1/D");
00139   t1->Branch("EtMCJ1",&EtMCJ1,"EtMCJ1/D");
00140   t1->Branch("EtZSP1",&EtZSP1,"EtZSP1/D");
00141   t1->Branch("EtJPT1",&EtJPT1,"EtJPT1/D");
00142   t1->Branch("DRMAXgjet1",&DRMAXgjet1,"DRMAXgjet1/D");
00143 
00144   t1->Branch("EtaGen2",&EtaGen2,"EtaGen2/D");
00145   t1->Branch("PhiGen2",&PhiGen2,"PhiGen2/D");
00146   t1->Branch("EtaRaw2",&EtaRaw2,"EtaRaw2/D");
00147   t1->Branch("PhiRaw2",&PhiRaw2,"PhiRaw2/D");
00148   t1->Branch("EtGen2",&EtGen2,"EtGen2/D");
00149   t1->Branch("EtRaw2",&EtRaw2,"EtRaw2/D");
00150   t1->Branch("EtMCJ2",&EtMCJ2,"EtMCJ2/D");
00151   t1->Branch("EtZSP2",&EtZSP2,"EtZSP2/D");
00152   t1->Branch("EtJPT2",&EtJPT2,"EtJPT2/D");
00153   t1->Branch("DRMAXgjet2",&DRMAXgjet2,"DRMAXgjet2/D");
00154 
00155   return ;
00156 }

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

Reimplemented from edm::EDAnalyzer.

Definition at line 160 of file JPTAnalyzer.cc.

References hOutputFile.

00160                     {
00161 
00162   hOutputFile->Write() ;
00163   hOutputFile->Close() ;
00164   
00165   return ;
00166 }


Member Data Documentation

string JPTAnalyzer::calojetsSrc [private]

Definition at line 94 of file JPTAnalyzer.cc.

Referenced by analyze(), and JPTAnalyzer().

double JPTAnalyzer::DRMAXgjet1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::DRMAXgjet2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtaGen1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtaGen2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtaRaw1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtaRaw2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtGen1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtGen2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtJPT1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtJPT2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtMCJ1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtMCJ2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtRaw1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtRaw2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtZSP1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::EtZSP2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

string JPTAnalyzer::fOutputFileName [private]

Definition at line 91 of file JPTAnalyzer.cc.

Referenced by beginJob(), and JPTAnalyzer().

string JPTAnalyzer::genjetsSrc [private]

Definition at line 98 of file JPTAnalyzer.cc.

Referenced by analyze(), and JPTAnalyzer().

TFile* JPTAnalyzer::hOutputFile [private]

Definition at line 110 of file JPTAnalyzer.cc.

Referenced by beginJob(), and endJob().

string JPTAnalyzer::JetCorrectionJPT [private]

Definition at line 105 of file JPTAnalyzer.cc.

Referenced by analyze(), and JPTAnalyzer().

double JPTAnalyzer::PhiGen1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::PhiGen2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::PhiRaw1 [private]

Definition at line 107 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

double JPTAnalyzer::PhiRaw2 [private]

Definition at line 108 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

TTree* JPTAnalyzer::t1 [private]

Definition at line 111 of file JPTAnalyzer.cc.

Referenced by analyze(), and beginJob().

string JPTAnalyzer::zspjetsSrc [private]

Definition at line 96 of file JPTAnalyzer.cc.

Referenced by analyze(), and JPTAnalyzer().


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