CMS 3D CMS Logo

EcalClusterLocal.cc
Go to the documentation of this file.
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  auto cell=geom->getGeometry(crystal);
58  const TruncatedPyramid* cpyr = dynamic_cast<const TruncatedPyramid*>(cell.get());
59  GlobalPoint center_pos = cpyr->getPosition(depth);
60  double EtaCentr = center_pos.eta();
61  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
62 
63  float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
64  if (dr<drmin) {
65  drmin = dr;
66  crystalseed = crystal;
67  }
68 
69  }
70 
71  ieta = crystalseed.ieta();
72  iphi = crystalseed.iphi();
73 
74  // Get center cell position from shower depth
75  auto cell=geom->getGeometry(crystalseed);
76  const TruncatedPyramid* cpyr = dynamic_cast<const TruncatedPyramid*>(cell.get());
77 
78  thetatilt = cpyr->getThetaAxis();
79  phitilt = cpyr->getPhiAxis();
80 
81  GlobalPoint center_pos = cpyr->getPosition(depth);
82 
83  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
84  double PhiWidth = (TMath::Pi()/180.);
85  phicry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
86  //Some flips to take into account ECAL barrel symmetries:
87  if (ieta<0) phicry *= -1.;
88 
89  double ThetaCentr = -center_pos.theta()+0.5*TMath::Pi();
90  double ThetaWidth = (TMath::Pi()/180.)*TMath::Cos(ThetaCentr);
91  etacry = (Theta-ThetaCentr)/ThetaWidth;
92  //flip to take into account ECAL barrel symmetries:
93  if (ieta<0) etacry *= -1.;
94 
95  return;
96 
97 }
98 
99 void EcalClusterLocal::localCoordsEE( const reco::CaloCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
100 {
102  es.get<CaloGeometryRecord>().get(pG);
103  localCoordsEE( bclus, *pG, xcry, ycry, ix, iy, thetatilt, phitilt);
104 }
105 
106 void EcalClusterLocal::localCoordsEE( const reco::CaloCluster &bclus, const CaloGeometry & caloGeometry, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
107 {
108 
109  assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalEndcap);
110 
111  const CaloSubdetectorGeometry* geom=caloGeometry.getSubdetectorGeometry(DetId::Ecal,EcalEndcap);//EcalBarrel = 1
112 
113  const math::XYZPoint& position_ = bclus.position();
114  //double Theta = -position_.theta()+0.5*TMath::Pi();
115  double Eta = position_.eta();
116  double Phi = TVector2::Phi_mpi_pi(position_.phi());
117  double X = position_.x();
118  double Y = position_.y();
119 
120  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
121  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
122  const float X0 = 0.89; float T0 = 1.2;
123  //different T0 value if outside of preshower coverage
124  if (TMath::Abs(bclus.eta())<1.653) T0 = 3.1;
125 
126  double depth = X0 * (T0 + log(bclus.energy()));
127 
128  //find max energy crystal
129  std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
130  float drmin = 999.;
131  EEDetId crystalseed;
132  //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
133  for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {
134 
135  EEDetId crystal(crystals_vector[icry].first);
136 
137  auto cell=geom->getGeometry(crystal);
138  const TruncatedPyramid* cpyr = dynamic_cast<const TruncatedPyramid*>(cell.get());
139  GlobalPoint center_pos = cpyr->getPosition(depth);
140  double EtaCentr = center_pos.eta();
141  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
142 
143  float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
144  if (dr<drmin) {
145  drmin = dr;
146  crystalseed = crystal;
147  }
148 
149  }
150 
151  ix = crystalseed.ix();
152  iy = crystalseed.iy();
153 
154  // Get center cell position from shower depth
155  auto cell=geom->getGeometry(crystalseed);
156  const TruncatedPyramid* cpyr = dynamic_cast<const TruncatedPyramid*>(cell.get());
157 
158  thetatilt = cpyr->getThetaAxis();
159  phitilt = cpyr->getPhiAxis();
160 
161  GlobalPoint center_pos = cpyr->getPosition(depth);
162 
163  double XCentr = center_pos.x();
164  double XWidth = 2.59;
165  xcry = (X-XCentr)/XWidth;
166 
167 
168  double YCentr = center_pos.y();
169  double YWidth = 2.59;
170  ycry = (Y-YCentr)/YWidth;
171 
172  return;
173 
174 }
const double Pi
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:131
GlobalPoint getPosition(CCGFloat depth) const override
int ix() const
Definition: EEDetId.h:76
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:197
CCGFloat getPhiAxis() const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
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:126
int iy() const
Definition: EEDetId.h:82
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
CCGFloat getThetaAxis() const
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
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
T get() const
Definition: EventSetup.h:68
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