CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoParticleFlow/PFClusterProducer/src/PFClusterAlgo.cc

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