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