CMS 3D CMS Logo

Multi5x5SuperClusterProducer.cc

Go to the documentation of this file.
00001 // C/C++ headers
00002 #include <iostream>
00003 #include <vector>
00004 #include <memory>
00005 #include <sstream>
00006 
00007 // Framework
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/Framework/interface/EventSetup.h"
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "FWCore/Utilities/interface/Exception.h"
00013 
00014 // Reconstruction Classes
00015 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00016 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00017 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00018 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00019 
00020 // Class header file
00021 #include "RecoEcal/EgammaClusterProducers/interface/Multi5x5SuperClusterProducer.h"
00022 
00023 
00024 Multi5x5SuperClusterProducer::Multi5x5SuperClusterProducer(const edm::ParameterSet& ps)
00025 {
00026   // The verbosity level
00027   std::string verbosityString = ps.getParameter<std::string>("VerbosityLevel");
00028   if      (verbosityString == "DEBUG")   verbosity = Multi5x5BremRecoveryClusterAlgo::pDEBUG;
00029   else if (verbosityString == "WARNING") verbosity = Multi5x5BremRecoveryClusterAlgo::pWARNING;
00030   else if (verbosityString == "INFO")    verbosity = Multi5x5BremRecoveryClusterAlgo::pINFO;
00031   else                                   verbosity = Multi5x5BremRecoveryClusterAlgo::pERROR;
00032 
00033   endcapClusterProducer_ = ps.getParameter<std::string>("endcapClusterProducer");
00034   barrelClusterProducer_ = ps.getParameter<std::string>("barrelClusterProducer");
00035 
00036   endcapClusterCollection_ = ps.getParameter<std::string>("endcapClusterCollection");
00037   barrelClusterCollection_ = ps.getParameter<std::string>("barrelClusterCollection");
00038 
00039   endcapSuperclusterCollection_ = ps.getParameter<std::string>("endcapSuperclusterCollection");
00040   barrelSuperclusterCollection_ = ps.getParameter<std::string>("barrelSuperclusterCollection");
00041 
00042   doBarrel_ = ps.getParameter<bool>("doBarrel");
00043   doEndcaps_ = ps.getParameter<bool>("doEndcaps");
00044 
00045   barrelEtaSearchRoad_ = ps.getParameter<double>("barrelEtaSearchRoad");
00046   barrelPhiSearchRoad_ = ps.getParameter<double>("barrelPhiSearchRoad");
00047   endcapEtaSearchRoad_ = ps.getParameter<double>("endcapEtaSearchRoad");
00048   endcapPhiSearchRoad_ = ps.getParameter<double>("endcapPhiSearchRoad");
00049   seedTransverseEnergyThreshold_ = ps.getParameter<double>("seedTransverseEnergyThreshold");
00050 
00051   const edm::ParameterSet bremRecoveryPset = ps.getParameter<edm::ParameterSet>("bremRecoveryPset");
00052   bool dynamicPhiRoad = ps.getParameter<bool>("dynamicPhiRoad");
00053 
00054   bremAlgo_p = new Multi5x5BremRecoveryClusterAlgo(bremRecoveryPset, barrelEtaSearchRoad_, barrelPhiSearchRoad_, 
00055                                          endcapEtaSearchRoad_, endcapPhiSearchRoad_, 
00056                                          dynamicPhiRoad, seedTransverseEnergyThreshold_, verbosity);
00057 
00058   produces< reco::SuperClusterCollection >(endcapSuperclusterCollection_);
00059   produces< reco::SuperClusterCollection >(barrelSuperclusterCollection_);
00060 
00061   totalE = 0;
00062   noSuperClusters = 0;
00063   nEvt_ = 0;
00064 }
00065 
00066 
00067 Multi5x5SuperClusterProducer::~Multi5x5SuperClusterProducer()
00068 {
00069   delete bremAlgo_p;
00070 }
00071 
00072 void
00073 Multi5x5SuperClusterProducer::endJob() {
00074   double averEnergy = 0.;
00075   std::ostringstream str;
00076   str << "Multi5x5SuperClusterProducer::endJob()\n"
00077       << "  total # reconstructed super clusters: " << noSuperClusters << "\n"
00078       << "  total energy of all clusters: " << totalE << "\n";
00079   if(noSuperClusters>0) { 
00080     averEnergy = totalE / noSuperClusters;
00081     str << "  average SuperCluster energy = " << averEnergy << "\n";
00082   }
00083   edm::LogInfo("Multi5x5SuperClusterProducerInfo") << str.str() << "\n";
00084 
00085 }
00086 
00087 
00088 void Multi5x5SuperClusterProducer::produce(edm::Event& evt, const edm::EventSetup& es)
00089 {
00090   if(doEndcaps_)
00091     produceSuperclustersForECALPart(evt, endcapClusterProducer_, endcapClusterCollection_, endcapSuperclusterCollection_);
00092 
00093   if(doBarrel_)
00094     produceSuperclustersForECALPart(evt, barrelClusterProducer_, barrelClusterCollection_, barrelSuperclusterCollection_);
00095 
00096   nEvt_++;
00097 }
00098 
00099 
00100 void Multi5x5SuperClusterProducer::produceSuperclustersForECALPart(edm::Event& evt, 
00101                                                            std::string clusterProducer, 
00102                                                            std::string clusterCollection,
00103                                                            std::string superclusterCollection)
00104 {
00105   // get the cluster collection out and turn it to a BasicClusterRefVector:
00106   reco::BasicClusterRefVector *clusterRefVector_p = new reco::BasicClusterRefVector;
00107   getClusterRefVector(evt, clusterProducer, clusterCollection, clusterRefVector_p);
00108 
00109   // run the brem recovery and get the SC collection
00110   std::auto_ptr<reco::SuperClusterCollection> 
00111     superclusters_ap(new reco::SuperClusterCollection(bremAlgo_p->makeSuperClusters(*clusterRefVector_p)));
00112 
00113   // count the total energy and the number of superclusters
00114   reco::SuperClusterCollection::iterator it;
00115   for (it = superclusters_ap->begin(); it != superclusters_ap->end(); it++)
00116     {
00117       totalE += it->energy();
00118       noSuperClusters++;
00119     }
00120 
00121   // put the SC collection in the event
00122   evt.put(superclusters_ap, superclusterCollection);
00123 
00124   delete clusterRefVector_p;
00125 }
00126 
00127 
00128 void Multi5x5SuperClusterProducer::getClusterRefVector(edm::Event& evt, std::string clusterProducer_, std::string clusterCollection_, reco::BasicClusterRefVector *clusterRefVector_p)
00129 {  
00130   edm::Handle<reco::BasicClusterCollection> bccHandle;
00131   try
00132     {
00133       evt.getByLabel(clusterProducer_, clusterCollection_, bccHandle);
00134       if (!(bccHandle.isValid()))
00135         {
00136           edm::LogError("Multi5x5SuperClusterProducerError") << "could not get a handle on the BasicCluster Collection!";
00137           clusterRefVector_p = 0;
00138         }
00139     } 
00140   catch ( cms::Exception& ex )
00141     {
00142       edm::LogError("Multi5x5SuperClusterProducerError") << "Error! can't get the product " << clusterCollection_.c_str(); 
00143       clusterRefVector_p = 0;
00144     }
00145 
00146   const reco::BasicClusterCollection *clusterCollection_p = bccHandle.product();
00147   for (unsigned int i = 0; i < clusterCollection_p->size(); i++)
00148     {
00149       clusterRefVector_p->push_back(reco::BasicClusterRef(bccHandle, i));
00150     }
00151 }                               
00152 
00153 
00154 
00155 

Generated on Tue Jun 9 17:43:15 2009 for CMSSW by  doxygen 1.5.4