CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoEcal/EgammaClusterProducers/src/CleanAndMergeProducer.cc

Go to the documentation of this file.
00001 // C/C++ headers
00002 #include <iostream>
00003 #include <vector>
00004 #include <memory>
00005 
00006 // Framework
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 // Reconstruction Classes
00013 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
00014 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00015 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00016 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00017 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00018 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00019 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00020 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
00021 #include "CondFormats/EcalObjects/interface/EcalChannelStatus.h"
00022 
00023 // Geometry
00024 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00025 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00026 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00027 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00028 #include "Geometry/CaloTopology/interface/EcalBarrelTopology.h"
00029 #include "Geometry/CaloTopology/interface/EcalEndcapTopology.h"
00030 #include "Geometry/CaloTopology/interface/EcalPreshowerTopology.h"
00031 
00032 // Class header file
00033 #include "RecoEcal/EgammaClusterProducers/interface/CleanAndMergeProducer.h"
00034 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00035 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00036 
00037 
00038 /*
00039 CleanAndMergeProducer:
00040 ^^^^^^^^^^^^^^^^^^^^^^
00041 
00042 Takes  as  input the cleaned  and the  uncleaned collection of SC
00043 and produces an uncleaned collection of SC without the duplicates
00044 and  a collection of  references  to the cleaned collection of SC
00045 
00046 25 June 2010
00047 Nikolaos Rompotis and Chris Seez  - Imperial College London
00048 many thanks to Shahram Rahatlou and Federico Ferri
00049 
00050 
00051 */
00052 
00053 
00054 CleanAndMergeProducer::CleanAndMergeProducer(const edm::ParameterSet& ps)
00055 {
00056 
00057         // get the parameters
00058         // the cleaned collection:
00059         cleanScInputTag_   = ps.getParameter<edm::InputTag>("cleanScInputTag");
00060         uncleanScInputTag_ = ps.getParameter<edm::InputTag>("uncleanScInputTag");          
00061 
00062         // the names of the products to be produced:
00063         bcCollection_ = ps.getParameter<std::string>("bcCollection");
00064         scCollection_ = ps.getParameter<std::string>("scCollection");
00065         refScCollection_ = ps.getParameter<std::string>("refScCollection");
00066 
00067         std::map<std::string,double> providedParameters;  
00068         providedParameters.insert(std::make_pair("LogWeighted",ps.getParameter<bool>("posCalc_logweight")));
00069         providedParameters.insert(std::make_pair("T0_barl",ps.getParameter<double>("posCalc_t0")));
00070         providedParameters.insert(std::make_pair("W0",ps.getParameter<double>("posCalc_w0")));
00071         providedParameters.insert(std::make_pair("X0",ps.getParameter<double>("posCalc_x0")));
00072 
00073         hitproducer_ = ps.getParameter<std::string>("ecalhitproducer");
00074         hitcollection_ =ps.getParameter<std::string>("ecalhitcollection");
00075 
00076         // the products:
00077         produces< reco::BasicClusterCollection >(bcCollection_);
00078         produces< reco::SuperClusterCollection >(scCollection_);
00079         produces< reco::SuperClusterRefVector >(refScCollection_);
00080 
00081 }
00082 
00083 CleanAndMergeProducer::~CleanAndMergeProducer() {;}
00084 
00085 void CleanAndMergeProducer::produce(edm::Event& evt, 
00086                                     const edm::EventSetup& es)
00087 {
00088         // get the input collections
00089         // ______________________________________________________________________________________
00090         //
00091         // cluster collections:
00092 
00093         edm::Handle<reco::SuperClusterCollection> pCleanSC;
00094         edm::Handle<reco::SuperClusterCollection> pUncleanSC;
00095 
00096         evt.getByLabel(cleanScInputTag_, pCleanSC);
00097         if (!(pCleanSC.isValid())) 
00098         {
00099           edm::LogWarning("MissingInput")<< "could not get a handle on the clean Super Clusters!";
00100           return;
00101         }
00102 
00103         evt.getByLabel(uncleanScInputTag_, pUncleanSC);
00104         if (!(pUncleanSC.isValid())) 
00105         {
00106           
00107           edm::LogWarning("MissingInput")<< "could not get a handle on the unclean Super Clusters!";
00108           return;
00109         }
00110 
00111         //
00112         // collections are all taken now _____________________________________________________
00113         //
00114         //
00115         // the collections to be produced:
00116         reco::BasicClusterCollection basicClusters;
00117         reco::SuperClusterCollection superClusters;
00118         reco::SuperClusterRefVector * scRefs = new reco::SuperClusterRefVector;
00119 
00120         //
00121         // run over the uncleaned SC and check how many of them are matched to the cleaned ones
00122         // if you find a matched one, create a reference to the cleaned collection and store it
00123         // if you find an unmatched one, then keep all its basic clusters in the basic Clusters 
00124         // vector
00125         int uncleanSize =  pUncleanSC->size();
00126         int cleanSize =    pCleanSC->size();
00127 
00128         LogDebug("EcalCleaning") << "Size of Clean Collection: " << cleanSize 
00129                 << ", uncleanSize: " << uncleanSize  << std::endl;
00130         //
00131         // keep whether the SC in unique in the uncleaned collection
00132         std::vector<int> isUncleanOnly;     // 1 if unique in the uncleaned
00133         std::vector<int> basicClusterOwner; // contains the index of the SC that owns that BS
00134         std::vector<int> isSeed; // if this basic cluster is a seed it is 1
00135         for (int isc =0; isc< uncleanSize; ++isc) {
00136           reco::SuperClusterRef unscRef( pUncleanSC, isc );    
00137           const std::vector< std::pair<DetId, float> > & uhits = unscRef->hitsAndFractions();
00138           int uhitsSize = uhits.size();
00139           bool foundTheSame = false;
00140           for (int jsc=0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
00141             reco::SuperClusterRef cscRef( pCleanSC, jsc );
00142             const std::vector< std::pair<DetId, float> > & chits = cscRef->hitsAndFractions();
00143             int chitsSize =  chits.size();
00144             foundTheSame = true;
00145             if (unscRef->seed()->seed() == cscRef->seed()->seed() && chitsSize == uhitsSize) { 
00146               // if the clusters are exactly the same then because the clustering
00147               // algorithm works in a deterministic way, the order of the rechits
00148               // will be the same
00149               for (int i=0; i< chitsSize; ++i) {
00150                 if (uhits[i].first != chits[i].first ) { foundTheSame=false;  break;}
00151               }
00152               if (foundTheSame) { // ok you have found it!
00153                 // make the reference
00154                 //scRefs->push_back( edm::Ref<reco::SuperClusterCollection>(pCleanSC, jsc) );
00155                 scRefs->push_back( cscRef );
00156                 isUncleanOnly.push_back(0);
00157                 break;
00158               }
00159             }
00160           }
00161           if (not foundTheSame) {
00162             // mark it as unique in the uncleaned
00163             isUncleanOnly.push_back(1);
00164             // keep all its basic clusters
00165             for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00166               // the basic clusters
00167               basicClusters.push_back(**bciter);
00168               basicClusterOwner.push_back(isc);
00169             }
00170           }
00171         }
00172         int bcSize =  basicClusters.size();
00173         
00174         LogDebug("EcalCleaning") << "Found cleaned SC: " << cleanSize <<  " uncleaned SC: " 
00175                 << uncleanSize << " from which " << scRefs->size() 
00176                 << " will become refs to the cleaned collection" ;
00177      
00178 
00179         // now you have the collection of basic clusters of the SC to be remain in the
00180         // in the clean collection, export them to the event
00181         // you will need to reread them later in order to make correctly the refs to the SC
00182         std::auto_ptr< reco::BasicClusterCollection> basicClusters_p(new reco::BasicClusterCollection);
00183         basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
00184         edm::OrphanHandle<reco::BasicClusterCollection> bccHandle =  
00185                 evt.put(basicClusters_p, bcCollection_);
00186         if (!(bccHandle.isValid())) {
00187           edm::LogWarning("MissingInput")<<"could not get a handle on the BasicClusterCollection!" << std::endl;
00188           return;
00189         }
00190 
00191         LogDebug("EcalCleaning")<< "Got the BasicClusters from the event again";
00192         // now you have to create again your superclusters
00193         // you run over the uncleaned SC, but now you know which of them are
00194         // the ones that are needed and which are their basic clusters
00195         for (int isc=0; isc< uncleanSize; ++isc) {
00196           if (isUncleanOnly[isc]==1) { // look for sc that are unique in the unclean collection
00197             // make the basic cluster collection
00198             reco::CaloClusterPtrVector clusterPtrVector;
00199             reco::CaloClusterPtr seed; // the seed is the basic cluster with the highest energy
00200             double energy = -1;
00201             for (int jbc=0; jbc< bcSize; ++jbc) {
00202               if (basicClusterOwner[jbc]==isc) {
00203                 reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandle, jbc);
00204                 clusterPtrVector.push_back(currentClu);
00205                 if (energy< currentClu->energy()) {
00206                   energy = currentClu->energy(); seed = currentClu;
00207                 }
00208               }
00209             }
00210             reco::SuperClusterRef unscRef( pUncleanSC, isc ); 
00211             reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), seed, clusterPtrVector );
00212             superClusters.push_back(newSC);
00213           }
00214           
00215         }
00216         // export the collection of references to the clean collection
00217         std::auto_ptr< reco::SuperClusterRefVector >  scRefs_p( scRefs );
00218         evt.put(scRefs_p, refScCollection_);
00219 
00220         // the collection of basic clusters is already in the event
00221         // the collection of uncleaned SC
00222         std::auto_ptr< reco::SuperClusterCollection > superClusters_p(new reco::SuperClusterCollection);
00223         superClusters_p->assign(superClusters.begin(), superClusters.end());
00224         evt.put(superClusters_p, scCollection_);
00225         LogDebug("EcalCleaning")<< "Hybrid Clusters (Basic/Super) added to the Event! :-)";
00226 }