00001 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterTools.h"
00002
00003 #include "DataFormats/DetId/interface/DetId.h"
00004 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00005 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00006 #include "RecoCaloTools/Navigation/interface/CaloNavigator.h"
00007
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009
00010
00011 #include "CLHEP/Geometry/Transform3D.h"
00012
00013 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00014 #include "Geometry/CaloTopology/interface/CaloTopology.h"
00015 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00016 #include "Geometry/EcalAlgo/interface/EcalBarrelGeometry.h"
00017
00018 std::pair<DetId, float> EcalClusterTools::getMaximum( const std::vector< std::pair<DetId, float> > &v_id, const EcalRecHitCollection *recHits)
00019 {
00020 float max = 0;
00021 DetId id(0);
00022 for ( size_t i = 0; i < v_id.size(); ++i ) {
00023 float energy = recHitEnergy( v_id[i].first, recHits ) * v_id[i].second;
00024 if ( energy > max ) {
00025 max = energy;
00026 id = v_id[i].first;
00027 }
00028 }
00029 return std::pair<DetId, float>(id, max);
00030 }
00031
00032
00033
00034 std::pair<DetId, float> EcalClusterTools::getMaximum( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
00035 {
00036 return getMaximum( cluster.hitsAndFractions(), recHits );
00037 }
00038
00039
00040
00041 float EcalClusterTools::recHitEnergy(DetId id, const EcalRecHitCollection *recHits)
00042 {
00043 if ( id == DetId(0) ) {
00044 return 0;
00045 } else {
00046 EcalRecHitCollection::const_iterator it = recHits->find( id );
00047 if ( it != recHits->end() ) {
00048 return (*it).energy();
00049 } else {
00050
00051
00052 return 0;
00053 }
00054 }
00055 return 0;
00056 }
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071 float EcalClusterTools::matrixEnergy( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, DetId id, int ixMin, int ixMax, int iyMin, int iyMax )
00072 {
00073
00074 CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00075 float energy = 0;
00076 for ( int i = ixMin; i <= ixMax; ++i ) {
00077 for ( int j = iyMin; j <= iyMax; ++j ) {
00078 cursor.home();
00079 cursor.offsetBy( i, j );
00080 energy += recHitEnergy( *cursor, recHits );
00081 }
00082 }
00083
00084
00085
00086
00087
00088
00089 return energy;
00090 }
00091
00092
00093
00094 std::vector<DetId> EcalClusterTools::matrixDetId( const CaloTopology* topology, DetId id, int ixMin, int ixMax, int iyMin, int iyMax )
00095 {
00096 CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00097 std::vector<DetId> v;
00098 for ( int i = ixMin; i <= ixMax; ++i ) {
00099 for ( int j = iyMin; j <= iyMax; ++j ) {
00100 cursor.home();
00101 cursor.offsetBy( i, j );
00102 if ( *cursor != DetId(0) ) v.push_back( *cursor );
00103 }
00104 }
00105 return v;
00106 }
00107
00108
00109
00110 float EcalClusterTools::e2x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00111 {
00112 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00113 std::list<float> energies;
00114 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 0 ) );
00115 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, 0, 1 ) );
00116 energies.push_back( matrixEnergy( cluster, recHits, topology, id, 0, 1, 0, 1 ) );
00117 energies.push_back( matrixEnergy( cluster, recHits, topology, id, 0, 1, -1, 0 ) );
00118
00119
00120 return *std::max_element(energies.begin(),energies.end());
00121
00122 }
00123
00124
00125
00126 float EcalClusterTools::e3x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00127 {
00128 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00129 std::list<float> energies;
00130 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 0 ) );
00131 energies.push_back( matrixEnergy( cluster, recHits, topology, id, 0, 1, -1, 1 ) );
00132 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1, 0, 1 ) );
00133 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 1 ) );
00134 return *std::max_element(energies.begin(),energies.end());
00135 }
00136
00137
00138
00139 float EcalClusterTools::e3x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00140 {
00141 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00142 return matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 1 );
00143 }
00144
00145
00146
00147 float EcalClusterTools::e4x4( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00148 {
00149 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00150 std::list<float> energies;
00151 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -2, 1 ) );
00152 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -2, 1 ) );
00153 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -1, 2 ) );
00154 energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -1, 2 ) );
00155 return *std::max_element(energies.begin(),energies.end());
00156 }
00157
00158
00159
00160 float EcalClusterTools::e5x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00161 {
00162 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00163 return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, 2 );
00164 }
00165
00166
00167
00168 float EcalClusterTools::eMax( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00169 {
00170 return getMaximum( cluster.hitsAndFractions(), recHits ).second;
00171 }
00172
00173
00174
00175 float EcalClusterTools::e2nd( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00176 {
00177 std::list<float> energies;
00178 std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00179 if ( v_id.size() < 2 ) return 0;
00180 for ( size_t i = 0; i < v_id.size(); ++i ) {
00181 energies.push_back( recHitEnergy( v_id[i].first, recHits ) * v_id[i].second );
00182 }
00183 energies.sort();
00184 return *--(--energies.end());
00185
00186
00187 }
00188
00189
00190
00191 float EcalClusterTools::e2x5Right( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00192 {
00193 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00194 return matrixEnergy( cluster, recHits, topology, id, 1, 2, -2, 2 );
00195 }
00196
00197
00198
00199 float EcalClusterTools::e2x5Left( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00200 {
00201 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00202 return matrixEnergy( cluster, recHits, topology, id, -2, -1, -2, 2 );
00203 }
00204
00205
00206
00207 float EcalClusterTools::e2x5Top( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00208 {
00209 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00210 return matrixEnergy( cluster, recHits, topology, id, -2, 2, 1, 2 );
00211 }
00212
00213
00214
00215 float EcalClusterTools::e2x5Bottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00216 {
00217 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00218 return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, -1 );
00219 }
00220
00221
00222
00223 float EcalClusterTools::e2x5Max( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00224 {
00225 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00226
00227
00228 float left = matrixEnergy( cluster, recHits, topology, id, -1, -1, -2, 2 );
00229
00230 float right = matrixEnergy( cluster, recHits, topology, id, 1, 1, -2, 2 );
00231
00232 float centre = matrixEnergy( cluster, recHits, topology, id, 0, 0, -2, 2 );
00233
00234
00235 return left > right ? left+centre : right+centre;
00236 }
00237
00238
00239 float EcalClusterTools::e1x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00240 {
00241 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00242 return matrixEnergy( cluster, recHits, topology, id, 0, 0, -2, 2 );
00243 }
00244
00245
00246
00247 float EcalClusterTools::e5x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00248 {
00249 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00250 return matrixEnergy( cluster, recHits, topology, id, -2, 2, 0, 0 );
00251 }
00252
00253
00254
00255 float EcalClusterTools::e1x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00256 {
00257 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00258 return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, 1 );
00259 }
00260
00261
00262
00263 float EcalClusterTools::e3x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00264 {
00265 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00266 return matrixEnergy( cluster, recHits, topology, id, -1, 1, 0, 0 );
00267 }
00268
00269
00270
00271 float EcalClusterTools::eLeft( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00272 {
00273 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00274 return matrixEnergy( cluster, recHits, topology, id, -1, -1, 0, 0 );
00275 }
00276
00277
00278
00279 float EcalClusterTools::eRight( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00280 {
00281 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00282 return matrixEnergy( cluster, recHits, topology, id, 1, 1, 0, 0 );
00283 }
00284
00285
00286
00287 float EcalClusterTools::eTop( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00288 {
00289 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00290 return matrixEnergy( cluster, recHits, topology, id, 0, 0, 1, 1 );
00291 }
00292
00293
00294
00295 float EcalClusterTools::eBottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00296 {
00297 DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00298 return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, -1 );
00299 }
00300
00301
00302
00303 std::vector<float> EcalClusterTools::energyBasketFractionEta( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00304 {
00305 std::vector<float> basketFraction( 2 * EBDetId::kModulesPerSM );
00306 float clusterEnergy = cluster.energy();
00307 std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00308 if ( v_id[0].first.subdetId() != EcalBarrel ) {
00309 edm::LogWarning("EcalClusterTools::energyBasketFractionEta") << "Trying to get basket fraction for endcap basic-clusters. Basket fractions can be obtained ONLY for barrel basic-clusters. Returning empty vector.";
00310 return basketFraction;
00311 }
00312 for ( size_t i = 0; i < v_id.size(); ++i ) {
00313 basketFraction[ EBDetId(v_id[i].first).im()-1 + EBDetId(v_id[i].first).positiveZ()*EBDetId::kModulesPerSM ] += recHitEnergy( v_id[i].first, recHits ) * v_id[i].second / clusterEnergy;
00314 }
00315 std::sort( basketFraction.rbegin(), basketFraction.rend() );
00316 return basketFraction;
00317 }
00318
00319
00320
00321 std::vector<float> EcalClusterTools::energyBasketFractionPhi( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00322 {
00323 std::vector<float> basketFraction( 2 * (EBDetId::MAX_IPHI / EBDetId::kCrystalsInPhi) );
00324 float clusterEnergy = cluster.energy();
00325 std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00326 if ( v_id[0].first.subdetId() != EcalBarrel ) {
00327 edm::LogWarning("EcalClusterTools::energyBasketFractionPhi") << "Trying to get basket fraction for endcap basic-clusters. Basket fractions can be obtained ONLY for barrel basic-clusters. Returning empty vector.";
00328 return basketFraction;
00329 }
00330 for ( size_t i = 0; i < v_id.size(); ++i ) {
00331 basketFraction[ (EBDetId(v_id[i].first).iphi()-1)/EBDetId::kCrystalsInPhi + EBDetId(v_id[i].first).positiveZ()*EBDetId::kTowersInPhi] += recHitEnergy( v_id[i].first, recHits ) * v_id[i].second / clusterEnergy;
00332 }
00333 std::sort( basketFraction.rbegin(), basketFraction.rend() );
00334 return basketFraction;
00335 }
00336
00337
00338
00339 std::vector<EcalClusterTools::EcalClusterEnergyDeposition> EcalClusterTools::getEnergyDepTopology( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW, float w0 )
00340 {
00341 std::vector<EcalClusterTools::EcalClusterEnergyDeposition> energyDistribution;
00342
00343
00344 CLHEP::Hep3Vector clVect(cluster.position().x(), cluster.position().y(), cluster.position().z());
00345 CLHEP::Hep3Vector clDir(clVect);
00346 clDir*=1.0/clDir.mag();
00347
00348 CLHEP::Hep3Vector theta_axis(clDir.y(),-clDir.x(),0.0);
00349 theta_axis *= 1.0/theta_axis.mag();
00350 CLHEP::Hep3Vector phi_axis = theta_axis.cross(clDir);
00351
00352 std::vector< std::pair<DetId, float> > clusterDetIds = cluster.hitsAndFractions();
00353
00354 EcalClusterEnergyDeposition clEdep;
00355 EcalRecHit testEcalRecHit;
00356 std::vector< std::pair<DetId, float> >::iterator posCurrent;
00357
00358 for(posCurrent=clusterDetIds.begin(); posCurrent!=clusterDetIds.end(); ++posCurrent) {
00359 EcalRecHitCollection::const_iterator itt = recHits->find( (*posCurrent).first );
00360 testEcalRecHit=*itt;
00361
00362 if(( (*posCurrent).first != DetId(0)) && (recHits->find( (*posCurrent).first ) != recHits->end())) {
00363 clEdep.deposited_energy = testEcalRecHit.energy() * (*posCurrent).second;
00364
00365 if(logW) {
00366
00367
00368 double weight = std::max(0.0, w0 + log(fabs(clEdep.deposited_energy)/cluster.energy()) );
00369 if(weight==0) {
00370 LogDebug("ClusterShapeAlgo") << "Crystal has insufficient energy: E = "
00371 << clEdep.deposited_energy << " GeV; skipping... ";
00372 continue;
00373 }
00374 else LogDebug("ClusterShapeAlgo") << "===> got crystal. Energy = " << clEdep.deposited_energy << " GeV. ";
00375 }
00376 DetId id_ = (*posCurrent).first;
00377 const CaloCellGeometry *this_cell = geometry->getSubdetectorGeometry(id_)->getGeometry(id_);
00378 GlobalPoint cellPos = this_cell->getPosition();
00379 CLHEP::Hep3Vector gblPos (cellPos.x(),cellPos.y(),cellPos.z());
00380
00381 CLHEP::Hep3Vector diff = gblPos - clVect;
00382
00383
00384
00385 CLHEP::Hep3Vector DigiVect = diff - diff.dot(clDir)*clDir;
00386 clEdep.r = DigiVect.mag();
00387 LogDebug("ClusterShapeAlgo") << "E = " << clEdep.deposited_energy
00388 << "\tdiff = " << diff.mag()
00389 << "\tr = " << clEdep.r;
00390 clEdep.phi = DigiVect.angle(theta_axis);
00391 if(DigiVect.dot(phi_axis)<0) clEdep.phi = 2 * M_PI - clEdep.phi;
00392 energyDistribution.push_back(clEdep);
00393 }
00394 }
00395 return energyDistribution;
00396 }
00397
00398
00399
00400 std::vector<float> EcalClusterTools::lat( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW, float w0 )
00401 {
00402 std::vector<EcalClusterTools::EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
00403
00404 std::vector<float> lat;
00405 double r, redmoment=0;
00406 double phiRedmoment = 0 ;
00407 double etaRedmoment = 0 ;
00408 int n,n1,n2,tmp;
00409 int clusterSize=energyDistribution.size();
00410 float etaLat_, phiLat_, lat_;
00411 if (clusterSize<3) {
00412 etaLat_ = 0.0 ;
00413 lat_ = 0.0;
00414 lat.push_back(0.);
00415 lat.push_back(0.);
00416 lat.push_back(0.);
00417 return lat;
00418 }
00419
00420 n1=0; n2=1;
00421 if (energyDistribution[1].deposited_energy >
00422 energyDistribution[0].deposited_energy)
00423 {
00424 tmp=n2; n2=n1; n1=tmp;
00425 }
00426 for (int i=2; i<clusterSize; i++) {
00427 n=i;
00428 if (energyDistribution[i].deposited_energy >
00429 energyDistribution[n1].deposited_energy)
00430 {
00431 tmp = n2;
00432 n2 = n1; n1 = i; n=tmp;
00433 } else {
00434 if (energyDistribution[i].deposited_energy >
00435 energyDistribution[n2].deposited_energy)
00436 {
00437 tmp=n2; n2=i; n=tmp;
00438 }
00439 }
00440
00441 r = energyDistribution[n].r;
00442 redmoment += r*r* energyDistribution[n].deposited_energy;
00443 double rphi = r * cos (energyDistribution[n].phi) ;
00444 phiRedmoment += rphi * rphi * energyDistribution[n].deposited_energy;
00445 double reta = r * sin (energyDistribution[n].phi) ;
00446 etaRedmoment += reta * reta * energyDistribution[n].deposited_energy;
00447 }
00448 double e1 = energyDistribution[n1].deposited_energy;
00449 double e2 = energyDistribution[n2].deposited_energy;
00450
00451 lat_ = redmoment/(redmoment+2.19*2.19*(e1+e2));
00452 phiLat_ = phiRedmoment/(phiRedmoment+2.19*2.19*(e1+e2));
00453 etaLat_ = etaRedmoment/(etaRedmoment+2.19*2.19*(e1+e2));
00454
00455 lat.push_back(etaLat_);
00456 lat.push_back(phiLat_);
00457 lat.push_back(lat_);
00458 return lat;
00459 }
00460
00461
00462
00463 math::XYZVector EcalClusterTools::meanClusterPosition( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, const CaloGeometry *geometry )
00464 {
00465
00466 math::XYZVector meanPosition(0.0, 0.0, 0.0);
00467 std::vector<DetId> v_id = matrixDetId( topology, getMaximum( cluster, recHits ).first, -2, 2, -2, 2 );
00468 for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00469 GlobalPoint positionGP = geometry->getSubdetectorGeometry( *it )->getGeometry( *it )->getPosition();
00470 math::XYZVector position(positionGP.x(),positionGP.y(),positionGP.z());
00471 meanPosition = meanPosition + recHitEnergy( *it, recHits ) * position;
00472 }
00473 return meanPosition / e5x5( cluster, recHits, topology );
00474 }
00475
00476
00477
00478
00479
00480
00481
00482
00483
00484
00485 std::pair<float,float> EcalClusterTools::mean5x5PositionInLocalCrysCoord(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology)
00486 {
00487 DetId seedId = getMaximum( cluster, recHits ).first;
00488 float meanDEta=0.;
00489 float meanDPhi=0.;
00490 float energySum=0.;
00491
00492 std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00493 for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00494 float energy = recHitEnergy(*it,recHits);
00495 if(energy<0.) continue;
00496 meanDEta += energy * getNrCrysDiffInEta(*it,seedId);
00497 meanDPhi += energy * getNrCrysDiffInPhi(*it,seedId);
00498 energySum +=energy;
00499 }
00500 meanDEta /=energySum;
00501 meanDPhi /=energySum;
00502 return std::pair<float,float>(meanDEta,meanDPhi);
00503 }
00504
00505
00506
00507
00508
00509
00510
00511 std::pair<float,float> EcalClusterTools::mean5x5PositionInXY(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology)
00512 {
00513 DetId seedId = getMaximum( cluster, recHits ).first;
00514
00515 std::pair<float,float> meanXY(0.,0.);
00516 if(seedId.subdetId()==EcalBarrel) return meanXY;
00517
00518 float energySum=0.;
00519
00520 std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00521 for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00522 float energy = recHitEnergy(*it,recHits);
00523 if(energy<0.) continue;
00524 meanXY.first += energy * getNormedIX(*it);
00525 meanXY.second += energy * getNormedIY(*it);
00526 energySum +=energy;
00527 }
00528 meanXY.first/=energySum;
00529 meanXY.second/=energySum;
00530 return meanXY;
00531 }
00532
00533
00534 std::vector<float> EcalClusterTools::covariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits, const CaloTopology *topology, const CaloGeometry* geometry, float w0)
00535 {
00536 float e_5x5 = e5x5( cluster, recHits, topology );
00537 float covEtaEta, covEtaPhi, covPhiPhi;
00538 if (e_5x5 >= 0.) {
00539
00540 std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00541 math::XYZVector meanPosition = meanClusterPosition( cluster, recHits, topology, geometry );
00542
00543
00544 double numeratorEtaEta = 0;
00545 double numeratorEtaPhi = 0;
00546 double numeratorPhiPhi = 0;
00547 double denominator = 0;
00548
00549 DetId id = getMaximum( v_id, recHits ).first;
00550 CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00551 for ( int i = -2; i <= 2; ++i ) {
00552 for ( int j = -2; j <= 2; ++j ) {
00553 cursor.home();
00554 cursor.offsetBy( i, j );
00555 float energy = recHitEnergy( *cursor, recHits );
00556
00557 if ( energy <= 0 ) continue;
00558
00559 GlobalPoint position = geometry->getSubdetectorGeometry(*cursor)->getGeometry(*cursor)->getPosition();
00560
00561 double dPhi = position.phi() - meanPosition.phi();
00562 if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00563 if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() + dPhi; }
00564
00565 double dEta = position.eta() - meanPosition.eta();
00566 double w = 0.;
00567 w = std::max(0.0, w0 + log( energy / e_5x5 ));
00568
00569 denominator += w;
00570 numeratorEtaEta += w * dEta * dEta;
00571 numeratorEtaPhi += w * dEta * dPhi;
00572 numeratorPhiPhi += w * dPhi * dPhi;
00573 }
00574 }
00575
00576 if (denominator != 0.0) {
00577 covEtaEta = numeratorEtaEta / denominator;
00578 covEtaPhi = numeratorEtaPhi / denominator;
00579 covPhiPhi = numeratorPhiPhi / denominator;
00580 } else {
00581 covEtaEta = 999.9;
00582 covEtaPhi = 999.9;
00583 covPhiPhi = 999.9;
00584 }
00585
00586 } else {
00587
00588
00589 covEtaEta = 0;
00590 covEtaPhi = 0;
00591 covPhiPhi = 0;
00592 }
00593 std::vector<float> v;
00594 v.push_back( covEtaEta );
00595 v.push_back( covEtaPhi );
00596 v.push_back( covPhiPhi );
00597 return v;
00598 }
00599
00600
00601
00602
00603
00604
00605
00606
00607
00608
00609 std::vector<float> EcalClusterTools::localCovariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology,float w0)
00610 {
00611
00612 float e_5x5 = e5x5( cluster, recHits, topology );
00613 float covEtaEta, covEtaPhi, covPhiPhi;
00614
00615 if (e_5x5 >= 0.) {
00616
00617 std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00618 std::pair<float,float> mean5x5PosInNrCrysFromSeed = mean5x5PositionInLocalCrysCoord( cluster, recHits, topology );
00619 std::pair<float,float> mean5x5XYPos = mean5x5PositionInXY(cluster,recHits,topology);
00620
00621
00622 double numeratorEtaEta = 0;
00623 double numeratorEtaPhi = 0;
00624 double numeratorPhiPhi = 0;
00625 double denominator = 0;
00626
00627
00628
00629 const double barrelCrysSize = 0.01745;
00630 const double endcapCrysSize = 0.0447;
00631
00632 DetId seedId = getMaximum( v_id, recHits ).first;
00633
00634 bool isBarrel=seedId.subdetId()==EcalBarrel;
00635 const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
00636
00637 CaloNavigator<DetId> cursor = CaloNavigator<DetId>( seedId, topology->getSubdetectorTopology( seedId ) );
00638
00639 for ( int eastNr = -2; eastNr <= 2; ++eastNr ) {
00640 for ( int northNr = -2; northNr <= 2; ++northNr ) {
00641 cursor.home();
00642 cursor.offsetBy( eastNr, northNr);
00643 float energy = recHitEnergy( *cursor, recHits );
00644 if ( energy <= 0 ) continue;
00645
00646 float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
00647 float dPhi = 0;
00648
00649 if(isBarrel) dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
00650 else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
00651
00652
00653 double w = std::max(0.0, w0 + log( energy / e_5x5 ));
00654
00655 denominator += w;
00656 numeratorEtaEta += w * dEta * dEta;
00657 numeratorEtaPhi += w * dEta * dPhi;
00658 numeratorPhiPhi += w * dPhi * dPhi;
00659 }
00660 }
00661
00662
00663
00664 if (denominator != 0.0) {
00665 covEtaEta = crysSize*crysSize* numeratorEtaEta / denominator;
00666 covEtaPhi = crysSize*crysSize* numeratorEtaPhi / denominator;
00667 covPhiPhi = crysSize*crysSize* numeratorPhiPhi / denominator;
00668 } else {
00669 covEtaEta = 999.9;
00670 covEtaPhi = 999.9;
00671 covPhiPhi = 999.9;
00672 }
00673
00674
00675 } else {
00676
00677
00678 covEtaEta = 0;
00679 covEtaPhi = 0;
00680 covPhiPhi = 0;
00681 }
00682 std::vector<float> v;
00683 v.push_back( covEtaEta );
00684 v.push_back( covEtaPhi );
00685 v.push_back( covPhiPhi );
00686 return v;
00687 }
00688
00689
00690
00691
00692 double EcalClusterTools::zernike20( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0, bool logW, float w0 )
00693 {
00694 return absZernikeMoment( cluster, recHits, geometry, 2, 0, R0, logW, w0 );
00695 }
00696
00697
00698
00699 double EcalClusterTools::zernike42( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0, bool logW, float w0 )
00700 {
00701 return absZernikeMoment( cluster, recHits, geometry, 4, 2, R0, logW, w0 );
00702 }
00703
00704
00705
00706 double EcalClusterTools::absZernikeMoment( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
00707 {
00708
00709 if ((m>n) || ((n-m)%2 != 0) || (n<0) || (m<0)) return -1;
00710
00711
00712
00713 if ((n>20) || (R0<=2.19)) return -1;
00714 if (n<=5) return fast_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0 );
00715 else return calc_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0 );
00716 }
00717
00718
00719 double EcalClusterTools::fast_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
00720 {
00721 double r,ph,e,Re=0,Im=0;
00722 double TotalEnergy = cluster.energy();
00723 int index = (n/2)*(n/2)+(n/2)+m;
00724 std::vector<EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
00725 int clusterSize = energyDistribution.size();
00726 if(clusterSize < 3) return 0.0;
00727
00728 for (int i=0; i<clusterSize; i++)
00729 {
00730 r = energyDistribution[i].r / R0;
00731 if (r<1) {
00732 std::vector<double> pol;
00733 pol.push_back( f00(r) );
00734 pol.push_back( f11(r) );
00735 pol.push_back( f20(r) );
00736 pol.push_back( f22(r) );
00737 pol.push_back( f31(r) );
00738 pol.push_back( f33(r) );
00739 pol.push_back( f40(r) );
00740 pol.push_back( f42(r) );
00741 pol.push_back( f44(r) );
00742 pol.push_back( f51(r) );
00743 pol.push_back( f53(r) );
00744 pol.push_back( f55(r) );
00745 ph = (energyDistribution[i]).phi;
00746 e = energyDistribution[i].deposited_energy;
00747 Re = Re + e/TotalEnergy * pol[index] * cos( (double) m * ph);
00748 Im = Im - e/TotalEnergy * pol[index] * sin( (double) m * ph);
00749 }
00750 }
00751 return sqrt(Re*Re+Im*Im);
00752 }
00753
00754
00755
00756 double EcalClusterTools::calc_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
00757 {
00758 double r, ph, e, Re=0, Im=0, f_nm;
00759 double TotalEnergy = cluster.energy();
00760 std::vector<EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
00761 int clusterSize=energyDistribution.size();
00762 if(clusterSize<3) return 0.0;
00763
00764 for (int i = 0; i < clusterSize; ++i)
00765 {
00766 r = energyDistribution[i].r / R0;
00767 if (r < 1) {
00768 ph = energyDistribution[i].phi;
00769 e = energyDistribution[i].deposited_energy;
00770 f_nm = 0;
00771 for (int s=0; s<=(n-m)/2; s++) {
00772 if (s%2==0) {
00773 f_nm = f_nm + factorial(n-s)*pow(r,(double) (n-2*s))/(factorial(s)*factorial((n+m)/2-s)*factorial((n-m)/2-s));
00774 } else {
00775 f_nm = f_nm - factorial(n-s)*pow(r,(double) (n-2*s))/(factorial(s)*factorial((n+m)/2-s)*factorial((n-m)/2-s));
00776 }
00777 }
00778 Re = Re + e/TotalEnergy * f_nm * cos( (double) m*ph);
00779 Im = Im - e/TotalEnergy * f_nm * sin( (double) m*ph);
00780 }
00781 }
00782 return sqrt(Re*Re+Im*Im);
00783 }
00784
00785
00786
00787
00788
00789 float EcalClusterTools::getIEta(const DetId& id)
00790 {
00791 if(id.det()==DetId::Ecal){
00792 if(id.subdetId()==EcalBarrel){
00793 EBDetId ebId(id);
00794 return ebId.ieta();
00795 }else if(id.subdetId()==EcalEndcap){
00796 float iXNorm = getNormedIX(id);
00797 float iYNorm = getNormedIY(id);
00798
00799 return std::sqrt(iXNorm*iXNorm+iYNorm*iYNorm);
00800 }
00801 }
00802 return 0.;
00803 }
00804
00805
00806
00807
00808
00809
00810 float EcalClusterTools::getIPhi(const DetId& id)
00811 {
00812 if(id.det()==DetId::Ecal){
00813 if(id.subdetId()==EcalBarrel){
00814 EBDetId ebId(id);
00815 return ebId.iphi();
00816 }
00817 }
00818 return 0.;
00819 }
00820
00821
00822 float EcalClusterTools::getNormedIX(const DetId& id)
00823 {
00824 if(id.det()==DetId::Ecal && id.subdetId()==EcalEndcap){
00825 EEDetId eeId(id);
00826 int iXNorm = eeId.ix()-50;
00827 if(iXNorm<=0) iXNorm--;
00828 return iXNorm;
00829 }
00830 return 0;
00831 }
00832
00833
00834 float EcalClusterTools::getNormedIY(const DetId& id)
00835 {
00836 if(id.det()==DetId::Ecal && id.subdetId()==EcalEndcap){
00837 EEDetId eeId(id);
00838 int iYNorm = eeId.iy()-50;
00839 if(iYNorm<=0) iYNorm--;
00840 return iYNorm;
00841 }
00842 return 0;
00843 }
00844
00845
00846 float EcalClusterTools::getNrCrysDiffInEta(const DetId& crysId,const DetId& orginId)
00847 {
00848 float crysIEta = getIEta(crysId);
00849 float orginIEta = getIEta(orginId);
00850 bool isBarrel = orginId.subdetId()==EcalBarrel;
00851
00852 float nrCrysDiff = crysIEta-orginIEta;
00853
00854
00855
00856 if(isBarrel){
00857 if(crysIEta*orginIEta<0){
00858 if(crysIEta>0) nrCrysDiff--;
00859 else nrCrysDiff++;
00860 }
00861 }
00862 return nrCrysDiff;
00863 }
00864
00865
00866 float EcalClusterTools::getNrCrysDiffInPhi(const DetId& crysId,const DetId& orginId)
00867 {
00868 float crysIPhi = getIPhi(crysId);
00869 float orginIPhi = getIPhi(orginId);
00870 bool isBarrel = orginId.subdetId()==EcalBarrel;
00871
00872 float nrCrysDiff = crysIPhi-orginIPhi;
00873
00874 if(isBarrel){
00875 if (nrCrysDiff > + 180) { nrCrysDiff = nrCrysDiff - 360; }
00876 if (nrCrysDiff < - 180) { nrCrysDiff = nrCrysDiff + 360; }
00877 }
00878 return nrCrysDiff;
00879 }
00880
00881
00882 float EcalClusterTools::getDPhiEndcap(const DetId& crysId,float meanX,float meanY)
00883 {
00884 float iXNorm = getNormedIX(crysId);
00885 float iYNorm = getNormedIY(crysId);
00886
00887 float hitLocalR2 = (iXNorm-meanX)*(iXNorm-meanX)+(iYNorm-meanY)*(iYNorm-meanY);
00888 float hitR2 = iXNorm*iXNorm+iYNorm*iYNorm;
00889 float meanR2 = meanX*meanX+meanY*meanY;
00890 float hitR = sqrt(hitR2);
00891 float meanR = sqrt(meanR2);
00892
00893 float phi = acos((hitR2+meanR2-hitLocalR2)/(2*hitR*meanR));
00894 float dPhi = hitR*phi;
00895
00896 return dPhi;
00897 }
00898
00899 std::vector<float> EcalClusterTools::scLocalCovariances(const reco::SuperCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology, float w0)
00900 {
00901 const reco::BasicCluster bcluster = *(cluster.seed());
00902
00903 float e_5x5 = e5x5(bcluster, recHits, topology);
00904 float covEtaEta, covEtaPhi, covPhiPhi;
00905
00906 if (e_5x5 >= 0.) {
00907 std::vector<std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00908 std::pair<float,float> mean5x5PosInNrCrysFromSeed = mean5x5PositionInLocalCrysCoord(bcluster, recHits, topology);
00909 std::pair<float,float> mean5x5XYPos = mean5x5PositionInXY(cluster,recHits,topology);
00910
00911 double numeratorEtaEta = 0;
00912 double numeratorEtaPhi = 0;
00913 double numeratorPhiPhi = 0;
00914 double denominator = 0;
00915
00916 const double barrelCrysSize = 0.01745;
00917 const double endcapCrysSize = 0.0447;
00918
00919 DetId seedId = getMaximum(v_id, recHits).first;
00920 bool isBarrel=seedId.subdetId()==EcalBarrel;
00921
00922 const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
00923
00924 for (size_t i = 0; i < v_id.size(); ++i) {
00925 CaloNavigator<DetId> cursor = CaloNavigator<DetId>(v_id[i].first, topology->getSubdetectorTopology(v_id[i].first));
00926 float energy = recHitEnergy(*cursor, recHits);
00927
00928 if (energy <= 0) continue;
00929
00930 float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
00931 float dPhi = 0;
00932 if(isBarrel) dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
00933 else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
00934
00935
00936
00937 double w = 0.;
00938 w = std::max(0.0, w0 + log( energy / e_5x5 ));
00939
00940 denominator += w;
00941 numeratorEtaEta += w * dEta * dEta;
00942 numeratorEtaPhi += w * dEta * dPhi;
00943 numeratorPhiPhi += w * dPhi * dPhi;
00944 }
00945
00946
00947 if (denominator != 0.0) {
00948 covEtaEta = crysSize*crysSize* numeratorEtaEta / denominator;
00949 covEtaPhi = crysSize*crysSize* numeratorEtaPhi / denominator;
00950 covPhiPhi = crysSize*crysSize* numeratorPhiPhi / denominator;
00951 } else {
00952 covEtaEta = 999.9;
00953 covEtaPhi = 999.9;
00954 covPhiPhi = 999.9;
00955 }
00956
00957 } else {
00958
00959
00960 covEtaEta = 0;
00961 covEtaPhi = 0;
00962 covPhiPhi = 0;
00963 }
00964
00965 std::vector<float> v;
00966 v.push_back( covEtaEta );
00967 v.push_back( covEtaPhi );
00968 v.push_back( covPhiPhi );
00969
00970 return v;
00971 }
00972
00973
00974
00975
00976
00977
00978
00979
00980 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( const reco::BasicCluster &basicCluster, const EcalRecHitCollection &recHits, double phiCorrectionFactor, double w0, bool useLogWeights) {
00981
00982 Cluster2ndMoments returnMoments;
00983 returnMoments.sMaj = -1.;
00984 returnMoments.sMin = -1.;
00985 returnMoments.alpha = 0.;
00986
00987
00988 if( fabs( basicCluster.eta() ) < 1.479 ) {
00989
00990 std::vector<const EcalRecHit*> RH_ptrs;
00991
00992 std::vector< std::pair<DetId, float> > myHitsPair = basicCluster.hitsAndFractions();
00993 std::vector<DetId> usedCrystals;
00994 for(unsigned int i=0; i< myHitsPair.size(); i++){
00995 usedCrystals.push_back(myHitsPair[i].first);
00996 }
00997
00998 for(unsigned int i=0; i<usedCrystals.size(); i++){
00999
01000 EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01001 RH_ptrs.push_back( &(*myRH) );
01002 }
01003
01004 returnMoments = EcalClusterTools::cluster2ndMoments(RH_ptrs, phiCorrectionFactor, w0, useLogWeights);
01005
01006 }
01007
01008 return returnMoments;
01009
01010 }
01011
01012
01013 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( const reco::SuperCluster &superCluster, const EcalRecHitCollection &recHits, double phiCorrectionFactor, double w0, bool useLogWeights) {
01014
01015
01016 Cluster2ndMoments returnMoments;
01017 returnMoments.sMaj = -1.;
01018 returnMoments.sMin = -1.;
01019 returnMoments.alpha = 0.;
01020
01021
01022 if( fabs( superCluster.eta() ) < 1.479 ) {
01023 returnMoments = EcalClusterTools::cluster2ndMoments( *(superCluster.seed()), recHits, phiCorrectionFactor, w0, useLogWeights);
01024 }
01025
01026 return returnMoments;
01027
01028 }
01029
01030
01031 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( std::vector<const EcalRecHit*> RH_ptrs, double phiCorrectionFactor, double w0, bool useLogWeights) {
01032
01033 double mid_eta,mid_phi;
01034 mid_eta=mid_phi=0.;
01035
01036 double Etot = EcalClusterTools::getSumEnergy( RH_ptrs );
01037
01038 double max_phi=-10.;
01039 double min_phi=100.;
01040
01041 std::vector<double> etaDetId;
01042 std::vector<double> phiDetId;
01043 std::vector<double> wiDetId;
01044
01045 int nCry=0;
01046 double denominator=0.;
01047
01048
01049
01050 for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01051
01052
01053 EBDetId temp_EBDetId( (*rh_ptr)->detid() );
01054 double temp_eta=(temp_EBDetId.ieta() > 0. ? temp_EBDetId.ieta() + 84.5 : temp_EBDetId.ieta() + 85.5);
01055 double temp_phi=temp_EBDetId.iphi() - 0.5;
01056 double temp_ene=(*rh_ptr)->energy();
01057
01058 double temp_wi=((useLogWeights) ?
01059 std::max(0., w0 + log( fabs(temp_ene)/Etot ))
01060 : temp_ene);
01061
01062 if(temp_phi>max_phi) max_phi=temp_phi;
01063 if(temp_phi<min_phi) min_phi=temp_phi;
01064 etaDetId.push_back(temp_eta);
01065 phiDetId.push_back(temp_phi);
01066 wiDetId.push_back(temp_wi);
01067 denominator+=temp_wi;
01068 nCry++;
01069 }
01070
01071
01072
01073 if(max_phi==359.5 && min_phi==0.5){
01074 for(int i=0; i<nCry; i++){
01075 if(phiDetId[i] - 179. > 0.) phiDetId[i]-=360.;
01076 mid_phi+=phiDetId[i]*wiDetId[i];
01077 mid_eta+=etaDetId[i]*wiDetId[i];
01078 }
01079 }
01080
01081 else{
01082 for(int i=0; i<nCry; i++){
01083 mid_phi+=phiDetId[i]*wiDetId[i];
01084 mid_eta+=etaDetId[i]*wiDetId[i];
01085 }
01086 }
01087
01088 mid_eta/=denominator;
01089 mid_phi/=denominator;
01090
01091
01092
01093
01094 double See=0.;
01095 double Spp=0.;
01096 double Sep=0.;
01097
01098
01099 for(int i=0; i<nCry; i++) {
01100 See += (wiDetId[i]*(etaDetId[i]-mid_eta)*(etaDetId[i]-mid_eta)) / denominator;
01101 Spp += phiCorrectionFactor*(wiDetId[i]*(phiDetId[i]-mid_phi)*(phiDetId[i]-mid_phi)) / denominator;
01102 Sep += sqrt(phiCorrectionFactor)*(wiDetId[i]*(etaDetId[i]-mid_eta)*(phiDetId[i]-mid_phi)) / denominator;
01103 }
01104
01105 Cluster2ndMoments returnMoments;
01106
01107
01108 returnMoments.sMaj = ((See + Spp) + sqrt((See - Spp)*(See - Spp) + 4.*Sep*Sep)) / 2.;
01109 returnMoments.sMin = ((See + Spp) - sqrt((See - Spp)*(See - Spp) + 4.*Sep*Sep)) / 2.;
01110
01111 returnMoments.alpha = atan( (See - Spp + sqrt( (Spp - See)*(Spp - See) + 4.*Sep*Sep )) / (2.*Sep));
01112
01113 return returnMoments;
01114
01115 }
01116
01117
01118
01119
01120
01121
01122
01123
01124
01125
01126 std::vector<float> EcalClusterTools::roundnessBarrelSuperClusters( const reco::SuperCluster &superCluster ,const EcalRecHitCollection &recHits, int weightedPositionMethod, float energyThreshold){
01127 std::vector<const EcalRecHit*>RH_ptrs;
01128 std::vector< std::pair<DetId, float> > myHitsPair = superCluster.hitsAndFractions();
01129 std::vector<DetId> usedCrystals;
01130 for(unsigned int i=0; i< myHitsPair.size(); i++){
01131 usedCrystals.push_back(myHitsPair[i].first);
01132 }
01133 for(unsigned int i=0; i<usedCrystals.size(); i++){
01134
01135 EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01136 if( myRH != recHits.end() && myRH->energy() > energyThreshold){
01137 RH_ptrs.push_back( &(*myRH) );
01138 }
01139 }
01140 std::vector<float> temp = EcalClusterTools::roundnessSelectedBarrelRecHits(RH_ptrs,weightedPositionMethod);
01141 return temp;
01142 }
01143
01144
01145
01146
01147
01148
01149 std::vector<float> EcalClusterTools::roundnessBarrelSuperClustersUserExtended( const reco::SuperCluster &superCluster ,const EcalRecHitCollection &recHits, int ieta_delta, int iphi_delta, float energyRHThresh, int weightedPositionMethod){
01150
01151 std::vector<const EcalRecHit*>RH_ptrs;
01152 std::vector< std::pair<DetId, float> > myHitsPair = superCluster.hitsAndFractions();
01153 std::vector<DetId> usedCrystals;
01154 for(unsigned int i=0; i< myHitsPair.size(); i++){
01155 usedCrystals.push_back(myHitsPair[i].first);
01156 }
01157
01158 for(unsigned int i=0; i<usedCrystals.size(); i++){
01159
01160 EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01161 if(myRH != recHits.end() && myRH->energy() > energyRHThresh)
01162 RH_ptrs.push_back( &(*myRH) );
01163 }
01164
01165
01166 std::vector<int> seedPosition = EcalClusterTools::getSeedPosition( RH_ptrs );
01167
01168 for(EcalRecHitCollection::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++){
01169 EBDetId EBdetIdi( rh->detid() );
01170
01171 bool inEtaWindow = ( abs( deltaIEta(seedPosition[0],EBdetIdi.ieta()) ) <= ieta_delta );
01172 bool inPhiWindow = ( abs( deltaIPhi(seedPosition[1],EBdetIdi.iphi()) ) <= iphi_delta );
01173 bool passEThresh = ( rh->energy() > energyRHThresh );
01174 bool alreadyCounted = false;
01175
01176
01177 bool is_SCrh_inside_recHits = false;
01178 for(unsigned int i=0; i<usedCrystals.size(); i++){
01179 EcalRecHitCollection::const_iterator SCrh = recHits.find(usedCrystals[i]);
01180 if(SCrh != recHits.end()){
01181 is_SCrh_inside_recHits = true;
01182 if( rh->detid() == SCrh->detid() ) alreadyCounted = true;
01183 }
01184 }
01185
01186 if( is_SCrh_inside_recHits && !alreadyCounted && passEThresh && inEtaWindow && inPhiWindow){
01187 RH_ptrs.push_back( &(*rh) );
01188 }
01189
01190 }
01191 return EcalClusterTools::roundnessSelectedBarrelRecHits(RH_ptrs,weightedPositionMethod);
01192 }
01193
01194
01195
01196
01197 std::vector<float> EcalClusterTools::roundnessSelectedBarrelRecHits( std::vector<const EcalRecHit*> RH_ptrs, int weightedPositionMethod){
01198
01199
01200 std::vector<float> shapes;
01201
01202
01203 if(RH_ptrs.size()<2){
01204 shapes.push_back( -3 );
01205 shapes.push_back( -3 );
01206 return shapes;
01207 }
01208
01209
01210 std::vector<int> seedPosition = EcalClusterTools::getSeedPosition( RH_ptrs );
01211 int tempInt = seedPosition[0];
01212 if(tempInt <0) tempInt++;
01213 float energyTotal = EcalClusterTools::getSumEnergy( RH_ptrs );
01214
01215
01216 float centerIEta = 0.;
01217 float centerIPhi = 0.;
01218 float denominator = 0.;
01219
01220 for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01221
01222 EBDetId EBdetIdi( (*rh_ptr)->detid() );
01223 if(fabs(energyTotal) < 0.0001){
01224
01225 shapes.push_back( -2 );
01226 shapes.push_back( -2 );
01227 return shapes;
01228 }
01229 float weight = 0;
01230 if(fabs(weightedPositionMethod)<0.0001){
01231 weight = (*rh_ptr)->energy()/energyTotal;
01232 }else{
01233 weight = std::max(0.0, 4.2 + log((*rh_ptr)->energy()/energyTotal));
01234 }
01235 denominator += weight;
01236 centerIEta += weight*deltaIEta(seedPosition[0],EBdetIdi.ieta());
01237 centerIPhi += weight*deltaIPhi(seedPosition[1],EBdetIdi.iphi());
01238 }
01239 if(fabs(denominator) < 0.0001){
01240
01241 shapes.push_back( -2 );
01242 shapes.push_back( -2 );
01243 return shapes;
01244 }
01245 centerIEta = centerIEta / denominator;
01246 centerIPhi = centerIPhi / denominator;
01247
01248
01249
01250 TMatrixDSym inertia(2);
01251 double inertia00 = 0.;
01252 double inertia01 = 0.;
01253 double inertia11 = 0.;
01254 int i = 0;
01255 for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01256
01257 EBDetId EBdetIdi( (*rh_ptr)->detid() );
01258
01259 if(fabs(energyTotal) < 0.0001){
01260
01261 shapes.push_back( -2 );
01262 shapes.push_back( -2 );
01263 return shapes;
01264 }
01265 float weight = 0;
01266 if(fabs(weightedPositionMethod) < 0.0001){
01267 weight = (*rh_ptr)->energy()/energyTotal;
01268 }else{
01269 weight = std::max(0.0, 4.2 + log((*rh_ptr)->energy()/energyTotal));
01270 }
01271
01272 float ieta_rh_to_center = deltaIEta(seedPosition[0],EBdetIdi.ieta()) - centerIEta;
01273 float iphi_rh_to_center = deltaIPhi(seedPosition[1],EBdetIdi.iphi()) - centerIPhi;
01274
01275 inertia00 += weight*iphi_rh_to_center*iphi_rh_to_center;
01276 inertia01 -= weight*iphi_rh_to_center*ieta_rh_to_center;
01277 inertia11 += weight*ieta_rh_to_center*ieta_rh_to_center;
01278 i++;
01279 }
01280
01281 inertia[0][0] = inertia00;
01282 inertia[0][1] = inertia01;
01283 inertia[1][0] = inertia01;
01284 inertia[1][1] = inertia11;
01285
01286
01287
01288 TMatrixD eVectors(2,2);
01289 TVectorD eValues(2);
01290
01291 eVectors=inertia.EigenVectors(eValues);
01292
01293
01294
01295
01296
01297
01298 TVectorD smallerAxis(2);
01299 smallerAxis[0]=eVectors[0][1];
01300 smallerAxis[1]=eVectors[1][1];
01301
01302
01303 Double_t temp = fabs(smallerAxis[0]);
01304 if(fabs(eValues[0]) < 0.0001){
01305
01306 shapes.push_back( -2 );
01307 shapes.push_back( -2 );
01308 return shapes;
01309 }
01310
01311 float Roundness = eValues[1]/eValues[0];
01312 float Angle=acos(temp);
01313
01314 if( -0.00001 < Roundness && Roundness < 0) Roundness = 0.;
01315 if( -0.00001 < Angle && Angle < 0 ) Angle = 0.;
01316
01317 shapes.push_back( Roundness );
01318 shapes.push_back( Angle );
01319 return shapes;
01320
01321 }
01322
01323
01324
01325
01326 int EcalClusterTools::deltaIPhi(int seed_iphi, int rh_iphi){
01327 int rel_iphi = rh_iphi - seed_iphi;
01328
01329 if(rel_iphi > 180) rel_iphi = rel_iphi - 360;
01330 if(rel_iphi < -180) rel_iphi = rel_iphi + 360;
01331 return rel_iphi;
01332 }
01333
01334
01335
01336
01337 int EcalClusterTools::deltaIEta(int seed_ieta, int rh_ieta){
01338
01339 if(seed_ieta < 0) seed_ieta++;
01340 if(rh_ieta < 0) rh_ieta++;
01341 int rel_ieta = rh_ieta - seed_ieta;
01342 return rel_ieta;
01343 }
01344
01345
01346 std::vector<int> EcalClusterTools::getSeedPosition(std::vector<const EcalRecHit*> RH_ptrs){
01347 std::vector<int> seedPosition;
01348 float eSeedRH = 0;
01349 int iEtaSeedRH = 0;
01350 int iPhiSeedRH = 0;
01351
01352 for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01353
01354
01355 EBDetId EBdetIdi( (*rh_ptr)->detid() );
01356
01357 if(eSeedRH < (*rh_ptr)->energy()){
01358 eSeedRH = (*rh_ptr)->energy();
01359 iEtaSeedRH = EBdetIdi.ieta();
01360 iPhiSeedRH = EBdetIdi.iphi();
01361 }
01362
01363 }
01364
01365 seedPosition.push_back(iEtaSeedRH);
01366 seedPosition.push_back(iPhiSeedRH);
01367 return seedPosition;
01368 }
01369
01370
01371 float EcalClusterTools::getSumEnergy(std::vector<const EcalRecHit*> RH_ptrs){
01372
01373 float sumE = 0.;
01374
01375 for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01376 sumE += (*rh_ptr)->energy();
01377 }
01378
01379 return sumE;
01380 }