CMS 3D CMS Logo

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