CMS 3D CMS Logo

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