Go to the documentation of this file.00001
00002 #include <iostream>
00003 #include <vector>
00004 #include <memory>
00005
00006
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
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
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
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053 UnifiedSCCollectionProducer::UnifiedSCCollectionProducer(const edm::ParameterSet& ps)
00054 {
00055
00056
00057
00058 cleanBcCollection_ = ps.getParameter<edm::InputTag>("cleanBcCollection");
00059 cleanScCollection_ = ps.getParameter<edm::InputTag>("cleanScCollection");
00060
00061 uncleanBcCollection_ = ps.getParameter<edm::InputTag>("uncleanBcCollection");
00062 uncleanScCollection_ = ps.getParameter<edm::InputTag>("uncleanScCollection");
00063
00064
00065
00066 bcCollection_ = ps.getParameter<std::string>("bcCollection");
00067 scCollection_ = ps.getParameter<std::string>("scCollection");
00068
00069 bcCollectionUncleanOnly_ = ps.getParameter<std::string>("bcCollectionUncleanOnly");
00070 scCollectionUncleanOnly_ = ps.getParameter<std::string>("scCollectionUncleanOnly");
00071
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
00089
00090
00091
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
00115 reco::BasicClusterCollection basicClusters;
00116 reco::SuperClusterCollection superClusters;
00117
00118 reco::BasicClusterCollection basicClustersUncleanOnly;
00119 reco::SuperClusterCollection superClustersUncleanOnly;
00120
00121
00122
00123
00124
00125
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
00136 std::vector<int> inUncleanOnlyInd;
00137 std::vector<int> inCleanInd;
00138 std::vector<int> inCleanOnlyInd;
00139 std::vector<DetId> scUncleanSeedDetId;
00140 std::vector<DetId> scCleanSeedDetId;
00141
00142
00143 std::vector< std::pair<int, int> > basicClusterOwner;
00144 std::vector< std::pair<int, int> > basicClusterOwnerUncleanOnly;
00145
00146 std::vector<int> uncleanBasicClusterIsSeed;
00147
00148
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) {
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
00161
00162
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) {
00169
00170 inUncleanOnlyInd.push_back(0);
00171 inCleanInd.push_back(jsc);
00172 scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00173
00174
00175 for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00176
00177 basicClusters.push_back(**bciter);
00178
00179 basicClusterOwner.push_back( std::make_pair(isc,0) );
00180 }
00181 break;
00182 }
00183 }
00184 if (not foundTheSame) {
00185
00186 inUncleanOnlyInd.push_back(1);
00187 scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00188
00189 for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00190
00191 basicClustersUncleanOnly.push_back(**bciter);
00192 basicClusterOwnerUncleanOnly.push_back( std::make_pair(isc,0) );
00193 }
00194 }
00195 }
00196
00197 int inCleanSize = inCleanInd.size();
00198
00199
00200
00201 for (int jsc =0; jsc< cleanSize; ++jsc) {
00202
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
00217 basicClusters.push_back(**bciter);
00218 basicClusterOwner.push_back( std::make_pair(jsc,1) );
00219 }
00220 }
00221
00222
00223
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
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
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
00263
00264
00265
00266
00267 for (int isc=0; isc< uncleanSize; ++isc) {
00268 reco::CaloClusterPtrVector clusterPtrVector;
00269
00270 reco::CaloClusterPtr seed;
00271 if (inUncleanOnlyInd[isc] == 1) {
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 {
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
00297 reco::SuperClusterRef unscRef( pUncleanSC, isc );
00298 reco::SuperCluster newSC(unscRef->energy(), unscRef->position(),
00299 seed, clusterPtrVector );
00300
00301 if (inUncleanOnlyInd[isc] == 1) {
00302
00303 newSC.setFlags(reco::CaloCluster::uncleanOnly);
00304 superClustersUncleanOnly.push_back(newSC);
00305 }
00306 else {
00307
00308 newSC.setFlags(reco::CaloCluster::common);
00309 superClusters.push_back(newSC);
00310 }
00311
00312
00313 }
00314
00315
00316
00317
00318
00319
00320
00321 for (int jsc=0; jsc< cleanSize; ++jsc) {
00322
00323
00324 if (inCleanOnlyInd[jsc] == 0) continue;
00325 reco::CaloClusterPtrVector clusterPtrVector;
00326
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
00344 superClusters.push_back(newSC);
00345
00346 }
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
00367
00368
00369
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
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
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 }