CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFPhotonClusters.cc
Go to the documentation of this file.
5 
7 
8 #include <TMath.h>
9 #include <TVector2.h>
10 using namespace reco;
12  PFClusterRef_(PFClusterRef)
13 {
14  if(PFClusterRef_->layer()==PFLayer:: ECAL_BARREL )isEB_=true;
15  else isEB_=false;
16  SetSeed();
17  PFCrystalCoor();
18  for(int i=0; i<5; ++i)
19  for(int j=0; j<5; ++j)e5x5_[i][j]=0;
22 }
23 
25  double PFSeedE=0;
26  math::XYZVector axis;
28  DetId idseed;
29  const std::vector< reco::PFRecHitFraction >& PFRecHits=
30  PFClusterRef_->recHitFractions();
31 
32  for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin();
33  it != PFRecHits.end(); ++it){
34  const PFRecHitRef& RefPFRecHit = it->recHitRef();
35  double frac=it->fraction();
36  float E= RefPFRecHit->energy()* frac;
37  if(E>PFSeedE){
38  PFSeedE=E;
39  axis=RefPFRecHit->getAxisXYZ();
40  position=RefPFRecHit->position();
41  idseed = RefPFRecHit->detId();
42  }
43  }
44  idseed_=idseed;
46  seedAxis_=axis;
47 }
48 
50  if(PFClusterRef_->layer()==PFLayer:: ECAL_BARREL ){//is Barrel
51  isEB_=true;
52  EBDetId EBidSeed=EBDetId(idseed_.rawId());
53  CrysIEta_=EBidSeed.ieta();
54  CrysIPhi_=EBidSeed.iphi();
55  double depth = PFClusterRef_->getDepthCorrection(PFClusterRef_->energy(), false, false);
56  math::XYZVector center_pos = seedPosition_+depth*seedAxis_;
57  //Crystal Coordinates:
58  double Pi=TMath::Pi();
59  float Phi=PFClusterRef_->position().phi();
60  double Theta = -(PFClusterRef_->position().theta())+0.5* Pi;
61  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
62  double PhiWidth = (Pi/180.);
63  double PhiCry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
64  double ThetaCentr = -center_pos.theta()+0.5*Pi;
65  double ThetaWidth = (Pi/180.)*cos(ThetaCentr);
66 
67  double EtaCry = (Theta-ThetaCentr)/ThetaWidth;
68  CrysEta_=EtaCry;
69  CrysPhi_=PhiCry;
70 
71  if(abs(CrysIEta_)==1 || abs(CrysIEta_)==2 )
72  CrysIEtaCrack_=abs(CrysIEta_);
73  if(abs(CrysIEta_)>2 && abs(CrysIEta_)<24)
75  if(abs(CrysIEta_)==24)
77  if(abs(CrysIEta_)==25)
79  if(abs(CrysIEta_)==26)
81  if(abs(CrysIEta_)==27)
83  if(abs(CrysIEta_)>27 && abs(CrysIEta_)<44)
85  if(abs(CrysIEta_)==44)
87  if(abs(CrysIEta_)==45)
88  CrysIEtaCrack_=10;
89  if(abs(CrysIEta_)==46)
90  CrysIEtaCrack_=11;
91  if(abs(CrysIEta_)==47)
92  CrysIEtaCrack_=12;
93  if(abs(CrysIEta_)>47 && abs(CrysIEta_)<64)
94  CrysIEtaCrack_=13;
95  if(abs(CrysIEta_)==64)
96  CrysIEtaCrack_=14;
97  if(abs(CrysIEta_)==65)
98  CrysIEtaCrack_=15;
99  if(abs(CrysIEta_)==66)
100  CrysIEtaCrack_=16;
101  if(abs(CrysIEta_)==67)
102  CrysIEtaCrack_=17;
103  if(abs(CrysIEta_)>67 && abs(CrysIEta_)<84)
104  CrysIEtaCrack_=18;
105  if(abs(CrysIEta_)==84)
106  CrysIEtaCrack_=19;
107  if(abs(CrysIEta_)==85)
108  CrysIEtaCrack_=20;
109  }
110  else{
111  isEB_=false;
112  EEDetId EEidSeed=EEDetId(idseed_.rawId());
113  CrysIX_=EEidSeed.ix();
114  CrysIY_=EEidSeed.iy();
115  float X0 = 0.89; float T0 = 1.2;
116  if(fabs(PFClusterRef_->eta())>1.653)T0=3.1;
117  double depth = X0 * (T0 + log(PFClusterRef_->energy()));
118  math::XYZVector center_pos=(seedPosition_)+depth*seedAxis_;
119  double XCentr = center_pos.x();
120  double YCentr = center_pos.y();
121  double XWidth = 2.59;
122  double YWidth = 2.59;
123 
124  CrysX_=(PFClusterRef_->x()-XCentr)/XWidth;
125  CrysY_=(PFClusterRef_->y()-YCentr)/YWidth;
126  }
127 }
128 
130  const std::vector< reco::PFRecHitFraction >& PFRecHits=PFClusterRef_->recHitFractions();
131  for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin(); it != PFRecHits.end(); ++it){
132  const PFRecHitRef& RefPFRecHit = it->recHitRef();
133  DetId id=RefPFRecHit->detId();
134  double frac=it->fraction();
135  float E=RefPFRecHit->energy()*frac;
136  if(isEB_){
137  int deta=EBDetId::distanceEta(id,idseed_);
138  int dphi=EBDetId::distancePhi(id,idseed_);
139  if(abs(deta)>2 ||abs(dphi)>2)continue;
140  //f(abs(dphi)<=2 && abs(deta)<=2){
141  EBDetId EBidSeed=EBDetId(idseed_.rawId());
142  EBDetId EBid=EBDetId(id.rawId());
143  int ind1=EBidSeed.ieta()-EBid.ieta();
144  int ind2=EBid.iphi()-EBidSeed.iphi();
145  if(EBidSeed.ieta() * EBid.ieta() > 0){
146  ind1=EBid.ieta()-EBidSeed.ieta();
147  }
148  else{ //near EB+ EB-
149  ind1=(1-(EBidSeed.ieta()-EBid.ieta()));
150  }
151  int iEta=ind1+2;
152  int iPhi=ind2+2;
153  //std::cout<<"IEta, IPhi "<<iEta<<", "<<iPhi<<std::endl;
154  if (iEta*5 + iPhi > 24 || iEta*5 + iPhi < 0){ //really check just writing outside the array
155  edm::LogInfo("OutOfBounds")<<"iEta = "<<iEta<<" iPhi = "<<iPhi;
156  } else {
157  e5x5_[iEta][iPhi]=E;
158  }
159  }
160  else{
161  int dx=EEDetId::distanceX(id,idseed_);
162  int dy=EEDetId::distanceY(id,idseed_);
163  if(abs(dx)>2 ||abs(dy>2))continue;
164  EEDetId EEidSeed=EEDetId(idseed_.rawId());
165  EEDetId EEid=EEDetId(id.rawId());
166  int ind1=EEid.ix()-EEidSeed.ix();
167  int ind2=EEid.iy()-EEidSeed.iy();
168  int ix=ind1+2;
169  int iy=ind2+2;
170  //std::cout<<"IX, IY "<<ix<<", "<<iy<<std::endl;
171  if (ix*5 + iy > 24 || ix*5 + iy < 0 ){ //no indication there was a problem here, just sanity-protection as above
172  edm::LogInfo("OutOfBounds")<<"ix = "<<ix<<" iy = "<<iy;
173  } else {
174  e5x5_[ix][iy]=E;
175  }
176  }
177  }
178 }
179 
181  double numeratorEtaWidth = 0.;
182  double numeratorPhiWidth = 0.;
183  double numeratorEtaPhiWidth = 0.;
184  double ClustEta=PFClusterRef_->eta();
185  double ClustPhi=PFClusterRef_->phi();
186  const std::vector< reco::PFRecHitFraction >& PFRecHits=PFClusterRef_->recHitFractions();
187  for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin(); it != PFRecHits.end(); ++it){
188  const PFRecHitRef& RefPFRecHit = it->recHitRef();
189  float E=RefPFRecHit->energy() * it->fraction();
190  double dEta = RefPFRecHit->position().eta() - ClustEta;
191  double dPhi = RefPFRecHit->position().phi() - ClustPhi;
192  if (dPhi > + TMath::Pi()) { dPhi = TMath::TwoPi() - dPhi; }
193  if (dPhi < - TMath::Pi()) { dPhi = TMath::TwoPi() + dPhi; }
194  numeratorEtaWidth += E * dEta * dEta;
195  numeratorPhiWidth += E * dPhi * dPhi;
196  numeratorEtaPhiWidth += E * fabs(dPhi) * fabs(dEta);
197  }
198  double denominator=PFClusterRef_->energy();
199  sigetaeta_ = sqrt(numeratorEtaWidth / denominator);
200  sigphiphi_ = sqrt(numeratorPhiWidth / denominator);
201  sigetaphi_ = sqrt(numeratorEtaPhiWidth / denominator);
202 }
const double TwoPi
const double Pi
int i
Definition: DBlmapReader.cc:9
int ix() const
Definition: EEDetId.h:71
static int distanceX(const EEDetId &a, const EEDetId &b)
Definition: EEDetId.cc:602
#define abs(x)
Definition: mlp_lapack.h:159
double Phi_mpi_pi(double x)
Definition: JetUtil.h:24
int iphi() const
get the crystal iphi
Definition: EBDetId.h:46
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
static int distanceEta(const EBDetId &a, const EBDetId &b)
Definition: EBDetId.cc:185
math::XYZVector seedAxis_
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
T sqrt(T t)
Definition: SSEVec.h:46
static int distancePhi(const EBDetId &a, const EBDetId &b)
Definition: EBDetId.cc:193
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int j
Definition: DBlmapReader.cc:9
int iy() const
Definition: EEDetId.h:77
int ieta() const
get the crystal ieta
Definition: EBDetId.h:44
static int distanceY(const EEDetId &a, const EEDetId &b)
Definition: EEDetId.cc:607
Definition: DetId.h:20
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
double e5x5_[5][5]
static int position[264][3]
Definition: ReadPGInfo.cc:509
PFPhotonClusters(PFClusterRef PFClusterRef)
static const double X0
PFClusterRef PFClusterRef_
math::XYZVector seedPosition_