CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoEcal/EgammaClusterProducers/src/PreshowerAnalyzer.cc

Go to the documentation of this file.
00001 
00008 //
00009 // Original Author:  Shahram Rahatlou
00010 //         Created:  10 May 200
00011 // $Id: PreshowerAnalyzer.cc,v 1.6 2009/12/18 20:45:01 wmtan Exp $
00012 //
00013 
00014 #include "RecoEcal/EgammaClusterProducers/interface/PreshowerAnalyzer.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Utilities/interface/Exception.h"
00018 #include "DataFormats/Common/interface/Handle.h"
00019 
00020 #include "TFile.h"
00021 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00022 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00023 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00024 
00025 #include "DataFormats/EgammaReco/interface/PreshowerCluster.h"
00026 #include "DataFormats/EgammaReco/interface/PreshowerClusterFwd.h"
00027 #include "RecoEcal/EgammaClusterProducers/interface/PreshowerClusterProducer.h"
00028 
00029 //========================================================================
00030 PreshowerAnalyzer::PreshowerAnalyzer( const edm::ParameterSet& ps )
00031 //========================================================================
00032 {
00033  
00034   EminDE_ = ps.getParameter<double>("EminDE");
00035   EmaxDE_ = ps.getParameter<double>("EmaxDE");
00036   nBinDE_ = ps.getParameter<int>("nBinDE");
00037 
00038   EminSC_ = ps.getParameter<double>("EminSC");
00039   EmaxSC_ = ps.getParameter<double>("EmaxSC");
00040   nBinSC_ = ps.getParameter<int>("nBinSC");
00041 
00042   preshClusterCollectionX_ = ps.getParameter<std::string>("preshClusterCollectionX");
00043   preshClusterCollectionY_ = ps.getParameter<std::string>("preshClusterCollectionY");
00044   preshClusterProducer_   = ps.getParameter<std::string>("preshClusterProducer");
00045 
00046   islandEndcapSuperClusterCollection1_ = ps.getParameter<std::string>("islandEndcapSuperClusterCollection1");
00047   islandEndcapSuperClusterProducer1_   = ps.getParameter<std::string>("islandEndcapSuperClusterProducer1");
00048 
00049   islandEndcapSuperClusterCollection2_ = ps.getParameter<std::string>("islandEndcapSuperClusterCollection2");
00050   islandEndcapSuperClusterProducer2_   = ps.getParameter<std::string>("islandEndcapSuperClusterProducer2");
00051 
00052   outputFile_   = ps.getParameter<std::string>("outputFile");
00053   rootFile_ = TFile::Open(outputFile_.c_str(),"RECREATE"); // open output file to store histograms
00054 
00055   // calibration parameters:
00056   calib_planeX_ = ps.getParameter<double>("preshCalibPlaneX");
00057   calib_planeY_ = ps.getParameter<double>("preshCalibPlaneY");
00058   gamma_        = ps.getParameter<double>("preshCalibGamma");
00059   mip_          = ps.getParameter<double>("preshCalibMIP");
00060 
00061   nEvt_ = 0; 
00062 }
00063 
00064 //========================================================================
00065 PreshowerAnalyzer::~PreshowerAnalyzer()
00066 //========================================================================
00067 {
00068   delete rootFile_;
00069 }
00070 
00071 //========================================================================
00072 void PreshowerAnalyzer::beginJob() {
00073 //========================================================================
00074 
00075   rootFile_->cd();
00076 
00077   h1_esE_x = new TH1F("esE_x"," ES cluster Energy in  X-plane",20, 0, 0.03);
00078   h1_esE_y = new TH1F("esE_y"," ES cluster Energy in  Y-plane",20, 0, 0.03);
00079   h1_esEta_x = new TH1F("esEta_x"," ES cluster Eta in X-plane",12, 1.5, 2.7);
00080   h1_esEta_y = new TH1F("esEta_y"," ES cluster Eta in Y-plane",12, 1.5, 2.7);
00081   h1_esPhi_x = new TH1F("esPhi_x"," ES cluster Phi in X-plane",20, 0, 6.28);
00082   h1_esPhi_y = new TH1F("esPhi_y"," ES cluster Phi in Y-plane",20, 0, 6.28);
00083   h1_esNhits_x = new TH1F("esNhits_x"," ES cluster Nhits in  X-plane",10, 0, 10);
00084   h1_esNhits_y = new TH1F("esNhits_y"," ES cluster Nhits in  Y-plane",10, 0, 10);
00085   h1_esDeltaE = new TH1F("esDeltaE"," DeltaE", nBinDE_, EminDE_, EmaxDE_); 
00086   h1_nclu_x = new TH1F("esNclu_x"," number of ES clusters (for one SC) in X-plane",20, 0, 80);
00087   h1_nclu_y = new TH1F("esNclu_y"," number of ES clusters (for one SC) in Y-plane",20, 0, 80);
00088 
00089   h1_islandEESCEnergy1 = new TH1F("islandEESCEnergy1","Energy of super clusters with island algo - endcap1",nBinSC_,EminSC_,EmaxSC_);
00090   h1_islandEESCEnergy2 = new TH1F("islandEESCEnergy2","Energy of super clusters with island algo - endcap2",nBinSC_,EminSC_,EmaxSC_);
00091 }
00092 
00093 
00094 //========================================================================
00095 void
00096 PreshowerAnalyzer::analyze( const edm::Event& evt, const edm::EventSetup& es ) {
00097 //========================================================================
00098 
00099   using namespace edm; // needed for all fwk related classe
00100 
00101   //std::cout << "\n .......  Event # " << nEvt_+1 << " is analyzing ....... " << std::endl << std::endl;
00102 
00103   // Get island super clusters
00104   Handle<reco::SuperClusterCollection> pIslandEndcapSuperClusters1;
00105   evt.getByLabel(islandEndcapSuperClusterProducer1_, islandEndcapSuperClusterCollection1_, pIslandEndcapSuperClusters1);
00106   const reco::SuperClusterCollection* islandEndcapSuperClusters1 = pIslandEndcapSuperClusters1.product();
00107   //std::cout << "\n islandEndcapSuperClusters1->size() = " << islandEndcapSuperClusters1->size() << std::endl;
00108 
00109   // loop over the super clusters and fill the histogram
00110   for(reco::SuperClusterCollection::const_iterator aClus = islandEndcapSuperClusters1->begin();
00111                                                     aClus != islandEndcapSuperClusters1->end(); aClus++) {
00112     h1_islandEESCEnergy1->Fill( aClus->energy() );
00113   }
00114 
00115   // Get island super clusters
00116   Handle<reco::SuperClusterCollection> pIslandEndcapSuperClusters2;
00117   evt.getByLabel(islandEndcapSuperClusterProducer2_, islandEndcapSuperClusterCollection2_, pIslandEndcapSuperClusters2);
00118   const reco::SuperClusterCollection* islandEndcapSuperClusters2 = pIslandEndcapSuperClusters2.product();
00119   //std::cout << "\n islandEndcapSuperClusters2->size() = " << islandEndcapSuperClusters2->size() << std::endl;
00120 
00121   // loop over the super clusters and fill the histogram
00122   for(reco::SuperClusterCollection::const_iterator aClus = islandEndcapSuperClusters2->begin();
00123                                                     aClus != islandEndcapSuperClusters2->end(); aClus++) {
00124     h1_islandEESCEnergy2->Fill( aClus->energy() );
00125   }
00126 
00127 
00128   // Get ES clusters in X plane
00129   Handle<reco::PreshowerClusterCollection> pPreshowerClustersX;
00130   evt.getByLabel(preshClusterProducer_, preshClusterCollectionX_, pPreshowerClustersX);
00131   const reco::PreshowerClusterCollection *clustersX = pPreshowerClustersX.product();
00132   h1_nclu_x->Fill( clustersX->size() );
00133   //std::cout << "\n pPreshowerClustersX->size() = " << clustersX->size() << std::endl;
00134 
00135   Handle<reco::PreshowerClusterCollection> pPreshowerClustersY;
00136   evt.getByLabel(preshClusterProducer_, preshClusterCollectionY_, pPreshowerClustersY);
00137   const reco::PreshowerClusterCollection *clustersY = pPreshowerClustersY.product();
00138   h1_nclu_y->Fill( clustersY->size() );
00139   //std::cout << "\n pPreshowerClustersY->size() = " << clustersY->size() << std::endl;
00140 
00141 
00142   // loop over the ES clusters and fill the histogram
00143   float e1 = 0;
00144   for(reco::PreshowerClusterCollection::const_iterator esClus = clustersX->begin();
00145                                                        esClus !=clustersX->end(); esClus++) {
00146       e1 += esClus->energy();
00147       h1_esE_x->Fill( esClus->energy() );  
00148       h1_esEta_x->Fill( esClus->eta() );
00149       h1_esPhi_x->Fill( esClus->phi() );
00150       h1_esNhits_x->Fill( esClus->nhits() );     
00151   }
00152 
00153   float e2 = 0;
00154   for(reco::PreshowerClusterCollection::const_iterator esClus = clustersY->begin();
00155                                                        esClus !=clustersY->end(); esClus++) {
00156       e2 += esClus->energy();
00157       h1_esE_y->Fill( esClus->energy() );  
00158       h1_esEta_y->Fill( esClus->eta() );
00159       h1_esPhi_y->Fill( esClus->phi() );
00160       h1_esNhits_y->Fill( esClus->nhits() );     
00161   }
00162   
00163   float deltaE = 0;       
00164   if(e1+e2 > 1.0e-10) {
00165       // GeV to #MIPs
00166       e1 = e1 / mip_;
00167       e2 = e2 / mip_;
00168       deltaE = gamma_*(calib_planeX_*e1+calib_planeY_*e2);       
00169    }
00170 
00171   h1_esDeltaE->Fill(deltaE);
00172 
00173    nEvt_++;
00174 
00175 }
00176 
00177 //========================================================================
00178 void PreshowerAnalyzer::endJob() {
00179 //========================================================================
00180 
00181    rootFile_->cd();
00182 
00183    h1_esE_x->Write();     
00184    h1_esE_y->Write();
00185    h1_esEta_x->Write();
00186    h1_esEta_y->Write();
00187    h1_esPhi_x->Write();
00188    h1_esPhi_y->Write();
00189    h1_esNhits_x->Write();
00190    h1_esNhits_y->Write();          
00191    h1_esDeltaE->Write();
00192    h1_nclu_x->Write();
00193    h1_nclu_y->Write();
00194 
00195    h1_islandEESCEnergy1->Write();
00196    h1_islandEESCEnergy2->Write();
00197 
00198    rootFile_->Close();
00199 }