CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/Calibration/EcalAlCaRecoProducers/src/AlCaECALRecHitReducer.cc

Go to the documentation of this file.
00001 #include "Calibration/EcalAlCaRecoProducers/interface/AlCaECALRecHitReducer.h"
00002 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00003 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00004 #include "DataFormats/EgammaCandidates/interface/SiStripElectron.h"
00005 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00006 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00009 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00010 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00011 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00012 #include "DataFormats/EcalDetId/interface/ESDetId.h"
00013 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00014 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00015 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00016 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00017 #include "Geometry/EcalAlgo/interface/EcalPreshowerGeometry.h"
00018 #include "Geometry/CaloTopology/interface/EcalPreshowerTopology.h"
00019 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00020 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 #include "Geometry/CaloTopology/interface/CaloTopology.h"
00023 #include "Geometry/CaloTopology/interface/EcalBarrelTopology.h"
00024 
00025 //#define ALLrecHits
00026 //#define DEBUG
00027 
00028 //#define QUICK -> if commented loop over the recHits of the SC and add them to the list of recHits to be saved
00029 //                 comment it if you want a faster module but be sure the window is large enough
00030 
00031 AlCaECALRecHitReducer::AlCaECALRecHitReducer(const edm::ParameterSet& iConfig) 
00032 {
00033 
00034   ebRecHitsLabel_ = iConfig.getParameter< edm::InputTag > ("ebRecHitsLabel");
00035   eeRecHitsLabel_ = iConfig.getParameter< edm::InputTag > ("eeRecHitsLabel");
00036   //  esRecHitsLabel_ = iConfig.getParameter< edm::InputTag > ("esRecHitsLabel");
00037   electronLabel_ = iConfig.getParameter< edm::InputTag > ("electronLabel");
00038 
00039   alcaBarrelHitsCollection_ = iConfig.getParameter<std::string>("alcaBarrelHitCollection");
00040   alcaEndcapHitsCollection_ = iConfig.getParameter<std::string>("alcaEndcapHitCollection");
00041   //  alcaPreshowerHitsCollection_ = iConfig.getParameter<std::string>("alcaPreshowerHitCollection");
00042   
00043   etaSize_ = iConfig.getParameter<int> ("etaSize");
00044   phiSize_ = iConfig.getParameter<int> ("phiSize");
00045   // FIXME: minimum size of etaSize_ and phiSize_
00046   if ( phiSize_ % 2 == 0 ||  etaSize_ % 2 == 0)
00047     edm::LogError("AlCaECALRecHitReducerError") << "Size of eta/phi should be odd numbers";
00048  
00049   weight_= iConfig.getParameter<double> ("eventWeight");
00050  
00051   //  esNstrips_  = iConfig.getParameter<int> ("esNstrips");
00052   //  esNcolumns_ = iConfig.getParameter<int> ("esNcolumns");
00053   
00054   //register your products
00055   produces< EBRecHitCollection > (alcaBarrelHitsCollection_) ;
00056   produces< EERecHitCollection > (alcaEndcapHitsCollection_) ;
00057   //  produces< ESRecHitCollection > (alcaPreshowerHitsCollection_) ;
00058   produces< double > ("weight") ;
00059 }
00060 
00061 
00062 AlCaECALRecHitReducer::~AlCaECALRecHitReducer()
00063 {}
00064 
00065 
00066 // ------------ method called to produce the data  ------------
00067 void
00068 AlCaECALRecHitReducer::produce (edm::Event& iEvent, 
00069                                 const edm::EventSetup& iSetup)
00070 {
00071   using namespace edm;
00072   using namespace std;
00073   
00074   EcalRecHitCollection::const_iterator recHit_itr;
00075 
00076   // get the ECAL geometry:
00077   ESHandle<CaloGeometry> geoHandle;
00078   iSetup.get<CaloGeometryRecord>().get(geoHandle);
00079 
00080   edm::ESHandle<CaloTopology> theCaloTopology;
00081   iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
00082   const CaloTopology *caloTopology = theCaloTopology.product();
00083   
00084   // Get GSFElectrons
00085   Handle<reco::GsfElectronCollection> pElectrons;
00086   iEvent.getByLabel(electronLabel_, pElectrons);
00087   if (!pElectrons.isValid()) {
00088     edm::LogError ("reading") << electronLabel_ << " not found" ; 
00089     //      std::cerr << "[AlCaECALRecHitReducer]" << electronLabel_ << " not found" ; 
00090     return ;
00091   }
00092   
00093   const reco::GsfElectronCollection * electronCollection = 
00094     pElectrons.product();
00095   
00096   // get RecHits
00097   Handle<EBRecHitCollection> barrelRecHitsHandle;
00098   bool barrelIsFull = true ;
00099   
00100   iEvent.getByLabel(ebRecHitsLabel_,barrelRecHitsHandle);
00101   if (!barrelRecHitsHandle.isValid()) {
00102     edm::LogError ("reading") << ebRecHitsLabel_ << " not found" ; 
00103     barrelIsFull = false ;
00104   }
00105   
00106   const EBRecHitCollection * barrelHitsCollection = 0 ;
00107   if (barrelIsFull)  
00108     barrelHitsCollection = barrelRecHitsHandle.product () ;
00109   
00110   // get RecHits
00111   Handle<EERecHitCollection> endcapRecHitsHandle;
00112   bool endcapIsFull = true ;
00113   
00114   iEvent.getByLabel(eeRecHitsLabel_,endcapRecHitsHandle);
00115   if (!endcapRecHitsHandle.isValid()) {
00116     edm::LogError ("reading") << eeRecHitsLabel_ << " not found" ; 
00117     endcapIsFull = false ;
00118   }
00119   
00120   const EERecHitCollection * endcapHitsCollection = 0 ;
00121   if (endcapIsFull)  
00122     endcapHitsCollection = endcapRecHitsHandle.product () ;
00123   //  const EERecHitCollection * endcapHitsCollection = endcapRecHitsHandle.product();
00124   
00125   //   // get ES RecHits
00126   //   Handle<ESRecHitCollection> preshowerRecHitsHandle;
00127   //   bool preshowerIsFull = true ;
00128   
00129   //   iEvent.getByLabel(esRecHitsLabel_,preshowerRecHitsHandle);
00130   //   if (!preshowerRecHitsHandle.isValid()) {
00131   //     edm::LogError ("reading") << esRecHitsLabel_ << " not found" ; 
00132   //     preshowerIsFull = false ;
00133   //   }
00134   
00135   //   const ESRecHitCollection * preshowerHitsCollection = 0 ;
00136   //   if (preshowerIsFull)  
00137   //     preshowerHitsCollection = preshowerRecHitsHandle.product () ;
00138 
00139   //   // make a vector to store the used ES rechits:
00140   //   set<ESDetId> used_strips;
00141   //   used_strips.clear();
00142  
00143   //Create empty output collections
00144   std::auto_ptr< EBRecHitCollection > miniEBRecHitCollection (new EBRecHitCollection) ;
00145   std::auto_ptr< EERecHitCollection > miniEERecHitCollection (new EERecHitCollection) ;  
00146   //  std::auto_ptr< ESRecHitCollection > miniESRecHitCollection (new ESRecHitCollection) ;  
00147   std::auto_ptr< double > weight (new double(1));
00148   (*weight) = weight_;
00149   
00150   //  loop on SiStrip Electrons
00151   
00152   reco::GsfElectronCollection::const_iterator eleIt;
00153   int nEle_EB=0;
00154   int nEle_EE=0;
00155 
00156   std::set<DetId> reducedRecHit_EBmap;
00157   std::set<DetId> reducedRecHit_EEmap;
00158 
00159   for (eleIt=electronCollection->begin(); eleIt!=electronCollection->end(); eleIt++) {
00160     // barrel
00161     const reco::SuperCluster& sc = *(eleIt->superCluster()) ;
00162 
00163     if (eleIt->isEB()) {
00164       nEle_EB++;
00165 
00166       // find the seed 
00167       EBDetId seed=(sc.seed()->seed());
00168 
00169       std::vector<DetId> recHit_window = caloTopology->getWindow(seed, phiSize_, etaSize_);
00170       for(unsigned int i =0; i < recHit_window.size(); i++){
00171 #ifdef DEBUG2
00172         std::cout << i << "/" << recHit_window.size() << "\t" << recHit_window[i]() << std::endl;
00173 #endif
00174         reducedRecHit_EBmap.insert(recHit_window[i]);
00175 #ifdef DEBUG2
00176         EBDetId ebrechit(recHit_window[i]);
00177         std::cout << ebrechit.ieta() << "\t" << ebrechit.iphi() << std::endl;
00178 #endif
00179       }
00180       
00181 #ifdef DEBUG      
00182       // find the most energetic crystal 
00183       float energy_recHit_max=-999;
00184 
00185       if(reducedRecHit_EBmap.size() < sc.size())
00186         std::cerr << "[WARNING] number of recHit in selected window < RECO SC recHits!" << std::endl;
00187 #endif
00188 
00189 #ifndef QUICK
00190       const std::vector< std::pair<DetId, float> > & scHits = sc.hitsAndFractions();
00191 
00192 #ifdef DEBUG
00193       std::vector< std::pair<DetId, float> >::const_iterator scHit_max_itr = scHits.end();
00194 #endif
00195       for(std::vector< std::pair<DetId, float> >::const_iterator scHit_itr = scHits.begin(); 
00196           scHit_itr != scHits.end(); scHit_itr++){
00197         // the map fills just one time (avoiding double insert of recHits)
00198         reducedRecHit_EBmap.insert(scHit_itr->first);
00199 
00200 #ifdef DEBUG2
00201         const EcalRecHit ecalRecHit = *(barrelHitsCollection->find( (*scHit_itr).first ));
00202         if(energy_recHit_max <  ecalRecHit.energy()){
00203           scHit_max_itr = scHit_itr;
00204           energy_recHit_max=ecalRecHit.energy();
00205         }
00206 #endif
00207       }
00208 #endif
00209 
00210 #ifdef DEBUG2
00211       // cross check, the seed should be the highest energetic crystal in the SC
00212       if(EBDetId(scHit_max_itr->first) != seed)
00213         std::cerr << "[ERROR] highest energetic crystal is not the seed of the SC" << std::endl;
00214 
00215       else{
00216         
00217         std::cout << "[DEBUG] highest energetic crystal = " << EBDetId(scHit_max_itr->first) << std::endl;
00218         std::cout << "[DEBUG] seed of the SC = " << seed << std::endl;
00219       }
00220 #endif
00221       //                                                                   (id, phi, eta)
00222 
00223       if(reducedRecHit_EBmap.size() < sc.size()){
00224         if(eleIt->ecalDrivenSeed())
00225           edm::LogError("AlCaSavedRecHitsEB") << "[ERROR] ecalDrivenSeed: number of saved recHits < RECO SC recHits!: " << reducedRecHit_EBmap.size() << " < " << sc.size() << std::endl;
00226         else
00227           edm::LogWarning("AlCaSavedRecHitsEB") << "[WARNING] trackerDrivenSeed: number of saved recHits < RECO SC recHits!: " << reducedRecHit_EBmap.size() << " < " << sc.size() << std::endl;
00228 
00229       }
00230 
00231     } else { // endcap
00232       nEle_EE++;
00233 
00234       // find the seed 
00235       EEDetId seed=(sc.seed()->seed());
00236 
00237       // get detId for a window around the seed of the SC
00238       int sideSize = std::max(phiSize_,etaSize_);
00239       std::vector<DetId> recHit_window = caloTopology->getWindow(seed, sideSize, sideSize);
00240 
00241       // fill the recHit map with the DetId of the window
00242       for(std::vector<DetId>::const_iterator window_itr = recHit_window.begin(); 
00243           window_itr != recHit_window.end(); window_itr++){
00244         reducedRecHit_EEmap.insert(*window_itr);
00245       }
00246 #ifdef DEBUG
00247       if(reducedRecHit_EEmap.size() < sc.size())
00248         std::cerr << "[WARNING] number of recHit in selected window < RECO SC recHits!" << std::endl;
00249 #endif
00250 
00251       const std::vector< std::pair<DetId, float> > & scHits = sc.hitsAndFractions();
00252 
00253 #ifndef QUICK
00254       // fill the recHit map with the DetId of the SC recHits
00255 
00256       for(std::vector< std::pair<DetId, float> >::const_iterator scHit_itr = scHits.begin(); 
00257           scHit_itr != scHits.end(); scHit_itr++){
00258         // the map fills just one time (avoiding double insert of recHits)
00259         reducedRecHit_EEmap.insert(scHit_itr->first);   
00260 
00261       }
00262 #endif
00263 
00264       if(reducedRecHit_EEmap.size() < sc.size()){
00265         if(eleIt->ecalDrivenSeed())
00266           edm::LogError("AlCaSavedRecHitsEE") << "[ERROR] ecalDrivenSeed: number of saved recHits < RECO SC recHits!: " << reducedRecHit_EEmap.size() << " < " << sc.size() << std::endl;
00267         else
00268           edm::LogWarning("AlCaSavedRecHitsEE") << "[WARNING] trackerDrivenSeed: number of saved recHits < RECO SC recHits!: " << reducedRecHit_EEmap.size() << " < " << sc.size() << std::endl;
00269 
00270       }
00271 
00272     } // end of endcap
00273   }
00274 
00275 #ifndef ALLrecHits
00276   for(std::set<DetId>::const_iterator itr = reducedRecHit_EBmap.begin();
00277       itr != reducedRecHit_EBmap.end(); itr++){
00278     if (barrelHitsCollection->find(*itr) != barrelHitsCollection->end())
00279       miniEBRecHitCollection->push_back(*(barrelHitsCollection->find(*itr)));
00280   }
00281 #else
00282   for(EcalRecHitCollection::const_iterator recHits_itr = barrelHitsCollection->begin();
00283       recHits_itr != barrelHitsCollection->end();
00284       recHits_itr++){
00285     miniEBRecHitCollection->push_back(*recHits_itr);
00286   }
00287 #endif
00288 
00289   // fill the alcareco reduced recHit collection
00290   for(std::set<DetId>::const_iterator itr = reducedRecHit_EEmap.begin();
00291       itr != reducedRecHit_EEmap.end(); itr++){
00292     if (endcapHitsCollection->find(*itr) != endcapHitsCollection->end())
00293       miniEERecHitCollection->push_back(*(endcapHitsCollection->find(*itr)));
00294   }
00295 
00296 
00297 
00298 #ifdef DEBUG
00299   std::cout << "nEle_EB= " << nEle_EB << "\tnEle_EE = " << nEle_EE << std::endl;
00300   if(nEle_EB > 0 && miniEBRecHitCollection->size() < (unsigned int) phiSize_*etaSize_)
00301     edm::LogError("AlCaECALRecHitReducerError") << "Size EBRecHitCollection < " << phiSize_*etaSize_ << ": "  << miniEBRecHitCollection->size() ;
00302 
00303   int side = phiSize_ ;
00304   if (phiSize_ < etaSize_) side = etaSize_ ;
00305 
00306   if(nEle_EE > 0 && miniEERecHitCollection->size() < (unsigned int )side*side)
00307     edm::LogError("AlCaECALRecHitReducerError") << "Size EERecHitCollection < " << side*side << ": "  << miniEERecHitCollection->size() ;
00308 #endif
00309   
00310   //Put selected information in the event
00311   iEvent.put( miniEBRecHitCollection,alcaBarrelHitsCollection_ );
00312   iEvent.put( miniEERecHitCollection,alcaEndcapHitsCollection_ );     
00313   //  iEvent.put( miniESRecHitCollection,alcaPreshowerHitsCollection_ );     
00314   iEvent.put( weight, "weight");     
00315 }
00316 
00317 
00318 DEFINE_FWK_MODULE(AlCaECALRecHitReducer);
00319