CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/SLHCUpgradeSimulations/L1CaloTrigger/plugins/CaloTriggerAnalyzer.cc

Go to the documentation of this file.
00001 
00002 #include "SLHCUpgradeSimulations/L1CaloTrigger/plugins/CaloTriggerAnalyzer.h"
00003 #include "Math/GenVector/VectorUtil.h"
00004 #include <iostream>
00005 #include <iomanip>
00006 
00007 
00008 CaloTriggerAnalyzer::CaloTriggerAnalyzer(const edm::ParameterSet& iConfig):
00009   src_(iConfig.getParameter<edm::InputTag>("src")),
00010   ref_(iConfig.getParameter<edm::InputTag>("ref")),
00011   DR_(iConfig.getParameter<double>("deltaR")),
00012   threshold_(iConfig.getParameter<double>("threshold")),
00013   maxEta_(iConfig.getUntrackedParameter<double>("maxEta",2.5))
00014 {
00015    //now do what ever initialization is needed
00016 
00017   edm::Service<TFileService> fs;
00018 
00019   ptNum    = fs->make<TH1F>( "ptNum"   , "ptNum"   , 20  ,  0., 100. );
00020   ptDenom  = fs->make<TH1F>( "ptDenom" , "ptDenom" , 20  ,  0., 100. );
00021   etaNum   = fs->make<TH1F>( "etaNum"  , "etaNum"  , 20  ,  -2.5, 2.5 );
00022   etaDenom = fs->make<TH1F>( "etaDenom", "etaDenom", 20  ,  -2.5, 2.5 );
00023   pt       = fs->make<TH1F>( "pt"      , "pt", 20  ,  0. , 100. );
00024   highestPt= fs->make<TH1F>( "highestPt"      , "highestPt", 50  ,  0. , 100. );
00025   secondPt = fs->make<TH1F>( "secondHighestPt", "secondHighestPt", 50  ,  0. , 100. );
00026   highestPtGen=fs->make<TH1F>("highestPtGen","highestPtGen",100, 0., 200.);
00027   secondPtGen=fs->make<TH1F>("secondPtGen","secondPtGen",100, 0., 200.);
00028   dPt      = fs->make<TH1F>( "dPt"      , "dPt", 50  , -1  , 1 );
00029   dEta     = fs->make<TH1F>( "dEta"      , "dEta", 50  , -0.5  , 0.5 );
00030   dPhi     = fs->make<TH1F>( "dPhi"      , "dPhi", 50  , -0.5  , 0.5 );
00031   RPt = fs->make<TH1F>("RPt", "Pt_Ratio", 10, 0, 2.);
00032   RPtEta = fs->make<TProfile>("RPtEta","Pt_Ratio as fcn of abs(eta)",13,0.0,2.6,0,2);
00033   RPtEtaFull = fs->make<TProfile>("RPtEtaFull","Pt_Ratio as fcn of eta",26,-2.6,2.6,0,2);
00034   absEta = fs->make<TH1F>("abseta","abs(eta)",13,0.,2.6); 
00035   dR = fs->make<TH1F>("dR","delta(R)",50,0,2);
00036 }
00037 
00038 
00039 CaloTriggerAnalyzer::~CaloTriggerAnalyzer()
00040 {
00041 
00042    // do anything here that needs to be done at desctruction time
00043    // (e.g. close files, deallocate resources etc.)
00044 
00045 }
00046 
00047 
00048 //
00049 // member functions
00050 //
00051 
00052 // ------------ method called to for each event  ------------
00053 void
00054 CaloTriggerAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00055 {
00056   //  std::cout << src_ << std::endl;
00057   edm::Handle<edm::View<reco::Candidate> > ref;
00058   edm::Handle<edm::View<reco::Candidate> > src;
00059 
00060   bool gotRef = iEvent.getByLabel(ref_,ref);
00061   bool gotSrc = iEvent.getByLabel(src_,src);
00062 
00063   if(gotSrc) {
00064     highPt=0;
00065     secondPtf=0;
00066     for(edm::View<reco::Candidate>::const_iterator i = src->begin(); i!= src->end();++i)
00067       {
00068         pt->Fill(i->pt());
00069         if (i->pt()>highPt){
00070           secondPtf=highPt;
00071           highPt=i->pt();
00072         } else if (i->pt()>secondPtf){
00073           secondPtf=i->pt();
00074         }
00075       }
00076 
00077     if(src->size()>0)
00078       highestPt->Fill(highPt);
00079     else
00080       highestPt->Fill(0.0);
00081 
00082 
00083     if(src->size()>1)
00084       secondPt->Fill(secondPtf);
00085     else
00086       secondPt->Fill(0.0);
00087   }
00088 
00089   if(gotRef) {
00090     //get highest Pt gen object--loop over all to make sure it's the highest!
00091     highestGenPt=-1.0;
00092     secondGenPt=-1.0;
00093     for(edm::View<reco::Candidate>::const_iterator i = ref->begin(); i!= ref->end();++i){
00094 
00095       if (i->pt()>highestGenPt){
00096         highestGenPt=i->pt();
00097         highestPtGen->Fill(i->pt());
00098       } else if (i->pt()>secondGenPt){
00099         secondGenPt=i->pt();
00100         secondPtGen->Fill(i->pt());
00101       }
00102       
00103       if(fabs(i->eta())<maxEta_&&i->pt()>threshold_/2.)
00104         {
00105           ptDenom->Fill(i->pt());
00106           etaDenom->Fill(i->eta());
00107           
00108           //      printf("ref pt = %f  eta = %f phi = %f\n",i->pt(),i->eta(),i->phi());
00109           
00110           if(gotSrc)
00111             {
00112               bool matched=false;
00113               math::XYZTLorentzVector highestV(0.0001,0.,0.,0.0002);
00114               for(edm::View<reco::Candidate>::const_iterator j = src->begin(); j!= src->end();++j)
00115                 {
00116                   dR->Fill(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4()));
00117                   if(j->pt()>threshold_)
00118                     {
00119                       if(ROOT::Math::VectorUtil::DeltaR(i->p4(),j->p4())<DR_) {
00120                         //      printf("matched pt = %f  eta = %f phi = %f\n",j->pt(),j->eta(),j->phi());
00121                         
00122                         if(j->pt()>highestV.pt())
00123                           highestV = j->p4();
00124                         
00125                         matched=true;
00126                       }
00127                       
00128                     }
00129                 }
00130             if(matched)
00131               {
00132                 RPt->Fill(i->pt()/highestV.pt());
00133                 RPtEtaFull->Fill(highestV.eta(), i->pt()/highestV.pt());
00134                 RPtEta->Fill(fabs(highestV.eta()), i->pt()/highestV.pt());
00135 
00136                 //              printf("matched abs(eta) = %f \n", fabs(highestV.eta()) );
00137                 absEta->Fill(fabs(highestV.eta()));
00138                 dPt->Fill((highestV.pt()-i->pt())/i->pt());
00139                 dEta->Fill(highestV.eta()-i->eta());
00140                 dPhi->Fill(highestV.phi()-i->phi());
00141                 
00142                 ptNum->Fill(i->pt());
00143                 etaNum->Fill(i->eta());
00144               }
00145             }
00146         }
00147     }
00148   }
00149 }
00150 
00151 DEFINE_FWK_MODULE(CaloTriggerAnalyzer);