CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoEcal/EgammaCoreTools/src/EcalClusterTools.cc

Go to the documentation of this file.
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 float EcalClusterTools::getFraction( const std::vector< std::pair<DetId, float> > &v_id, DetId id
00019                           ){
00020   float frac = 1.0;
00021   for ( size_t i = 0; i < v_id.size(); ++i ) {
00022     if(v_id[i].first.rawId()==id.rawId()){
00023       frac=v_id[i].second;
00024     }
00025   }
00026   return frac;
00027 }
00028 
00029 
00030 std::pair<DetId, float> EcalClusterTools::getMaximum( const std::vector< std::pair<DetId, float> > &v_id, const EcalRecHitCollection *recHits)
00031 {
00032     float max = 0;
00033     DetId id(0);
00034     for ( size_t i = 0; i < v_id.size(); ++i ) {
00035         float energy = recHitEnergy( v_id[i].first, recHits ) * v_id[i].second;
00036         if ( energy > max ) {
00037             max = energy;
00038             id = v_id[i].first;
00039         }
00040     }
00041     return std::pair<DetId, float>(id, max);
00042 }
00043 
00044 std::pair<DetId, float> EcalClusterTools::getMaximum( const std::vector< std::pair<DetId, float> > &v_id, const EcalRecHitCollection *recHits,std::vector<int> flagsexcl,  std::vector<int> severitiesexcl, const  EcalSeverityLevelAlgo *sevLv)
00045 {
00046     float max = 0;
00047     DetId id(0);
00048     for ( size_t i = 0; i < v_id.size(); ++i ) {
00049         float energy = recHitEnergy( v_id[i].first, recHits,flagsexcl, severitiesexcl, sevLv ) * v_id[i].second;
00050         if ( energy > max ) {
00051             max = energy;
00052             id = v_id[i].first;
00053         }
00054     }
00055     return std::pair<DetId, float>(id, max);
00056 }
00057 
00058 
00059 std::pair<DetId, float> EcalClusterTools::getMaximum( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits)
00060 {
00061     return getMaximum( cluster.hitsAndFractions(), recHits );
00062 }
00063 
00064 std::pair<DetId, float> EcalClusterTools::getMaximum( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits,std::vector<int> flagsexcl,  std::vector<int> severitiesexcl, const  EcalSeverityLevelAlgo *sevLv )
00065 {
00066     return getMaximum( cluster.hitsAndFractions(), recHits, flagsexcl, severitiesexcl, sevLv );
00067 }
00068 
00069 
00070 float EcalClusterTools::recHitEnergy(DetId id, const EcalRecHitCollection *recHits)
00071 {
00072     if ( id == DetId(0) ) {
00073         return 0;
00074     } else {
00075         EcalRecHitCollection::const_iterator it = recHits->find( id );
00076         if ( it != recHits->end() ) {
00077             return (*it).energy();
00078         } else {
00079             //throw cms::Exception("EcalRecHitNotFound") << "The recHit corresponding to the DetId" << id.rawId() << " not found in the EcalRecHitCollection";
00080             // the recHit is not in the collection (hopefully zero suppressed)
00081             return 0;
00082         }
00083     }
00084     return 0;
00085 }
00086 
00087 float EcalClusterTools::recHitEnergy(DetId id, const EcalRecHitCollection *recHits, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const  EcalSeverityLevelAlgo *sevLv)
00088 {
00089     if ( id == DetId(0) ) {
00090         return 0;
00091     } else {
00092         EcalRecHitCollection::const_iterator it = recHits->find( id );
00093         if ( it != recHits->end() ) {
00094           // avoid anomalous channels (recoFlag based)
00095           uint32_t rhFlag = (*it).recoFlag();
00096           std::vector<int>::const_iterator vit = std::find( flagsexcl.begin(), flagsexcl.end(), rhFlag );
00097           //if your flag was found to be one which is excluded, zero out
00098           //this energy.
00099           if ( vit != flagsexcl.end() ) return 0;
00100             
00101           int severityFlag =  sevLv->severityLevel( it->id(), *recHits);
00102           std::vector<int>::const_iterator sit = std::find(severitiesexcl.begin(), severitiesexcl.end(), severityFlag);
00103           //if you were flagged by some condition (kWeird etc.)
00104           //zero out this energy.
00105           if (sit!= severitiesexcl.end())
00106             return 0; 
00107           //If we make it here, you're a found, clean hit.
00108           return (*it).energy();
00109         } else {
00110           //throw cms::Exception("EcalRecHitNotFound") << "The recHit corresponding to the DetId" << id.rawId() << " not found in the EcalRecHitCollection";
00111           // the recHit is not in the collection (hopefully zero suppressed)
00112           return 0;
00113         }
00114     }
00115     return 0;
00116 }
00117 
00118 
00119 // Returns the energy in a rectangle of crystals
00120 // specified in eta by ixMin and ixMax
00121 //       and in phi by iyMin and iyMax
00122 //
00123 // Reference picture (X=seed crystal)
00124 //    iy ___________
00125 //     2 |_|_|_|_|_|
00126 //     1 |_|_|_|_|_|
00127 //     0 |_|_|X|_|_|
00128 //    -1 |_|_|_|_|_|
00129 //    -2 |_|_|_|_|_|
00130 //      -2 -1 0 1 2 ix
00131 float EcalClusterTools::matrixEnergy( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, DetId id, int ixMin, int ixMax, int iyMin, int iyMax )
00132 {
00133   //take into account fractions
00134     // fast version
00135     CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00136     float energy = 0;
00137     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00138     for ( int i = ixMin; i <= ixMax; ++i ) {
00139         for ( int j = iyMin; j <= iyMax; ++j ) {
00140           cursor.home();
00141           cursor.offsetBy( i, j );
00142           float frac=getFraction(v_id,*cursor);
00143           energy += recHitEnergy( *cursor, recHits )*frac;
00144         }
00145     }
00146     // slow elegant version
00147     //float energy = 0;
00148     //std::vector<DetId> v_id = matrixDetId( topology, id, ixMin, ixMax, iyMin, iyMax );
00149     //for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00150     //        energy += recHitEnergy( *it, recHits );
00151     //}
00152     return energy;
00153 }
00154 
00155 float EcalClusterTools::matrixEnergy( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, DetId id, int ixMin, int ixMax, int iyMin, int iyMax, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00156 {
00157     // fast version
00158     CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00159     float energy = 0;
00160     for ( int i = ixMin; i <= ixMax; ++i ) {
00161         for ( int j = iyMin; j <= iyMax; ++j ) {
00162             cursor.home();
00163             cursor.offsetBy( i, j );
00164             energy += recHitEnergy( *cursor, recHits, flagsexcl, severitiesexcl, sevLv );
00165         }
00166     }
00167     return energy;
00168 }
00169 
00170 
00171 
00172 std::vector<DetId> EcalClusterTools::matrixDetId( const CaloTopology* topology, DetId id, int ixMin, int ixMax, int iyMin, int iyMax )
00173 {
00174     CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00175     std::vector<DetId> v;
00176     for ( int i = ixMin; i <= ixMax; ++i ) {
00177         for ( int j = iyMin; j <= iyMax; ++j ) {
00178             cursor.home();
00179             cursor.offsetBy( i, j );
00180             if ( *cursor != DetId(0) ) v.push_back( *cursor );
00181         }
00182     }
00183     return v;
00184 }
00185 
00186 
00187 
00188 float EcalClusterTools::e2x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00189 {
00190     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00191     std::list<float> energies;
00192     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 0 ) );
00193     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0,  0, 1 ) );
00194     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1,  0, 1 ) );
00195     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1, -1, 0 ) );
00196 
00197 
00198     return *std::max_element(energies.begin(),energies.end());
00199 
00200 }
00201 
00202 
00203 float EcalClusterTools::e2x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00204 {
00205     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00206     std::list<float> energies;
00207     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 0,flagsexcl, severitiesexcl, sevLv ) );
00208     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0,  0, 1,flagsexcl, severitiesexcl, sevLv ) );
00209     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1,  0, 1,flagsexcl, severitiesexcl, sevLv ) );
00210     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1, -1, 0,flagsexcl, severitiesexcl, sevLv ) );
00211 
00212 
00213     return *std::max_element(energies.begin(),energies.end());
00214 
00215 }
00216 
00217 
00218 float EcalClusterTools::e3x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00219 {
00220     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00221     std::list<float> energies;
00222     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 0 ) );
00223     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1, -1, 1 ) );
00224     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1,  0, 1 ) );
00225     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 1 ) );
00226     return *std::max_element(energies.begin(),energies.end());
00227 }
00228 
00229 float EcalClusterTools::e3x2( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00230 {
00231     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00232     std::list<float> energies;
00233     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 0,flagsexcl, severitiesexcl, sevLv ) );
00234     energies.push_back( matrixEnergy( cluster, recHits, topology, id,  0, 1, -1, 1,flagsexcl, severitiesexcl, sevLv ) );
00235     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 1,  0, 1,flagsexcl, severitiesexcl, sevLv ) );
00236     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 0, -1, 1,flagsexcl, severitiesexcl, sevLv ) );
00237     return *std::max_element(energies.begin(),energies.end());
00238 }
00239 
00240 float EcalClusterTools::e3x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00241 {
00242     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00243     return matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 1 );
00244 }
00245 
00246 
00247 float EcalClusterTools::e3x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00248 {
00249     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00250     return matrixEnergy( cluster, recHits, topology, id, -1, 1, -1, 1,flagsexcl, severitiesexcl, sevLv );
00251 }
00252 
00253 
00254 float EcalClusterTools::e4x4( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00255 {
00256     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00257     std::list<float> energies;
00258     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -2, 1 ) );
00259     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -2, 1 ) );
00260     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -1, 2 ) );
00261     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -1, 2 ) );
00262     return *std::max_element(energies.begin(),energies.end());
00263 }
00264 
00265 float EcalClusterTools::e4x4( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00266 {
00267     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00268     std::list<float> energies;
00269     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -2, 1,flagsexcl, severitiesexcl, sevLv ) );
00270     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -2, 1,flagsexcl, severitiesexcl, sevLv ) );
00271     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -2, 1, -1, 2,flagsexcl, severitiesexcl, sevLv ) );
00272     energies.push_back( matrixEnergy( cluster, recHits, topology, id, -1, 2, -1, 2,flagsexcl, severitiesexcl, sevLv ) );
00273     return *std::max_element(energies.begin(),energies.end());
00274 }
00275 
00276 
00277 
00278 float EcalClusterTools::e5x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00279 {
00280     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00281     return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, 2 );
00282 }
00283 
00284 float EcalClusterTools::e5x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00285 {
00286     DetId id = getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).first;
00287     return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, 2,flagsexcl, severitiesexcl, sevLv );
00288 }
00289 
00290 float EcalClusterTools::eMax( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00291 {
00292     return getMaximum( cluster.hitsAndFractions(), recHits ).second;
00293 }
00294 
00295 float EcalClusterTools::eMax( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv  )
00296 {
00297     return getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).second;
00298 }
00299 
00300 
00301 float EcalClusterTools::e2nd( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00302 {
00303     std::list<float> energies;
00304     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00305     if ( v_id.size() < 2 ) return 0;
00306     for ( size_t i = 0; i < v_id.size(); ++i ) {
00307         energies.push_back( recHitEnergy( v_id[i].first, recHits ) * v_id[i].second );
00308     }
00309     energies.sort();             
00310     return *--(--energies.end());
00311 
00312 
00313 }
00314 
00315 float EcalClusterTools::e2nd( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00316 {
00317     std::list<float> energies;
00318     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00319     if ( v_id.size() < 2 ) return 0;
00320     for ( size_t i = 0; i < v_id.size(); ++i ) {
00321         energies.push_back( recHitEnergy( v_id[i].first, recHits,flagsexcl, severitiesexcl, sevLv ) * v_id[i].second );
00322     }
00323     energies.sort();             
00324     return *--(--energies.end());
00325 
00326 
00327 }
00328 
00329 
00330 
00331 float EcalClusterTools::e2x5Right( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00332 {
00333     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00334     return matrixEnergy( cluster, recHits, topology, id, 1, 2, -2, 2 );
00335 }
00336 
00337 
00338 float EcalClusterTools::e2x5Right( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00339 {
00340     DetId id = getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).first;
00341     return matrixEnergy( cluster, recHits, topology, id, 1, 2, -2, 2,flagsexcl, severitiesexcl, sevLv );
00342 }
00343 
00344 
00345 float EcalClusterTools::e2x5Left( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00346 {
00347     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00348     return matrixEnergy( cluster, recHits, topology, id, -2, -1, -2, 2 );
00349 }
00350 
00351 float EcalClusterTools::e2x5Left( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00352 {
00353     DetId id = getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).first;
00354     return matrixEnergy( cluster, recHits, topology, id, -2, -1, -2, 2,flagsexcl, severitiesexcl, sevLv );
00355 }
00356 
00357 
00358 // 
00359 float EcalClusterTools::e2x5Top( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00360 {
00361     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00362     return matrixEnergy( cluster, recHits, topology, id, -2, 2, 1, 2 );
00363 }
00364 
00365 float EcalClusterTools::e2x5Top( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00366 {
00367     DetId id = getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).first;
00368     return matrixEnergy( cluster, recHits, topology, id, -2, 2, 1, 2,flagsexcl, severitiesexcl, sevLv );
00369 }
00370 
00371 
00372 float EcalClusterTools::e2x5Bottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00373 {
00374     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00375     return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, -1 );
00376 }
00377 
00378 float EcalClusterTools::e2x5Bottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00379 {
00380     DetId id = getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv  ).first;
00381     return matrixEnergy( cluster, recHits, topology, id, -2, 2, -2, -1,flagsexcl, severitiesexcl, sevLv );
00382 }
00383 
00384 // Energy in 2x5 strip containing the max crystal.
00385 // Adapted from code by Sam Harper
00386 float EcalClusterTools::e2x5Max( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00387 {
00388     DetId id =      getMaximum( cluster.hitsAndFractions(), recHits ).first;
00389 
00390     // 1x5 strip left of seed
00391     float left   = matrixEnergy( cluster, recHits, topology, id, -1, -1, -2, 2 );
00392     // 1x5 strip right of seed
00393     float right  = matrixEnergy( cluster, recHits, topology, id,  1,  1, -2, 2 );
00394     // 1x5 strip containing seed
00395     float centre = matrixEnergy( cluster, recHits, topology, id,  0,  0, -2, 2 );
00396 
00397     // Return the maximum of (left+center) or (right+center) strip
00398     return left > right ? left+centre : right+centre;
00399 }
00400 
00401 float EcalClusterTools::e2x5Max( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00402 {
00403     DetId id =      getMaximum( cluster.hitsAndFractions(), recHits,flagsexcl, severitiesexcl, sevLv ).first;
00404 
00405     // 1x5 strip left of seed
00406     float left   = matrixEnergy( cluster, recHits, topology, id, -1, -1, -2, 2,flagsexcl, severitiesexcl, sevLv );
00407     // 1x5 strip right of seed
00408     float right  = matrixEnergy( cluster, recHits, topology, id,  1,  1, -2, 2,flagsexcl, severitiesexcl, sevLv );
00409     // 1x5 strip containing seed
00410     float centre = matrixEnergy( cluster, recHits, topology, id,  0,  0, -2, 2,flagsexcl, severitiesexcl, sevLv );
00411 
00412     // Return the maximum of (left+center) or (right+center) strip
00413     return left > right ? left+centre : right+centre;
00414 }
00415 
00416 
00417 float EcalClusterTools::e1x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00418 {
00419     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00420     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -2, 2 );
00421 }
00422 
00423 float EcalClusterTools::e1x5( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv)
00424 {
00425     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00426     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -2, 2 ,
00427                          flagsexcl, severitiesexcl, sevLv);
00428 }
00429 
00430 
00431 
00432 float EcalClusterTools::e5x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00433 {
00434     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00435     return matrixEnergy( cluster, recHits, topology, id, -2, 2, 0, 0 );
00436 }
00437 
00438 float EcalClusterTools::e5x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00439 {
00440     DetId id = getMaximum( cluster.hitsAndFractions(), recHits).first;
00441     return matrixEnergy( cluster, recHits, topology, id, -2, 2, 0, 0,flagsexcl, severitiesexcl, sevLv );
00442 }
00443 
00444 
00445 float EcalClusterTools::e1x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00446 {
00447     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00448     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, 1 );
00449 }
00450 
00451 float EcalClusterTools::e1x3( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00452 {
00453     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00454     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, 1,flagsexcl, severitiesexcl, sevLv );
00455 }
00456 
00457 
00458 float EcalClusterTools::e3x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00459 {
00460     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00461     return matrixEnergy( cluster, recHits, topology, id, -1, 1, 0, 0 );
00462 }
00463 
00464 float EcalClusterTools::e3x1( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00465 {
00466     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00467     return matrixEnergy( cluster, recHits, topology, id, -1, 1, 0, 0,flagsexcl, severitiesexcl, sevLv );
00468 }
00469 
00470 
00471 float EcalClusterTools::eLeft( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00472 {
00473     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00474     return matrixEnergy( cluster, recHits, topology, id, -1, -1, 0, 0 );
00475 }
00476 
00477 float EcalClusterTools::eLeft( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00478 {
00479     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00480     return matrixEnergy( cluster, recHits, topology, id, -1, -1, 0, 0,flagsexcl, severitiesexcl, sevLv );
00481 }
00482 
00483 
00484 float EcalClusterTools::eRight( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00485 {
00486     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00487     return matrixEnergy( cluster, recHits, topology, id, 1, 1, 0, 0 );
00488 }
00489 
00490 float EcalClusterTools::eRight( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00491 {
00492     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00493     return matrixEnergy( cluster, recHits, topology, id, 1, 1, 0, 0,flagsexcl, severitiesexcl, sevLv );
00494 }
00495 
00496 
00497 float EcalClusterTools::eTop( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00498 {
00499     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00500     return matrixEnergy( cluster, recHits, topology, id, 0, 0, 1, 1 );
00501 }
00502 
00503 float EcalClusterTools::eTop( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00504 {
00505     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00506     return matrixEnergy( cluster, recHits, topology, id, 0, 0, 1, 1,flagsexcl, severitiesexcl, sevLv );
00507 }
00508 
00509 
00510 
00511 float EcalClusterTools::eBottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology )
00512 {
00513     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00514     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, -1 );
00515 }
00516 
00517 float EcalClusterTools::eBottom( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology* topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00518 {
00519     DetId id = getMaximum( cluster.hitsAndFractions(), recHits ).first;
00520     return matrixEnergy( cluster, recHits, topology, id, 0, 0, -1, -1,flagsexcl, severitiesexcl, sevLv  );
00521 }
00522 
00523 std::vector<float> EcalClusterTools::energyBasketFractionEta( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00524 {
00525     std::vector<float> basketFraction( 2 * EBDetId::kModulesPerSM );
00526     float clusterEnergy = cluster.energy();
00527     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00528     if ( v_id[0].first.subdetId() != EcalBarrel ) {
00529         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.";
00530         return basketFraction;
00531     }
00532     for ( size_t i = 0; i < v_id.size(); ++i ) {
00533         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;
00534     }
00535     std::sort( basketFraction.rbegin(), basketFraction.rend() );
00536     return basketFraction;
00537 }
00538 
00539 
00540 std::vector<float> EcalClusterTools::energyBasketFractionEta( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00541 {
00542     std::vector<float> basketFraction( 2 * EBDetId::kModulesPerSM );
00543     float clusterEnergy = cluster.energy();
00544     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00545     if ( v_id[0].first.subdetId() != EcalBarrel ) {
00546         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.";
00547         return basketFraction;
00548     }
00549     for ( size_t i = 0; i < v_id.size(); ++i ) {
00550         basketFraction[ EBDetId(v_id[i].first).im()-1 + EBDetId(v_id[i].first).positiveZ()*EBDetId::kModulesPerSM ] += recHitEnergy( v_id[i].first, recHits,flagsexcl, severitiesexcl, sevLv ) * v_id[i].second / clusterEnergy;
00551     }
00552     std::sort( basketFraction.rbegin(), basketFraction.rend() );
00553     return basketFraction;
00554 }
00555 
00556 std::vector<float> EcalClusterTools::energyBasketFractionPhi( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits )
00557 {
00558     std::vector<float> basketFraction( 2 * (EBDetId::MAX_IPHI / EBDetId::kCrystalsInPhi) );
00559     float clusterEnergy = cluster.energy();
00560     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00561     if ( v_id[0].first.subdetId() != EcalBarrel ) {
00562         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.";
00563         return basketFraction;
00564     }
00565     for ( size_t i = 0; i < v_id.size(); ++i ) {
00566         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;
00567     }
00568     std::sort( basketFraction.rbegin(), basketFraction.rend() );
00569     return basketFraction;
00570 }
00571 
00572 
00573 std::vector<float> EcalClusterTools::energyBasketFractionPhi( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00574 {
00575     std::vector<float> basketFraction( 2 * (EBDetId::MAX_IPHI / EBDetId::kCrystalsInPhi) );
00576     float clusterEnergy = cluster.energy();
00577     std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00578     if ( v_id[0].first.subdetId() != EcalBarrel ) {
00579         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.";
00580         return basketFraction;
00581     }
00582     for ( size_t i = 0; i < v_id.size(); ++i ) {
00583         basketFraction[ (EBDetId(v_id[i].first).iphi()-1)/EBDetId::kCrystalsInPhi + EBDetId(v_id[i].first).positiveZ()*EBDetId::kTowersInPhi] += recHitEnergy( v_id[i].first, recHits,flagsexcl, severitiesexcl, sevLv ) * v_id[i].second / clusterEnergy;
00584     }
00585     std::sort( basketFraction.rbegin(), basketFraction.rend() );
00586     return basketFraction;
00587 }
00588 
00589 
00590 std::vector<EcalClusterTools::EcalClusterEnergyDeposition> EcalClusterTools::getEnergyDepTopology( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW, float w0 )
00591 {
00592     std::vector<EcalClusterTools::EcalClusterEnergyDeposition> energyDistribution;
00593     // init a map of the energy deposition centered on the
00594     // cluster centroid. This is for momenta calculation only.
00595     CLHEP::Hep3Vector clVect(cluster.position().x(), cluster.position().y(), cluster.position().z());
00596     CLHEP::Hep3Vector clDir(clVect);
00597     clDir*=1.0/clDir.mag();
00598     // in the transverse plane, axis perpendicular to clusterDir
00599     CLHEP::Hep3Vector theta_axis(clDir.y(),-clDir.x(),0.0);
00600     theta_axis *= 1.0/theta_axis.mag();
00601     CLHEP::Hep3Vector phi_axis = theta_axis.cross(clDir);
00602 
00603     std::vector< std::pair<DetId, float> > clusterDetIds = cluster.hitsAndFractions();
00604 
00605     EcalClusterEnergyDeposition clEdep;
00606     EcalRecHit testEcalRecHit;
00607     std::vector< std::pair<DetId, float> >::iterator posCurrent;
00608     // loop over crystals
00609     for(posCurrent=clusterDetIds.begin(); posCurrent!=clusterDetIds.end(); ++posCurrent) {
00610         EcalRecHitCollection::const_iterator itt = recHits->find( (*posCurrent).first );
00611         testEcalRecHit=*itt;
00612 
00613         if(( (*posCurrent).first != DetId(0)) && (recHits->find( (*posCurrent).first ) != recHits->end())) {
00614             clEdep.deposited_energy = testEcalRecHit.energy() * (*posCurrent).second;
00615             // if logarithmic weight is requested, apply cut on minimum energy of the recHit
00616             if(logW) {
00617                 //double w0 = parameterMap_.find("W0")->second;
00618 
00619                 double weight = std::max(0.0, w0 + log(fabs(clEdep.deposited_energy)/cluster.energy()) );
00620                 if(weight==0) {
00621                     LogDebug("ClusterShapeAlgo") << "Crystal has insufficient energy: E = " 
00622                         << clEdep.deposited_energy << " GeV; skipping... ";
00623                     continue;
00624                 }
00625                 else LogDebug("ClusterShapeAlgo") << "===> got crystal. Energy = " << clEdep.deposited_energy << " GeV. ";
00626             }
00627             DetId id_ = (*posCurrent).first;
00628             const CaloCellGeometry *this_cell = geometry->getSubdetectorGeometry(id_)->getGeometry(id_);
00629             GlobalPoint cellPos = this_cell->getPosition();
00630             CLHEP::Hep3Vector gblPos (cellPos.x(),cellPos.y(),cellPos.z()); //surface position?
00631             // Evaluate the distance from the cluster centroid
00632             CLHEP::Hep3Vector diff = gblPos - clVect;
00633             // Important: for the moment calculation, only the "lateral distance" is important
00634             // "lateral distance" r_i = distance of the digi position from the axis Origin-Cluster Center
00635             // ---> subtract the projection on clDir
00636             CLHEP::Hep3Vector DigiVect = diff - diff.dot(clDir)*clDir;
00637             clEdep.r = DigiVect.mag();
00638             LogDebug("ClusterShapeAlgo") << "E = " << clEdep.deposited_energy
00639                 << "\tdiff = " << diff.mag()
00640                 << "\tr = " << clEdep.r;
00641             clEdep.phi = DigiVect.angle(theta_axis);
00642             if(DigiVect.dot(phi_axis)<0) clEdep.phi = 2 * M_PI - clEdep.phi;
00643             energyDistribution.push_back(clEdep);
00644         }
00645     } 
00646     return energyDistribution;
00647 }
00648 
00649 
00650 
00651 std::vector<float> EcalClusterTools::lat( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, bool logW, float w0 )
00652 {
00653     std::vector<EcalClusterTools::EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
00654 
00655     std::vector<float> lat;
00656     double r, redmoment=0;
00657     double phiRedmoment = 0 ;
00658     double etaRedmoment = 0 ;
00659     int n,n1,n2,tmp;
00660     int clusterSize=energyDistribution.size();
00661     float etaLat_, phiLat_, lat_;
00662     if (clusterSize<3) {
00663         etaLat_ = 0.0 ; 
00664         lat_ = 0.0;
00665         lat.push_back(0.);
00666         lat.push_back(0.);
00667         lat.push_back(0.);
00668         return lat; 
00669     }
00670 
00671     n1=0; n2=1;
00672     if (energyDistribution[1].deposited_energy > 
00673             energyDistribution[0].deposited_energy) 
00674     {
00675         tmp=n2; n2=n1; n1=tmp;
00676     }
00677     for (int i=2; i<clusterSize; i++) {
00678         n=i;
00679         if (energyDistribution[i].deposited_energy > 
00680                 energyDistribution[n1].deposited_energy) 
00681         {
00682             tmp = n2;
00683             n2 = n1; n1 = i; n=tmp;
00684         } else {
00685             if (energyDistribution[i].deposited_energy > 
00686                     energyDistribution[n2].deposited_energy) 
00687             {
00688                 tmp=n2; n2=i; n=tmp;
00689             }
00690         }
00691 
00692         r = energyDistribution[n].r;
00693         redmoment += r*r* energyDistribution[n].deposited_energy;
00694         double rphi = r * cos (energyDistribution[n].phi) ;
00695         phiRedmoment += rphi * rphi * energyDistribution[n].deposited_energy;
00696         double reta = r * sin (energyDistribution[n].phi) ;
00697         etaRedmoment += reta * reta * energyDistribution[n].deposited_energy;
00698     } 
00699     double e1 = energyDistribution[n1].deposited_energy;
00700     double e2 = energyDistribution[n2].deposited_energy;
00701 
00702     lat_ = redmoment/(redmoment+2.19*2.19*(e1+e2));
00703     phiLat_ = phiRedmoment/(phiRedmoment+2.19*2.19*(e1+e2));
00704     etaLat_ = etaRedmoment/(etaRedmoment+2.19*2.19*(e1+e2));
00705 
00706     lat.push_back(etaLat_);
00707     lat.push_back(phiLat_);
00708     lat.push_back(lat_);
00709     return lat;
00710 }
00711 
00712 
00713 
00714 math::XYZVector EcalClusterTools::meanClusterPosition( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, const CaloGeometry *geometry )
00715 {
00716     // find mean energy position of a 5x5 cluster around the maximum
00717     math::XYZVector meanPosition(0.0, 0.0, 0.0);
00718     std::vector<DetId> v_id = matrixDetId( topology, getMaximum( cluster, recHits ).first, -2, 2, -2, 2 );
00719     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00720         GlobalPoint positionGP = geometry->getSubdetectorGeometry( *it )->getGeometry( *it )->getPosition();
00721         math::XYZVector position(positionGP.x(),positionGP.y(),positionGP.z());
00722         meanPosition = meanPosition + recHitEnergy( *it, recHits ) * position;
00723     }
00724     return meanPosition / e5x5( cluster, recHits, topology );
00725 }
00726 
00727 //================================================= meanClusterPosition===================================================================================
00728 
00729 math::XYZVector EcalClusterTools::meanClusterPosition( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloTopology *topology, const CaloGeometry *geometry, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv )
00730 {
00731     // find mean energy position of a 5x5 cluster around the maximum
00732     math::XYZVector meanPosition(0.0, 0.0, 0.0);
00733     std::vector<DetId> v_id = matrixDetId( topology, getMaximum( cluster, recHits ).first, -2, 2, -2, 2 );
00734     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {
00735         GlobalPoint positionGP = geometry->getSubdetectorGeometry( *it )->getGeometry( *it )->getPosition();
00736         math::XYZVector position(positionGP.x(),positionGP.y(),positionGP.z());
00737         meanPosition = meanPosition + recHitEnergy( *it, recHits,flagsexcl, severitiesexcl, sevLv ) * position;
00738     }
00739     return meanPosition / e5x5( cluster, recHits, topology,flagsexcl, severitiesexcl, sevLv );
00740 }
00741 
00742 
00743 //returns mean energy weighted eta/phi in crystals from the seed
00744 //iPhi is not defined for endcap and is returned as zero
00745 //return <eta,phi>
00746 //we have an issue in working out what to do for negative energies
00747 //I (Sam Harper) think it makes sense to ignore crystals with E<0 in the calculation as they are ignored
00748 //in the sigmaIEtaIEta calculation (well they arent fully ignored, they do still contribute to the e5x5 sum
00749 //in the sigmaIEtaIEta calculation but not here)
00750 std::pair<float,float>  EcalClusterTools::mean5x5PositionInLocalCrysCoord(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology)
00751 {
00752     DetId seedId =  getMaximum( cluster, recHits ).first;
00753     float meanDEta=0.;
00754     float meanDPhi=0.;
00755     float energySum=0.;
00756 
00757     std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00758     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {  
00759         float energy = recHitEnergy(*it,recHits);
00760         if(energy<0.) continue;//skipping negative energy crystals
00761         meanDEta += energy * getNrCrysDiffInEta(*it,seedId);
00762         meanDPhi += energy * getNrCrysDiffInPhi(*it,seedId);    
00763         energySum +=energy;
00764     }
00765     meanDEta /=energySum;
00766     meanDPhi /=energySum;
00767     return std::pair<float,float>(meanDEta,meanDPhi);
00768 }
00769 
00770 std::pair<float,float>  EcalClusterTools::mean5x5PositionInLocalCrysCoord(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv)
00771 {
00772     DetId seedId =  getMaximum( cluster, recHits ).first;
00773     float meanDEta=0.;
00774     float meanDPhi=0.;
00775     float energySum=0.;
00776 
00777     std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00778     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {  
00779         float energy = recHitEnergy(*it,recHits,flagsexcl, severitiesexcl, sevLv);
00780         if(energy<0.) continue;//skipping negative energy crystals
00781         meanDEta += energy * getNrCrysDiffInEta(*it,seedId);
00782         meanDPhi += energy * getNrCrysDiffInPhi(*it,seedId);    
00783         energySum +=energy;
00784     }
00785     meanDEta /=energySum;
00786     meanDPhi /=energySum;
00787     return std::pair<float,float>(meanDEta,meanDPhi);
00788 }
00789 
00790 //returns mean energy weighted x/y in normalised crystal coordinates
00791 //only valid for endcap, returns 0,0 for barrel
00792 //we have an issue in working out what to do for negative energies
00793 //I (Sam Harper) think it makes sense to ignore crystals with E<0 in the calculation as they are ignored
00794 //in the sigmaIEtaIEta calculation (well they arent fully ignored, they do still contribute to the e5x5 sum
00795 //in the sigmaIEtaIEta calculation but not here)
00796 std::pair<float,float> EcalClusterTools::mean5x5PositionInXY(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology)
00797 {
00798     DetId seedId =  getMaximum( cluster, recHits ).first;
00799 
00800     std::pair<float,float> meanXY(0.,0.);
00801     if(seedId.subdetId()==EcalBarrel) return meanXY;
00802 
00803     float energySum=0.;
00804 
00805     std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00806     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {  
00807         float energy = recHitEnergy(*it,recHits);
00808         if(energy<0.) continue;//skipping negative energy crystals
00809         meanXY.first += energy * getNormedIX(*it);
00810         meanXY.second += energy * getNormedIY(*it);
00811         energySum +=energy;
00812     }
00813     meanXY.first/=energySum;
00814     meanXY.second/=energySum;
00815     return meanXY;
00816 }
00817 
00818 std::pair<float,float> EcalClusterTools::mean5x5PositionInXY(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology, std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv)
00819 {
00820     DetId seedId =  getMaximum( cluster, recHits ).first;
00821 
00822     std::pair<float,float> meanXY(0.,0.);
00823     if(seedId.subdetId()==EcalBarrel) return meanXY;
00824 
00825     float energySum=0.;
00826 
00827     std::vector<DetId> v_id = matrixDetId( topology,seedId, -2, 2, -2, 2 );
00828     for ( std::vector<DetId>::const_iterator it = v_id.begin(); it != v_id.end(); ++it ) {  
00829         float energy = recHitEnergy(*it,recHits,flagsexcl, severitiesexcl, sevLv);
00830         if(energy<0.) continue;//skipping negative energy crystals
00831         meanXY.first += energy * getNormedIX(*it);
00832         meanXY.second += energy * getNormedIY(*it);
00833         energySum +=energy;
00834     }
00835     meanXY.first/=energySum;
00836     meanXY.second/=energySum;
00837     return meanXY;
00838 }
00839 
00840 
00841 std::vector<float> EcalClusterTools::covariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits, const CaloTopology *topology, const CaloGeometry* geometry, float w0)
00842 {
00843     float e_5x5 = e5x5( cluster, recHits, topology );
00844     float covEtaEta, covEtaPhi, covPhiPhi;
00845     if (e_5x5 >= 0.) {
00846         //double w0_ = parameterMap_.find("W0")->second;
00847         std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00848         math::XYZVector meanPosition = meanClusterPosition( cluster, recHits, topology, geometry );
00849 
00850         // now we can calculate the covariances
00851         double numeratorEtaEta = 0;
00852         double numeratorEtaPhi = 0;
00853         double numeratorPhiPhi = 0;
00854         double denominator     = 0;
00855 
00856         DetId id = getMaximum( v_id, recHits ).first;
00857         CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00858         for ( int i = -2; i <= 2; ++i ) {
00859             for ( int j = -2; j <= 2; ++j ) {
00860                 cursor.home();
00861                 cursor.offsetBy( i, j );
00862                 float energy = recHitEnergy( *cursor, recHits );
00863 
00864                 if ( energy <= 0 ) continue;
00865 
00866                 GlobalPoint position = geometry->getSubdetectorGeometry(*cursor)->getGeometry(*cursor)->getPosition();
00867 
00868                 double dPhi = position.phi() - meanPosition.phi();
00869                 if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00870                 if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() + dPhi; }
00871 
00872                 double dEta = position.eta() - meanPosition.eta();
00873                 double w = 0.;
00874                 w = std::max(0.0, w0 + log( energy / e_5x5 ));
00875 
00876                 denominator += w;
00877                 numeratorEtaEta += w * dEta * dEta;
00878                 numeratorEtaPhi += w * dEta * dPhi;
00879                 numeratorPhiPhi += w * dPhi * dPhi;
00880             }
00881         }
00882 
00883         if (denominator != 0.0) {
00884             covEtaEta =  numeratorEtaEta / denominator;
00885             covEtaPhi =  numeratorEtaPhi / denominator;
00886             covPhiPhi =  numeratorPhiPhi / denominator;
00887         } else {
00888             covEtaEta = 999.9;
00889             covEtaPhi = 999.9;
00890             covPhiPhi = 999.9;
00891         }
00892 
00893     } else {
00894         // Warn the user if there was no energy in the cells and return zeroes.
00895         //       std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
00896         covEtaEta = 0;
00897         covEtaPhi = 0;
00898         covPhiPhi = 0;
00899     }
00900     std::vector<float> v;
00901     v.push_back( covEtaEta );
00902     v.push_back( covEtaPhi );
00903     v.push_back( covPhiPhi );
00904     return v;
00905 }
00906 
00907 //==================================================== Covariances===========================================================================
00908 
00909 std::vector<float> EcalClusterTools::covariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits, const CaloTopology *topology, const CaloGeometry* geometry,std::vector<int> flagsexcl,std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0)
00910 {
00911     float e_5x5 = e5x5( cluster, recHits, topology,flagsexcl, severitiesexcl, sevLv );
00912     float covEtaEta, covEtaPhi, covPhiPhi;
00913     if (e_5x5 >= 0.) {
00914         //double w0_ = parameterMap_.find("W0")->second;
00915         std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00916         math::XYZVector meanPosition = meanClusterPosition( cluster, recHits, topology, geometry,flagsexcl, severitiesexcl, sevLv );
00917 
00918         // now we can calculate the covariances
00919         double numeratorEtaEta = 0;
00920         double numeratorEtaPhi = 0;
00921         double numeratorPhiPhi = 0;
00922         double denominator     = 0;
00923 
00924         DetId id = getMaximum( v_id, recHits ).first;
00925         CaloNavigator<DetId> cursor = CaloNavigator<DetId>( id, topology->getSubdetectorTopology( id ) );
00926         for ( int i = -2; i <= 2; ++i ) {
00927             for ( int j = -2; j <= 2; ++j ) {
00928                 cursor.home();
00929                 cursor.offsetBy( i, j );
00930                 float energy = recHitEnergy( *cursor, recHits,flagsexcl, severitiesexcl, sevLv );
00931 
00932                 if ( energy <= 0 ) continue;
00933 
00934                 GlobalPoint position = geometry->getSubdetectorGeometry(*cursor)->getGeometry(*cursor)->getPosition();
00935 
00936                 double dPhi = position.phi() - meanPosition.phi();
00937                 if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00938                 if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() + dPhi; }
00939 
00940                 double dEta = position.eta() - meanPosition.eta();
00941                 double w = 0.;
00942                 w = std::max(0.0, w0 + log( energy / e_5x5 ));
00943 
00944                 denominator += w;
00945                 numeratorEtaEta += w * dEta * dEta;
00946                 numeratorEtaPhi += w * dEta * dPhi;
00947                 numeratorPhiPhi += w * dPhi * dPhi;
00948             }
00949         }
00950 
00951         if (denominator != 0.0) {
00952             covEtaEta =  numeratorEtaEta / denominator;
00953             covEtaPhi =  numeratorEtaPhi / denominator;
00954             covPhiPhi =  numeratorPhiPhi / denominator;
00955         } else {
00956             covEtaEta = 999.9;
00957             covEtaPhi = 999.9;
00958             covPhiPhi = 999.9;
00959         }
00960 
00961     } else {
00962         // Warn the user if there was no energy in the cells and return zeroes.
00963         //       std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
00964         covEtaEta = 0;
00965         covEtaPhi = 0;
00966         covPhiPhi = 0;
00967     }
00968     std::vector<float> v;
00969     v.push_back( covEtaEta );
00970     v.push_back( covEtaPhi );
00971     v.push_back( covPhiPhi );
00972     return v;
00973 }
00974 
00975 
00976 
00977 
00978 
00979 //for covIEtaIEta,covIEtaIPhi and covIPhiIPhi are defined but only covIEtaIEta has been actively studied
00980 //instead of using absolute eta/phi it counts crystals normalised so that it gives identical results to normal covariances except near the cracks where of course its better 
00981 //it also does not require any eta correction function in the endcap
00982 //it is multipled by an approprate crystal size to ensure it gives similar values to covariances(...)
00983 std::vector<float> EcalClusterTools::localCovariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology,float w0)
00984 {
00985 
00986     float e_5x5 = e5x5( cluster, recHits, topology );
00987     float covEtaEta, covEtaPhi, covPhiPhi;
00988 
00989     if (e_5x5 >= 0.) {
00990         //double w0_ = parameterMap_.find("W0")->second;
00991         std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
00992         std::pair<float,float> mean5x5PosInNrCrysFromSeed =  mean5x5PositionInLocalCrysCoord( cluster, recHits, topology );
00993         std::pair<float,float> mean5x5XYPos =  mean5x5PositionInXY(cluster,recHits,topology);
00994 
00995         // now we can calculate the covariances
00996         double numeratorEtaEta = 0;
00997         double numeratorEtaPhi = 0;
00998         double numeratorPhiPhi = 0;
00999         double denominator     = 0;
01000 
01001         //these allow us to scale the localCov by the crystal size 
01002         //so that the localCovs have the same average value as the normal covs
01003         const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
01004         const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
01005 
01006         DetId seedId = getMaximum( v_id, recHits ).first;
01007 
01008         bool isBarrel=seedId.subdetId()==EcalBarrel;
01009         const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
01010 
01011         CaloNavigator<DetId> cursor = CaloNavigator<DetId>( seedId, topology->getSubdetectorTopology( seedId ) );
01012 
01013         for ( int eastNr = -2; eastNr <= 2; ++eastNr ) { //east is eta in barrel
01014             for ( int northNr = -2; northNr <= 2; ++northNr ) { //north is phi in barrel
01015                 cursor.home();
01016                 cursor.offsetBy( eastNr, northNr);
01017                 float energy = recHitEnergy( *cursor, recHits );
01018                 if ( energy <= 0 ) continue;
01019 
01020                 float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
01021                 float dPhi = 0;
01022 
01023                 if(isBarrel)  dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
01024                 else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
01025 
01026 
01027                 double w = std::max(0.0, w0 + log( energy / e_5x5 ));
01028 
01029                 denominator += w;
01030                 numeratorEtaEta += w * dEta * dEta;
01031                 numeratorEtaPhi += w * dEta * dPhi;
01032                 numeratorPhiPhi += w * dPhi * dPhi;
01033             } //end east loop
01034         }//end north loop
01035 
01036 
01037         //multiplying by crysSize to make the values compariable to normal covariances
01038         if (denominator != 0.0) {
01039             covEtaEta =  crysSize*crysSize* numeratorEtaEta / denominator;
01040             covEtaPhi =  crysSize*crysSize* numeratorEtaPhi / denominator;
01041             covPhiPhi =  crysSize*crysSize* numeratorPhiPhi / denominator;
01042         } else {
01043             covEtaEta = 999.9;
01044             covEtaPhi = 999.9;
01045             covPhiPhi = 999.9;
01046         }
01047 
01048 
01049     } else {
01050         // Warn the user if there was no energy in the cells and return zeroes.
01051         //       std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
01052         covEtaEta = 0;
01053         covEtaPhi = 0;
01054         covPhiPhi = 0;
01055     }
01056     std::vector<float> v;
01057     v.push_back( covEtaEta );
01058     v.push_back( covEtaPhi );
01059     v.push_back( covPhiPhi );
01060     return v;
01061 }
01062 
01063 //==================================================================localCovariances======================================================================
01064 
01065 std::vector<float> EcalClusterTools::localCovariances(const reco::BasicCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology,std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv,float w0)
01066 {
01067 
01068     float e_5x5 = e5x5( cluster, recHits, topology,flagsexcl, severitiesexcl, sevLv );
01069     float covEtaEta, covEtaPhi, covPhiPhi;
01070 
01071     if (e_5x5 >= 0.) {
01072         //double w0_ = parameterMap_.find("W0")->second;
01073         std::vector< std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
01074         std::pair<float,float> mean5x5PosInNrCrysFromSeed =  mean5x5PositionInLocalCrysCoord( cluster, recHits, topology,flagsexcl, severitiesexcl, sevLv );
01075         std::pair<float,float> mean5x5XYPos =  mean5x5PositionInXY(cluster,recHits,topology,flagsexcl, severitiesexcl, sevLv);
01076 
01077         // now we can calculate the covariances
01078         double numeratorEtaEta = 0;
01079         double numeratorEtaPhi = 0;
01080         double numeratorPhiPhi = 0;
01081         double denominator     = 0;
01082 
01083         //these allow us to scale the localCov by the crystal size 
01084         //so that the localCovs have the same average value as the normal covs
01085         const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
01086         const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
01087 
01088         DetId seedId = getMaximum( v_id, recHits ).first;
01089 
01090         bool isBarrel=seedId.subdetId()==EcalBarrel;
01091         const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
01092 
01093         CaloNavigator<DetId> cursor = CaloNavigator<DetId>( seedId, topology->getSubdetectorTopology( seedId ) );
01094 
01095         for ( int eastNr = -2; eastNr <= 2; ++eastNr ) { //east is eta in barrel
01096             for ( int northNr = -2; northNr <= 2; ++northNr ) { //north is phi in barrel
01097                 cursor.home();
01098                 cursor.offsetBy( eastNr, northNr);
01099                 float energy = recHitEnergy( *cursor, recHits,flagsexcl, severitiesexcl, sevLv);
01100                 if ( energy <= 0 ) continue;
01101 
01102                 float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
01103                 float dPhi = 0;
01104 
01105                 if(isBarrel)  dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
01106                 else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
01107 
01108 
01109                 double w = std::max(0.0, w0 + log( energy / e_5x5 ));
01110 
01111                 denominator += w;
01112                 numeratorEtaEta += w * dEta * dEta;
01113                 numeratorEtaPhi += w * dEta * dPhi;
01114                 numeratorPhiPhi += w * dPhi * dPhi;
01115             } //end east loop
01116         }//end north loop
01117 
01118 
01119         //multiplying by crysSize to make the values compariable to normal covariances
01120         if (denominator != 0.0) {
01121             covEtaEta =  crysSize*crysSize* numeratorEtaEta / denominator;
01122             covEtaPhi =  crysSize*crysSize* numeratorEtaPhi / denominator;
01123             covPhiPhi =  crysSize*crysSize* numeratorPhiPhi / denominator;
01124         } else {
01125             covEtaEta = 999.9;
01126             covEtaPhi = 999.9;
01127             covPhiPhi = 999.9;
01128         }
01129 
01130 
01131     } else {
01132         // Warn the user if there was no energy in the cells and return zeroes.
01133         //       std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
01134         covEtaEta = 0;
01135         covEtaPhi = 0;
01136         covPhiPhi = 0;
01137     }
01138     std::vector<float> v;
01139     v.push_back( covEtaEta );
01140     v.push_back( covEtaPhi );
01141     v.push_back( covPhiPhi );
01142     return v;
01143 }
01144 
01145 
01146 double EcalClusterTools::zernike20( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0, bool logW, float w0 )
01147 {
01148     return absZernikeMoment( cluster, recHits, geometry, 2, 0, R0, logW, w0 );
01149 }
01150 
01151 
01152 
01153 double EcalClusterTools::zernike42( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, double R0, bool logW, float w0 )
01154 {
01155     return absZernikeMoment( cluster, recHits, geometry, 4, 2, R0, logW, w0 );
01156 }
01157 
01158 
01159 
01160 double EcalClusterTools::absZernikeMoment( const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
01161 {
01162     // 1. Check if n,m are correctly
01163     if ((m>n) || ((n-m)%2 != 0) || (n<0) || (m<0)) return -1;
01164 
01165     // 2. Check if n,R0 are within validity Range :
01166     // n>20 or R0<2.19cm  just makes no sense !
01167     if ((n>20) || (R0<=2.19)) return -1;
01168     if (n<=5) return fast_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0 );
01169     else return calc_AbsZernikeMoment(cluster, recHits, geometry, n, m, R0, logW, w0 );
01170 }
01171 
01172 
01173 double EcalClusterTools::fast_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
01174 {
01175     double r,ph,e,Re=0,Im=0;
01176     double TotalEnergy = cluster.energy();
01177     int index = (n/2)*(n/2)+(n/2)+m;
01178     std::vector<EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
01179     int clusterSize = energyDistribution.size();
01180     if(clusterSize < 3) return 0.0;
01181 
01182     for (int i=0; i<clusterSize; i++)
01183     { 
01184         r = energyDistribution[i].r / R0;
01185         if (r<1) {
01186             std::vector<double> pol;
01187             pol.push_back( f00(r) );
01188             pol.push_back( f11(r) );
01189             pol.push_back( f20(r) );
01190             pol.push_back( f22(r) );
01191             pol.push_back( f31(r) );
01192             pol.push_back( f33(r) );
01193             pol.push_back( f40(r) );
01194             pol.push_back( f42(r) );
01195             pol.push_back( f44(r) );
01196             pol.push_back( f51(r) );
01197             pol.push_back( f53(r) );
01198             pol.push_back( f55(r) );
01199             ph = (energyDistribution[i]).phi;
01200             e = energyDistribution[i].deposited_energy;
01201             Re = Re + e/TotalEnergy * pol[index] * cos( (double) m * ph);
01202             Im = Im - e/TotalEnergy * pol[index] * sin( (double) m * ph);
01203         }
01204     }
01205     return sqrt(Re*Re+Im*Im);
01206 }
01207 
01208 
01209 
01210 double EcalClusterTools::calc_AbsZernikeMoment(const reco::BasicCluster &cluster, const EcalRecHitCollection *recHits, const CaloGeometry *geometry, int n, int m, double R0, bool logW, float w0 )
01211 {
01212     double r, ph, e, Re=0, Im=0, f_nm;
01213     double TotalEnergy = cluster.energy();
01214     std::vector<EcalClusterEnergyDeposition> energyDistribution = getEnergyDepTopology( cluster, recHits, geometry, logW, w0 );
01215     int clusterSize=energyDistribution.size();
01216     if(clusterSize<3) return 0.0;
01217 
01218     for (int i = 0; i < clusterSize; ++i)
01219     { 
01220         r = energyDistribution[i].r / R0;
01221         if (r < 1) {
01222             ph = energyDistribution[i].phi;
01223             e = energyDistribution[i].deposited_energy;
01224             f_nm = 0;
01225             for (int s=0; s<=(n-m)/2; s++) {
01226                 if (s%2==0) { 
01227                     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));
01228                 } else {
01229                     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));
01230                 }
01231             }
01232             Re = Re + e/TotalEnergy * f_nm * cos( (double) m*ph);
01233             Im = Im - e/TotalEnergy * f_nm * sin( (double) m*ph);
01234         }
01235     }
01236     return sqrt(Re*Re+Im*Im);
01237 }
01238 
01239 //returns the crystal 'eta' from the det id
01240 //it is defined as the number of crystals from the centre in the eta direction
01241 //for the barrel with its eta/phi geometry it is always integer
01242 //for the endcap it is fractional due to the x/y geometry
01243 float  EcalClusterTools::getIEta(const DetId& id)
01244 {
01245     if(id.det()==DetId::Ecal){
01246         if(id.subdetId()==EcalBarrel){
01247             EBDetId ebId(id);
01248             return ebId.ieta();
01249         }else if(id.subdetId()==EcalEndcap){
01250             float iXNorm = getNormedIX(id);
01251             float iYNorm = getNormedIY(id);
01252 
01253             return std::sqrt(iXNorm*iXNorm+iYNorm*iYNorm);
01254         }
01255     }
01256     return 0.;    
01257 }
01258 
01259 
01260 //returns the crystal 'phi' from the det id
01261 //it is defined as the number of crystals from the centre in the phi direction
01262 //for the barrel with its eta/phi geometry it is always integer
01263 //for the endcap it is not defined 
01264 float  EcalClusterTools::getIPhi(const DetId& id)
01265 {
01266     if(id.det()==DetId::Ecal){
01267         if(id.subdetId()==EcalBarrel){
01268             EBDetId ebId(id);
01269             return ebId.iphi();
01270         }
01271     }
01272     return 0.;    
01273 }
01274 
01275 //want to map 1=-50,50=-1,51=1 and 100 to 50 so sub off one if zero or neg
01276 float EcalClusterTools::getNormedIX(const DetId& id)
01277 {
01278     if(id.det()==DetId::Ecal && id.subdetId()==EcalEndcap){
01279         EEDetId eeId(id);      
01280         int iXNorm  = eeId.ix()-50;
01281         if(iXNorm<=0) iXNorm--;
01282         return iXNorm;
01283     }
01284     return 0;
01285 }
01286 
01287 //want to map 1=-50,50=-1,51=1 and 100 to 50 so sub off one if zero or neg
01288 float EcalClusterTools::getNormedIY(const DetId& id)
01289 {
01290     if(id.det()==DetId::Ecal && id.subdetId()==EcalEndcap){
01291         EEDetId eeId(id);      
01292         int iYNorm  = eeId.iy()-50;
01293         if(iYNorm<=0) iYNorm--;
01294         return iYNorm;
01295     }
01296     return 0;
01297 }
01298 
01299 //nr crystals crysId is away from orgin id in eta
01300 float EcalClusterTools::getNrCrysDiffInEta(const DetId& crysId,const DetId& orginId)
01301 {
01302     float crysIEta = getIEta(crysId);
01303     float orginIEta = getIEta(orginId);
01304     bool isBarrel = orginId.subdetId()==EcalBarrel;
01305 
01306     float nrCrysDiff = crysIEta-orginIEta;
01307 
01308     //no iEta=0 in barrel, so if go from positive to negative
01309     //need to reduce abs(detEta) by 1
01310     if(isBarrel){ 
01311         if(crysIEta*orginIEta<0){ // -1 to 1 transition
01312             if(crysIEta>0) nrCrysDiff--;
01313             else nrCrysDiff++;
01314         }
01315     }
01316     return nrCrysDiff;
01317 }
01318 
01319 //nr crystals crysId is away from orgin id in phi
01320 float EcalClusterTools::getNrCrysDiffInPhi(const DetId& crysId,const DetId& orginId)
01321 {
01322     float crysIPhi = getIPhi(crysId);
01323     float orginIPhi = getIPhi(orginId);
01324     bool isBarrel = orginId.subdetId()==EcalBarrel;
01325 
01326     float nrCrysDiff = crysIPhi-orginIPhi;
01327 
01328     if(isBarrel){ //if barrel, need to map into 0-180 
01329         if (nrCrysDiff > + 180) { nrCrysDiff = nrCrysDiff - 360; }
01330         if (nrCrysDiff < - 180) { nrCrysDiff = nrCrysDiff + 360; }
01331     }
01332     return nrCrysDiff;
01333 }
01334 
01335 //nr crystals crysId is away from mean phi in 5x5 in phi
01336 float EcalClusterTools::getDPhiEndcap(const DetId& crysId,float meanX,float meanY)
01337 {
01338     float iXNorm  = getNormedIX(crysId);
01339     float iYNorm  = getNormedIY(crysId);
01340 
01341     float hitLocalR2 = (iXNorm-meanX)*(iXNorm-meanX)+(iYNorm-meanY)*(iYNorm-meanY);
01342     float hitR2 = iXNorm*iXNorm+iYNorm*iYNorm;
01343     float meanR2 = meanX*meanX+meanY*meanY;
01344     float hitR = sqrt(hitR2);
01345     float meanR = sqrt(meanR2);
01346 
01347     float tmp = (hitR2+meanR2-hitLocalR2)/(2*hitR*meanR);
01348     if (tmp<-1) tmp =-1;
01349     if (tmp>1)  tmp=1;
01350     float phi = acos(tmp);
01351     float dPhi = hitR*phi;
01352 
01353     return dPhi;
01354 }
01355 
01356 std::vector<float> EcalClusterTools::scLocalCovariances(const reco::SuperCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology, float w0)
01357 {
01358     const reco::BasicCluster bcluster = *(cluster.seed());
01359 
01360     float e_5x5 = e5x5(bcluster, recHits, topology);
01361     float covEtaEta, covEtaPhi, covPhiPhi;
01362 
01363     if (e_5x5 >= 0.) {
01364         std::vector<std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
01365         std::pair<float,float> mean5x5PosInNrCrysFromSeed =  mean5x5PositionInLocalCrysCoord(bcluster, recHits, topology);
01366         std::pair<float,float> mean5x5XYPos =  mean5x5PositionInXY(cluster,recHits,topology);
01367         // now we can calculate the covariances
01368         double numeratorEtaEta = 0;
01369         double numeratorEtaPhi = 0;
01370         double numeratorPhiPhi = 0;
01371         double denominator     = 0;
01372 
01373         const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
01374         const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
01375 
01376         DetId seedId = getMaximum(v_id, recHits).first;  
01377         bool isBarrel=seedId.subdetId()==EcalBarrel;
01378 
01379         const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
01380 
01381         for (size_t i = 0; i < v_id.size(); ++i) {
01382             CaloNavigator<DetId> cursor = CaloNavigator<DetId>(v_id[i].first, topology->getSubdetectorTopology(v_id[i].first));
01383             float energy = recHitEnergy(*cursor, recHits);
01384 
01385             if (energy <= 0) continue;
01386 
01387             float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
01388             float dPhi = 0;
01389             if(isBarrel)  dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
01390             else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
01391 
01392 
01393 
01394             double w = 0.;
01395             w = std::max(0.0, w0 + log( energy / e_5x5 ));
01396 
01397             denominator += w;
01398             numeratorEtaEta += w * dEta * dEta;
01399             numeratorEtaPhi += w * dEta * dPhi;
01400             numeratorPhiPhi += w * dPhi * dPhi;
01401         }
01402 
01403         //multiplying by crysSize to make the values compariable to normal covariances
01404         if (denominator != 0.0) {
01405             covEtaEta =  crysSize*crysSize* numeratorEtaEta / denominator;
01406             covEtaPhi =  crysSize*crysSize* numeratorEtaPhi / denominator;
01407             covPhiPhi =  crysSize*crysSize* numeratorPhiPhi / denominator;
01408         } else {
01409             covEtaEta = 999.9;
01410             covEtaPhi = 999.9;
01411             covPhiPhi = 999.9;
01412         }
01413 
01414     } else {
01415         // Warn the user if there was no energy in the cells and return zeroes.
01416         // std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
01417         covEtaEta = 0;
01418         covEtaPhi = 0;
01419         covPhiPhi = 0;
01420     }
01421 
01422     std::vector<float> v;
01423     v.push_back( covEtaEta );
01424     v.push_back( covEtaPhi );
01425     v.push_back( covPhiPhi );
01426 
01427     return v;
01428 }
01429 
01430 
01431 //================================================================== scLocalCovariances==============================================================
01432 
01433 std::vector<float> EcalClusterTools::scLocalCovariances(const reco::SuperCluster &cluster, const EcalRecHitCollection* recHits,const CaloTopology *topology,std::vector<int> flagsexcl, std::vector<int> severitiesexcl, const EcalSeverityLevelAlgo *sevLv, float w0)
01434 {
01435     const reco::BasicCluster bcluster = *(cluster.seed());
01436 
01437     float e_5x5 = e5x5(bcluster, recHits, topology);
01438     float covEtaEta, covEtaPhi, covPhiPhi;
01439 
01440     if (e_5x5 >= 0.) {
01441         std::vector<std::pair<DetId, float> > v_id = cluster.hitsAndFractions();
01442         std::pair<float,float> mean5x5PosInNrCrysFromSeed =  mean5x5PositionInLocalCrysCoord(bcluster, recHits, topology,flagsexcl, severitiesexcl, sevLv);
01443         std::pair<float,float> mean5x5XYPos =  mean5x5PositionInXY(cluster,recHits,topology,flagsexcl, severitiesexcl, sevLv);
01444         // now we can calculate the covariances
01445         double numeratorEtaEta = 0;
01446         double numeratorEtaPhi = 0;
01447         double numeratorPhiPhi = 0;
01448         double denominator     = 0;
01449 
01450         const double barrelCrysSize = 0.01745; //approximate size of crystal in eta,phi in barrel
01451         const double endcapCrysSize = 0.0447; //the approximate crystal size sigmaEtaEta was corrected to in the endcap
01452 
01453         DetId seedId = getMaximum(v_id, recHits).first;  
01454         bool isBarrel=seedId.subdetId()==EcalBarrel;
01455 
01456         const double crysSize = isBarrel ? barrelCrysSize : endcapCrysSize;
01457 
01458         for (size_t i = 0; i < v_id.size(); ++i) {
01459             CaloNavigator<DetId> cursor = CaloNavigator<DetId>(v_id[i].first, topology->getSubdetectorTopology(v_id[i].first));
01460             float energy = recHitEnergy(*cursor, recHits,flagsexcl, severitiesexcl, sevLv);
01461 
01462             if (energy <= 0) continue;
01463 
01464             float dEta = getNrCrysDiffInEta(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.first;
01465             float dPhi = 0;
01466             if(isBarrel)  dPhi = getNrCrysDiffInPhi(*cursor,seedId) - mean5x5PosInNrCrysFromSeed.second;
01467             else dPhi = getDPhiEndcap(*cursor,mean5x5XYPos.first,mean5x5XYPos.second);
01468 
01469 
01470 
01471             double w = 0.;
01472             w = std::max(0.0, w0 + log( energy / e_5x5 ));
01473 
01474             denominator += w;
01475             numeratorEtaEta += w * dEta * dEta;
01476             numeratorEtaPhi += w * dEta * dPhi;
01477             numeratorPhiPhi += w * dPhi * dPhi;
01478         }
01479 
01480         //multiplying by crysSize to make the values compariable to normal covariances
01481         if (denominator != 0.0) {
01482             covEtaEta =  crysSize*crysSize* numeratorEtaEta / denominator;
01483             covEtaPhi =  crysSize*crysSize* numeratorEtaPhi / denominator;
01484             covPhiPhi =  crysSize*crysSize* numeratorPhiPhi / denominator;
01485         } else {
01486             covEtaEta = 999.9;
01487             covEtaPhi = 999.9;
01488             covPhiPhi = 999.9;
01489         }
01490 
01491     } else {
01492         // Warn the user if there was no energy in the cells and return zeroes.
01493         // std::cout << "\ClusterShapeAlgo::Calculate_Covariances:  no energy in supplied cells.\n";
01494         covEtaEta = 0;
01495         covEtaPhi = 0;
01496         covPhiPhi = 0;
01497     }
01498 
01499     std::vector<float> v;
01500     v.push_back( covEtaEta );
01501     v.push_back( covEtaPhi );
01502     v.push_back( covPhiPhi );
01503 
01504     return v;
01505 }
01506 
01507 // compute cluster second moments with respect to principal axes (eigenvectors of sEtaEta, sPhiPhi, sEtaPhi matrix)
01508 // store also angle alpha between major axis and phi.
01509 // takes into account shower elongation in phi direction due to magnetic field effect: 
01510 // default value of 0.8 ensures sMaj = sMin for unconverted photons 
01511 // (if phiCorrectionFactor=1 sMaj > sMin and alpha=0 also for unconverted photons)
01512 
01513 
01514 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( const reco::BasicCluster &basicCluster, const EcalRecHitCollection &recHits, double phiCorrectionFactor, double w0, bool useLogWeights) {
01515 
01516   Cluster2ndMoments returnMoments;
01517   returnMoments.sMaj = -1.;
01518   returnMoments.sMin = -1.;
01519   returnMoments.alpha = 0.;
01520 
01521   // for now implemented only for EB:
01522   //  if( fabs( basicCluster.eta() ) < 1.479 ) { 
01523 
01524     std::vector<const EcalRecHit*> RH_ptrs;
01525     
01526     std::vector< std::pair<DetId, float> > myHitsPair = basicCluster.hitsAndFractions();
01527     std::vector<DetId> usedCrystals;
01528     for(unsigned int i=0; i< myHitsPair.size(); i++){
01529       usedCrystals.push_back(myHitsPair[i].first);
01530     }
01531     
01532     for(unsigned int i=0; i<usedCrystals.size(); i++){
01533       //get pointer to recHit object
01534       EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01535       RH_ptrs.push_back(  &(*myRH)  );
01536     }
01537 
01538       returnMoments = EcalClusterTools::cluster2ndMoments(RH_ptrs, phiCorrectionFactor, w0, useLogWeights);
01539 
01540       //  }
01541 
01542   return returnMoments;
01543 
01544 }
01545 
01546 
01547 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( const reco::SuperCluster &superCluster, const EcalRecHitCollection &recHits, double phiCorrectionFactor, double w0, bool useLogWeights) {
01548 
01549   // for now returns second moments of supercluster seed cluster:
01550   Cluster2ndMoments returnMoments;
01551   returnMoments.sMaj = -1.;
01552   returnMoments.sMin = -1.;
01553   returnMoments.alpha = 0.;
01554 
01555   // for now implemented only for EB:
01556   //  if( fabs( superCluster.eta() ) < 1.479 ) { 
01557     returnMoments = EcalClusterTools::cluster2ndMoments( *(superCluster.seed()), recHits, phiCorrectionFactor, w0, useLogWeights);
01558     //  }
01559 
01560   return returnMoments;
01561 
01562 }
01563 
01564 
01565 Cluster2ndMoments EcalClusterTools::cluster2ndMoments( std::vector<const EcalRecHit*> RH_ptrs, double phiCorrectionFactor, double w0, bool useLogWeights) {
01566 
01567   double mid_eta(0),mid_phi(0),mid_x(0),mid_y(0);
01568   
01569   double Etot = EcalClusterTools::getSumEnergy(  RH_ptrs  );
01570  
01571   double max_phi=-10.;
01572   double min_phi=100.;
01573   
01574   
01575   std::vector<double> etaDetId;
01576   std::vector<double> phiDetId;
01577   std::vector<double> xDetId;
01578   std::vector<double> yDetId;
01579   std::vector<double> wiDetId;
01580  
01581   unsigned int nCry=0;
01582   double denominator=0.;
01583   bool isBarrel(1);
01584 
01585   // loop over rechits and compute weights:
01586   for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01587 
01588     //get iEta, iPhi
01589     double temp_eta(0),temp_phi(0),temp_x(0),temp_y(0);
01590     isBarrel = (*rh_ptr)->detid().subdetId()==EcalBarrel;
01591     
01592     if(isBarrel) {
01593       temp_eta = (getIEta((*rh_ptr)->detid()) > 0. ? getIEta((*rh_ptr)->detid()) + 84.5 : getIEta((*rh_ptr)->detid()) + 85.5);
01594       temp_phi= getIPhi((*rh_ptr)->detid()) - 0.5;
01595     }
01596     else {
01597       temp_eta = getIEta((*rh_ptr)->detid());  
01598       temp_x =  getNormedIX((*rh_ptr)->detid());
01599       temp_y =  getNormedIY((*rh_ptr)->detid());
01600     }     
01601 
01602     double temp_ene=(*rh_ptr)->energy();
01603     
01604     double temp_wi=((useLogWeights) ?
01605                     std::max(0., w0 + log( fabs(temp_ene)/Etot ))
01606                     :  temp_ene);
01607 
01608 
01609     if(temp_phi>max_phi) max_phi=temp_phi;
01610     if(temp_phi<min_phi) min_phi=temp_phi;
01611     etaDetId.push_back(temp_eta);
01612     phiDetId.push_back(temp_phi);
01613     xDetId.push_back(temp_x);
01614     yDetId.push_back(temp_y);
01615     wiDetId.push_back(temp_wi);
01616     denominator+=temp_wi;
01617     nCry++;
01618   }
01619 
01620   if(isBarrel){
01621     // correct phi wrap-around:
01622     if(max_phi==359.5 && min_phi==0.5){ 
01623       for(unsigned int i=0; i<nCry; i++){
01624         if(phiDetId[i] - 179. > 0.) phiDetId[i]-=360.; 
01625         mid_phi+=phiDetId[i]*wiDetId[i];
01626         mid_eta+=etaDetId[i]*wiDetId[i];
01627       }
01628     } else{
01629       for(unsigned int i=0; i<nCry; i++){
01630         mid_phi+=phiDetId[i]*wiDetId[i];
01631         mid_eta+=etaDetId[i]*wiDetId[i];
01632       }
01633     }
01634   }else{
01635     for(unsigned int i=0; i<nCry; i++){
01636       mid_eta+=etaDetId[i]*wiDetId[i];      
01637       mid_x+=xDetId[i]*wiDetId[i];
01638       mid_y+=yDetId[i]*wiDetId[i];
01639     }
01640   }
01641   
01642   mid_eta/=denominator;
01643   mid_phi/=denominator;
01644   mid_x/=denominator;
01645   mid_y/=denominator;
01646 
01647 
01648   // See = sigma eta eta
01649   // Spp = (B field corrected) sigma phi phi
01650   // See = (B field corrected) sigma eta phi
01651   double See=0.;
01652   double Spp=0.;
01653   double Sep=0.;
01654   double deta(0),dphi(0);
01655   // compute (phi-corrected) covariance matrix:
01656   for(unsigned int i=0; i<nCry; i++) {
01657     if(isBarrel) {
01658       deta = etaDetId[i]-mid_eta;
01659       dphi = phiDetId[i]-mid_phi;
01660     } else {
01661       deta = etaDetId[i]-mid_eta;
01662       float hitLocalR2 = (xDetId[i]-mid_x)*(xDetId[i]-mid_x)+(yDetId[i]-mid_y)*(yDetId[i]-mid_y);
01663       float hitR2 = xDetId[i]*xDetId[i]+yDetId[i]*yDetId[i];
01664       float meanR2 = mid_x*mid_x+mid_y*mid_y;
01665       float hitR = sqrt(hitR2);
01666       float meanR = sqrt(meanR2);
01667       float phi = acos((hitR2+meanR2-hitLocalR2)/(2*hitR*meanR));
01668       dphi = hitR*phi;
01669 
01670     }
01671     See += (wiDetId[i]* deta * deta) / denominator;
01672     Spp += phiCorrectionFactor*(wiDetId[i]* dphi * dphi) / denominator;
01673     Sep += sqrt(phiCorrectionFactor)*(wiDetId[i]*deta*dphi) / denominator;
01674   }
01675 
01676   Cluster2ndMoments returnMoments;
01677 
01678   // compute matrix eigenvalues:
01679   returnMoments.sMaj = ((See + Spp) + sqrt((See - Spp)*(See - Spp) + 4.*Sep*Sep)) / 2.;
01680   returnMoments.sMin = ((See + Spp) - sqrt((See - Spp)*(See - Spp) + 4.*Sep*Sep)) / 2.;
01681 
01682   returnMoments.alpha = atan( (See - Spp + sqrt( (Spp - See)*(Spp - See) + 4.*Sep*Sep )) / (2.*Sep));
01683 
01684   return returnMoments;
01685 
01686 }
01687 
01688 
01689 
01690 
01691 //compute shower shapes: roundness and angle in a vector. Roundness is 0th element, Angle is 1st element.
01692 //description: uses classical mechanics inertia tensor.
01693 //             roundness is smaller_eValue/larger_eValue
01694 //             angle is the angle from the iEta axis to the smallest eVector (a.k.a. the shower's elongated axis)
01695 // this function uses only recHits belonging to a SC above energyThreshold (default 0)
01696 // you can select linear weighting = energy_recHit/total_energy         (weightedPositionMethod=0) default
01697 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
01698 std::vector<float> EcalClusterTools::roundnessBarrelSuperClusters( const reco::SuperCluster &superCluster ,const EcalRecHitCollection &recHits, int weightedPositionMethod, float energyThreshold){//int positionWeightingMethod=0){
01699     std::vector<const EcalRecHit*>RH_ptrs;
01700     std::vector< std::pair<DetId, float> > myHitsPair = superCluster.hitsAndFractions();
01701     std::vector<DetId> usedCrystals;
01702     for(unsigned int i=0; i< myHitsPair.size(); i++){
01703         usedCrystals.push_back(myHitsPair[i].first);
01704     }
01705     for(unsigned int i=0; i<usedCrystals.size(); i++){
01706         //get pointer to recHit object
01707         EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01708         if( myRH != recHits.end() && myRH->energy() > energyThreshold){ //require rec hit to have positive energy
01709             RH_ptrs.push_back(  &(*myRH)  );
01710         }
01711     }
01712     std::vector<float> temp = EcalClusterTools::roundnessSelectedBarrelRecHits(RH_ptrs,weightedPositionMethod); 
01713     return temp;
01714 }
01715 
01716 // this function uses all recHits within specified window ( with default values ieta_delta=2, iphi_delta=5) around SC's highest recHit.
01717 // recHits must pass an energy threshold "energyRHThresh" (default 0)
01718 // you can select linear weighting = energy_recHit/total_energy         (weightedPositionMethod=0)
01719 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
01720 
01721 std::vector<float> EcalClusterTools::roundnessBarrelSuperClustersUserExtended( const reco::SuperCluster &superCluster ,const EcalRecHitCollection &recHits, int ieta_delta, int iphi_delta, float energyRHThresh, int weightedPositionMethod){
01722 
01723     std::vector<const EcalRecHit*>RH_ptrs;
01724     std::vector< std::pair<DetId, float> > myHitsPair = superCluster.hitsAndFractions();
01725     std::vector<DetId> usedCrystals;
01726     for(unsigned int i=0; i< myHitsPair.size(); i++){
01727         usedCrystals.push_back(myHitsPair[i].first);
01728     }
01729 
01730     for(unsigned int i=0; i<usedCrystals.size(); i++){
01731         //get pointer to recHit object
01732         EcalRecHitCollection::const_iterator myRH = recHits.find(usedCrystals[i]);
01733         if(myRH != recHits.end() && myRH->energy() > energyRHThresh)
01734             RH_ptrs.push_back(  &(*myRH)  );
01735     }
01736 
01737 
01738     std::vector<int> seedPosition = EcalClusterTools::getSeedPosition(  RH_ptrs  );
01739 
01740     for(EcalRecHitCollection::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++){
01741         EBDetId EBdetIdi( rh->detid() );
01742         //if(rh != recHits.end())
01743         bool inEtaWindow = (   abs(  deltaIEta(seedPosition[0],EBdetIdi.ieta())  ) <= ieta_delta   );
01744         bool inPhiWindow = (   abs(  deltaIPhi(seedPosition[1],EBdetIdi.iphi())  ) <= iphi_delta   );
01745         bool passEThresh = (  rh->energy() > energyRHThresh  );
01746         bool alreadyCounted = false;
01747 
01748         // figure out if the rechit considered now is already inside the SC
01749         bool is_SCrh_inside_recHits = false;
01750         for(unsigned int i=0; i<usedCrystals.size(); i++){
01751             EcalRecHitCollection::const_iterator SCrh = recHits.find(usedCrystals[i]);
01752             if(SCrh != recHits.end()){
01753                 is_SCrh_inside_recHits = true;
01754                 if( rh->detid() == SCrh->detid()  ) alreadyCounted = true;
01755             }
01756         }//for loop over SC's recHits
01757 
01758         if( is_SCrh_inside_recHits && !alreadyCounted && passEThresh && inEtaWindow && inPhiWindow){
01759             RH_ptrs.push_back( &(*rh) );
01760         }
01761 
01762     }//for loop over rh
01763     return EcalClusterTools::roundnessSelectedBarrelRecHits(RH_ptrs,weightedPositionMethod);
01764 }
01765 
01766 // this function computes the roundness and angle variables for vector of pointers to recHits you pass it
01767 // you can select linear weighting = energy_recHit/total_energy         (weightedPositionMethod=0)
01768 // or log weighting = max( 0.0, 4.2 + log(energy_recHit/total_energy) ) (weightedPositionMethod=1)
01769 std::vector<float> EcalClusterTools::roundnessSelectedBarrelRecHits( std::vector<const EcalRecHit*> RH_ptrs, int weightedPositionMethod){//int weightedPositionMethod = 0){
01770     //positionWeightingMethod = 0 linear weighting, 1 log energy weighting
01771 
01772     std::vector<float> shapes; // this is the returning vector
01773 
01774     //make sure photon has more than one crystal; else roundness and angle suck
01775     if(RH_ptrs.size()<2){
01776         shapes.push_back( -3 );
01777         shapes.push_back( -3 );
01778         return shapes;
01779     }
01780 
01781     //Find highest E RH (Seed) and save info, compute sum total energy used
01782     std::vector<int> seedPosition = EcalClusterTools::getSeedPosition(  RH_ptrs  );// *recHits);
01783     int tempInt = seedPosition[0];
01784     if(tempInt <0) tempInt++;
01785     float energyTotal = EcalClusterTools::getSumEnergy(  RH_ptrs  );
01786 
01787     //1st loop over rechits: compute new weighted center position in coordinates centered on seed
01788     float centerIEta = 0.;
01789     float centerIPhi = 0.;
01790     float denominator = 0.;
01791 
01792     for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01793         //get iEta, iPhi
01794         EBDetId EBdetIdi( (*rh_ptr)->detid() );
01795         if(fabs(energyTotal) < 0.0001){
01796             // don't /0, bad!
01797             shapes.push_back( -2 );
01798             shapes.push_back( -2 );
01799             return shapes;
01800         }
01801         float weight = 0;
01802         if(fabs(weightedPositionMethod)<0.0001){ //linear
01803             weight = (*rh_ptr)->energy()/energyTotal;
01804         }else{ //logrithmic
01805             weight = std::max(0.0, 4.2 + log((*rh_ptr)->energy()/energyTotal));
01806         }
01807         denominator += weight;
01808         centerIEta += weight*deltaIEta(seedPosition[0],EBdetIdi.ieta()); 
01809         centerIPhi += weight*deltaIPhi(seedPosition[1],EBdetIdi.iphi());
01810     }
01811     if(fabs(denominator) < 0.0001){
01812         // don't /0, bad!
01813         shapes.push_back( -2 );
01814         shapes.push_back( -2 );
01815         return shapes;
01816     }
01817     centerIEta = centerIEta / denominator;
01818     centerIPhi = centerIPhi / denominator;
01819 
01820 
01821     //2nd loop over rechits: compute inertia tensor
01822     TMatrixDSym inertia(2); //initialize 2d inertia tensor
01823     double inertia00 = 0.;
01824     double inertia01 = 0.;// = inertia10 b/c matrix should be symmetric
01825     double inertia11 = 0.;
01826     int i = 0;
01827     for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01828         //get iEta, iPhi
01829         EBDetId EBdetIdi( (*rh_ptr)->detid() );
01830 
01831         if(fabs(energyTotal) < 0.0001){
01832             // don't /0, bad!
01833             shapes.push_back( -2 );
01834             shapes.push_back( -2 );
01835             return shapes;
01836         }
01837         float weight = 0;
01838         if(fabs(weightedPositionMethod) < 0.0001){ //linear
01839             weight = (*rh_ptr)->energy()/energyTotal;
01840         }else{ //logrithmic
01841             weight = std::max(0.0, 4.2 + log((*rh_ptr)->energy()/energyTotal));
01842         }
01843 
01844         float ieta_rh_to_center = deltaIEta(seedPosition[0],EBdetIdi.ieta()) - centerIEta;
01845         float iphi_rh_to_center = deltaIPhi(seedPosition[1],EBdetIdi.iphi()) - centerIPhi;
01846 
01847         inertia00 += weight*iphi_rh_to_center*iphi_rh_to_center;
01848         inertia01 -= weight*iphi_rh_to_center*ieta_rh_to_center;
01849         inertia11 += weight*ieta_rh_to_center*ieta_rh_to_center;
01850         i++;
01851     }
01852 
01853     inertia[0][0] = inertia00;
01854     inertia[0][1] = inertia01; // use same number here
01855     inertia[1][0] = inertia01; // and here to insure symmetry
01856     inertia[1][1] = inertia11;
01857 
01858 
01859     //step 1 find principal axes of inertia
01860     TMatrixD eVectors(2,2);
01861     TVectorD eValues(2);
01862     //std::cout<<"EcalClusterTools::showerRoundness- about to compute eVectors"<<std::endl;
01863     eVectors=inertia.EigenVectors(eValues); //ordered highest eV to lowest eV (I checked!)
01864     //and eVectors are in columns of matrix! I checked!
01865     //and they are normalized to 1
01866 
01867 
01868 
01869     //step 2 select eta component of smaller eVal's eVector
01870     TVectorD smallerAxis(2);//easiest to spin SC on this axis (smallest eVal)
01871     smallerAxis[0]=eVectors[0][1];//row,col  //eta component
01872     smallerAxis[1]=eVectors[1][1];           //phi component
01873 
01874     //step 3 compute interesting quatities
01875     Double_t temp = fabs(smallerAxis[0]);// closer to 1 ->beamhalo, closer to 0 something else
01876     if(fabs(eValues[0]) < 0.0001){
01877         // don't /0, bad!
01878         shapes.push_back( -2 );
01879         shapes.push_back( -2 );
01880         return shapes;
01881     }
01882 
01883     float Roundness = eValues[1]/eValues[0];
01884     float Angle=acos(temp);
01885 
01886     if( -0.00001 < Roundness && Roundness < 0) Roundness = 0.;
01887     if( -0.00001 < Angle && Angle < 0 ) Angle = 0.;
01888 
01889     shapes.push_back( Roundness );
01890     shapes.push_back( Angle );
01891     return shapes;
01892 
01893 }
01894 //private functions useful for roundnessBarrelSuperClusters etc.
01895 //compute delta iphi between a seed and a particular recHit
01896 //iphi [1,360]
01897 //safe gaurds are put in to ensure the difference is between [-180,180]
01898 int EcalClusterTools::deltaIPhi(int seed_iphi, int rh_iphi){
01899     int rel_iphi = rh_iphi - seed_iphi;
01900     // take care of cyclic variable iphi [1,360]
01901     if(rel_iphi >  180) rel_iphi = rel_iphi - 360;
01902     if(rel_iphi < -180) rel_iphi = rel_iphi + 360;
01903     return rel_iphi;
01904 }
01905 
01906 //compute delta ieta between a seed and a particular recHit
01907 //ieta [-85,-1] and [1,85]
01908 //safe gaurds are put in to shift the negative ieta +1 to make an ieta=0 so differences are computed correctly
01909 int EcalClusterTools::deltaIEta(int seed_ieta, int rh_ieta){
01910     // get rid of the fact that there is no ieta=0
01911     if(seed_ieta < 0) seed_ieta++;
01912     if(rh_ieta   < 0) rh_ieta++;
01913     int rel_ieta = rh_ieta - seed_ieta;
01914     return rel_ieta;
01915 }
01916 
01917 //return (ieta,iphi) of highest energy recHit of the recHits passed to this function
01918 std::vector<int> EcalClusterTools::getSeedPosition(std::vector<const EcalRecHit*> RH_ptrs){
01919     std::vector<int> seedPosition;
01920     float eSeedRH = 0;
01921     int iEtaSeedRH = 0;
01922     int iPhiSeedRH = 0;
01923 
01924     for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01925 
01926         //get iEta, iPhi
01927         EBDetId EBdetIdi( (*rh_ptr)->detid() );
01928 
01929         if(eSeedRH < (*rh_ptr)->energy()){
01930             eSeedRH = (*rh_ptr)->energy();
01931             iEtaSeedRH = EBdetIdi.ieta();
01932             iPhiSeedRH = EBdetIdi.iphi();
01933         }
01934 
01935     }// for loop
01936 
01937     seedPosition.push_back(iEtaSeedRH);
01938     seedPosition.push_back(iPhiSeedRH);
01939     return seedPosition;
01940 }
01941 
01942 // return the total energy of the recHits passed to this function
01943 float EcalClusterTools::getSumEnergy(std::vector<const EcalRecHit*> RH_ptrs){
01944 
01945     float sumE = 0.;
01946 
01947     for(std::vector<const EcalRecHit*>::const_iterator rh_ptr = RH_ptrs.begin(); rh_ptr != RH_ptrs.end(); rh_ptr++){
01948         sumE += (*rh_ptr)->energy();
01949     }// for loop
01950 
01951     return sumE;
01952 }