CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoEcal/EgammaClusterProducers/src/EcalDigiSelector.cc

Go to the documentation of this file.
00001 #include <vector>
00002 #include "RecoEcal/EgammaClusterProducers/interface/EcalDigiSelector.h"
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "DataFormats/DetId/interface/DetId.h"
00006 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00007 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00008 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00009 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00010 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00011 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00012 #include "DataFormats/EcalDigi/interface/EcalDataFrame.h"
00013 #include "DataFormats/EcalDigi/interface/EEDataFrame.h"
00014 #include "DataFormats/EcalDigi/interface/EBDataFrame.h"
00015 #include "Geometry/CaloTopology/interface/CaloTopology.h"
00016 #include "Geometry/CaloEventSetup/interface/CaloTopologyRecord.h"
00017 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00018 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
00019 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
00020 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include <TMath.h>
00023 #include <iostream>
00024 #include <set>
00025 
00026 
00027 using namespace reco;
00028 EcalDigiSelector::EcalDigiSelector(const edm::ParameterSet& ps)
00029 {
00030   selectedEcalEBDigiCollection_ = ps.getParameter<std::string>("selectedEcalEBDigiCollection");
00031   selectedEcalEEDigiCollection_ = ps.getParameter<std::string>("selectedEcalEEDigiCollection");
00032 
00033   barrelSuperClusterProducer_ = ps.getParameter<edm::InputTag>("barrelSuperClusterProducer");
00034   endcapSuperClusterProducer_ = ps.getParameter<edm::InputTag>("endcapSuperClusterProducer");
00035 
00036   EcalEBDigiTag_ = ps.getParameter<edm::InputTag>("EcalEBDigiTag");
00037   EcalEEDigiTag_ = ps.getParameter<edm::InputTag>("EcalEEDigiTag");
00038 
00039   EcalEBRecHitTag_ = ps.getParameter<edm::InputTag>("EcalEBRecHitTag");
00040   EcalEERecHitTag_ = ps.getParameter<edm::InputTag>("EcalEERecHitTag");
00041   
00042   cluster_pt_thresh_ = ps.getParameter<double>("cluster_pt_thresh");
00043   single_cluster_thresh_ = ps.getParameter<double>("single_cluster_thresh");
00044 
00045   nclus_sel_ = ps.getParameter<int>("nclus_sel");
00046   produces<EBDigiCollection>(selectedEcalEBDigiCollection_);
00047   produces<EEDigiCollection>(selectedEcalEEDigiCollection_);
00048 
00049 }
00050 
00051 
00052 EcalDigiSelector::~EcalDigiSelector()
00053 {
00054 }
00055 
00056 void EcalDigiSelector::produce(edm::Event& evt, const edm::EventSetup& es)
00057 {
00058   
00059   //Get BarrelSuperClusters to start.
00060   edm::Handle<reco::SuperClusterCollection> pBarrelSuperClusters;
00061   
00062   evt.getByLabel(barrelSuperClusterProducer_, pBarrelSuperClusters);
00063   if (!pBarrelSuperClusters.isValid()){
00064     edm::LogError("EcalDigiSelector")
00065       << "can't get collection with label " << barrelSuperClusterProducer_ ;
00066   }
00067   const reco::SuperClusterCollection & BarrelSuperClusters = *pBarrelSuperClusters;
00068   //Got BarrelSuperClusters
00069 
00070   //Get BarrelSuperClusters to start.
00071   edm::Handle<reco::SuperClusterCollection> pEndcapSuperClusters;
00072   
00073   evt.getByLabel(endcapSuperClusterProducer_, pEndcapSuperClusters);
00074   if (!pEndcapSuperClusters.isValid()){
00075     edm::LogError("EcalDigiSelector")
00076       << "Error! can't get collection with label " << endcapSuperClusterProducer_;
00077   }
00078   const reco::SuperClusterCollection & EndcapSuperClusters = *pEndcapSuperClusters;
00079   //Got EndcapSuperClusters
00080 
00081   reco::SuperClusterCollection saveBarrelSuperClusters;
00082   reco::SuperClusterCollection saveEndcapSuperClusters;
00083   bool meet_single_thresh = false;
00084   //Loop over barrel superclusters, and apply threshold
00085   for (int loop=0;loop<int(BarrelSuperClusters.size());loop++){
00086     SuperCluster clus1 = BarrelSuperClusters[loop];
00087     float eta1 = clus1.eta();
00088     float energy1 = clus1.energy();
00089     float theta1 = 2*atan(exp(-1.*eta1));
00090     float cluspt1 = energy1 * sin(theta1);
00091     if (cluspt1 > cluster_pt_thresh_){
00092       saveBarrelSuperClusters.push_back(clus1);
00093       if (cluspt1 > single_cluster_thresh_)
00094         meet_single_thresh = true;
00095       
00096     }
00097   }
00098 
00099   //Loop over endcap superclusters, and apply threshold
00100   for (int loop=0;loop<int(EndcapSuperClusters.size());loop++){
00101     SuperCluster clus1 = EndcapSuperClusters[loop];
00102     float eta1 = clus1.eta();
00103     float energy1 = clus1.energy();
00104     float theta1 = 2*atan(exp(-1.*eta1));
00105     float cluspt1 = energy1 * sin(theta1);
00106     if (cluspt1 > cluster_pt_thresh_){
00107       saveEndcapSuperClusters.push_back(clus1);
00108       if (cluspt1 > single_cluster_thresh_)
00109         meet_single_thresh = true;
00110     }
00111   }
00112   
00113   std::auto_ptr<EBDigiCollection> SEBDigiCol(new EBDigiCollection);
00114   std::auto_ptr<EEDigiCollection> SEEDigiCol(new EEDigiCollection);
00115   int TotClus = saveBarrelSuperClusters.size() + saveEndcapSuperClusters.size();
00116   //  std::cout << "Barrel Clusters: " << saveBarrelSuperClusters.size();
00117   //  std::cout << " Endcap Clusters: " << saveEndcapSuperClusters.size();
00118   //  std::cout << " Total: " << TotClus << std::endl;
00119   if (TotClus >= nclus_sel_ || meet_single_thresh){
00120     EcalClusterLazyTools tooly(evt, es, EcalEBRecHitTag_, EcalEERecHitTag_);
00121     if (saveBarrelSuperClusters.size() > 0){
00122       
00123       //Get barrel rec hit collection
00124 //       edm::Handle<EcalRecHitCollection> ecalhitsCollH;
00125 //       evt.getByLabel(EcalEBRecHitTag_, ecalhitsCollH);
00126 //       const EcalRecHitCollection* rechitsCollection = ecalhitsCollH.product();
00127 //      std::set<DetId> saveTheseDetIds;
00128       
00129       //get barrel digi collection
00130       edm::Handle<EBDigiCollection> pdigis;
00131       const EBDigiCollection* digis=0;
00132       evt.getByLabel(EcalEBDigiTag_,pdigis);
00133       if (!pdigis.isValid()){
00134               edm::LogError("EcalDigiSelector")
00135                       << "can't get collection with label " << EcalEBDigiTag_ ;
00136       } else {
00137               digis = pdigis.product(); // get a ptr to the product
00138       }
00139       if ( digis ) {
00140          std::vector<DetId> saveTheseDetIds;
00141          //pick out the detids for the 3x3 in each of the selected superclusters
00142          for (int loop = 0;loop < int(saveBarrelSuperClusters.size());loop++){
00143            SuperCluster clus1 = saveBarrelSuperClusters[loop];
00144            CaloClusterPtr bcref = clus1.seed();
00145            const BasicCluster *bc = bcref.get();
00146            //Get the maximum detid
00147            std::pair<DetId, float> EDetty = tooly.getMaximum(*bc);
00148            //get the 3x3 array centered on maximum detid.
00149            std::vector<DetId> detvec = tooly.matrixDetId(EDetty.first, -1, 1, -1, 1);
00150            //Loop over the 3x3
00151            for (int ik = 0;ik<int(detvec.size());++ik)
00152              saveTheseDetIds.push_back(detvec[ik]);
00153            //saveTheseDetIds.insert(detvec[ik]);
00154          }
00155          for (int detloop=0; detloop < int(saveTheseDetIds.size());++detloop){
00156                 EBDetId detL = EBDetId(saveTheseDetIds[detloop]);
00157            //      int ebcounter=0;
00158            for (EBDigiCollection::const_iterator blah = digis->begin();
00159                 blah!=digis->end();blah++){
00160              
00161              if (detL == blah->id()){
00162                EBDataFrame myDigi = (*blah);
00163                SEBDigiCol->push_back(detL);
00164                //SEBDigiCol->push_back(detL, myDigi);
00165                EBDataFrame df( SEBDigiCol->back());
00166                for (int iq =0;iq<myDigi.size();++iq){
00167                  df.setSample(iq, myDigi.sample(iq).raw());
00168                }
00169                //ebcounter++;
00170              } 
00171            }
00172            //if (ebcounter >= int(saveTheseDetIds.size())) break;
00173          }//loop over dets
00174     
00175 
00176          //      std::cout << "size of new digi container contents: " << SEBDigiCol->size() << std::endl;
00177 
00178       }
00179 
00180     }//If barrel superclusters need saving.
00181     
00182     
00183     if (saveEndcapSuperClusters.size() > 0){
00184       //Get endcap rec hit collection
00185       //get endcap digi collection
00186       edm::Handle<EEDigiCollection> pdigis;
00187       const EEDigiCollection* digis=0;
00188       evt.getByLabel(EcalEEDigiTag_,pdigis);
00189       if (!pdigis.isValid()){
00190               edm::LogError("EcalDigiSelector")
00191                       << "can't get collection with label " << EcalEEDigiTag_ ;
00192       } else {
00193               digis = pdigis.product(); // get a ptr to the product
00194       }
00195       if ( digis ) {
00196          //std::vector<DetId> saveTheseDetIds;
00197          std::set<DetId> saveTheseDetIds;
00198          //pick out the digis for the 3x3 in each of the selected superclusters
00199          for (int loop = 0;loop < int(saveEndcapSuperClusters.size());loop++){
00200            SuperCluster clus1 = saveEndcapSuperClusters[loop];
00201            CaloClusterPtr bcref = clus1.seed();
00202            const BasicCluster *bc = bcref.get();
00203            //Get the maximum detid
00204            std::pair<DetId, float> EDetty = tooly.getMaximum(*bc);
00205            //get the 3x3 array centered on maximum detid.
00206            std::vector<DetId> detvec = tooly.matrixDetId(EDetty.first, -1, 1, -1, 1);
00207            //Loop over the 3x3
00208            for (int ik = 0;ik<int(detvec.size());++ik)
00209              //saveTheseDetIds.push_back(detvec[ik]);
00210              saveTheseDetIds.insert(detvec[ik]);
00211          }
00212          int eecounter=0;
00213          for (EEDigiCollection::const_iterator blah = digis->begin();
00214               blah!=digis->end();blah++){
00215            std::set<DetId>::const_iterator finder = saveTheseDetIds.find(blah->id());
00216            if (finder!=saveTheseDetIds.end()){
00217              EEDetId detL = EEDetId(*finder);
00218 
00219              if (detL == blah->id()){
00220                EEDataFrame myDigi = (*blah);
00221                SEEDigiCol->push_back(detL);
00222                EEDataFrame df( SEEDigiCol->back());
00223                for (int iq =0;iq<myDigi.size();++iq){
00224                  df.setSample(iq, myDigi.sample(iq).raw());           
00225                }
00226                eecounter++;
00227                
00228              }
00229            }
00230            if (eecounter >= int(saveTheseDetIds.size())) break;
00231          }//loop over digis
00232          //      std::cout << "Current new digi container contents: " << SEEDigiCol->size() << std::endl;
00233       }
00234     }//If endcap superclusters need saving.
00235     
00236   }//If we're actually saving stuff
00237   
00238   //Okay, either my collections have been filled with the requisite Digis, or they haven't.
00239   
00240   //std::cout << "Saving: " << SEBDigiCol->size() << " barrel digis and " << SEEDigiCol->size() << " endcap digis." << std::endl;
00241 
00242   //Empty collection, or full, still put in event.
00243   SEBDigiCol->sort();
00244   SEEDigiCol->sort();
00245   evt.put(SEBDigiCol, selectedEcalEBDigiCollection_);
00246   evt.put(SEEDigiCol, selectedEcalEEDigiCollection_);
00247   
00248 }