00001 #include "RecoParticleFlow/PFClusterProducer/interface/PFClusterAlgo.h"
00002 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
00003 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00004 #include "Math/GenVector/VectorUtil.h"
00005 #include "TFile.h"
00006 #include "TH2F.h"
00007 #include "TROOT.h"
00008
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010
00011 #include <stdexcept>
00012 #include <string>
00013 #include <sstream>
00014
00015 using namespace std;
00016
00017 unsigned PFClusterAlgo::prodNum_ = 1;
00018
00019
00020
00021
00022 PFClusterAlgo::PFClusterAlgo() :
00023 pfClusters_( new vector<reco::PFCluster> ),
00024 pfRecHitsCleaned_( new vector<reco::PFRecHit> ),
00025 threshBarrel_(0.),
00026 threshPtBarrel_(0.),
00027 threshSeedBarrel_(0.2),
00028 threshPtSeedBarrel_(0.),
00029 threshEndcap_(0.),
00030 threshPtEndcap_(0.),
00031 threshSeedEndcap_(0.6),
00032 threshPtSeedEndcap_(0.),
00033 threshCleanBarrel_(1E5),
00034 minS4S1Barrel_(0.),
00035 threshDoubleSpikeBarrel_(1E9),
00036 minS6S2DoubleSpikeBarrel_(-1.),
00037 threshCleanEndcap_(1E5),
00038 minS4S1Endcap_(0.),
00039 threshDoubleSpikeEndcap_(1E9),
00040 minS6S2DoubleSpikeEndcap_(-1.),
00041 nNeighbours_(4),
00042 posCalcNCrystal_(-1),
00043 posCalcP1_(-1),
00044 showerSigma_(5),
00045 useCornerCells_(false),
00046 cleanRBXandHPDs_(false),
00047 debug_(false)
00048 {
00049 file_ = 0;
00050 hBNeighbour = 0;
00051 hENeighbour = 0;
00052 }
00053
00054 void
00055 PFClusterAlgo::write() {
00056
00057 if ( file_ ) {
00058 file_->Write();
00059 cout << "Benchmark output written to file " << file_->GetName() << endl;
00060 file_->Close();
00061 }
00062
00063 }
00064
00065
00066 void PFClusterAlgo::doClustering( const PFRecHitHandle& rechitsHandle ) {
00067 rechitsHandle_ = rechitsHandle;
00068 doClustering( *rechitsHandle );
00069 }
00070
00071 void PFClusterAlgo::doClustering( const reco::PFRecHitCollection& rechits ) {
00072
00073
00074 if(pfClusters_.get() ) pfClusters_->clear();
00075 else
00076 pfClusters_.reset( new std::vector<reco::PFCluster> );
00077
00078 if(pfRecHitsCleaned_.get() ) pfRecHitsCleaned_->clear();
00079 else
00080 pfRecHitsCleaned_.reset( new std::vector<reco::PFRecHit> );
00081
00082
00083 eRecHits_.clear();
00084
00085 bool initMask = false;
00086 if( mask_.size() != rechits.size() ) {
00087 initMask = true;
00088 mask_.clear();
00089 mask_.reserve( rechits.size() );
00090
00091 if( ! mask_.empty() )
00092 cerr<<"PClusterAlgo::doClustering: map size should be "<<mask_.size()
00093 <<". Will be reinitialized."<<endl;
00094 }
00095
00096 color_.clear();
00097 color_.reserve( rechits.size() );
00098 seedStates_.clear();
00099 seedStates_.reserve( rechits.size() );
00100 usedInTopo_.clear();
00101 usedInTopo_.reserve( rechits.size() );
00102
00103 for ( unsigned i = 0; i < rechits.size(); i++ ) {
00104 eRecHits_.insert( make_pair( rechit(i, rechits).energy(), i) );
00105 if(initMask) mask_.push_back( true );
00106 color_.push_back( 0 );
00107 seedStates_.push_back( UNKNOWN );
00108 usedInTopo_.push_back( false );
00109 }
00110
00111 if ( cleanRBXandHPDs_ ) cleanRBXAndHPD( rechits);
00112
00113
00114 findSeeds( rechits );
00115
00116
00117 buildTopoClusters( rechits );
00118
00119
00120 for(unsigned i=0; i<topoClusters_.size(); i++) {
00121
00122 const std::vector< unsigned >& topocluster = topoClusters_[i];
00123 buildPFClusters( topocluster, rechits );
00124 }
00125 }
00126
00127
00128 void PFClusterAlgo::setMask( const std::vector<bool>& mask ) {
00129 mask_ = mask;
00130 }
00131
00132
00133
00134
00135 double PFClusterAlgo::parameter( Parameter paramtype,
00136 PFLayer::Layer layer,
00137 unsigned iCoeff ) const {
00138
00139
00140 double value = 0;
00141
00142 switch( layer ) {
00143 case PFLayer::ECAL_BARREL:
00144 case PFLayer::HCAL_BARREL1:
00145 case PFLayer::HCAL_BARREL2:
00146
00147 switch(paramtype) {
00148 case THRESH:
00149 value = threshBarrel_;
00150 break;
00151 case SEED_THRESH:
00152 value = threshSeedBarrel_;
00153 break;
00154 case PT_THRESH:
00155 value = threshPtBarrel_;
00156 break;
00157 case SEED_PT_THRESH:
00158 value = threshPtSeedBarrel_;
00159 break;
00160 case CLEAN_THRESH:
00161 value = threshCleanBarrel_;
00162 break;
00163 case CLEAN_S4S1:
00164 value = minS4S1Barrel_[iCoeff];
00165 break;
00166 case DOUBLESPIKE_THRESH:
00167 value = threshDoubleSpikeBarrel_;
00168 break;
00169 case DOUBLESPIKE_S6S2:
00170 value = minS6S2DoubleSpikeBarrel_;
00171 break;
00172 default:
00173 cerr<<"PFClusterAlgo::parameter : unknown parameter type "
00174 <<paramtype<<endl;
00175 assert(0);
00176 }
00177 break;
00178 case PFLayer::ECAL_ENDCAP:
00179 case PFLayer::HCAL_ENDCAP:
00180 case PFLayer::PS1:
00181 case PFLayer::PS2:
00182 case PFLayer::HF_EM:
00183 case PFLayer::HF_HAD:
00184
00185 switch(paramtype) {
00186 case THRESH:
00187 value = threshEndcap_;
00188 break;
00189 case SEED_THRESH:
00190 value = threshSeedEndcap_;
00191 break;
00192 case PT_THRESH:
00193 value = threshPtEndcap_;
00194 break;
00195 case SEED_PT_THRESH:
00196 value = threshPtSeedEndcap_;
00197 break;
00198 case CLEAN_THRESH:
00199 value = threshCleanEndcap_;
00200 break;
00201 case CLEAN_S4S1:
00202 value = minS4S1Endcap_[iCoeff];
00203 break;
00204 case DOUBLESPIKE_THRESH:
00205 value = threshDoubleSpikeEndcap_;
00206 break;
00207 case DOUBLESPIKE_S6S2:
00208 value = minS6S2DoubleSpikeEndcap_;
00209 break;
00210 default:
00211 cerr<<"PFClusterAlgo::parameter : unknown parameter type "
00212 <<paramtype<<endl;
00213 assert(0);
00214 }
00215 break;
00216 default:
00217 cerr<<"PFClusterAlgo::parameter : unknown layer "<<layer<<endl;
00218 assert(0);
00219 break;
00220 }
00221
00222 return value;
00223 }
00224
00225
00226 void
00227 PFClusterAlgo::cleanRBXAndHPD( const reco::PFRecHitCollection& rechits ) {
00228
00229 std::map< int, std::vector<unsigned> > hpds;
00230 std::map< int, std::vector<unsigned> > rbxs;
00231
00232 for(EH ih = eRecHits_.begin(); ih != eRecHits_.end(); ih++ ) {
00233
00234 unsigned rhi = ih->second;
00235
00236 if(! masked(rhi) ) continue;
00237
00238 const reco::PFRecHit& rhit = rechit( rhi, rechits);
00239
00240 int layer = rhit.layer();
00241 if ( layer != PFLayer::HCAL_BARREL1 &&
00242 layer != PFLayer::HCAL_ENDCAP ) break;
00243 HcalDetId theHcalDetId = HcalDetId(rhit.detId());
00244 int ieta = theHcalDetId.ieta();
00245 int iphi = theHcalDetId.iphi();
00246 int ihpd = ieta < 0 ?
00247 ( layer == PFLayer::HCAL_ENDCAP ? -(iphi+1)/2-100 : -iphi ) :
00248 ( layer == PFLayer::HCAL_ENDCAP ? (iphi+1)/2+100 : iphi ) ;
00249 hpds[ihpd].push_back(rhi);
00250 int irbx = ieta < 0 ?
00251 ( layer == PFLayer::HCAL_ENDCAP ? -(iphi+5)/4 - 20 : -(iphi+5)/4 ) :
00252 ( layer == PFLayer::HCAL_ENDCAP ? (iphi+5)/4 + 20 : (iphi+5)/4 ) ;
00253 if ( irbx == 19 ) irbx = 1;
00254 else if ( irbx == -19 ) irbx = -1;
00255 else if ( irbx == 39 ) irbx = 21;
00256 else if ( irbx == -39 ) irbx = -21;
00257 rbxs[irbx].push_back(rhi);
00258 }
00259
00260
00261 for ( std::map<int, std::vector<unsigned> >::iterator itrbx = rbxs.begin();
00262 itrbx != rbxs.end(); ++itrbx ) {
00263 if ( ( abs(itrbx->first)<20 && itrbx->second.size() > 30 ) ||
00264 ( abs(itrbx->first)>20 && itrbx->second.size() > 30 ) ) {
00265 const std::vector<unsigned>& rhits = itrbx->second;
00266 double totalEta = 0.;
00267 double totalEtaW = 0.;
00268 double totalPhi = 0.;
00269 double totalPhiW = 0.;
00270 double totalEta2 = 1E-9;
00271 double totalEta2W = 1E-9;
00272 double totalPhi2 = 1E-9;
00273 double totalPhi2W = 1E-9;
00274 double totalEnergy = 0.;
00275 double totalEnergy2 = 1E-9;
00276 unsigned nSeeds = rhits.size();
00277 unsigned nSeeds0 = rhits.size();
00278 std::map< int,std::vector<unsigned> > theHPDs;
00279 std::multimap< double,unsigned > theEnergies;
00280 for ( unsigned jh=0; jh < rhits.size(); ++jh ) {
00281 const reco::PFRecHit& hit = rechit(rhits[jh], rechits);
00282
00283 unsigned nN = 0;
00284 bool isASeed = true;
00285 const vector<unsigned>& neighbours4 = *(& hit.neighbours4());
00286 for(unsigned in=0; in<neighbours4.size(); in++) {
00287 const reco::PFRecHit& neighbour = rechit( neighbours4[in], rechits );
00288
00289 if( neighbour.energy() > hit.energy() ) {
00290 --nSeeds;
00291 --nSeeds0;
00292 isASeed = false;
00293 break;
00294 } else {
00295 if ( neighbour.energy() > 0.4 ) ++nN;
00296 }
00297 }
00298 if ( isASeed && !nN ) --nSeeds0;
00299
00300 HcalDetId theHcalDetId = HcalDetId(hit.detId());
00301
00302 int iphi = theHcalDetId.iphi();
00303
00304 if ( hit.layer() == PFLayer::HCAL_BARREL1 )
00305 theHPDs[iphi].push_back(rhits[jh]);
00306 else
00307 theHPDs[(iphi-1)/2].push_back(rhits[jh]);
00308 theEnergies.insert(std::pair<double,unsigned>(hit.energy(),rhits[jh]));
00309 totalEnergy += hit.energy();
00310 totalPhi += fabs(hit.position().phi());
00311 totalPhiW += hit.energy()*fabs(hit.position().phi());
00312 totalEta += hit.position().eta();
00313 totalEtaW += hit.energy()*hit.position().eta();
00314 totalEnergy2 += hit.energy()*hit.energy();
00315 totalPhi2 += hit.position().phi()*hit.position().phi();
00316 totalPhi2W += hit.energy()*hit.position().phi()*hit.position().phi();
00317 totalEta2 += hit.position().eta()*hit.position().eta();
00318 totalEta2W += hit.energy()*hit.position().eta()*hit.position().eta();
00319 }
00320
00321 totalPhi /= rhits.size();
00322 totalEta /= rhits.size();
00323 totalPhiW /= totalEnergy;
00324 totalEtaW /= totalEnergy;
00325 totalPhi2 /= rhits.size();
00326 totalEta2 /= rhits.size();
00327 totalPhi2W /= totalEnergy;
00328 totalEta2W /= totalEnergy;
00329 totalPhi2 = std::sqrt(totalPhi2 - totalPhi*totalPhi);
00330 totalEta2 = std::sqrt(totalEta2 - totalEta*totalEta);
00331 totalPhi2W = std::sqrt(totalPhi2W - totalPhiW*totalPhiW);
00332 totalEta2W = std::sqrt(totalEta2W - totalEtaW*totalEtaW);
00333 totalEnergy /= rhits.size();
00334 totalEnergy2 /= rhits.size();
00335 totalEnergy2 = std::sqrt(totalEnergy2 - totalEnergy*totalEnergy);
00336
00337 if ( nSeeds0 > 6 ) {
00338 unsigned nHPD15 = 0;
00339 for ( std::map<int, std::vector<unsigned> >::iterator itHPD = theHPDs.begin();
00340 itHPD != theHPDs.end(); ++itHPD ) {
00341 int hpdN = itHPD->first;
00342 const std::vector<unsigned>& hpdHits = itHPD->second;
00343 if ( ( abs(hpdN) < 100 && hpdHits.size() > 14 ) ||
00344 ( abs(hpdN) > 100 && hpdHits.size() > 14 ) ) ++nHPD15;
00345 }
00346 if ( nHPD15 > 1 ) {
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362 std::multimap<double, unsigned >::iterator ntEn = theEnergies.end();
00363 --ntEn;
00364 unsigned nn = 0;
00365 double threshold = 1.;
00366 for ( std::multimap<double, unsigned >::iterator itEn = theEnergies.begin();
00367 itEn != theEnergies.end(); ++itEn ) {
00368 ++nn;
00369 if ( nn < 5 ) {
00370 mask_[itEn->second] = false;
00371 } else if ( nn == 5 ) {
00372 threshold = itEn->first * 5.;
00373 mask_[itEn->second] = false;
00374 } else {
00375 if ( itEn->first < threshold ) mask_[itEn->second] = false;
00376 }
00377 if ( !masked(itEn->second) ) {
00378 reco::PFRecHit theCleanedHit(rechit(itEn->second, rechits));
00379
00380 pfRecHitsCleaned_->push_back(theCleanedHit);
00381 }
00382
00383
00384
00385
00386
00387
00388
00389
00390 }
00391 }
00392 }
00393 }
00394 }
00395
00396
00397 std::map<int, std::vector<unsigned> >::iterator neighbour1;
00398 std::map<int, std::vector<unsigned> >::iterator neighbour2;
00399 std::map<int, std::vector<unsigned> >::iterator neighbour0;
00400 std::map<int, std::vector<unsigned> >::iterator neighbour3;
00401 unsigned size1 = 0;
00402 unsigned size2 = 0;
00403 for ( std::map<int, std::vector<unsigned> >::iterator ithpd = hpds.begin();
00404 ithpd != hpds.end(); ++ithpd ) {
00405
00406 const std::vector<unsigned>& rhits = ithpd->second;
00407 std::multimap< double,unsigned > theEnergies;
00408 double totalEnergy = 0.;
00409 double totalEnergy2 = 1E-9;
00410 for ( unsigned jh=0; jh < rhits.size(); ++jh ) {
00411 const reco::PFRecHit& hit = rechit(rhits[jh], rechits);
00412 totalEnergy += hit.energy();
00413 totalEnergy2 += hit.energy()*hit.energy();
00414 theEnergies.insert(std::pair<double,unsigned>(hit.energy(),rhits[jh]));
00415 }
00416 totalEnergy /= rhits.size();
00417 totalEnergy2 /= rhits.size();
00418 totalEnergy2 = std::sqrt(totalEnergy2 - totalEnergy*totalEnergy);
00419
00420 if ( ithpd->first == 1 ) neighbour1 = hpds.find(72);
00421 else if ( ithpd->first == -1 ) neighbour1 = hpds.find(-72);
00422 else if ( ithpd->first == 101 ) neighbour1 = hpds.find(136);
00423 else if ( ithpd->first == -101 ) neighbour1 = hpds.find(-136);
00424 else neighbour1 = ithpd->first > 0 ? hpds.find(ithpd->first-1) : hpds.find(ithpd->first+1) ;
00425
00426 if ( ithpd->first == 72 ) neighbour2 = hpds.find(1);
00427 else if ( ithpd->first == -72 ) neighbour2 = hpds.find(-1);
00428 else if ( ithpd->first == 136 ) neighbour2 = hpds.find(101);
00429 else if ( ithpd->first == -136 ) neighbour2 = hpds.find(-101);
00430 else neighbour2 = ithpd->first > 0 ? hpds.find(ithpd->first+1) : hpds.find(ithpd->first-1) ;
00431
00432 if ( neighbour1 != hpds.end() ) {
00433 if ( neighbour1->first == 1 ) neighbour0 = hpds.find(72);
00434 else if ( neighbour1->first == -1 ) neighbour0 = hpds.find(-72);
00435 else if ( neighbour1->first == 101 ) neighbour0 = hpds.find(136);
00436 else if ( neighbour1->first == -101 ) neighbour0 = hpds.find(-136);
00437 else neighbour0 = neighbour1->first > 0 ? hpds.find(neighbour1->first-1) : hpds.find(neighbour1->first+1) ;
00438 }
00439
00440 if ( neighbour2 != hpds.end() ) {
00441 if ( neighbour2->first == 72 ) neighbour3 = hpds.find(1);
00442 else if ( neighbour2->first == -72 ) neighbour3 = hpds.find(-1);
00443 else if ( neighbour2->first == 136 ) neighbour3 = hpds.find(101);
00444 else if ( neighbour2->first == -136 ) neighbour3 = hpds.find(-101);
00445 else neighbour3 = neighbour2->first > 0 ? hpds.find(neighbour2->first+1) : hpds.find(neighbour2->first-1) ;
00446 }
00447
00448 size1 = neighbour1 != hpds.end() ? neighbour1->second.size() : 0;
00449 size2 = neighbour2 != hpds.end() ? neighbour2->second.size() : 0;
00450
00451
00452 if ( size1 > 10 ) {
00453 if ( ( abs(neighbour1->first) > 100 && neighbour1->second.size() > 15 ) ||
00454 ( abs(neighbour1->first) < 100 && neighbour1->second.size() > 12 ) )
00455 size1 = neighbour0 != hpds.end() ? neighbour0->second.size() : 0;
00456 }
00457 if ( size2 > 10 ) {
00458 if ( ( abs(neighbour2->first) > 100 && neighbour2->second.size() > 15 ) ||
00459 ( abs(neighbour2->first) < 100 && neighbour2->second.size() > 12 ) )
00460 size2 = neighbour3 != hpds.end() ? neighbour3->second.size() : 0;
00461 }
00462
00463 if ( ( abs(ithpd->first) > 100 && ithpd->second.size() > 15 ) ||
00464 ( abs(ithpd->first) < 100 && ithpd->second.size() > 12 ) )
00465 if ( (float)(size1 + size2)/(float)ithpd->second.size() < 1.0 ) {
00466
00467
00468
00469
00470
00471
00472 std::multimap<double, unsigned >::iterator ntEn = theEnergies.end();
00473 --ntEn;
00474 unsigned nn = 0;
00475 double threshold = 1.;
00476 for ( std::multimap<double, unsigned >::iterator itEn = theEnergies.begin();
00477 itEn != theEnergies.end(); ++itEn ) {
00478 ++nn;
00479 if ( nn < 5 ) {
00480 mask_[itEn->second] = false;
00481 } else if ( nn == 5 ) {
00482 threshold = itEn->first * 2.5;
00483 mask_[itEn->second] = false;
00484 } else {
00485 if ( itEn->first < threshold ) mask_[itEn->second] = false;
00486 }
00487 if ( !masked(itEn->second) ) {
00488 reco::PFRecHit theCleanedHit(rechit(itEn->second, rechits));
00489
00490 pfRecHitsCleaned_->push_back(theCleanedHit);
00491 }
00492
00493
00494
00495
00496
00497
00498
00499
00500 }
00501 }
00502 }
00503
00504 }
00505
00506
00507 void PFClusterAlgo::findSeeds( const reco::PFRecHitCollection& rechits ) {
00508
00509 seeds_.clear();
00510
00511
00512 #ifdef PFLOW_DEBUG
00513 if(debug_)
00514 cout<<"PFClusterAlgo::findSeeds : start"<<endl;
00515 #endif
00516
00517
00518
00519 const vector<unsigned> noNeighbours(0, static_cast<unsigned>(0));
00520
00521
00522 for(EH ih = eRecHits_.begin(); ih != eRecHits_.end(); ih++ ) {
00523
00524 unsigned rhi = ih->second;
00525
00526 if(! masked(rhi) ) continue;
00527
00528
00529 double rhenergy = ih->first;
00530 const reco::PFRecHit& wannaBeSeed = rechit(rhi, rechits);
00531
00532 if( seedStates_[rhi] == NO ) continue;
00533
00534
00535
00536 int layer = wannaBeSeed.layer();
00537 double seedThresh = parameter( SEED_THRESH,
00538 static_cast<PFLayer::Layer>(layer) );
00539 double seedPtThresh = parameter( SEED_PT_THRESH,
00540 static_cast<PFLayer::Layer>(layer) );
00541 double cleanThresh = parameter( CLEAN_THRESH,
00542 static_cast<PFLayer::Layer>(layer) );
00543 double minS4S1_a = parameter( CLEAN_S4S1,
00544 static_cast<PFLayer::Layer>(layer), 0 );
00545 double minS4S1_b = parameter( CLEAN_S4S1,
00546 static_cast<PFLayer::Layer>(layer), 1 );
00547 double doubleSpikeThresh = parameter( DOUBLESPIKE_THRESH,
00548 static_cast<PFLayer::Layer>(layer) );
00549 double doubleSpikeS6S2 = parameter( DOUBLESPIKE_S6S2,
00550 static_cast<PFLayer::Layer>(layer) );
00551
00552
00553 #ifdef PFLOW_DEBUG
00554 if(debug_)
00555 cout<<"layer:"<<layer<<" seedThresh:"<<seedThresh<<endl;
00556 #endif
00557
00558
00559 if( rhenergy < seedThresh || (seedPtThresh>0. && wannaBeSeed.pt2() < seedPtThresh*seedPtThresh )) {
00560 seedStates_[rhi] = NO;
00561 continue;
00562 }
00563
00564
00565
00566 const vector<unsigned>* nbp;
00567 double tighterE = 1.0;
00568 double tighterF = 1.0;
00569
00570 switch ( layer ) {
00571 case PFLayer::ECAL_BARREL:
00572 case PFLayer::ECAL_ENDCAP:
00573 case PFLayer::HCAL_BARREL1:
00574 case PFLayer::HCAL_BARREL2:
00575 case PFLayer::HCAL_ENDCAP:
00576 tighterE = 2.0;
00577 tighterF = 3.0;
00578 case PFLayer::HF_EM:
00579 case PFLayer::HF_HAD:
00580 if( nNeighbours_ == 4 ) {
00581 nbp = & wannaBeSeed.neighbours4();
00582 }
00583 else if( nNeighbours_ == 8 ) {
00584 nbp = & wannaBeSeed.neighbours8();
00585 }
00586 else if( nNeighbours_ == 0 ) {
00587 nbp = & noNeighbours;
00588
00589
00590 }
00591 else {
00592 cerr<<"you're not allowed to set n neighbours to "
00593 <<nNeighbours_<<endl;
00594 assert(0);
00595 }
00596 break;
00597 case PFLayer::PS1:
00598 case PFLayer::PS2:
00599 nbp = & wannaBeSeed.neighbours4();
00600 break;
00601
00602 default:
00603 cerr<<"CellsEF::PhotonSeeds : unknown layer "<<layer<<endl;
00604 assert(0);
00605 }
00606
00607 const vector<unsigned>& neighbours = *nbp;
00608
00609
00610
00611
00612 seedStates_[rhi] = YES;
00613 for(unsigned in=0; in<neighbours.size(); in++) {
00614
00615 unsigned rhj = neighbours[in];
00616
00617 if ( !masked(rhj) ) continue;
00618 const reco::PFRecHit& neighbour = rechit( rhj, rechits );
00619
00620
00621 if( neighbour.energy() > wannaBeSeed.energy() ) {
00622 seedStates_[rhi] = NO;
00623 break;
00624 }
00625 }
00626
00627
00628
00629 if ( file_ || wannaBeSeed.energy() > cleanThresh ) {
00630
00631 const vector<unsigned>& neighbours4 = *(& wannaBeSeed.neighbours4());
00632
00633 double surroundingEnergy = wannaBeSeed.energyUp();
00634 double neighbourEnergy = 0.;
00635 double layerEnergy = 0.;
00636 for(unsigned in4=0; in4<neighbours4.size(); in4++) {
00637 unsigned rhj = neighbours4[in4];
00638
00639 if ( !masked(rhj) ) continue;
00640 const reco::PFRecHit& neighbour = rechit( rhj, rechits );
00641 surroundingEnergy += neighbour.energy() + neighbour.energyUp();
00642 neighbourEnergy += neighbour.energy() + neighbour.energyUp();
00643 layerEnergy += neighbour.energy();
00644 }
00645
00646
00647
00648
00649 double fraction1 = surroundingEnergy/wannaBeSeed.energy();
00650
00651
00652
00653
00654
00655
00656
00657
00658 if ( file_ ) {
00659 if ( layer == PFLayer::ECAL_BARREL || layer == PFLayer::HCAL_BARREL1 ) {
00660
00661
00662
00663
00664
00665
00666
00667 hBNeighbour->Fill(1./wannaBeSeed.energy(), fraction1);
00668 } else if ( fabs(wannaBeSeed.position().eta()) < 5.0 ) {
00669 if ( layer == PFLayer::ECAL_ENDCAP || layer == PFLayer::HCAL_ENDCAP ) {
00670 if ( wannaBeSeed.energy() > 1000 ) {
00671 if ( fabs(wannaBeSeed.position().eta()) < 2.8 )
00672 hENeighbour->Fill(1./wannaBeSeed.energy(), fraction1);
00673 }
00674 } else {
00675 hENeighbour->Fill(1./wannaBeSeed.energy(), fraction1);
00676 }
00677 }
00678 }
00679
00680
00681 if ( wannaBeSeed.energy() > cleanThresh ) {
00682 double f1Cut = minS4S1_a * log10(wannaBeSeed.energy()) + minS4S1_b;
00683 if ( fraction1 < f1Cut ) {
00684
00685 double eta = wannaBeSeed.position().eta();
00686 double phi = wannaBeSeed.position().phi();
00687 std::pair<double,double> dcr = dCrack(phi,eta);
00688 double dcrmin = std::min(dcr.first, dcr.second);
00689 eta = fabs(eta);
00690 if ( eta < 5.0 &&
00691 ( ( eta < 2.85 && dcrmin > 1. ) ||
00692 ( rhenergy > tighterE*cleanThresh &&
00693 fraction1 < f1Cut/tighterF ) )
00694 ) {
00695 seedStates_[rhi] = CLEAN;
00696 mask_[rhi] = false;
00697 reco::PFRecHit theCleanedHit(wannaBeSeed);
00698
00699 pfRecHitsCleaned_->push_back(theCleanedHit);
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710 }
00711 }
00712 }
00713 }
00714
00715
00716 if ( mask_[rhi] && wannaBeSeed.energy() > doubleSpikeThresh ) {
00717
00718 double surroundingEnergyi = 0.;
00719 double enmax = -999.;
00720 unsigned mostEnergeticNeighbour = 0;
00721 const vector<unsigned>& neighbours4i = *(& wannaBeSeed.neighbours4());
00722 for(unsigned in4=0; in4<neighbours4i.size(); in4++) {
00723 unsigned rhj = neighbours4i[in4];
00724 if ( !masked(rhj) ) continue;
00725 const reco::PFRecHit& neighbouri = rechit( rhj, rechits );
00726 surroundingEnergyi += neighbouri.energy();
00727 if ( neighbouri.energy() > enmax ) {
00728 enmax = neighbouri.energy();
00729 mostEnergeticNeighbour = rhj;
00730 }
00731 }
00732
00733 if ( enmax > 0. ) {
00734 unsigned rhj = mostEnergeticNeighbour;
00735 const reco::PFRecHit& neighbouri = rechit( rhj, rechits );
00736 double surroundingEnergyj = 0.;
00737
00738
00739 const vector<unsigned>& neighbours4j = *(& neighbouri.neighbours4());
00740 for(unsigned jn4=0; jn4<neighbours4j.size(); jn4++) {
00741 unsigned rhk = neighbours4j[jn4];
00742 const reco::PFRecHit& neighbourj = rechit( rhk, rechits );
00743 surroundingEnergyj += neighbourj.energy();
00744 }
00745
00746 double surroundingEnergyFraction =
00747 (surroundingEnergyi+surroundingEnergyj) / (wannaBeSeed.energy()+neighbouri.energy()) - 1.;
00748 if ( surroundingEnergyFraction < doubleSpikeS6S2 ) {
00749 double eta = wannaBeSeed.position().eta();
00750 double phi = wannaBeSeed.position().phi();
00751 std::pair<double,double> dcr = dCrack(phi,eta);
00752 double dcrmin = std::min(dcr.first, dcr.second);
00753 eta = fabs(eta);
00754 if ( ( eta < 5.0 && dcrmin > 1. ) ||
00755 ( wannaBeSeed.energy() > tighterE*doubleSpikeThresh &&
00756 surroundingEnergyFraction < doubleSpikeS6S2/tighterF ) ) {
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766 seedStates_[rhi] = CLEAN;
00767 mask_[rhi] = false;
00768 reco::PFRecHit theCleanedSeed(wannaBeSeed);
00769 pfRecHitsCleaned_->push_back(theCleanedSeed);
00770
00771 seedStates_[rhj] = CLEAN;
00772 mask_[rhj] = false;
00773 reco::PFRecHit theCleanedNeighbour(wannaBeSeed);
00774 pfRecHitsCleaned_->push_back(theCleanedNeighbour);
00775 }
00776 }
00777 } else {
00778
00779
00780
00781
00782 }
00783 }
00784
00785 if ( seedStates_[rhi] == YES ) {
00786
00787
00788 seeds_.push_back( rhi );
00789
00790
00791 paint(rhi, SEED);
00792
00793
00794 for(unsigned in=0; in<neighbours.size(); in++) {
00795 seedStates_[ neighbours[in] ] = NO;
00796 }
00797 }
00798 }
00799
00800 #ifdef PFLOW_DEBUG
00801 if(debug_)
00802 cout<<"PFClusterAlgo::findSeeds : done"<<endl;
00803 #endif
00804 }
00805
00806
00807
00808
00809 void PFClusterAlgo::buildTopoClusters( const reco::PFRecHitCollection& rechits ){
00810
00811 topoClusters_.clear();
00812
00813 #ifdef PFLOW_DEBUG
00814 if(debug_)
00815 cout<<"PFClusterAlgo::buildTopoClusters start"<<endl;
00816 #endif
00817
00818 for(unsigned is = 0; is<seeds_.size(); is++) {
00819
00820 unsigned rhi = seeds_[is];
00821
00822 if( !masked(rhi) ) continue;
00823
00824
00825
00826 if( usedInTopo_[rhi] ) {
00827 #ifdef PFLOW_DEBUG
00828 if(debug_)
00829 cout<<rhi<<" used"<<endl;
00830 #endif
00831 continue;
00832 }
00833
00834 vector< unsigned > topocluster;
00835 buildTopoCluster( topocluster, rhi, rechits );
00836
00837 if(topocluster.empty() ) continue;
00838
00839 topoClusters_.push_back( topocluster );
00840
00841 }
00842
00843 #ifdef PFLOW_DEBUG
00844 if(debug_)
00845 cout<<"PFClusterAlgo::buildTopoClusters done"<<endl;
00846 #endif
00847
00848 return;
00849 }
00850
00851
00852 void
00853 PFClusterAlgo::buildTopoCluster( vector< unsigned >& cluster,
00854 unsigned rhi,
00855 const reco::PFRecHitCollection& rechits ){
00856
00857
00858 #ifdef PFLOW_DEBUG
00859 if(debug_)
00860 cout<<"PFClusterAlgo::buildTopoCluster in"<<endl;
00861 #endif
00862
00863 const reco::PFRecHit& rh = rechit( rhi, rechits);
00864
00865 double e = rh.energy();
00866 int layer = rh.layer();
00867
00868 double thresh = parameter( THRESH,
00869 static_cast<PFLayer::Layer>(layer) );
00870 double ptThresh = parameter( PT_THRESH,
00871 static_cast<PFLayer::Layer>(layer) );
00872
00873
00874 if( e < thresh || (ptThresh > 0. && rh.pt2() < ptThresh*ptThresh) ) {
00875 #ifdef PFLOW_DEBUG
00876 if(debug_)
00877 cout<<"return : "<<e<<"<"<<thresh<<endl;
00878 #endif
00879 return;
00880 }
00881
00882
00883
00884 cluster.push_back( rhi );
00885
00886
00887 usedInTopo_[ rhi ] = true;
00888
00889
00890
00891
00892
00893 const std::vector< unsigned >& nbs
00894 = useCornerCells_ ? rh.neighbours8() : rh.neighbours4();
00895
00896 for(unsigned i=0; i<nbs.size(); i++) {
00897
00898
00899
00900
00901
00902
00903
00904
00905 if( usedInTopo_[ nbs[i] ] ) {
00906 #ifdef PFLOW_DEBUG
00907 if(debug_)
00908 cout<<rhi<<" used"<<endl;
00909 #endif
00910 continue;
00911 }
00912
00913 if( !masked(nbs[i]) ) continue;
00914 buildTopoCluster( cluster, nbs[i], rechits );
00915 }
00916 #ifdef PFLOW_DEBUG
00917 if(debug_)
00918 cout<<"PFClusterAlgo::buildTopoCluster out"<<endl;
00919 #endif
00920
00921 }
00922
00923
00924 void
00925 PFClusterAlgo::buildPFClusters( const std::vector< unsigned >& topocluster,
00926 const reco::PFRecHitCollection& rechits )
00927 {
00928
00929
00930
00931
00932
00933
00934
00935 vector<reco::PFCluster> curpfclusters;
00936 vector<reco::PFCluster> curpfclusterswodepthcor;
00937 vector< unsigned > seedsintopocluster;
00938
00939
00940 for(unsigned i=0; i<topocluster.size(); i++ ) {
00941
00942 unsigned rhi = topocluster[i];
00943
00944 if( seedStates_[rhi] == YES ) {
00945
00946 reco::PFCluster cluster;
00947 reco::PFCluster clusterwodepthcor;
00948
00949 double fraction = 1.0;
00950
00951 reco::PFRecHitRef recHitRef = createRecHitRef( rechits, rhi );
00952
00953 cluster.addRecHitFraction( reco::PFRecHitFraction( recHitRef,
00954 fraction ) );
00955
00956
00957
00958 calculateClusterPosition( cluster,
00959 clusterwodepthcor,
00960 true );
00961
00962
00963
00964 curpfclusters.push_back( cluster );
00965 curpfclusterswodepthcor.push_back( clusterwodepthcor );
00966 #ifdef PFLOW_DEBUG
00967 if(debug_) {
00968 cout << "PFClusterAlgo::buildPFClusters: seed "
00969 << rechit( rhi, rechits) <<endl;
00970 cout << "PFClusterAlgo::buildPFClusters: pfcluster initialized : "
00971 << cluster <<endl;
00972 }
00973 #endif
00974
00975
00976 seedsintopocluster.push_back( rhi );
00977
00978 }
00979 }
00980
00981
00982
00983
00984 int posCalcNCrystal = seedsintopocluster.size()>1 ? posCalcNCrystal_:-1;
00985 double ns2 = std::max(1.,(double)(seedsintopocluster.size())-1.);
00986 ns2 *= ns2;
00987
00988
00989
00990 unsigned iter = 0;
00991 unsigned niter = 50;
00992 double diff = ns2;
00993
00994
00995 vector<double> ener;
00996 vector<double> dist;
00997 vector<double> frac;
00998 vector<math::XYZVector> tmp;
00999 while ( iter++ < niter && diff > 1E-8*ns2 ) {
01000
01001
01002 ener.clear();
01003 tmp.clear();
01004
01005 for ( unsigned ic=0; ic<curpfclusters.size(); ic++ ) {
01006 ener.push_back( curpfclusters[ic].energy() );
01007
01008 math::XYZVector v;
01009 v = curpfclusters[ic].position();
01010
01011 tmp.push_back( v );
01012
01013 #ifdef PFLOW_DEBUG
01014 if(debug_) {
01015 cout<<"saving photon pos "<<ic<<" "<<curpfclusters[ic]<<endl;
01016 cout<<tmp[ic].X()<<" "<<tmp[ic].Y()<<" "<<tmp[ic].Z()<<endl;
01017 }
01018 #endif
01019
01020 curpfclusters[ic].reset();
01021 }
01022
01023
01024
01025 for( unsigned irh=0; irh<topocluster.size(); irh++ ) {
01026
01027 unsigned rhindex = topocluster[irh];
01028
01029 const reco::PFRecHit& rh = rechit( rhindex, rechits);
01030
01031
01032
01033 dist.clear();
01034 frac.clear();
01035 double fractot = 0.;
01036
01037 bool isaseed = isSeed(rhindex);
01038
01039 math::XYZVector cposxyzcell;
01040 cposxyzcell = rh.position();
01041
01042 #ifdef PFLOW_DEBUG
01043 if(debug_) {
01044 cout<<rh<<endl;
01045 cout<<"start loop on curpfclusters"<<endl;
01046 }
01047 #endif
01048
01049
01050 for ( unsigned ic=0; ic<tmp.size(); ic++) {
01051
01052 #ifdef PFLOW_DEBUG
01053 if(debug_) cout<<"pfcluster "<<ic<<endl;
01054 #endif
01055
01056 double frc=0.;
01057 bool seedexclusion=true;
01058
01059
01060
01061
01062 math::XYZVector cposxyzclust;
01063 cposxyzclust = curpfclusterswodepthcor[ic].position();
01064
01065 #ifdef PFLOW_DEBUG
01066 if(debug_) {
01067
01068 cout<<"CLUSTER "<<cposxyzclust.X()<<","
01069 <<cposxyzclust.Y()<<","
01070 <<cposxyzclust.Z()<<"\t\t"
01071 <<"CELL "<<cposxyzcell.X()<<","
01072 <<cposxyzcell.Y()<<","
01073 <<cposxyzcell.Z()<<endl;
01074 }
01075 #endif
01076
01077
01078
01079
01080 math::XYZVector deltav = cposxyzclust;
01081 deltav -= cposxyzcell;
01082 double d = deltav.R() / showerSigma_;
01083
01084
01085
01086
01087
01088 #ifdef PFLOW_DEBUG
01089 if( d > 10. && debug_ ) {
01090 paint(rhindex, SPECIAL);
01091 cout<<"PFClusterAlgo Warning: distance too large"<<d<<endl;
01092 }
01093 #endif
01094 dist.push_back( d );
01095
01096
01097 if( rhindex == seedsintopocluster[ic] && seedexclusion ) {
01098 frc = 1.;
01099 #ifdef PFLOW_DEBUG
01100 if(debug_) cout<<"this cell is a seed for the current photon"<<endl;
01101 #endif
01102 }
01103 else if( isaseed && seedexclusion ) {
01104 frc = 0.;
01105 #ifdef PFLOW_DEBUG
01106 if(debug_) cout<<"this cell is a seed for another photon"<<endl;
01107 #endif
01108 }
01109 else {
01110
01111
01112 frc = ener[ic] * exp ( - dist[ic]*dist[ic] / 2. );
01113
01114 #ifdef PFLOW_DEBUG
01115 if(debug_) {
01116 cout<<"dist["<<ic<<"] "<<dist[ic]
01117 <<", sigma="<<sigma
01118 <<", frc="<<frc<<endl;
01119 }
01120 #endif
01121
01122 }
01123 fractot += frc;
01124 frac.push_back(frc);
01125 }
01126
01127
01128 #ifdef PFLOW_DEBUG
01129 if(debug_) cout<<"start add cell"<<endl;
01130 #endif
01131 for ( unsigned ic=0; ic<tmp.size(); ++ic ) {
01132 #ifdef PFLOW_DEBUG
01133 if(debug_)
01134 cout<<" frac["<<ic<<"] "<<frac[ic]<<" "<<fractot<<" "<<rh<<endl;
01135 #endif
01136
01137 if( fractot )
01138 frac[ic] /= fractot;
01139 else {
01140 #ifdef PFLOW_DEBUG
01141 if( debug_ ) {
01142 int layer = rh.layer();
01143 cerr<<"fractot = 0 ! "<<layer<<endl;
01144
01145 for( unsigned trh=0; trh<topocluster.size(); trh++ ) {
01146 unsigned tindex = topocluster[trh];
01147 const reco::PFRecHit& rh = rechit( tindex, rechits);
01148 cout<<rh<<endl;
01149 }
01150
01151
01152 }
01153 #endif
01154
01155 continue;
01156 }
01157
01158
01159
01160
01161
01162
01163
01164
01165
01166
01167
01168
01169 if ( dist[ic] < 10. || frac[ic] > 0.99999 ) {
01170
01171 reco::PFRecHitRef recHitRef = createRecHitRef( rechits, rhindex );
01172 reco::PFRecHitFraction rhf( recHitRef,frac[ic] );
01173 curpfclusters[ic].addRecHitFraction( rhf );
01174 }
01175 }
01176
01177 }
01178
01179
01180
01181 diff = 0.;
01182 for ( unsigned ic=0; ic<tmp.size(); ++ic ) {
01183 calculateClusterPosition( curpfclusters[ic], curpfclusterswodepthcor[ic],
01184 true, posCalcNCrystal );
01185 #ifdef PFLOW_DEBUG
01186 if(debug_) cout<<"new iter "<<ic<<endl;
01187 if(debug_) cout<<curpfclusters[ic]<<endl;
01188 #endif
01189
01190 double delta = ROOT::Math::VectorUtil::DeltaR(curpfclusters[ic].position(),tmp[ic]);
01191 if ( delta > diff ) diff = delta;
01192 }
01193 }
01194
01195
01196
01197 #ifdef PFLOW_DEBUG
01198 if ( iter >= 50 && debug_ )
01199 cout << "PFClusterAlgo Warning: "
01200 << "more than "<<niter<<" iterations in pfcluster finding: "
01201 << setprecision(10) << diff << endl;
01202 #endif
01203
01204
01205
01206 for(unsigned ic=0; ic<curpfclusters.size(); ic++) {
01207 calculateClusterPosition(curpfclusters[ic], curpfclusterswodepthcor[ic],
01208 true, posCalcNCrystal);
01209 pfClusters_->push_back(curpfclusters[ic]);
01210 }
01211 }
01212
01213
01214
01215 void
01216 PFClusterAlgo::calculateClusterPosition(reco::PFCluster& cluster,
01217 reco::PFCluster& clusterwodepthcor,
01218 bool depcor,
01219 int posCalcNCrystal) {
01220
01221 if( posCalcNCrystal_ != -1 &&
01222 posCalcNCrystal_ != 5 &&
01223 posCalcNCrystal_ != 9 ) {
01224 throw "PFCluster::calculatePosition : posCalcNCrystal_ must be -1, 5, or 9.";
01225 }
01226
01227 if(!posCalcNCrystal) posCalcNCrystal = posCalcNCrystal_;
01228
01229 cluster.position_.SetXYZ(0,0,0);
01230
01231 cluster.energy_ = 0;
01232
01233 double normalize = 0;
01234
01235
01236
01237
01238 map <PFLayer::Layer, double> layers;
01239 unsigned seedIndex = 0;
01240 bool seedIndexFound = false;
01241
01242
01243
01244
01245 for (unsigned ic=0; ic<cluster.rechits_.size(); ic++ ) {
01246
01247 unsigned rhi = cluster.rechits_[ic].recHitRef().index();
01248
01249
01250 const reco::PFRecHit& rh = *(cluster.rechits_[ic].recHitRef());
01251 double fraction = cluster.rechits_[ic].fraction();
01252
01253
01254
01255 if( isSeed(rhi) && fraction > 1e-9 ) {
01256 seedIndex = rhi;
01257 seedIndexFound = true;
01258 }
01259
01260 double recHitEnergy = rh.energy() * fraction;
01261
01262
01263 if( recHitEnergy!=recHitEnergy ) {
01264 ostringstream ostr;
01265 edm::LogError("PFClusterAlgo")<<"rechit "<<rh.detId()<<" has a NaN energy... The input of the particle flow clustering seems to be corrupted.";
01266 }
01267
01268 cluster.energy_ += recHitEnergy;
01269
01270
01271 PFLayer::Layer layer = rh.layer();
01272 map <PFLayer::Layer, double>:: iterator it = layers.find(layer);
01273 if (it != layers.end())
01274 it->second += recHitEnergy;
01275 else
01276 layers.insert(make_pair(layer, recHitEnergy));
01277 }
01278
01279 assert(seedIndexFound);
01280
01281
01282 double Emax = 0.;
01283 PFLayer::Layer layer = PFLayer::NONE;
01284 for (map<PFLayer::Layer, double>::iterator it = layers.begin();
01285 it != layers.end(); ++it) {
01286 double e = it->second;
01287 if(e > Emax){
01288 Emax = e;
01289 layer = it->first;
01290 }
01291 }
01292
01293
01294 cluster.setLayer( layer );
01295
01296
01297
01298
01299 double p1 = posCalcP1_;
01300 if( p1 < 0 ) {
01301
01302
01303
01304
01305
01306 switch(cluster.layer() ) {
01307 case PFLayer::ECAL_BARREL:
01308 case PFLayer::HCAL_BARREL1:
01309 case PFLayer::HCAL_BARREL2:
01310 p1 = threshBarrel_;
01311 break;
01312 case PFLayer::ECAL_ENDCAP:
01313 case PFLayer::HCAL_ENDCAP:
01314 case PFLayer::HF_EM:
01315 case PFLayer::HF_HAD:
01316 case PFLayer::PS1:
01317 case PFLayer::PS2:
01318 p1 = threshEndcap_;
01319 break;
01320
01321
01322
01323
01324
01325
01326
01327
01328
01329
01330
01331
01332
01333
01334
01335
01336
01337
01338 default:
01339 cerr<<"Clusters weight_p1 -1 not yet allowed for layer "<<layer
01340 <<". Chose a better value in the opt file"<<endl;
01341 assert(0);
01342 break;
01343 }
01344 }
01345 else if( p1< 1e-9 ) {
01346 p1 = 1e-9;
01347 }
01348
01349
01350
01351 reco::PFCluster::REPPoint clusterpos;
01352 math::XYZPoint clusterposxyz;
01353 math::XYZPoint firstrechitposxyz;
01354
01355 double maxe = -9999;
01356 double x = 0;
01357 double y = 0;
01358 double z = 0;
01359
01360 for (unsigned ic=0; ic<cluster.rechits_.size(); ic++ ) {
01361
01362 unsigned rhi = cluster.rechits_[ic].recHitRef().index();
01363
01364
01365 const reco::PFRecHit& rh = *(cluster.rechits_[ic].recHitRef());
01366
01367 if(rhi != seedIndex) {
01368 if( posCalcNCrystal == 5 ) {
01369 if(!rh.isNeighbour4(seedIndex) ) {
01370 continue;
01371 }
01372 }
01373 if( posCalcNCrystal == 9 ) {
01374 if(!rh.isNeighbour8(seedIndex) ) {
01375 continue;
01376 }
01377 }
01378 }
01379
01380 double fraction = cluster.rechits_[ic].fraction();
01381 double recHitEnergy = rh.energy() * fraction;
01382
01383 double norm = fraction < 1E-9 ? 0. : max(0., log(recHitEnergy/p1 ));
01384
01385 const math::XYZPoint& rechitposxyz = rh.position();
01386
01387 if( recHitEnergy > maxe ) {
01388 firstrechitposxyz = rechitposxyz;
01389 maxe = recHitEnergy;
01390 }
01391
01392 x += rechitposxyz.X() * norm;
01393 y += rechitposxyz.Y() * norm;
01394 z += rechitposxyz.Z() * norm;
01395
01396
01397 normalize += norm;
01398 }
01399
01400
01401
01402 if( normalize < 1e-9 ) {
01403
01404
01405 cout << "Watch out : cluster too far from its seeding cell, set to 0,0,0" << endl;
01406 clusterposxyz.SetXYZ(0,0,0);
01407 clusterpos.SetCoordinates(0,0,0);
01408 return;
01409 }
01410 else {
01411 x /= normalize;
01412 y /= normalize;
01413 z /= normalize;
01414 clusterposxyz.SetCoordinates( x, y, z);
01415 clusterpos.SetCoordinates( clusterposxyz.Rho(), clusterposxyz.Eta(), clusterposxyz.Phi() );
01416 }
01417
01418 cluster.posrep_ = clusterpos;
01419 cluster.position_ = clusterposxyz;
01420
01421 clusterwodepthcor = cluster;
01422
01423
01424
01425
01426
01427
01428 if( depcor &&
01429 ( cluster.layer() == PFLayer::ECAL_BARREL ||
01430 cluster.layer() == PFLayer::ECAL_ENDCAP ) ) {
01431
01432
01433 double corra = reco::PFCluster::depthCorA_;
01434 double corrb = reco::PFCluster::depthCorB_;
01435 if( abs(clusterpos.Eta() ) < 2.6 &&
01436 abs(clusterpos.Eta() ) > 1.65 ) {
01437
01438
01439 corra = reco::PFCluster::depthCorAp_;
01440 corrb = reco::PFCluster::depthCorBp_;
01441 }
01442
01443 double depth = 0;
01444
01445 switch( reco::PFCluster::depthCorMode_ ) {
01446 case 1:
01447 depth = corra * ( corrb + log(cluster.energy_) );
01448 break;
01449 case 2:
01450 depth = corra;
01451 break;
01452 default:
01453 cerr<<"PFClusterAlgo::calculateClusterPosition : unknown function for depth correction! "<<endl;
01454 assert(0);
01455 }
01456
01457
01458
01459
01460
01461
01462
01463
01464
01465
01466
01467
01468
01469
01470
01471
01472
01473 math::XYZVector depthv( clusterposxyz.X(),
01474 clusterposxyz.Y(),
01475 clusterposxyz.Z() );
01476 depthv /= sqrt(depthv.Mag2() );
01477 depthv *= depth;
01478
01479
01480
01481 math::XYZPoint clusterposxyzcor;
01482
01483 maxe = -9999;
01484 x = 0;
01485 y = 0;
01486 z = 0;
01487 cluster.posrep_.SetXYZ(0,0,0);
01488 normalize = 0;
01489 for (unsigned ic=0; ic<cluster.rechits_.size(); ic++ ) {
01490
01491 unsigned rhi = cluster.rechits_[ic].recHitRef().index();
01492
01493
01494 const reco::PFRecHit& rh = *(cluster.rechits_[ic].recHitRef());
01495
01496 if(rhi != seedIndex) {
01497 if( posCalcNCrystal == 5 ) {
01498 if(!rh.isNeighbour4(seedIndex) ) {
01499 continue;
01500 }
01501 }
01502 if( posCalcNCrystal == 9 ) {
01503 if(!rh.isNeighbour8(seedIndex) ) {
01504 continue;
01505 }
01506 }
01507 }
01508
01509 double fraction = cluster.rechits_[ic].fraction();
01510 double recHitEnergy = rh.energy() * fraction;
01511
01512 const math::XYZPoint& rechitposxyz = rh.position();
01513
01514
01515 math::XYZVector rechitaxis = rh.getAxisXYZ();
01516
01517
01518 math::XYZVector rechitaxisu( rechitaxis );
01519 rechitaxisu /= sqrt( rechitaxis.Mag2() );
01520
01521 math::XYZVector displacement( rechitaxisu );
01522
01523 displacement *= rechitaxisu.Dot( depthv );
01524
01525 math::XYZPoint rechitposxyzcor( rechitposxyz );
01526 rechitposxyzcor += displacement;
01527
01528 if( recHitEnergy > maxe ) {
01529 firstrechitposxyz = rechitposxyzcor;
01530 maxe = recHitEnergy;
01531 }
01532
01533 double norm = fraction < 1E-9 ? 0. : max(0., log(recHitEnergy/p1 ));
01534
01535 x += rechitposxyzcor.X() * norm;
01536 y += rechitposxyzcor.Y() * norm;
01537 z += rechitposxyzcor.Z() * norm;
01538
01539
01540 normalize += norm;
01541 }
01542
01543
01544 if(normalize < 1e-9) {
01545 cerr<<"--------------------"<<endl;
01546 cerr<< cluster <<endl;
01547 assert(0);
01548 }
01549 else {
01550 x /= normalize;
01551 y /= normalize;
01552 z /= normalize;
01553
01554
01555 clusterposxyzcor.SetCoordinates(x,y,z);
01556 cluster.posrep_.SetCoordinates( clusterposxyzcor.Rho(),
01557 clusterposxyzcor.Eta(),
01558 clusterposxyzcor.Phi() );
01559 cluster.position_ = clusterposxyzcor;
01560 clusterposxyz = clusterposxyzcor;
01561 }
01562 }
01563 }
01564
01565
01566
01567 const reco::PFRecHit&
01568 PFClusterAlgo::rechit(unsigned i,
01569 const reco::PFRecHitCollection& rechits ) {
01570
01571 if(i >= rechits.size() ) {
01572 string err = "PFClusterAlgo::rechit : out of range";
01573 throw std::out_of_range(err);
01574 }
01575
01576 return rechits[i];
01577 }
01578
01579
01580
01581 bool PFClusterAlgo::masked(unsigned rhi) const {
01582
01583 if(rhi>=mask_.size() ) {
01584 string err = "PFClusterAlgo::masked : out of range";
01585 throw std::out_of_range(err);
01586 }
01587
01588 return mask_[rhi];
01589 }
01590
01591
01592 unsigned PFClusterAlgo::color(unsigned rhi) const {
01593
01594 if(rhi>=color_.size() ) {
01595 string err = "PFClusterAlgo::color : out of range";
01596 throw std::out_of_range(err);
01597 }
01598
01599 return color_[rhi];
01600 }
01601
01602
01603
01604 bool PFClusterAlgo::isSeed(unsigned rhi) const {
01605
01606 if(rhi>=seedStates_.size() ) {
01607 string err = "PFClusterAlgo::isSeed : out of range";
01608 throw std::out_of_range(err);
01609 }
01610
01611 return seedStates_[rhi]>0 ? true : false;
01612 }
01613
01614
01615 void PFClusterAlgo::paint(unsigned rhi, unsigned color ) {
01616
01617 if(rhi>=color_.size() ) {
01618 string err = "PFClusterAlgo::color : out of range";
01619 throw std::out_of_range(err);
01620 }
01621
01622 color_[rhi] = color;
01623 }
01624
01625
01626 reco::PFRecHitRef
01627 PFClusterAlgo::createRecHitRef( const reco::PFRecHitCollection& rechits,
01628 unsigned rhi ) {
01629
01630 if( rechitsHandle_.isValid() ) {
01631 return reco::PFRecHitRef( rechitsHandle_, rhi );
01632 }
01633 else {
01634 return reco::PFRecHitRef( &rechits, rhi );
01635 }
01636 }
01637
01638
01639 ostream& operator<<(ostream& out,const PFClusterAlgo& algo) {
01640 if(!out) return out;
01641 out<<"PFClusterAlgo parameters : "<<endl;
01642 out<<"-----------------------------------------------------"<<endl;
01643 out<<"threshBarrel : "<<algo.threshBarrel_ <<endl;
01644 out<<"threshSeedBarrel : "<<algo.threshSeedBarrel_ <<endl;
01645 out<<"threshPtBarrel : "<<algo.threshPtBarrel_ <<endl;
01646 out<<"threshPtSeedBarrel : "<<algo.threshPtSeedBarrel_ <<endl;
01647 out<<"threshCleanBarrel : "<<algo.threshCleanBarrel_ <<endl;
01648 out<<"minS4S1Barrel : "<<algo.minS4S1Barrel_[0]<<" x log10(E) + "<<algo.minS4S1Barrel_[1]<<endl;
01649 out<<"threshEndcap : "<<algo.threshEndcap_ <<endl;
01650 out<<"threshSeedEndcap : "<<algo.threshSeedEndcap_ <<endl;
01651 out<<"threshPtEndcap : "<<algo.threshPtEndcap_ <<endl;
01652 out<<"threshPtSeedEndcap : "<<algo.threshPtSeedEndcap_ <<endl;
01653 out<<"threshEndcap : "<<algo.threshEndcap_ <<endl;
01654 out<<"threshCleanEndcap : "<<algo.threshCleanEndcap_ <<endl;
01655 out<<"minS4S1Endcap : "<<algo.minS4S1Endcap_[0]<<" x log10(E) + "<<algo.minS4S1Endcap_[1]<<endl;
01656 out<<"nNeighbours : "<<algo.nNeighbours_ <<endl;
01657 out<<"posCalcNCrystal : "<<algo.posCalcNCrystal_ <<endl;
01658 out<<"posCalcP1 : "<<algo.posCalcP1_ <<endl;
01659 out<<"showerSigma : "<<algo.showerSigma_ <<endl;
01660 out<<"useCornerCells : "<<algo.useCornerCells_ <<endl;
01661
01662 out<<endl;
01663 out<<algo.pfClusters_->size()<<" clusters:"<<endl;
01664
01665 for(unsigned i=0; i<algo.pfClusters_->size(); i++) {
01666 out<<(*algo.pfClusters_)[i]<<endl;
01667
01668 if(!out) return out;
01669 }
01670
01671 return out;
01672 }
01673
01674
01675
01676 std::pair<double,double>
01677 PFClusterAlgo::dCrack(double phi, double eta){
01678
01679 static double pi= M_PI;
01680
01681
01682 static std::vector<double> cPhi;
01683 if(cPhi.size()==0)
01684 {
01685 cPhi.resize(18,0);
01686 cPhi[0]=2.97025;
01687 for(unsigned i=1;i<=17;++i) cPhi[i]=cPhi[0]-2*i*pi/18;
01688 }
01689
01690
01691 static double delta_cPhi=0.00638;
01692
01693 double defi;
01694
01695
01696 if(eta<0) phi +=delta_cPhi;
01697
01698 if (phi>=-pi && phi<=pi){
01699
01700
01701 if (phi<cPhi[17] || phi>=cPhi[0]){
01702 if (phi<0) phi+= 2*pi;
01703 defi = std::min(fabs(phi -cPhi[0]),fabs(phi-cPhi[17]-2*pi));
01704 }
01705
01706
01707 else{
01708 bool OK = false;
01709 unsigned i=16;
01710 while(!OK){
01711 if (phi<cPhi[i]){
01712 defi=std::min(fabs(phi-cPhi[i+1]),fabs(phi-cPhi[i]));
01713 OK=true;
01714 }
01715 else i-=1;
01716 }
01717 }
01718 }
01719 else{
01720 defi=0.;
01721 std::cout<<"Problem in dminphi"<<std::endl;
01722 }
01723
01724
01725 static std::vector<double> cEta;
01726 if ( cEta.size() == 0 ) {
01727 cEta.push_back(0.0);
01728 cEta.push_back(4.44747e-01) ; cEta.push_back(-4.44747e-01) ;
01729 cEta.push_back(7.92824e-01) ; cEta.push_back(-7.92824e-01) ;
01730 cEta.push_back(1.14090e+00) ; cEta.push_back(-1.14090e+00) ;
01731 cEta.push_back(1.47464e+00) ; cEta.push_back(-1.47464e+00) ;
01732 }
01733 double deta = 999.;
01734
01735 for ( unsigned ieta=0; ieta<cEta.size(); ++ieta ) {
01736 deta = std::min(deta,fabs(eta-cEta[ieta]));
01737 }
01738
01739 defi /= 0.0175;
01740 deta /= 0.0175;
01741 return std::pair<double,double>(defi,deta);
01742
01743 }