CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoEcal/EgammaClusterProducers/src/UnifiedSCCollectionProducer.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 // Class header file
00024 #include "RecoEcal/EgammaClusterProducers/interface/UnifiedSCCollectionProducer.h"
00025 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00026 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00027 
00028 
00029 /*
00030 UnifiedSCCollectionProducer:
00031 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^
00032 
00033 Takes  as  input the cleaned  and the  uncleaned collection of SC
00034 and  produces two collections of SC: one with the clean SC, but flagged
00035 such that with the algoID value one can identify the SC that are also
00036 in the unclean collection and a collection with the unclean only SC.
00037 This collection has the algoID enumeration of the SC altered
00038 such that:
00039 flags = 0   (cleanedOnly)     cluster is only in the cleaned collection
00040 flags = 100 (common)          cluster is common in both collections
00041 flags = 200 (uncleanedOnly)   cluster is only in the uncleaned collection
00042 
00043 In that way the user can get hold of objects from the
00044 -  cleaned   collection only if they choose flags <  200
00045 -  uncleaned collection only if they choose flags >= 100
00046 
00047 18 Aug 2010
00048 Nikolaos Rompotis and Chris Seez  - Imperial College London
00049 many thanks to David Wardrope, Shahram Rahatlou and Federico Ferri
00050 */
00051 
00052 
00053 UnifiedSCCollectionProducer::UnifiedSCCollectionProducer(const edm::ParameterSet& ps)
00054 {
00055 
00056         // get the parameters
00057         // the cleaned collection:
00058         cleanBcCollection_ = ps.getParameter<edm::InputTag>("cleanBcCollection");
00059         cleanScCollection_ = ps.getParameter<edm::InputTag>("cleanScCollection");
00060         // the uncleaned collection
00061         uncleanBcCollection_ = ps.getParameter<edm::InputTag>("uncleanBcCollection");
00062         uncleanScCollection_ = ps.getParameter<edm::InputTag>("uncleanScCollection");
00063         // the names of the products to be produced:
00064         //
00065         // the clean collection: this is as it was before, but labeled
00066         bcCollection_ = ps.getParameter<std::string>("bcCollection");
00067         scCollection_ = ps.getParameter<std::string>("scCollection");
00068         // the unclean only collection: SC unique to the unclean collection
00069         bcCollectionUncleanOnly_ = ps.getParameter<std::string>("bcCollectionUncleanOnly");
00070         scCollectionUncleanOnly_ = ps.getParameter<std::string>("scCollectionUncleanOnly");
00071         // the products:
00072         produces< reco::BasicClusterCollection >(bcCollection_);
00073         produces< reco::SuperClusterCollection >(scCollection_);
00074         produces< reco::BasicClusterCollection >(bcCollectionUncleanOnly_);
00075         produces< reco::SuperClusterCollection >(scCollectionUncleanOnly_);
00076 
00077 }
00078 
00079 
00080 UnifiedSCCollectionProducer::~UnifiedSCCollectionProducer() {;}
00081 
00082 
00083 void UnifiedSCCollectionProducer::produce(edm::Event& evt, 
00084                                           const edm::EventSetup& es)
00085 {
00086 
00087   edm::LogInfo("UnifiedSC")<< ">>>>> Entering UnifiedSCCollectionProducer <<<<<";
00088   // get the input collections
00089   // __________________________________________________________________________
00090   //
00091   // cluster collections:
00092   edm::Handle<reco::BasicClusterCollection> pCleanBC;
00093   edm::Handle<reco::SuperClusterCollection> pCleanSC;
00094   //
00095   edm::Handle<reco::BasicClusterCollection> pUncleanBC;
00096   edm::Handle<reco::SuperClusterCollection> pUncleanSC;
00097   
00098   evt.getByLabel(cleanScCollection_, pCleanSC);
00099   if (!(pCleanSC.isValid())) 
00100     {
00101   
00102       edm::LogError("MissingInput") << "could not handle clean super clusters";
00103       return;
00104     }
00105   
00106         evt.getByLabel(uncleanScCollection_, pUncleanSC);
00107         if (!(pUncleanSC.isValid())) 
00108         {
00109 
00110           edm::LogError("MissingInput")<< "could not handle unclean super clusters!" ;
00111           return;
00112         }
00113 
00114         // the collections to be produced ___________________________________________
00115         reco::BasicClusterCollection basicClusters;
00116         reco::SuperClusterCollection superClusters;
00117         //
00118         reco::BasicClusterCollection basicClustersUncleanOnly;
00119         reco::SuperClusterCollection superClustersUncleanOnly;
00120         //
00121         // run over the uncleaned SC and check how many of them are matched to 
00122         // the cleaned ones
00123         // if you find a matched one, then keep the info that it is matched 
00124         //    along with which clean SC was matched + its basic clusters
00125         // if you find an unmatched one, keep the info and store its basic clusters
00126         //
00127         // 
00128         int uncleanSize = pUncleanSC->size();
00129         int cleanSize =   pCleanSC->size();
00130         
00131         LogDebug("UnifiedSC") << "Size of Clean Collection: " << cleanSize 
00132                         << ", uncleanSize: " << uncleanSize;
00133 
00134 
00135         // keep the indices
00136         std::vector<int> inUncleanOnlyInd;      // counting the unclean
00137         std::vector<int> inCleanInd;            // counting the unclean
00138         std::vector<int> inCleanOnlyInd;        // counting the clean
00139         std::vector<DetId> scUncleanSeedDetId;  // counting the unclean
00140         std::vector<DetId> scCleanSeedDetId;    // counting the clean
00141         // ontains the index of the SC that owns that BS
00142         // first basic cluster index, second: 0 for unclean and 1 for clean
00143         std::vector< std::pair<int, int> > basicClusterOwner; 
00144         std::vector< std::pair<int, int> > basicClusterOwnerUncleanOnly; 
00145         // if this basic cluster is a seed it is 1
00146         std::vector<int> uncleanBasicClusterIsSeed;
00147 
00148         // loop over unclean SC _____________________________________________________
00149         for (int isc =0; isc< uncleanSize; ++isc) {
00150                 reco::SuperClusterRef unscRef( pUncleanSC, isc);    
00151                 const std::vector< std::pair<DetId, float> > & uhits = unscRef->hitsAndFractions();
00152                 int uhitsSize = uhits.size();
00153                 bool foundTheSame = false;
00154                 for (int jsc=0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
00155                         reco::SuperClusterRef cscRef( pCleanSC, jsc );
00156                         const std::vector<std::pair<DetId,float> > & chits = cscRef->hitsAndFractions();
00157                         int chitsSize =  chits.size();
00158                         foundTheSame = false;
00159                         if (unscRef->seed()->seed()==cscRef->seed()->seed() && chitsSize == uhitsSize) { 
00160                                 // if the clusters are exactly the same then because the clustering
00161                                 // algorithm works in a deterministic way, the order of the rechits
00162                                 // will be the same
00163                                 for (int i=0; i< chitsSize; ++i) {
00164                                         if (uhits[i].first != chits[i].first ) { break;}
00165                                 }
00166                                 foundTheSame = true;
00167                         }
00168                         if (foundTheSame) { // ok you have found it:
00169                                 // this supercluster belongs to both collections
00170                                 inUncleanOnlyInd.push_back(0);
00171                                 inCleanInd.push_back(jsc); // keeps the index of the clean SC
00172                                 scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00173                                 //
00174                                 // keep its basic clusters:
00175                                 for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00176                                         // the basic clusters
00177                                         basicClusters.push_back(**bciter);
00178                                         // index of the unclean SC
00179                                         basicClusterOwner.push_back( std::make_pair(isc,0) ); 
00180                                 } 
00181                                 break; // break the loop over unclean sc
00182                         }
00183                 }
00184                 if (not foundTheSame) { // this SC is only in the unclean collection
00185                         // mark it as unique in the uncleaned
00186                         inUncleanOnlyInd.push_back(1);
00187                         scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00188                         // keep all its basic clusters
00189                         for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00190                                 // the basic clusters
00191                                 basicClustersUncleanOnly.push_back(**bciter);
00192                                 basicClusterOwnerUncleanOnly.push_back( std::make_pair(isc,0) );
00193                         }
00194                 }
00195         } // loop over the unclean SC _______________________________________________
00196         //
00197         int inCleanSize = inCleanInd.size();
00198         //
00199         // loop over the clean SC, check that are not in common with the unclean
00200         // ones and then store their SC as before ___________________________________
00201         for (int jsc =0; jsc< cleanSize; ++jsc) {
00202                 // check whether this index is already in the common collection
00203                 bool takenAlready = false;
00204                 for (int j=0; j< inCleanSize; ++j) {
00205                         if (jsc == inCleanInd[j]) { takenAlready = true ;break;}
00206                 }
00207                 if (takenAlready) {
00208                         inCleanOnlyInd.push_back(0);
00209                         scCleanSeedDetId.push_back(DetId(0));
00210                         continue;
00211                 }
00212                 inCleanOnlyInd.push_back(1);
00213                 reco::SuperClusterRef cscRef( pCleanSC, jsc );
00214                 scCleanSeedDetId.push_back(cscRef->seed()->seed());
00215                 for (reco::CaloCluster_iterator bciter = cscRef->clustersBegin(); bciter != cscRef->clustersEnd(); ++bciter) {
00216                         // the basic clusters
00217                         basicClusters.push_back(**bciter);
00218                         basicClusterOwner.push_back( std::make_pair(jsc,1) );
00219                 }
00220         } // end loop over clean SC _________________________________________________
00221         //
00222         //
00223         // at this point we have the basic cluster collection ready
00224         // 
00225         int bcSize = (int) basicClusters.size();
00226         int bcSizeUncleanOnly = (int) basicClustersUncleanOnly.size();
00227 
00228         LogDebug("UnifiedSC") << "Found cleaned SC: " << cleanSize 
00229                               <<  " uncleaned SC: "   << uncleanSize ;
00230         //
00231         // export the clusters to the event from the clean clusters
00232         std::auto_ptr< reco::BasicClusterCollection> 
00233                 basicClusters_p(new reco::BasicClusterCollection);
00234         basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
00235         edm::OrphanHandle<reco::BasicClusterCollection> bccHandle =  
00236                 evt.put(basicClusters_p, bcCollection_);
00237         if (!(bccHandle.isValid())) {
00238               
00239           edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters!";
00240           return;
00241         }
00242         reco::BasicClusterCollection basicClustersProd = *bccHandle;
00243 
00244         LogDebug("UnifiedSC")<< "Got the BasicClusters from the event again" ;
00245         //
00246         // export the clusters to the event: from the unclean only clusters
00247         std::auto_ptr< reco::BasicClusterCollection> 
00248                 basicClustersUncleanOnly_p(new reco::BasicClusterCollection);
00249         basicClustersUncleanOnly_p->assign(basicClustersUncleanOnly.begin(), 
00250                                            basicClustersUncleanOnly.end());
00251         edm::OrphanHandle<reco::BasicClusterCollection> bccHandleUncleanOnly =  
00252                 evt.put(basicClustersUncleanOnly_p, bcCollectionUncleanOnly_);
00253         if (!(bccHandleUncleanOnly.isValid())) {
00254 
00255           edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters (Unclean Only)!" ;
00256           return;
00257         }
00258         reco::BasicClusterCollection basicClustersUncleanOnlyProd = *bccHandleUncleanOnly;
00259         LogDebug("UnifiedSC")<< "Got the BasicClusters from the event again  (Unclean Only)" ;
00260         //
00261 
00262         // now we can build the SC collection
00263         //
00264         // start again from the unclean collection
00265         // all the unclean SC will become members of the new collection
00266         // with different algoIDs ___________________________________________________
00267         for (int isc=0; isc< uncleanSize; ++isc) {
00268                 reco::CaloClusterPtrVector clusterPtrVector;
00269                 // the seed is the basic cluster with the highest energy
00270                 reco::CaloClusterPtr seed; 
00271                 if (inUncleanOnlyInd[isc] == 1) { // unclean SC Unique in Unclean
00272                         for (int jbc=0; jbc< bcSizeUncleanOnly; ++jbc) {
00273                                 std::pair<int, int> theBcOwner = basicClusterOwnerUncleanOnly[jbc];
00274                                 if (theBcOwner.first == isc && theBcOwner.second == 0) {
00275                                         reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandleUncleanOnly,jbc);
00276                                         clusterPtrVector.push_back(currentClu);
00277                                         if (scUncleanSeedDetId[isc] == currentClu->seed()) {
00278                                                 seed = currentClu;
00279                                         }
00280                                 }
00281                         }
00282 
00283                 }
00284                 else { // unclean SC common in clean and unclean
00285                         for (int jbc=0; jbc< bcSize; ++jbc) {
00286                                 std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
00287                                 if (theBcOwner.first == isc && theBcOwner.second == 0) {
00288                                         reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00289                                         clusterPtrVector.push_back(currentClu);
00290                                         if (scUncleanSeedDetId[isc] == currentClu->seed()) {
00291                                                 seed = currentClu;
00292                                         }
00293                                 }
00294                         }
00295                 }
00296                 //std::cout << "before getting the uncl" << std::endl;
00297                 reco::SuperClusterRef unscRef( pUncleanSC, isc ); 
00298                 reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), 
00299                                          seed, clusterPtrVector );
00300                 // now set the algoID for this SC again
00301                 if (inUncleanOnlyInd[isc] == 1) {
00302                         // set up the quality to unclean only .............
00303                         newSC.setFlags(reco::CaloCluster::uncleanOnly);
00304                         superClustersUncleanOnly.push_back(newSC);
00305                 }
00306                 else {
00307                         // set up the quality to common  .............
00308                         newSC.setFlags(reco::CaloCluster::common);
00309                         superClusters.push_back(newSC);
00310                 }
00311                 // now you can store your SC
00312 
00313         } // end loop over unclean SC _______________________________________________
00314         //  flags numbering scheme
00315         //  flags =   0 = cleanedOnly     cluster is only in the cleaned collection
00316         //  flags = 100 = common          cluster is common in both collections
00317         //  flags = 200 = uncleanedOnly   cluster is only in the uncleaned collection
00318 
00319         // now loop over the clean SC and do the same but now you have to avoid the
00320         // the duplicated ones ______________________________________________________
00321         for (int jsc=0; jsc< cleanSize; ++jsc) {
00322                 //std::cout << "working in cl #" << jsc << std::endl;
00323                 // check that the SC is not in the unclean collection
00324                 if (inCleanOnlyInd[jsc] == 0) continue;
00325                 reco::CaloClusterPtrVector clusterPtrVector;
00326                 // the seed is the basic cluster with the highest energy
00327                 reco::CaloClusterPtr seed; 
00328                 for (int jbc=0; jbc< bcSize; ++jbc) {
00329                         std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
00330                         if (theBcOwner.first == jsc && theBcOwner.second == 1) {
00331                                 reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00332                                 clusterPtrVector.push_back(currentClu);
00333                                 if (scCleanSeedDetId[jsc] == currentClu->seed()) {
00334                                         seed = currentClu;
00335                                 }
00336                         }
00337                 }
00338                 reco::SuperClusterRef cscRef( pCleanSC, jsc ); 
00339                 reco::SuperCluster newSC(cscRef->energy(), cscRef->position(),
00340                                          seed, clusterPtrVector );
00341                 newSC.setFlags(reco::CaloCluster::cleanOnly);
00342 
00343                 // add it to the collection:
00344                 superClusters.push_back(newSC);
00345 
00346         } // end loop over clean SC _________________________________________________
00347 
00348 
00349         LogDebug("UnifiedSC")<< "New SC collection was created";
00350 
00351         std::auto_ptr< reco::SuperClusterCollection> 
00352                 superClusters_p(new reco::SuperClusterCollection);
00353         superClusters_p->assign(superClusters.begin(), superClusters.end());
00354 
00355         evt.put(superClusters_p, scCollection_);
00356         
00357         LogDebug("UnifiedSC") << "Clusters (Basic/Super) added to the Event! :-)";
00358 
00359         std::auto_ptr< reco::SuperClusterCollection> 
00360                 superClustersUncleanOnly_p(new reco::SuperClusterCollection);
00361         superClustersUncleanOnly_p->assign(superClustersUncleanOnly.begin(), 
00362                                            superClustersUncleanOnly.end());
00363 
00364         evt.put(superClustersUncleanOnly_p, scCollectionUncleanOnly_);
00365 
00366         // ----- debugging ----------
00367         // print the new collection SC quantities
00368 
00369         // print out the clean collection SC
00370         LogDebug("UnifiedSC") << "Clean Collection SC ";
00371         for (int i=0; i < cleanSize; ++i) {
00372           reco::SuperClusterRef cscRef( pCleanSC, i );
00373           LogDebug("UnifiedSC") << " >>> clean    #" << i << "; Energy: " << cscRef->energy()
00374                                 << " eta: " << cscRef->eta() 
00375                                 << " sc seed detid: " << cscRef->seed()->seed().rawId()
00376                                 << std::endl;
00377                 }
00378                 // the unclean SC
00379                 LogDebug("UnifiedSC") << "Unclean Collection SC ";
00380                 for (int i=0; i < uncleanSize; ++i) {
00381                         reco::SuperClusterRef uscRef( pUncleanSC, i );
00382                         LogDebug("UnifiedSC") << " >>> unclean  #" << i << "; Energy: " << uscRef->energy()
00383                                 << " eta: " << uscRef->eta() 
00384                                 << " sc seed detid: " << uscRef->seed()->seed().rawId();
00385                 }
00386                 // the new collection
00387                 LogDebug("UnifiedSC")<< "The new SC clean collection with size "<< superClusters.size()  << std::endl;
00388 
00389                 int new_unclean = 0, new_clean=0;
00390                 for (int i=0; i < (int) superClusters.size(); ++i) {
00391                         const reco::SuperCluster & nsc = superClusters[i];
00392                         LogDebug("UnifiedSC") << "SC was got" << std::endl
00393                                 << " ---> energy: " << nsc.energy() << std::endl
00394                                 << " ---> eta: " << nsc.eta() << std::endl
00395                                 << " ---> inClean: " << nsc.isInClean() << std::endl
00396                                 << " ---> id: " << nsc.seed()->seed().rawId() << std::endl
00397                                 << " >>> newSC    #" << i << "; Energy: " << nsc.energy()
00398                                 << " eta: " << nsc.eta()  << " isClean=" 
00399                                 << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
00400                                 << " sc seed detid: " << nsc.seed()->seed().rawId();
00401 
00402                         if (nsc.isInUnclean()) ++new_unclean;
00403                         if (nsc.isInClean()) ++new_clean;
00404                 }
00405                 LogDebug("UnifiedSC")<< "The new SC unclean only collection with size "<< superClustersUncleanOnly.size();
00406                 for (int i=0; i < (int) superClustersUncleanOnly.size(); ++i) {
00407                         const reco::SuperCluster nsc = superClustersUncleanOnly[i];
00408                         LogDebug ("UnifiedSC") << " >>> newSC    #" << i << "; Energy: " << nsc.energy()
00409                                 << " eta: " << nsc.eta()  << " isClean=" 
00410                                 << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
00411                                 << " sc seed detid: " << nsc.seed()->seed().rawId();
00412                         if (nsc.isInUnclean()) ++new_unclean;
00413                         if (nsc.isInClean()) ++new_clean;      
00414                 }
00415                 if ( (new_unclean != uncleanSize) || (new_clean != cleanSize) ) {
00416                         LogDebug("UnifiedSC") << ">>>>!!!!!! MISMATCH: new unclean/ old unclean= " 
00417                                 << new_unclean << " / " << uncleanSize 
00418                                 << ", new clean/ old clean" << new_clean << " / " << cleanSize;
00419                 }
00420         
00421 }