CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/RecoEcal/EgammaClusterProducers/src/UncleanSCRecoveryProducer.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 // Class header file
00021 #include "RecoEcal/EgammaClusterProducers/interface/UncleanSCRecoveryProducer.h"
00022 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00023 
00024 
00025 /*
00026 UncleanSCRecoveryProducer:
00027 ^^^^^^^^^^^^^^^^^^^^^^^^^^
00028 
00029 takes the collections of flagged clean and unclean only SC 
00030 (this is the output of UnifiedSCCollectionProducer) and
00031 recovers the original collection of unclean SC.
00032 
00033 18 Aug 2010
00034 Nikolaos Rompotis and Chris Seez  - Imperial College London
00035 many thanks to David Wardrope, Shahram Rahatlou and Federico Ferri
00036 */
00037 
00038 
00039 UncleanSCRecoveryProducer::UncleanSCRecoveryProducer(const edm::ParameterSet& ps)
00040 {
00041 
00042         // get the parameters
00043         // the cleaned collection:
00044         cleanBcCollection_ = ps.getParameter<edm::InputTag>("cleanBcCollection");
00045         cleanScCollection_ = ps.getParameter<edm::InputTag>("cleanScCollection");
00046         // the uncleaned collection
00047         uncleanBcCollection_ = ps.getParameter<edm::InputTag>("uncleanBcCollection");
00048         uncleanScCollection_ = ps.getParameter<edm::InputTag>("uncleanScCollection");
00049         // the names of the products to be produced:
00050         bcCollection_ = ps.getParameter<std::string>("bcCollection");
00051         scCollection_ = ps.getParameter<std::string>("scCollection");
00052         // the products:
00053         produces< reco::BasicClusterCollection >(bcCollection_);
00054         produces< reco::SuperClusterCollection >(scCollection_);
00055 
00056 
00057 }
00058 
00059 UncleanSCRecoveryProducer::~UncleanSCRecoveryProducer() {;}
00060 
00061 void UncleanSCRecoveryProducer::produce(edm::Event& evt, 
00062                                         const edm::EventSetup& es)
00063 {
00064         // __________________________________________________________________________
00065         //
00066         // cluster collections:
00067         edm::Handle<reco::BasicClusterCollection> pCleanBC;
00068         edm::Handle<reco::SuperClusterCollection> pCleanSC;
00069         //
00070         edm::Handle<reco::BasicClusterCollection> pUncleanBC;
00071         edm::Handle<reco::SuperClusterCollection> pUncleanSC;
00072         // clean collections ________________________________________________________
00073         evt.getByLabel(cleanScCollection_, pCleanSC);
00074         if (!(pCleanSC.isValid())) 
00075         {
00076                 edm::LogWarning("MissingInput") << "could not handle clean super clusters";
00077                 return;
00078         }
00079         const  reco::SuperClusterCollection cleanSC = *(pCleanSC.product());
00080 
00081         // unclean collections ______________________________________________________
00082         evt.getByLabel(uncleanBcCollection_, pUncleanBC);
00083         if (!(pUncleanBC.isValid())) 
00084         {
00085                 edm::LogWarning("MissingInput") << "could not handle unclean Basic Clusters!";
00086                 return;
00087         }
00088         const  reco::BasicClusterCollection uncleanBC = *(pUncleanBC.product());
00089         //
00090         evt.getByLabel(uncleanScCollection_, pUncleanSC);
00091         if (!(pUncleanSC.isValid())) 
00092         {
00093 
00094                 edm::LogWarning("MissingInput") << "could not handle unclean super clusters!";
00095                 return;
00096         }
00097         const  reco::SuperClusterCollection uncleanSC = *(pUncleanSC.product());
00098         int uncleanSize = pUncleanSC->size();
00099         int cleanSize =   pCleanSC->size();
00100 
00101         LogDebug("EcalCleaning")  << "Size of Clean Collection: " << cleanSize 
00102                 << ", uncleanSize: " << uncleanSize;
00103 
00104         // collections are all taken now ____________________________________________
00105         //
00106         // the collections to be produced ___________________________________________
00107         reco::BasicClusterCollection basicClusters;
00108         reco::SuperClusterCollection superClusters;
00109         //
00110         //
00111         // collect all the basic clusters of the SC that belong to the unclean
00112         // collection and put them into the basicClusters vector
00113         // keep the information of which SC they belong to
00114         //
00115         // loop over the unclean sc: all SC will join the new collection
00116         std::vector< std::pair<int, int> > basicClusterOwner; // counting all
00117 
00118         std::vector<DetId> scUncleanSeedDetId;  // counting the unclean
00119         for (int isc =0; isc< uncleanSize; ++isc) {
00120                 const reco::SuperCluster unsc = uncleanSC[isc];    
00121                 scUncleanSeedDetId.push_back(unsc.seed()->seed());
00122                 reco::CaloCluster_iterator bciter = unsc.clustersBegin();
00123                 for (; bciter != unsc.clustersEnd(); ++bciter) {
00124                         // the basic clusters
00125                         basicClusters.push_back(**bciter);
00126                         // index of the unclean SC
00127                         basicClusterOwner.push_back( std::make_pair(isc,0) ); 
00128                 }
00129         }
00130         // loop over the clean: only the ones which are in common with the unclean
00131         // are taken into account
00132 
00133         std::vector<DetId> scCleanSeedDetId;  // counting the clean
00134         std::vector<int> isToBeKept;
00135         for (int isc =0; isc< cleanSize; ++isc) {
00136                 reco::SuperClusterRef cscRef( pCleanSC, isc );    
00137                 scCleanSeedDetId.push_back(cscRef->seed()->seed());
00138                 for (reco::CaloCluster_iterator bciter = cscRef->clustersBegin(); bciter != cscRef->clustersEnd(); ++bciter) {
00139                         // the basic clusters
00140                         basicClusters.push_back(**bciter);
00141                         // index of the clean SC
00142                         basicClusterOwner.push_back( std::make_pair(isc,1) ); 
00143                 }
00144                 if (cscRef->isInUnclean()) isToBeKept.push_back(1);
00145                 else isToBeKept.push_back(0);
00146         }
00147         //
00148         // now export the basic clusters into the event and get them back
00149         std::auto_ptr< reco::BasicClusterCollection> basicClusters_p(new reco::BasicClusterCollection);
00150         basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
00151         edm::OrphanHandle<reco::BasicClusterCollection> bccHandle =  
00152                 evt.put(basicClusters_p, bcCollection_);
00153         if (!(bccHandle.isValid())) {
00154 
00155                 edm::LogWarning("MissingInput") << "could not handle the new BasicClusters!";
00156                 return;
00157         }
00158         reco::BasicClusterCollection basicClustersProd = *bccHandle;
00159 
00160         LogDebug("EcalCleaning") <<"Got the BasicClusters from the event again";
00161         int bcSize = bccHandle->size();
00162         //
00163         // now we can create the SC collection
00164         //
00165         // run over the unclean SC: all to be kept here
00166         for (int isc=0; isc< uncleanSize; ++isc) {
00167                 reco::CaloClusterPtrVector clusterPtrVector;
00168                 // the seed is the basic cluster with the highest energy
00169                 reco::CaloClusterPtr seed; 
00170                 for (int jbc=0; jbc< bcSize; ++jbc) {
00171                         std::pair<int,int> theBcOwner = basicClusterOwner[jbc];
00172                         if (theBcOwner.first == isc && theBcOwner.second == 0) {
00173                                 reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00174                                 clusterPtrVector.push_back(currentClu);
00175                                 if (scUncleanSeedDetId[isc] == currentClu->seed()) {
00176                                         seed = currentClu;
00177                                 }
00178                         }
00179                 }
00180                 const reco::SuperCluster unsc = uncleanSC[isc]; 
00181                 reco::SuperCluster newSC(unsc.energy(), unsc.position(), seed, clusterPtrVector );
00182                 newSC.setFlags(reco::CaloCluster::uncleanOnly);
00183                 superClusters.push_back(newSC);
00184         }
00185         // run over the clean SC: only those who are in common between the
00186         // clean and unclean collection are kept
00187         for (int isc=0; isc< cleanSize; ++isc) {
00188                 reco::SuperClusterRef cscRef( pCleanSC, isc); 
00189                 if (not cscRef->isInUnclean()) continue;
00190                 reco::CaloClusterPtrVector clusterPtrVector;
00191                 // the seed is the basic cluster with the highest energy
00192                 reco::CaloClusterPtr seed; 
00193                 for (int jbc=0; jbc< bcSize; ++jbc) {
00194                         std::pair<int,int> theBcOwner = basicClusterOwner[jbc];
00195                         if (theBcOwner.first == isc && theBcOwner.second == 1) {
00196                                 reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00197                                 clusterPtrVector.push_back(currentClu);
00198                                 if (scCleanSeedDetId[isc] == currentClu->seed()) {
00199                                         seed = currentClu;
00200                                 }
00201                         }
00202                 }
00203                 reco::SuperCluster newSC(cscRef->energy(), cscRef->position(), 
00204                                          seed, clusterPtrVector );
00205                 newSC.setFlags(reco::CaloCluster::common);
00206                 superClusters.push_back(newSC);
00207         }
00208 
00209         std::auto_ptr< reco::SuperClusterCollection> 
00210                 superClusters_p(new reco::SuperClusterCollection);
00211         superClusters_p->assign(superClusters.begin(), superClusters.end());
00212 
00213         evt.put(superClusters_p, scCollection_);
00214 
00215         LogDebug("EcalCleaning")<<"Clusters (Basic/Super) added to the Event! :-)";
00216 
00217         // ----- debugging ----------
00218         // print the new collection SC quantities
00219         // print out the clean collection SC
00220         LogDebug("EcalCleaning") << "Clean Collection SC ";
00221         for (int i=0; i < cleanSize; ++i) {
00222                 const reco::SuperCluster csc = cleanSC[i];
00223                 LogDebug("EcalCleaning") << " >>> clean    #" << i << "; Energy: " << csc.energy()
00224                         << " eta: " << csc.eta() 
00225                         << " sc seed detid: " << csc.seed()->seed().rawId();
00226         }
00227         // the unclean SC
00228         LogDebug("EcalCleaning") << "Unclean Collection SC ";
00229         for (int i=0; i < uncleanSize; ++i) {
00230                 const reco::SuperCluster usc = uncleanSC[i];
00231                 LogDebug("EcalCleaning") << " >>> unclean  #" << i << "; Energy: " << usc.energy()
00232                         << " eta: " << usc.eta() 
00233                         << " sc seed detid: " << usc.seed()->seed().rawId();
00234         }
00235         // the new collection
00236         LogDebug("EcalCleaning")<<"The new SC clean collection with size "<< superClusters.size();
00237         for (unsigned int i=0; i <  superClusters.size(); ++i) {
00238                 const reco::SuperCluster nsc = superClusters[i];
00239                 LogDebug("EcalCleaning")<< " >>> newSC    #" << i << "; Energy: " << nsc.energy()
00240                         << " eta: " << nsc.eta()  << " isClean=" 
00241                         << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
00242                         << " sc seed detid: " << nsc.seed()->seed().rawId();
00243         }
00244 }