CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoHI/HiEgammaAlgos/plugins/HiSpikeCleaner.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //111
00003 // Package:    HiSpikeCleaner
00004 // Class:      HiSpikeCleaner
00005 // 
00013 //
00014 // Original Author:  Yong Kim,32 4-A08,+41227673039,
00015 //         Created:  Mon Nov  1 18:22:21 CET 2010
00016 // $Id: HiSpikeCleaner.cc,v 1.9 2011/01/21 10:45:18 vlimant Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027 #include "FWCore/Framework/interface/EventSetup.h"
00028 #include "FWCore/Framework/interface/ESHandle.h"
00029 
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/Framework/interface/MakerMacros.h"
00032 
00033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00034 #include "FWCore/Utilities/interface/InputTag.h"
00035 
00036 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00037 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00038 
00039 #include "RecoEcal/EgammaCoreTools/interface/EcalTools.h"
00040 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgoRcd.h"
00041 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00042 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00043 
00044 
00045 
00046 //
00047 // class declaration
00048 //
00049 
00050 class HiSpikeCleaner : public edm::EDProducer {
00051    public:
00052       explicit HiSpikeCleaner(const edm::ParameterSet&);
00053       ~HiSpikeCleaner();
00054 
00055    private:
00056       virtual void beginJob() ;
00057       virtual void produce(edm::Event&, const edm::EventSetup&);
00058       virtual void endJob() ;
00059       
00060       // ----------member data ---------------------------
00061 
00062    edm::InputTag sCInputProducer_;
00063    edm::InputTag rHInputProducerB_;
00064    edm::InputTag rHInputProducerE_;
00065 
00066    std::string outputCollection_;
00067    double TimingCut_;
00068    double swissCutThr_;
00069    double etCut_;
00070    
00071 };
00072 
00073 //
00074 // constants, enums and typedefs
00075 //
00076 
00077 
00078 //
00079 // static data member definitions
00080 //
00081 
00082 //
00083 // constructors and destructor
00084 //
00085 HiSpikeCleaner::HiSpikeCleaner(const edm::ParameterSet& iConfig)
00086 {
00087    //register your products
00088 /* Examples
00089    produces<ExampleData2>();
00090 
00091    //if do put with a label
00092    produces<ExampleData2>("label");
00093 */
00094    //now do what ever other initialization is needed
00095    
00096    rHInputProducerB_  = iConfig.getParameter<edm::InputTag>("recHitProducerBarrel");
00097    rHInputProducerE_  = iConfig.getParameter<edm::InputTag>("recHitProducerEndcap");
00098 
00099    sCInputProducer_  = iConfig.getParameter<edm::InputTag>("originalSuperClusterProducer");
00100    TimingCut_      = iConfig.getUntrackedParameter<double>  ("TimingCut",4.0);
00101    swissCutThr_      = iConfig.getUntrackedParameter<double>("swissCutThr",0.95);
00102    etCut_            = iConfig.getParameter<double>("etCut");
00103    
00104    outputCollection_ = iConfig.getParameter<std::string>("outputColl");
00105    produces<reco::SuperClusterCollection>(outputCollection_);
00106    
00107    
00108    
00109 }
00110 
00111 
00112 HiSpikeCleaner::~HiSpikeCleaner()
00113 {
00114    // do anything here that needs to be done at desctruction time
00115    // (e.g. close files, deallocate resources etc.)
00116 }
00117 
00118 
00119 //
00120 // member functions
00121 //
00122 
00123 // ------------ method called to produce the data  ------------
00124 void
00125 HiSpikeCleaner::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00126 {
00127    using namespace edm;
00128 
00129 
00130    // Get raw SuperClusters from the event    
00131    Handle<reco::SuperClusterCollection> pRawSuperClusters;
00132    try { 
00133       iEvent.getByLabel(sCInputProducer_, pRawSuperClusters);
00134    } catch ( cms::Exception& ex ) {
00135       edm::LogError("EgammaSCCorrectionMakerError") 
00136          << "Error! can't get the rawSuperClusters " 
00137          << sCInputProducer_.label() ;
00138    }    
00139    
00140    // Get the RecHits from the event
00141    Handle<EcalRecHitCollection> pRecHitsB;
00142    try { 
00143       iEvent.getByLabel(rHInputProducerB_, pRecHitsB);
00144    } catch ( cms::Exception& ex ) {
00145       edm::LogError("EgammaSCCorrectionMakerError") 
00146          << "Error! can't get the RecHits " 
00147          << rHInputProducerB_.label();
00148    }    
00149    // Get the RecHits from the event                                                                                                            
00150    Handle<EcalRecHitCollection> pRecHitsE;
00151    try {
00152       iEvent.getByLabel(rHInputProducerE_, pRecHitsE);
00153    } catch ( cms::Exception& ex ) {
00154       edm::LogError("EgammaSCCorrectionMakerError")
00155          << "Error! can't get the RecHits "
00156          << rHInputProducerE_.label();
00157    }
00158 
00159    
00160    // get the channel status from the DB                                                                                                     
00161    //   edm::ESHandle<EcalChannelStatus> chStatus;
00162    //   iSetup.get<EcalChannelStatusRcd>().get(chStatus);
00163    
00164    edm::ESHandle<EcalSeverityLevelAlgo> ecalSevLvlAlgoHndl;
00165    iSetup.get<EcalSeverityLevelAlgoRcd>().get(ecalSevLvlAlgoHndl);
00166    
00167    
00168    // Create a pointer to the RecHits and raw SuperClusters
00169    const reco::SuperClusterCollection *rawClusters = pRawSuperClusters.product();
00170    
00171    
00172    EcalClusterLazyTools lazyTool(iEvent, iSetup, rHInputProducerB_,rHInputProducerE_);
00173 
00174    // Define a collection of corrected SuperClusters to put back into the event
00175    std::auto_ptr<reco::SuperClusterCollection> corrClusters(new reco::SuperClusterCollection);
00176    
00177    //  Loop over raw clusters and make corrected ones
00178    reco::SuperClusterCollection::const_iterator aClus;
00179    for(aClus = rawClusters->begin(); aClus != rawClusters->end(); aClus++)
00180       {
00181          double theEt = aClus->energy()/cosh( aClus->eta() ) ;
00182          //      std::cout << " et of SC = " << theEt << std::endl;
00183 
00184          if ( theEt < etCut_ )  continue;   // cut off low pT superclusters 
00185          
00186          bool flagS = true;
00187          int severity(-100);
00188          float swissCrx(0);
00189          
00190          const reco::CaloClusterPtr seed = aClus->seed();
00191          DetId id = lazyTool.getMaximum(*seed).first;
00192          const EcalRecHitCollection & rechits = *pRecHitsB;
00193          EcalRecHitCollection::const_iterator it = rechits.find( id );
00194          
00195          if( it != rechits.end() ) {
00196             severity = ecalSevLvlAlgoHndl->severityLevel(id, rechits);
00197             swissCrx = EcalTools::swissCross   (id, rechits, 0.,true);
00198             //      std::cout << "swissCross = " << swissCrx <<std::endl;
00199             // std::cout << " timing = " << it->time() << std::endl;
00200          }
00201          
00202          if ( fabs(it->time()) > TimingCut_ ) {
00203             flagS = false;
00204             //      std::cout << " timing = " << it->time() << std::endl;
00205             //   std::cout << " timing is bad........" << std::endl; 
00206          }
00207          if ( swissCrx > (float)swissCutThr_ ) {
00208             flagS = false ;     // swissCross cut
00209             //      std::cout << "swissCross = " << swissCrx <<std::endl;   
00210             //   std::cout << " removed by swiss cross cut" << std::endl;
00211          }
00212          // - kGood        --> good channel
00213          // - kProblematic --> problematic (e.g. noisy)
00214          // - kRecovered   --> recovered (e.g. an originally dead or saturated)
00215          // - kTime        --> the channel is out of time (e.g. spike)
00216          // - kWeird       --> weird (e.g. spike)
00217          // - kBad         --> bad, not suitable to be used in the reconstruction
00218          //   enum EcalSeverityLevel { kGood=0, kProblematic, kRecovered, kTime, kWeird, kBad };
00219             
00220          
00221          reco::SuperCluster newClus;
00222          if ( flagS == true)
00223             newClus=*aClus;
00224          else
00225             continue;
00226          corrClusters->push_back(newClus);
00227       }
00228    
00229    // Put collection of corrected SuperClusters into the event
00230    iEvent.put(corrClusters, outputCollection_);   
00231    
00232 }
00233 
00234 // ------------ method called once each job just before starting event loop  ------------
00235 void 
00236 HiSpikeCleaner::beginJob()
00237 {
00238 }
00239 
00240 // ------------ method called once each job just after ending the event loop  ------------
00241 void 
00242 HiSpikeCleaner::endJob() {
00243 }
00244 
00245 //define this as a plug-in
00246 DEFINE_FWK_MODULE(HiSpikeCleaner);