#include <HcalHaloAlgo.h>
Public Member Functions | |
reco::HcalHaloData | Calculate (const CaloGeometry &TheCaloGeometry, edm::Handle< HBHERecHitCollection > &TheHBHERecHits) |
float | GetHBRecHitEnergyThreshold () |
float | GetHERecHitEnergyThreshold () |
float | GetPhiWedgeEnergyThreshold () |
int | GetPhiWedgeNHitsThreshold () |
HcalHaloAlgo () | |
void | SetPhiWedgeEnergyThreshold (float SumE) |
void | SetPhiWedgeNHitsThreshold (int nhits) |
void | SetPhiWedgeThresholds (float SumE, int nhits) |
void | SetRecHitEnergyThresholds (float HB, float HE) |
~HcalHaloAlgo () | |
Private Attributes | |
float | HBRecHitEnergyThreshold |
float | HERecHitEnergyThreshold |
int | NHitsThreshold |
float | SumEnergyThreshold |
Definition at line 27 of file HcalHaloAlgo.h.
HcalHaloAlgo::HcalHaloAlgo | ( | ) |
Definition at line 17 of file HcalHaloAlgo.cc.
{ HBRecHitEnergyThreshold = 0.; HERecHitEnergyThreshold = 0.; SumEnergyThreshold = 0.; NHitsThreshold = 0; }
HcalHaloAlgo::~HcalHaloAlgo | ( | ) | [inline] |
Definition at line 32 of file HcalHaloAlgo.h.
{}
HcalHaloData HcalHaloAlgo::Calculate | ( | const CaloGeometry & | TheCaloGeometry, |
edm::Handle< HBHERecHitCollection > & | TheHBHERecHits | ||
) |
Definition at line 25 of file HcalHaloAlgo.cc.
References CompareTime(), reco::HcalHaloData::GetPhiWedges(), HcalBarrel, HcalEndcap, i, hit::id, HcalDetId::ieta(), j, reco::PhiWedge::SetPlusZOriginConfidence(), python::multivaluedict::sort(), and cond::rpcobgas::time.
Referenced by reco::HcalHaloDataProducer::produce().
{ HcalHaloData TheHcalHaloData; // Store Energy sum of rechits as a function of iPhi (iPhi goes from 1 to 72) float SumE[73]; // Store Number of rechits as a function of iPhi int NumHits[73]; // Store minimum time of rechit as a function of iPhi float MinTimeHits[73]; // Store maximum time of rechit as a function of iPhi float MaxTimeHits[73]; for(unsigned int i = 0 ; i < 73 ; i++ ) { SumE[i] = 0; NumHits[i]= 0; MinTimeHits[i] = 0.; MaxTimeHits[i] = 0.; } for( HBHERecHitCollection::const_iterator hit = TheHBHERecHits->begin() ; hit != TheHBHERecHits->end() ; hit++ ) { HcalDetId id = HcalDetId(hit->id()); switch ( id.subdet() ) { case HcalBarrel: if(hit->energy() < HBRecHitEnergyThreshold )continue; break; case HcalEndcap: if(hit->energy() < HERecHitEnergyThreshold ) continue; break; default: continue; } int iEta = id.ieta(); int iPhi = id.iphi(); if(iPhi < 73 && TMath::Abs(iEta) < 23 ) { SumE[iPhi]+= hit->energy(); NumHits[iPhi] ++; float time = hit->time(); MinTimeHits[iPhi] = time < MinTimeHits[iPhi] ? time : MinTimeHits[iPhi]; MaxTimeHits[iPhi] = time > MaxTimeHits[iPhi] ? time : MaxTimeHits[iPhi]; } } for( int iPhi = 1 ; iPhi < 73 ; iPhi++ ) { if( SumE[iPhi] >= SumEnergyThreshold && NumHits[iPhi] > NHitsThreshold ) { // Build PhiWedge and store to HcalHaloData if energy or #hits pass thresholds PhiWedge wedge(SumE[iPhi], iPhi, NumHits[iPhi], MinTimeHits[iPhi], MaxTimeHits[iPhi]); // Loop over rechits again to calculate direction based on timing info std::vector<const HBHERecHit*> Hits; for( HBHERecHitCollection::const_iterator hit = TheHBHERecHits->begin() ; hit != TheHBHERecHits->end() ; hit++ ) { HcalDetId id = HcalDetId(hit->id()); if( id.iphi() != iPhi ) continue; if( TMath::Abs(id.ieta() ) > 22 ) continue; // has to overlap geometrically w/ HB switch ( id.subdet() ) { case HcalBarrel: if(hit->energy() < HBRecHitEnergyThreshold )continue; break; case HcalEndcap: if(hit->energy() < HERecHitEnergyThreshold ) continue; break; default: continue; } Hits.push_back(&(*hit)); } std::sort( Hits.begin() , Hits.end() , CompareTime); float MinusToPlus = 0.; float PlusToMinus = 0.; for( unsigned int i = 0 ; i < Hits.size() ; i++ ) { HcalDetId id_i = HcalDetId(Hits[i]->id() ); int ieta_i = id_i.ieta(); for( unsigned int j = (i+1) ; j < Hits.size() ; j++ ) { HcalDetId id_j = HcalDetId(Hits[j]->id() ); int ieta_j = id_j.ieta(); if( ieta_i > ieta_j ) PlusToMinus += TMath::Abs(ieta_i - ieta_j ) ; else MinusToPlus += TMath::Abs(ieta_i - ieta_j); } } float PlusZOriginConfidence = (PlusToMinus + MinusToPlus )? PlusToMinus / ( PlusToMinus + MinusToPlus ) : -1. ; wedge.SetPlusZOriginConfidence( PlusZOriginConfidence ); TheHcalHaloData.GetPhiWedges().push_back( wedge ); } } return TheHcalHaloData; }
float HcalHaloAlgo::GetHBRecHitEnergyThreshold | ( | ) | [inline] |
Definition at line 46 of file HcalHaloAlgo.h.
References HBRecHitEnergyThreshold.
{ return HBRecHitEnergyThreshold;}
float HcalHaloAlgo::GetHERecHitEnergyThreshold | ( | ) | [inline] |
Definition at line 47 of file HcalHaloAlgo.h.
References HERecHitEnergyThreshold.
{ return HERecHitEnergyThreshold;}
float HcalHaloAlgo::GetPhiWedgeEnergyThreshold | ( | ) | [inline] |
Definition at line 50 of file HcalHaloAlgo.h.
References SumEnergyThreshold.
{ return SumEnergyThreshold;}
int HcalHaloAlgo::GetPhiWedgeNHitsThreshold | ( | ) | [inline] |
void HcalHaloAlgo::SetPhiWedgeEnergyThreshold | ( | float | SumE | ) | [inline] |
Definition at line 41 of file HcalHaloAlgo.h.
References SumEnergyThreshold.
{ SumEnergyThreshold = SumE ;}
void HcalHaloAlgo::SetPhiWedgeNHitsThreshold | ( | int | nhits | ) | [inline] |
Definition at line 42 of file HcalHaloAlgo.h.
References NHitsThreshold.
{ NHitsThreshold = nhits ; }
void HcalHaloAlgo::SetPhiWedgeThresholds | ( | float | SumE, |
int | nhits | ||
) | [inline] |
Definition at line 43 of file HcalHaloAlgo.h.
References NHitsThreshold, and SumEnergyThreshold.
Referenced by reco::HcalHaloDataProducer::produce().
{ SumEnergyThreshold = SumE ; NHitsThreshold = nhits ;}
void HcalHaloAlgo::SetRecHitEnergyThresholds | ( | float | HB, |
float | HE | ||
) | [inline] |
Definition at line 38 of file HcalHaloAlgo.h.
References HBRecHitEnergyThreshold, and HERecHitEnergyThreshold.
Referenced by reco::HcalHaloDataProducer::produce().
{ HBRecHitEnergyThreshold = HB; HERecHitEnergyThreshold = HE;}
float HcalHaloAlgo::HBRecHitEnergyThreshold [private] |
Definition at line 55 of file HcalHaloAlgo.h.
Referenced by GetHBRecHitEnergyThreshold(), and SetRecHitEnergyThresholds().
float HcalHaloAlgo::HERecHitEnergyThreshold [private] |
Definition at line 56 of file HcalHaloAlgo.h.
Referenced by GetHERecHitEnergyThreshold(), and SetRecHitEnergyThresholds().
int HcalHaloAlgo::NHitsThreshold [private] |
Definition at line 60 of file HcalHaloAlgo.h.
Referenced by GetPhiWedgeNHitsThreshold(), SetPhiWedgeNHitsThreshold(), and SetPhiWedgeThresholds().
float HcalHaloAlgo::SumEnergyThreshold [private] |
Definition at line 59 of file HcalHaloAlgo.h.
Referenced by GetPhiWedgeEnergyThreshold(), SetPhiWedgeEnergyThreshold(), and SetPhiWedgeThresholds().