Go to the documentation of this file.00001 #include "RecoMET/METAlgorithms/interface/EcalHaloAlgo.h"
00002 #include "DataFormats/Common/interface/ValueMap.h"
00003
00004
00005
00006
00007
00008
00009
00010
00011 using namespace std;
00012 using namespace reco;
00013 using namespace edm;
00014
00015 bool CompareTime(const EcalRecHit* x, const EcalRecHit* y){ return x->time() < y->time(); }
00016
00017 EcalHaloAlgo::EcalHaloAlgo()
00018 {
00019 RoundnessCut = 0 ;
00020 AngleCut = 0 ;
00021 EBRecHitEnergyThreshold = 0.;
00022 EERecHitEnergyThreshold = 0.;
00023 ESRecHitEnergyThreshold = 0.;
00024 SumEnergyThreshold = 0.;
00025 NHitsThreshold =0;
00026 }
00027
00028 EcalHaloData EcalHaloAlgo::Calculate(const CaloGeometry& TheCaloGeometry, edm::Handle<reco::PhotonCollection>& ThePhotons, edm::Handle<reco::SuperClusterCollection>& TheSuperClusters, edm::Handle<EBRecHitCollection>& TheEBRecHits, edm::Handle<EERecHitCollection>& TheEERecHits, edm::Handle<ESRecHitCollection>& TheESRecHits)
00029 {
00030 EcalHaloData TheEcalHaloData;
00031
00032
00033 float SumE[361];
00034
00035 int NumHits[361];
00036
00037 float MinTimeHits[361];
00038
00039 float MaxTimeHits[361];
00040
00041
00042 for(int i = 0 ; i < 361 ; i++ )
00043 {
00044 SumE[i] = 0.;
00045 NumHits[i] = 0 ;
00046 MinTimeHits[i] = 9999.;
00047 MaxTimeHits[i] = -9999.;
00048 }
00049
00050
00051 for(EBRecHitCollection::const_iterator hit = TheEBRecHits->begin() ; hit != TheEBRecHits->end() ; hit++ )
00052 {
00053
00054 if (hit->energy() < EBRecHitEnergyThreshold ) continue;
00055
00056
00057 DetId id = DetId(hit->id());
00058 const CaloSubdetectorGeometry* TheSubGeometry = 0;
00059 const CaloCellGeometry* cell = 0 ;
00060
00061
00062 TheSubGeometry = TheCaloGeometry.getSubdetectorGeometry(DetId::Ecal, 1);
00063 EBDetId EcalID(id.rawId());
00064 if( TheSubGeometry )
00065 cell = TheSubGeometry->getGeometry(id);
00066
00067 if(cell)
00068 {
00069 GlobalPoint globalpos = cell->getPosition();
00070
00071 int iPhi = EcalID.iphi();
00072
00073 if( iPhi < 361 )
00074 {
00075
00076 SumE[iPhi] += hit->energy();
00077 NumHits[iPhi] ++;
00078
00079 float time = hit->time();
00080 MinTimeHits[iPhi] = time < MinTimeHits[iPhi] ? time : MinTimeHits[iPhi];
00081 MaxTimeHits[iPhi] = time > MaxTimeHits[iPhi] ? time : MaxTimeHits[iPhi];
00082 }
00083 }
00084 }
00085
00086
00087 for( int iPhi = 1 ; iPhi < 361; iPhi++ )
00088 {
00089 if( SumE[iPhi] >= SumEnergyThreshold && NumHits[iPhi] > NHitsThreshold )
00090 {
00091
00092 PhiWedge wedge(SumE[iPhi], iPhi, NumHits[iPhi], MinTimeHits[iPhi], MaxTimeHits[iPhi]);
00093
00094
00095
00096
00097 std::vector<const EcalRecHit*> Hits;
00098 for(EBRecHitCollection::const_iterator hit = TheEBRecHits->begin() ; hit != TheEBRecHits->end() ; hit++ )
00099 {
00100 if (hit->energy() < EBRecHitEnergyThreshold ) continue;
00101
00102
00103 DetId id = DetId(hit->id());
00104 EBDetId EcalID(id.rawId());
00105 int Hit_iPhi = EcalID.iphi();
00106
00107 if( Hit_iPhi != iPhi ) continue;
00108 Hits.push_back( &(*hit) );
00109
00110 }
00111 std::sort( Hits.begin() , Hits.end(), CompareTime);
00112 float MinusToPlus = 0.;
00113 float PlusToMinus = 0.;
00114 for( unsigned int i = 0 ; i < Hits.size() ; i++ )
00115 {
00116 DetId id_i = DetId(Hits[i]->id());
00117 EBDetId EcalID_i(id_i.rawId());
00118 int ieta_i = EcalID_i.ieta();
00119 for( unsigned int j = (i+1) ; j < Hits.size() ; j++ )
00120 {
00121 DetId id_j = DetId(Hits[j]->id() );
00122 EBDetId EcalID_j(id_j.rawId());
00123 int ieta_j = EcalID_j.ieta();
00124 if( ieta_i > ieta_j ) PlusToMinus += TMath::Abs(ieta_j - ieta_i );
00125 else MinusToPlus += TMath::Abs(ieta_j -ieta_i) ;
00126 }
00127 }
00128
00129 float PlusZOriginConfidence = (PlusToMinus+MinusToPlus) ? PlusToMinus / (PlusToMinus+MinusToPlus) : -1.;
00130 wedge.SetPlusZOriginConfidence(PlusZOriginConfidence);
00131 TheEcalHaloData.GetPhiWedges().push_back(wedge);
00132 }
00133 }
00134
00135 std::vector<float> vShowerShapes_Roundness;
00136 std::vector<float> vShowerShapes_Angle ;
00137 for(reco::SuperClusterCollection::const_iterator cluster = TheSuperClusters->begin() ; cluster != TheSuperClusters->end() ; cluster++ )
00138 {
00139 if( abs(cluster->eta()) <= 1.48 )
00140 {
00141 vector<float> shapes = EcalClusterTools::roundnessBarrelSuperClusters( *cluster, (*TheEBRecHits.product()));
00142 float roundness = shapes[0];
00143 float angle = shapes[1];
00144
00145
00146 if( (roundness >=0 && roundness < GetRoundnessCut()) && angle >= 0 && angle < GetAngleCut() )
00147 {
00148 edm::Ref<SuperClusterCollection> TheClusterRef( TheSuperClusters, cluster - TheSuperClusters->begin());
00149 bool BelongsToPhoton = false;
00150 if( ThePhotons.isValid() )
00151 {
00152 for(reco::PhotonCollection::const_iterator iPhoton = ThePhotons->begin() ; iPhoton != ThePhotons->end() ; iPhoton++ )
00153 {
00154 if(iPhoton->isEB())
00155 if ( TheClusterRef == iPhoton->superCluster() )
00156 {
00157 BelongsToPhoton = true;
00158 break;
00159 }
00160 }
00161 }
00162
00163
00164 if( BelongsToPhoton )
00165 {
00166 TheEcalHaloData.GetSuperClusters().push_back( TheClusterRef ) ;
00167 }
00168 }
00169 vShowerShapes_Roundness.push_back(shapes[0]);
00170 vShowerShapes_Angle.push_back(shapes[1]);
00171 }
00172 else
00173 {
00174 vShowerShapes_Roundness.push_back(-1.);
00175 vShowerShapes_Angle.push_back(-1.);
00176 }
00177 }
00178
00179 edm::ValueMap<float>::Filler TheRoundnessFiller( TheEcalHaloData.GetShowerShapesRoundness() );
00180 TheRoundnessFiller.insert( TheSuperClusters, vShowerShapes_Roundness.begin(), vShowerShapes_Roundness.end() );
00181 TheRoundnessFiller.fill();
00182
00183 edm::ValueMap<float>::Filler TheAngleFiller( TheEcalHaloData.GetShowerShapesAngle() );
00184 TheAngleFiller.insert( TheSuperClusters, vShowerShapes_Angle.begin() , vShowerShapes_Angle.end() );
00185 TheAngleFiller.fill();
00186
00187 return TheEcalHaloData;
00188 }
00189
00190
00191