CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalClusterLocal.cc
Go to the documentation of this file.
1 #include "../interface/EcalClusterLocal.h"
2 
13 #include "TMath.h"
14 #include "TVector2.h"
15 
17 {}
18 
20 {}
21 
22 void EcalClusterLocal::localCoordsEB( const reco::CaloCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const {
24  es.get<CaloGeometryRecord>().get(pG);
25  localCoordsEB( bclus, *pG, etacry, phicry, ieta, iphi, thetatilt, phitilt);
26 }
27 
28 
29 
30 
31 void EcalClusterLocal::localCoordsEB( const reco::CaloCluster &bclus, const CaloGeometry & caloGeometry, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
32 {
33 
34  assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalBarrel);
35 
36  const CaloSubdetectorGeometry* geom=caloGeometry.getSubdetectorGeometry(DetId::Ecal,EcalBarrel);//EcalBarrel = 1
37 
38  const math::XYZPoint position_ = bclus.position();
39  double Theta = -position_.theta()+0.5*TMath::Pi();
40  double Eta = position_.eta();
41  double Phi = TVector2::Phi_mpi_pi(position_.phi());
42 
43  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
44  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
45  const float X0 = 0.89; const float T0 = 7.4;
46  double depth = X0 * (T0 + log(bclus.energy()));
47 
48  //find max energy crystal
49  std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
50  float drmin = 999.;
51  EBDetId crystalseed;
52  //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
53  for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {
54 
55  EBDetId crystal(crystals_vector[icry].first);
56 
57  const CaloCellGeometry* cell=geom->getGeometry(crystal);
58  GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
59  double EtaCentr = center_pos.eta();
60  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
61 
62  float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
63  if (dr<drmin) {
64  drmin = dr;
65  crystalseed = crystal;
66  }
67 
68  }
69 
70  ieta = crystalseed.ieta();
71  iphi = crystalseed.iphi();
72 
73  // Get center cell position from shower depth
74  const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
75  const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);
76 
77  thetatilt = cpyr->getThetaAxis();
78  phitilt = cpyr->getPhiAxis();
79 
80  GlobalPoint center_pos = cpyr->getPosition(depth);
81 
82  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
83  double PhiWidth = (TMath::Pi()/180.);
84  phicry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
85  //Some flips to take into account ECAL barrel symmetries:
86  if (ieta<0) phicry *= -1.;
87 
88  double ThetaCentr = -center_pos.theta()+0.5*TMath::Pi();
89  double ThetaWidth = (TMath::Pi()/180.)*TMath::Cos(ThetaCentr);
90  etacry = (Theta-ThetaCentr)/ThetaWidth;
91  //flip to take into account ECAL barrel symmetries:
92  if (ieta<0) etacry *= -1.;
93 
94  return;
95 
96 }
97 
98 void EcalClusterLocal::localCoordsEE( const reco::CaloCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
99 {
101  es.get<CaloGeometryRecord>().get(pG);
102  localCoordsEE( bclus, *pG, xcry, ycry, ix, iy, thetatilt, phitilt);
103 }
104 
105 void EcalClusterLocal::localCoordsEE( const reco::CaloCluster &bclus, const CaloGeometry & caloGeometry, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
106 {
107 
108  assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalEndcap);
109 
110  const CaloSubdetectorGeometry* geom=caloGeometry.getSubdetectorGeometry(DetId::Ecal,EcalEndcap);//EcalBarrel = 1
111 
112  const math::XYZPoint position_ = bclus.position();
113  //double Theta = -position_.theta()+0.5*TMath::Pi();
114  double Eta = position_.eta();
115  double Phi = TVector2::Phi_mpi_pi(position_.phi());
116  double X = position_.x();
117  double Y = position_.y();
118 
119  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
120  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
121  const float X0 = 0.89; float T0 = 1.2;
122  //different T0 value if outside of preshower coverage
123  if (TMath::Abs(bclus.eta())<1.653) T0 = 3.1;
124 
125  double depth = X0 * (T0 + log(bclus.energy()));
126 
127  //find max energy crystal
128  std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
129  float drmin = 999.;
130  EEDetId crystalseed;
131  //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
132  for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {
133 
134  EEDetId crystal(crystals_vector[icry].first);
135 
136  const CaloCellGeometry* cell=geom->getGeometry(crystal);
137  GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
138  double EtaCentr = center_pos.eta();
139  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
140 
141  float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
142  if (dr<drmin) {
143  drmin = dr;
144  crystalseed = crystal;
145  }
146 
147  }
148 
149  ix = crystalseed.ix();
150  iy = crystalseed.iy();
151 
152  // Get center cell position from shower depth
153  const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
154  const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);
155 
156  thetatilt = cpyr->getThetaAxis();
157  phitilt = cpyr->getPhiAxis();
158 
159  GlobalPoint center_pos = cpyr->getPosition(depth);
160 
161  double XCentr = center_pos.x();
162  double XWidth = 2.59;
163  xcry = (X-XCentr)/XWidth;
164 
165 
166  double YCentr = center_pos.y();
167  double YWidth = 2.59;
168  ycry = (Y-YCentr)/YWidth;
169 
170  return;
171 
172 }
const double Pi
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:126
int ix() const
Definition: EEDetId.h:76
assert(m_qm.get())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:192
double Phi_mpi_pi(double x)
Definition: JetUtil.h:24
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
CCGFloat getPhiAxis() const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:163
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
T Abs(T a)
Definition: MathUtil.h:49
double energy() const
cluster energy
Definition: CaloCluster.h:121
int iy() const
Definition: EEDetId.h:82
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
CCGFloat getThetaAxis() const
const GlobalPoint getPosition(CCGFloat depth) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const T & get() const
Definition: EventSetup.h:56
void localCoordsEB(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
A base class to handle the particular shape of Ecal Xtals. Taken from ORCA Calorimetry Code...
T eta() const
Definition: PV3DBase.h:76
void localCoordsEE(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
static const double X0
T x() const
Definition: PV3DBase.h:62
tuple log
Definition: cmsBatch.py:341