CMS 3D CMS Logo

EgammaBasicClusters.cc

Go to the documentation of this file.
00001 #include "Validation/EcalClusters/interface/EgammaBasicClusters.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "FWCore/Utilities/interface/Exception.h"
00005 #include "FWCore/ServiceRegistry/interface/Service.h"
00006 
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00009 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00010 #include "DQMServices/Core/interface/DQMStore.h"
00011 
00012 EgammaBasicClusters::EgammaBasicClusters( const edm::ParameterSet& ps )
00013 {
00014         outputFile_ = ps.getUntrackedParameter<std::string>("outputFile", "");
00015         CMSSW_Version_ = ps.getUntrackedParameter<std::string>("CMSSW_Version", "");
00016 
00017         verboseDBE_ = ps.getUntrackedParameter<bool>("verboseDBE", false);
00018 
00019         hist_min_Size_ = ps.getParameter<double>("hist_min_Size");
00020         hist_max_Size_ = ps.getParameter<double>("hist_max_Size");
00021         hist_bins_Size_ = ps.getParameter<int>   ("hist_bins_Size");
00022 
00023         hist_min_NumRecHits_ = ps.getParameter<double>("hist_min_NumRecHits");
00024         hist_max_NumRecHits_ = ps.getParameter<double>("hist_max_NumRecHits");
00025         hist_bins_NumRecHits_ = ps.getParameter<int>   ("hist_bins_NumRecHits");
00026 
00027         hist_min_ET_ = ps.getParameter<double>("hist_min_ET");
00028         hist_max_ET_ = ps.getParameter<double>("hist_max_ET");
00029         hist_bins_ET_ = ps.getParameter<int>   ("hist_bins_ET");
00030 
00031         hist_min_Eta_ = ps.getParameter<double>("hist_min_Eta");
00032         hist_max_Eta_ = ps.getParameter<double>("hist_max_Eta");
00033         hist_bins_Eta_ = ps.getParameter<int>   ("hist_bins_Eta");
00034 
00035         hist_min_Phi_ = ps.getParameter<double>("hist_min_Phi");
00036         hist_max_Phi_ = ps.getParameter<double>("hist_max_Phi");
00037         hist_bins_Phi_ = ps.getParameter<int>   ("hist_bins_Phi");
00038 
00039         barrelBasicClusterCollection_ = ps.getParameter<edm::InputTag>("barrelBasicClusterCollection");
00040         endcapBasicClusterCollection_ = ps.getParameter<edm::InputTag>("endcapBasicClusterCollection");
00041 }
00042 
00043 EgammaBasicClusters::~EgammaBasicClusters() {}
00044 
00045 void EgammaBasicClusters::beginJob(edm::EventSetup const&) 
00046 {
00047         dbe_ = edm::Service<DQMStore>().operator->();                   
00048 
00049         if ( verboseDBE_ )
00050         {
00051                 dbe_->setVerbose(1);
00052                 dbe_->showDirStructure();
00053         }
00054         else 
00055                 dbe_->setVerbose(0);
00056 
00057         dbe_->setCurrentFolder("EcalClustersV/CMSSW_"+CMSSW_Version_+"/EcalClusters/BasicClusters/");
00058 
00059         hist_EB_BC_Size_ 
00060                 = dbe_->book1D("hist_EB_BC_Size_","# Basic Clusters in Barrel",
00061                         hist_bins_Size_,hist_min_Size_,hist_max_Size_);
00062         hist_EE_BC_Size_ 
00063                 = dbe_->book1D("hist_EE_BC_Size_","# Basic Clusters in Endcap",
00064                         hist_bins_Size_,hist_min_Size_,hist_max_Size_);
00065 
00066         hist_EB_BC_NumRecHits_ 
00067                 = dbe_->book1D("hist_EB_BC_NumRecHits_","# of RecHits in Basic Clusters in Barrel",
00068                         hist_bins_NumRecHits_,hist_min_NumRecHits_,hist_max_NumRecHits_);
00069         hist_EE_BC_NumRecHits_ 
00070                 = dbe_->book1D("hist_EE_BC_NumRecHits_","# of RecHits in Basic Clusters in Endcap",
00071                         hist_bins_NumRecHits_,hist_min_NumRecHits_,hist_max_NumRecHits_);
00072 
00073         hist_EB_BC_ET_ 
00074                 = dbe_->book1D("hist_EB_BC_ET_","ET of Basic Clusters in Barrel",
00075                         hist_bins_ET_,hist_min_ET_,hist_max_ET_);
00076         hist_EE_BC_ET_ 
00077                 = dbe_->book1D("hist_EE_BC_ET_","ET of Basic Clusters in Endcap",
00078                         hist_bins_ET_,hist_min_ET_,hist_max_ET_);
00079 
00080         hist_EB_BC_Eta_ 
00081                 = dbe_->book1D("hist_EB_BC_Eta_","Eta of Basic Clusters in Barrel",
00082                         hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_);
00083         hist_EE_BC_Eta_ 
00084                 = dbe_->book1D("hist_EE_BC_Eta_","Eta of Basic Clusters in Endcap",
00085                         hist_bins_Eta_,hist_min_Eta_,hist_max_Eta_);
00086 
00087         hist_EB_BC_Phi_
00088                 = dbe_->book1D("hist_EB_BC_Phi_","Phi of Basic Clusters in Barrel",
00089                         hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_);
00090         hist_EE_BC_Phi_ 
00091                 = dbe_->book1D("hist_EE_BC_Phi_","Phi of Basic Clusters in Endcap",
00092                         hist_bins_Phi_,hist_min_Phi_,hist_max_Phi_);
00093 }
00094 
00095 void EgammaBasicClusters::analyze( const edm::Event& evt, const edm::EventSetup& es )
00096 {
00097         edm::Handle<reco::BasicClusterCollection> pBarrelBasicClusters;
00098         evt.getByLabel(barrelBasicClusterCollection_, pBarrelBasicClusters);
00099         if (!pBarrelBasicClusters.isValid()) {
00100           edm::LogError("EgammaBasicClusters") << "Error! can't get collection with label " 
00101                                                << barrelBasicClusterCollection_.label();
00102         }
00103 
00104         const reco::BasicClusterCollection* barrelBasicClusters = pBarrelBasicClusters.product();
00105         hist_EB_BC_Size_->Fill(barrelBasicClusters->size());
00106 
00107         for(reco::BasicClusterCollection::const_iterator aClus = barrelBasicClusters->begin(); 
00108                 aClus != barrelBasicClusters->end(); aClus++)
00109         {
00110                 hist_EB_BC_NumRecHits_->Fill(aClus->getHitsByDetId().size());
00111                 hist_EB_BC_ET_->Fill(aClus->energy()*aClus->position().theta());
00112                 hist_EB_BC_Eta_->Fill(aClus->position().eta());
00113                 hist_EB_BC_Phi_->Fill(aClus->position().phi());
00114         }
00115 
00116         edm::Handle<reco::BasicClusterCollection> pEndcapBasicClusters;
00117 
00118         evt.getByLabel(endcapBasicClusterCollection_, pEndcapBasicClusters);
00119         if (!pEndcapBasicClusters.isValid()) {
00120           edm::LogError("EgammaBasicClusters") << "Error! can't get collection with label " 
00121                                                << endcapBasicClusterCollection_.label();
00122         }
00123 
00124         const reco::BasicClusterCollection* endcapBasicClusters = pEndcapBasicClusters.product();
00125         hist_EE_BC_Size_->Fill(endcapBasicClusters->size());
00126 
00127         for(reco::BasicClusterCollection::const_iterator aClus = endcapBasicClusters->begin(); 
00128                 aClus != endcapBasicClusters->end(); aClus++)
00129         {
00130                 hist_EE_BC_NumRecHits_->Fill(aClus->getHitsByDetId().size());
00131                 hist_EE_BC_ET_->Fill(aClus->energy()*aClus->position().theta());
00132                 hist_EE_BC_Eta_->Fill(aClus->position().eta());
00133                 hist_EE_BC_Phi_->Fill(aClus->position().phi());
00134         }
00135 }
00136 
00137 void EgammaBasicClusters::endJob()
00138 {
00139         if (outputFile_.size() != 0 && dbe_) dbe_->save(outputFile_);
00140 }

Generated on Tue Jun 9 17:49:05 2009 for CMSSW by  doxygen 1.5.4