CMS 3D CMS Logo

Jet.cc
Go to the documentation of this file.
1 // Jet.cc
2 // Fedor Ratnikov, UMd
3 
4 #include <sstream>
5 #include <cmath>
6 
10 
11 //Own header file
13 
14 using namespace reco;
15 
16 namespace {
17  // approximate simple CALO geometry
18  // abstract baseclass for geometry.
19 
20  class CaloPoint {
21  public:
22  static const double depth; // one for all relative depth of the reference point between ECAL begin and HCAL end
23  static const double R_BARREL;
24  static const double R_BARREL2;
25  static const double Z_ENDCAP; // 1/2(EEz+HEz)
26  static const double R_FORWARD; // eta=3
27  static const double R_FORWARD2;
28  static const double Z_FORWARD;
29  static const double Z_BIG;
30  };
31 
32  // one for all relative depth of the reference point between ECAL begin and HCAL end
33  const double CaloPoint::depth = 0.1;
34  const double CaloPoint::R_BARREL = (1. - depth) * 143. + depth * 407.;
35  const double CaloPoint::R_BARREL2 = R_BARREL * R_BARREL;
36  const double CaloPoint::Z_ENDCAP = (1. - depth) * 320. + depth * 568.; // 1/2(EEz+HEz)
37  const double CaloPoint::R_FORWARD = Z_ENDCAP / std::sqrt(std::cosh(3.) * std::cosh(3.) - 1.); // eta=3
38  const double CaloPoint::R_FORWARD2 = R_FORWARD * R_FORWARD;
39  const double CaloPoint::Z_FORWARD = 1100. + depth * 165.;
40  const double CaloPoint::Z_BIG = 1.e5;
41 
42  //old zvertex only implementation:
43  class CaloPointZ : private CaloPoint {
44  public:
45  CaloPointZ(double fZ, double fEta) {
46  static const double ETA_MAX = 5.2;
47 
48  if (fZ > Z_ENDCAP)
49  fZ = Z_ENDCAP - 1.;
50  if (fZ < -Z_ENDCAP)
51  fZ = -Z_ENDCAP + 1; // sanity check
52 
53  double tanThetaAbs = std::sqrt(std::cosh(fEta) * std::cosh(fEta) - 1.);
54  double tanTheta = fEta >= 0 ? tanThetaAbs : -tanThetaAbs;
55 
56  double rEndcap = tanTheta == 0 ? 1.e10 : fEta > 0 ? (Z_ENDCAP - fZ) / tanTheta : (-Z_ENDCAP - fZ) / tanTheta;
57  if (rEndcap > R_BARREL) { // barrel
58  mR = R_BARREL;
59  mZ = fZ + R_BARREL * tanTheta;
60  } else {
61  double zRef = Z_BIG; // very forward;
62  if (rEndcap > R_FORWARD)
63  zRef = Z_ENDCAP; // endcap
64  else if (std::fabs(fEta) < ETA_MAX)
65  zRef = Z_FORWARD; // forward
66 
67  mZ = fEta > 0 ? zRef : -zRef;
68  mR = std::fabs((mZ - fZ) / tanTheta);
69  }
70  }
71  double etaReference(double fZ) {
72  Jet::Point p(r(), 0., z() - fZ);
73  return p.eta();
74  }
75 
76  double thetaReference(double fZ) {
77  Jet::Point p(r(), 0., z() - fZ);
78  return p.theta();
79  }
80 
81  double z() const { return mZ; }
82  double r() const { return mR; }
83 
84  private:
85  CaloPointZ(){};
86  double mZ;
87  double mR;
88  };
89 
90  //new implementation to derive CaloPoint for free 3d vertex.
91  //code provided thanks to Christophe Saout
92  template <typename Point>
93  class CaloPoint3D : private CaloPoint {
94  public:
95  template <typename Vector, typename Point2>
96  CaloPoint3D(const Point2 &vertex, const Vector &dir) {
97  // note: no sanity checks here, make sure vertex is inside the detector!
98 
99  // check if positive or negative (or none) endcap should be tested
100  int side = dir.z() < -1e-9 ? -1 : dir.z() > 1e-9 ? +1 : 0;
101 
102  double dirR = dir.Rho();
103 
104  // normalized direction in x-y plane
105  double dirUnit[2] = {dir.x() / dirR, dir.y() / dirR};
106 
107  // rotate the vertex into a coordinate system where dir lies along x
108 
109  // vtxLong is the longitudinal coordinate of the vertex wrt/ dir
110  double vtxLong = dirUnit[0] * vertex.x() + dirUnit[1] * vertex.y();
111 
112  // tIP is the (signed) transverse impact parameter
113  double tIP = dirUnit[0] * vertex.y() - dirUnit[1] * vertex.x();
114 
115  // r and z coordinate
116  double r, z;
117 
118  if (side) {
119  double slope = dirR / dir.z();
120 
121  // check extrapolation to endcap
122  r = vtxLong + slope * (side * Z_ENDCAP - vertex.z());
123  double r2 = sqr(r) + sqr(tIP);
124 
125  if (r2 < R_FORWARD2) {
126  // we are in the forward calorimeter, recompute
127  r = vtxLong + slope * (side * Z_FORWARD - vertex.z());
128  z = side * Z_FORWARD;
129  } else if (r2 < R_BARREL2) {
130  // we are in the endcap
131  z = side * Z_ENDCAP;
132  } else {
133  // we are in the barrel, do the intersection below
134  side = 0;
135  }
136  }
137 
138  if (!side) {
139  // we are in the barrel
140  double slope = dir.z() / dirR;
141  r = std::sqrt(R_BARREL2 - sqr(tIP));
142  z = vertex.z() + slope * (r - vtxLong);
143  }
144 
145  // rotate (r, tIP, z) back into original x-y coordinate system
146  point = Point(dirUnit[0] * r - dirUnit[1] * tIP, dirUnit[1] * r + dirUnit[0] * tIP, z);
147  }
148 
149  const Point &caloPoint() const { return point; }
150 
151  private:
152  template <typename T>
153  static inline T sqr(const T &value) {
154  return value * value;
155  }
156 
157  Point point;
158  };
159 
160 } // namespace
161 
162 Jet::Jet(const LorentzVector &fP4, const Point &fVertex, const Constituents &fConstituents)
163  : CompositePtrCandidate(0, fP4, fVertex), mJetArea(0), mPileupEnergy(0), mPassNumber(0), mIsWeighted(false) {
164  for (unsigned i = 0; i < fConstituents.size(); i++) {
165  addDaughter(fConstituents[i]);
166  }
167 }
168 
169 Jet::Jet(const LorentzVector &fP4, const Point &fVertex)
170  : CompositePtrCandidate(0, fP4, fVertex), mJetArea(0), mPileupEnergy(0), mPassNumber(0), mIsWeighted(false) {}
171 
173 Jet::EtaPhiMoments Jet::etaPhiStatistics() const {
174  std::vector<const Candidate *> towers = getJetConstituentsQuick();
175  double phiRef = phi();
176  double sumw = 0;
177  double sumEta = 0;
178  double sumEta2 = 0;
179  double sumPhi = 0;
180  double sumPhi2 = 0;
181  double sumEtaPhi = 0;
182  int i = towers.size();
183  while (--i >= 0) {
184  double eta = towers[i]->eta();
185  double phi = deltaPhi(towers[i]->phi(), phiRef);
186  double weight = towers[i]->et();
187  sumw += weight;
188  sumEta += eta * weight;
189  sumEta2 += eta * eta * weight;
190  sumPhi += phi * weight;
191  sumPhi2 += phi * phi * weight;
192  sumEtaPhi += eta * phi * weight;
193  }
194  Jet::EtaPhiMoments result;
195  if (sumw > 0) {
196  result.etaMean = sumEta / sumw;
197  result.phiMean = deltaPhi(phiRef + sumPhi, 0.);
198  result.etaEtaMoment = (sumEta2 - sumEta * sumEta / sumw) / sumw;
199  result.phiPhiMoment = (sumPhi2 - sumPhi * sumPhi / sumw) / sumw;
200  result.etaPhiMoment = (sumEtaPhi - sumEta * sumPhi / sumw) / sumw;
201  } else {
202  result.etaMean = 0;
203  result.phiMean = 0;
204  result.etaEtaMoment = 0;
205  result.phiPhiMoment = 0;
206  result.etaPhiMoment = 0;
207  }
208  return result;
209 }
210 
212 float Jet::etaetaMoment() const {
213  std::vector<const Candidate *> towers = getJetConstituentsQuick();
214  double sumw = 0;
215  double sum = 0;
216  double sum2 = 0;
217  int i = towers.size();
218  while (--i >= 0) {
219  double value = towers[i]->eta();
220  double weight = towers[i]->et();
221  sumw += weight;
222  sum += value * weight;
223  sum2 += value * value * weight;
224  }
225  return sumw > 0 ? (sum2 - sum * sum / sumw) / sumw : 0;
226 }
227 
229 float Jet::phiphiMoment() const {
230  double phiRef = phi();
231  std::vector<const Candidate *> towers = getJetConstituentsQuick();
232  double sumw = 0;
233  double sum = 0;
234  double sum2 = 0;
235  int i = towers.size();
236  while (--i >= 0) {
237  double value = deltaPhi(towers[i]->phi(), phiRef);
238  double weight = towers[i]->et();
239  sumw += weight;
240  sum += value * weight;
241  sum2 += value * value * weight;
242  }
243  return sumw > 0 ? (sum2 - sum * sum / sumw) / sumw : 0;
244 }
245 
247 float Jet::etaphiMoment() const {
248  double phiRef = phi();
249  std::vector<const Candidate *> towers = getJetConstituentsQuick();
250  double sumw = 0;
251  double sumA = 0;
252  double sumB = 0;
253  double sumAB = 0;
254  int i = towers.size();
255  while (--i >= 0) {
256  double valueA = towers[i]->eta();
257  double valueB = deltaPhi(towers[i]->phi(), phiRef);
258  double weight = towers[i]->et();
259  sumw += weight;
260  sumA += valueA * weight;
261  sumB += valueB * weight;
262  sumAB += valueA * valueB * weight;
263  }
264  return sumw > 0 ? (sumAB - sumA * sumB / sumw) / sumw : 0;
265 }
266 
268 float Jet::etInAnnulus(float fRmin, float fRmax) const {
269  float result = 0;
270  std::vector<const Candidate *> towers = getJetConstituentsQuick();
271  int i = towers.size();
272  while (--i >= 0) {
273  double r = deltaR(*this, *(towers[i]));
274  if (r >= fRmin && r < fRmax)
275  result += towers[i]->et();
276  }
277  return result;
278 }
279 
281 int Jet::nCarrying(float fFraction) const {
282  std::vector<const Candidate *> towers = getJetConstituentsQuick();
283  if (fFraction >= 1)
284  return towers.size();
285  double totalEt = 0;
286  for (unsigned i = 0; i < towers.size(); ++i)
287  totalEt += towers[i]->et();
288  double fractionEnergy = totalEt * fFraction;
289  unsigned result = 0;
290  for (; result < towers.size(); ++result) {
291  fractionEnergy -= towers[result]->et();
292  if (fractionEnergy <= 0)
293  return result + 1;
294  }
295  return 0;
296 }
297 
299 float Jet::maxDistance() const {
300  float result = 0;
301  std::vector<const Candidate *> towers = getJetConstituentsQuick();
302  for (unsigned i = 0; i < towers.size(); ++i) {
303  float dR = deltaR(*(towers[i]), *this);
304  if (dR > result)
305  result = dR;
306  }
307  return result;
308 }
309 
311 // kept for backwards compatibility, use detector/physicsP4 instead!
312 float Jet::physicsEta(float fZVertex, float fDetectorEta) {
313  CaloPointZ refPoint(0., fDetectorEta);
314  return refPoint.etaReference(fZVertex);
315 }
316 
318 // kept for backwards compatibility, use detector/physicsP4 instead!
319 float Jet::detectorEta(float fZVertex, float fPhysicsEta) {
320  CaloPointZ refPoint(fZVertex, fPhysicsEta);
321  return refPoint.etaReference(0.);
322 }
323 
325  const Candidate &inParticle,
326  const Candidate::Point &oldVertex) {
327  CaloPoint3D<Point> caloPoint(oldVertex, inParticle.momentum()); // Jet position in Calo.
328  Vector physicsDir = caloPoint.caloPoint() - newVertex;
329  double p = inParticle.momentum().r();
330  Vector p3 = p * physicsDir.unit();
331  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
332  return returnVector;
333 }
334 
336  CaloPoint3D<Point> caloPoint(vertex, inParticle.momentum()); // Jet position in Calo.
337  static const Point np(0, 0, 0);
338  Vector detectorDir = caloPoint.caloPoint() - np;
339  double p = inParticle.momentum().r();
340  Vector p3 = p * detectorDir.unit();
341  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
342  return returnVector;
343 }
344 
347  for (unsigned i = 0; i < numberOfDaughters(); i++) {
348  result.push_back(daughterPtr(i));
349  }
350  return result;
351 }
352 
353 std::vector<const Candidate *> Jet::getJetConstituentsQuick() const {
354  std::vector<const Candidate *> result;
355  int nDaughters = numberOfDaughters();
356  for (int i = 0; i < nDaughters; ++i) {
357  if (!(this->sourceCandidatePtr(i).isNonnull() && this->sourceCandidatePtr(i).isAvailable())) {
358  edm::LogInfo("JetConstituentPointer")
359  << "Bad pointer to jet constituent found. Check for possible dropped particle.";
360  continue;
361  }
362  result.push_back(daughter(i));
363  }
364  return result;
365 }
366 
367 float Jet::constituentPtDistribution() const {
368  Jet::Constituents constituents = this->getJetConstituents();
369 
370  float sum_pt2 = 0.;
371  float sum_pt = 0.;
372 
373  for (unsigned iConst = 0; iConst < constituents.size(); ++iConst) {
374  float pt = constituents[iConst]->p4().Pt();
375  float pt2 = pt * pt;
376 
377  sum_pt += pt;
378  sum_pt2 += pt2;
379 
380  } //for constituents
381 
382  float ptD_value = (sum_pt > 0.) ? sqrt(sum_pt2 / (sum_pt * sum_pt)) : 0.;
383 
384  return ptD_value;
385 
386 } //constituentPtDistribution
387 
388 float Jet::constituentEtaPhiSpread() const {
389  Jet::Constituents constituents = this->getJetConstituents();
390 
391  float sum_pt2 = 0.;
392  float sum_pt2deltaR2 = 0.;
393 
394  for (unsigned iConst = 0; iConst < constituents.size(); ++iConst) {
395  LorentzVector thisConstituent = constituents[iConst]->p4();
396 
397  float pt = thisConstituent.Pt();
398  float pt2 = pt * pt;
399  double dR = deltaR(*this, *(constituents[iConst]));
400  float pt2deltaR2 = pt * pt * dR * dR;
401 
402  sum_pt2 += pt2;
403  sum_pt2deltaR2 += pt2deltaR2;
404 
405  } //for constituents
406 
407  float rmsCand_value = (sum_pt2 > 0.) ? sum_pt2deltaR2 / sum_pt2 : 0.;
408 
409  return rmsCand_value;
410 
411 } //constituentEtaPhiSpread
412 
413 std::string Jet::print() const {
414  std::ostringstream out;
415  out << "Jet p/px/py/pz/pt: " << p() << '/' << px() << '/' << py() << '/' << pz() << '/' << pt() << std::endl
416  << " eta/phi: " << eta() << '/' << phi() << std::endl
417  << " # of constituents: " << nConstituents() << std::endl;
418  out << " Constituents:" << std::endl;
419  for (unsigned index = 0; index < numberOfDaughters(); index++) {
420  Constituent constituent = daughterPtr(index); // deref
421  if (constituent.isNonnull()) {
422  out << " #" << index << " p/pt/eta/phi: " << constituent->p() << '/' << constituent->pt() << '/'
423  << constituent->eta() << '/' << constituent->phi() << " productId/index: " << constituent.id() << '/'
424  << constituent.key() << std::endl;
425  } else {
426  out << " #" << index << " constituent is not available in the event" << std::endl;
427  }
428  }
429  return out.str();
430 }
431 
432 void Jet::scaleEnergy(double fScale) { setP4(p4() * fScale); }
433 
434 bool Jet::isJet() const { return true; }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
float etaetaMoment() const
eta-eta second moment, ET weighted
virtual CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
virtual std::vector< const reco::Candidate * > getJetConstituentsQuick() const
quick list of constituents
virtual double energy() const =0
energy
static float detectorEta(float fZVertex, float fPhysicsEta)
static function to convert physics eta to detector eta
double pz() const final
z coordinate of momentum vector
virtual void scaleEnergy(double fScale)
scale energy of the jet
double pt() const final
transverse momentum
float etInAnnulus(float fRmin, float fRmax) const
ET in annulus between rmin and rmax around jet direction.
int nCarrying(float fFraction) const
return # of constituent carrying fraction of energy
CandidatePtr sourceCandidatePtr(size_type i) const override
static const double slope[3]
std::vector< Constituent > Constituents
Definition: Jet.h:23
Definition: weight.py:1
const Point & vertex() const override
vertex position (overwritten by PF...)
float maxDistance() const
maximum distance from jet to constituent
float float float z
const LorentzVector & p4() const final
four-momentum Lorentz vector
float phiphiMoment() const
phi-phi second moment, ET weighted
virtual Constituents getJetConstituents() const
list of constituents
double px() const final
x coordinate of momentum vector
virtual Vector momentum() const =0
spatial momentum vector
double p() const final
magnitude of momentum vector
int np
Definition: AMPTWrapper.h:43
float constituentPtDistribution() const
fixed size vector
Definition: Vector.h:24
T sqrt(T t)
Definition: SSEVec.h:19
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
virtual std::string print() const
Print object.
bool isJet() const override
math::XYZPoint Point
math::XYZTLorentzVector LorentzVector
Definition: value.py:1
unsigned int index
Definition: LeafCandidate.h:31
double py() const final
y coordinate of momentum vector
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
edm::Ptr< Candidate > Constituent
Definition: Jet.h:22
Log< level::Info, false > LogInfo
size_t numberOfDaughters() const override
number of daughters
static Candidate::LorentzVector detectorP4(const Candidate::Point &vertex, const Candidate &inParticle)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
fixed size matrix
double et() const final
transverse energy
Structure Point Contains parameters of Gaussian fits to DMRs.
Square< F >::type sqr(const F &f)
Definition: Square.h:14
void addDaughter(const CandidatePtr &)
add a daughter via a reference
Jet()
Default constructor.
Definition: Jet.h:36
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
static Candidate::LorentzVector physicsP4(const Candidate::Point &newVertex, const Candidate &inParticle, const Candidate::Point &oldVertex=Candidate::Point(0, 0, 0))
long double T
double phi() const final
momentum azimuthal angle
EtaPhiMoments etaPhiStatistics() const
eta-phi statistics, ET weighted
static float physicsEta(float fZVertex, float fDetectorEta)
static function to convert detector eta to physics eta
float constituentEtaPhiSpread() const
void setP4(const LorentzVector &p4) final
set 4-momentum
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
math::XYZPoint Point
point in the space
Definition: LeafCandidate.h:27
float etaphiMoment() const
eta-phi second moment, ET weighted
double eta() const final
momentum pseudorapidity