test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
15 
16 using namespace reco;
17 
18 namespace {
19  // approximate simple CALO geometry
20  // abstract baseclass for geometry.
21 
22  class CaloPoint {
23  public:
24  static const double depth; // one for all relative depth of the reference point between ECAL begin and HCAL end
25  static const double R_BARREL;
26  static const double R_BARREL2;
27  static const double Z_ENDCAP; // 1/2(EEz+HEz)
28  static const double R_FORWARD; // eta=3
29  static const double R_FORWARD2;
30  static const double Z_FORWARD;
31  static const double Z_BIG;
32  };
33 
34  const double CaloPoint::depth = 0.1; // one for all relative depth of the reference point between ECAL begin and HCAL end
35  const double CaloPoint::R_BARREL = (1.-depth)*143.+depth*407.;
36  const double CaloPoint::R_BARREL2 = R_BARREL * R_BARREL;
37  const double CaloPoint::Z_ENDCAP = (1.-depth)*320.+depth*568.; // 1/2(EEz+HEz)
38  const double CaloPoint::R_FORWARD = Z_ENDCAP / std::sqrt (std::cosh(3.)*std::cosh(3.) -1.); // eta=3
39  const double CaloPoint::R_FORWARD2 = R_FORWARD * R_FORWARD;
40  const double CaloPoint::Z_FORWARD = 1100.+depth*165.;
41  const double CaloPoint::Z_BIG = 1.e5;
42 
43  //old zvertex only implementation:
44  class CaloPointZ: private CaloPoint{
45  public:
46  CaloPointZ (double fZ, double fEta){
47 
48  static const double ETA_MAX = 5.2;
49 
50  if (fZ > Z_ENDCAP) fZ = Z_ENDCAP-1.;
51  if (fZ < -Z_ENDCAP) 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 :
57  fEta > 0 ? (Z_ENDCAP - fZ) / tanTheta : (-Z_ENDCAP - fZ) / tanTheta;
58  if (rEndcap > R_BARREL) { // barrel
59  mR = R_BARREL;
60  mZ = fZ + R_BARREL * tanTheta;
61  }
62  else {
63  double zRef = Z_BIG; // very forward;
64  if (rEndcap > R_FORWARD) zRef = Z_ENDCAP; // endcap
65  else if (std::fabs (fEta) < ETA_MAX) 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  private:
84  CaloPointZ(){};
85  double mZ;
86  double mR;
87  };
88 
89  //new implementation to derive CaloPoint for free 3d vertex.
90  //code provided thanks to Christophe Saout
91  template<typename Point>
92  class CaloPoint3D : private CaloPoint {
93  public:
94  template<typename Vector, typename Point2>
95  CaloPoint3D(const Point2 &vertex, const Vector &dir)
96  {
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,
147  dirUnit[1] * r + dirUnit[0] * tIP,
148  z);
149  }
150 
151  const Point &caloPoint() const { return point; }
152 
153  private:
154  template<typename T>
155  static inline T sqr(const T &value) { return value * value; }
156 
157  Point point;
158  };
159 
160 }
161 
162 Jet::Jet (const LorentzVector& fP4,
163  const Point& fVertex,
164  const Constituents& fConstituents)
165  : CompositePtrCandidate (0, fP4, fVertex),
166  mJetArea (0),
167  mPileupEnergy (0),
168  mPassNumber (0)
169 {
170  for (unsigned i = 0; i < fConstituents.size (); i++) {
171  addDaughter (fConstituents [i]);
172  }
173 }
174 
175 Jet::Jet (const LorentzVector& fP4,
176  const Point& fVertex)
177  : CompositePtrCandidate (0, fP4, fVertex),
178  mJetArea (0),
179  mPileupEnergy (0),
180  mPassNumber (0)
181 {}
182 
184 Jet::EtaPhiMoments Jet::etaPhiStatistics () const {
185  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
186  double phiRef = phi();
187  double sumw = 0;
188  double sumEta = 0;
189  double sumEta2 = 0;
190  double sumPhi = 0;
191  double sumPhi2 = 0;
192  double sumEtaPhi = 0;
193  int i = towers.size();
194  while (--i >= 0) {
195  double eta = towers[i]->eta();
196  double phi = deltaPhi (towers[i]->phi(), phiRef);
197  double weight = towers[i]->et();
198  sumw += weight;
199  sumEta += eta * weight;
200  sumEta2 += eta * eta * weight;
201  sumPhi += phi * weight;
202  sumPhi2 += phi * phi * weight;
203  sumEtaPhi += eta * phi * weight;
204  }
205  Jet::EtaPhiMoments result;
206  if (sumw > 0) {
207  result.etaMean = sumEta / sumw;
208  result.phiMean = deltaPhi (phiRef + sumPhi, 0.);
209  result.etaEtaMoment = (sumEta2 - sumEta * sumEta / sumw) / sumw;
210  result.phiPhiMoment = (sumPhi2 - sumPhi * sumPhi / sumw) / sumw;
211  result.etaPhiMoment = (sumEtaPhi - sumEta * sumPhi / sumw) / sumw;
212  }
213  else {
214  result.etaMean = 0;
215  result.phiMean = 0;
216  result.etaEtaMoment = 0;
217  result.phiPhiMoment = 0;
218  result.etaPhiMoment = 0;
219  }
220  return result;
221 }
222 
224 float Jet::etaetaMoment () const {
225  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
226  double sumw = 0;
227  double sum = 0;
228  double sum2 = 0;
229  int i = towers.size();
230  while (--i >= 0) {
231  double value = towers[i]->eta();
232  double weight = towers[i]->et();
233  sumw += weight;
234  sum += value * weight;
235  sum2 += value * value * weight;
236  }
237  return sumw > 0 ? (sum2 - sum*sum/sumw ) / sumw : 0;
238 }
239 
241 float Jet::phiphiMoment () const {
242  double phiRef = phi();
243  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
244  double sumw = 0;
245  double sum = 0;
246  double sum2 = 0;
247  int i = towers.size();
248  while (--i >= 0) {
249  double value = deltaPhi (towers[i]->phi(), phiRef);
250  double weight = towers[i]->et();
251  sumw += weight;
252  sum += value * weight;
253  sum2 += value * value * weight;
254  }
255  return sumw > 0 ? (sum2 - sum*sum/sumw ) / sumw : 0;
256 }
257 
259 float Jet::etaphiMoment () const {
260  double phiRef = phi();
261  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
262  double sumw = 0;
263  double sumA = 0;
264  double sumB = 0;
265  double sumAB = 0;
266  int i = towers.size();
267  while (--i >= 0) {
268  double valueA = towers[i]->eta();
269  double valueB = deltaPhi (towers[i]->phi(), phiRef);
270  double weight = towers[i]->et();
271  sumw += weight;
272  sumA += valueA * weight;
273  sumB += valueB * weight;
274  sumAB += valueA * valueB * weight;
275  }
276  return sumw > 0 ? (sumAB - sumA*sumB/sumw ) / sumw : 0;
277 }
278 
280 float Jet::etInAnnulus (float fRmin, float fRmax) const {
281  float result = 0;
282  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
283  int i = towers.size ();
284  while (--i >= 0) {
285  double r = deltaR (*this, *(towers[i]));
286  if (r >= fRmin && r < fRmax) result += towers[i]->et ();
287  }
288  return result;
289 }
290 
292 int Jet::nCarrying (float fFraction) const {
293  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
294  if (fFraction >= 1) return towers.size();
295  double totalEt = 0;
296  for (unsigned i = 0; i < towers.size(); ++i) totalEt += towers[i]->et();
297  double fractionEnergy = totalEt * fFraction;
298  unsigned result = 0;
299  for (; result < towers.size(); ++result) {
300  fractionEnergy -= towers[result]->et();
301  if (fractionEnergy <= 0) return result+1;
302  }
303  return 0;
304 }
305 
307 float Jet::maxDistance () const {
308  float result = 0;
309  std::vector<const Candidate*> towers = getJetConstituentsQuick ();
310  for (unsigned i = 0; i < towers.size(); ++i) {
311  float dR = deltaR (*(towers[i]), *this);
312  if (dR > result) result = dR;
313  }
314  return result;
315 }
316 
318 // kept for backwards compatibility, use detector/physicsP4 instead!
319 float Jet::physicsEta (float fZVertex, float fDetectorEta) {
320  CaloPointZ refPoint (0., fDetectorEta);
321  return refPoint.etaReference (fZVertex);
322 }
323 
325 // kept for backwards compatibility, use detector/physicsP4 instead!
326 float Jet::detectorEta (float fZVertex, float fPhysicsEta) {
327  CaloPointZ refPoint (fZVertex, fPhysicsEta);
328  return refPoint.etaReference (0.);
329 }
330 
331 Candidate::LorentzVector Jet::physicsP4 (const Candidate::Point &newVertex, const Candidate &inParticle,const Candidate::Point &oldVertex) {
332  CaloPoint3D<Point> caloPoint(oldVertex,inParticle.momentum()); // Jet position in Calo.
333  Vector physicsDir = caloPoint.caloPoint() - newVertex;
334  double p = inParticle.momentum().r();
335  Vector p3 = p * physicsDir.unit();
336  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
337  return returnVector;
338 }
339 
340 Candidate::LorentzVector Jet::detectorP4 (const Candidate::Point &vertex, const Candidate &inParticle) {
341  CaloPoint3D<Point> caloPoint(vertex,inParticle.momentum()); // Jet position in Calo.
342  static const Point np(0,0,0);
343  Vector detectorDir = caloPoint.caloPoint() - np;
344  double p = inParticle.momentum().r();
345  Vector p3 = p * detectorDir.unit();
346  LorentzVector returnVector(p3.x(), p3.y(), p3.z(), inParticle.energy());
347  return returnVector;
348 }
349 
350 
353  for (unsigned i = 0; i < numberOfDaughters(); i++) {
354  result.push_back (daughterPtr (i));
355  }
356  return result;
357 }
358 
359 std::vector<const Candidate*> Jet::getJetConstituentsQuick () const {
360  std::vector<const Candidate*> result;
361  int nDaughters = numberOfDaughters();
362  for (int i = 0; i < nDaughters; ++i) {
363  if(!(this->sourceCandidatePtr(i).isNonnull() && this->sourceCandidatePtr(i).isAvailable()))
364  {
365  edm::LogInfo("JetConstituentPointer") << "Bad pointer to jet constituent found. Check for possible dropped particle.";
366  continue;
367  }
368  result.push_back (daughter (i));
369  }
370  return result;
371 }
372 
373 
374 
375 float Jet::constituentPtDistribution() const {
376 
377  Jet::Constituents constituents = this->getJetConstituents();
378 
379  float sum_pt2 = 0.;
380  float sum_pt = 0.;
381 
382  for( unsigned iConst=0; iConst<constituents.size(); ++iConst ) {
383 
384  float pt = constituents[iConst]->p4().Pt();
385  float pt2 = pt*pt;
386 
387  sum_pt += pt;
388  sum_pt2 += pt2;
389 
390  } //for constituents
391 
392  float ptD_value = (sum_pt>0.) ? sqrt( sum_pt2 / (sum_pt*sum_pt) ) : 0.;
393 
394  return ptD_value;
395 
396 } //constituentPtDistribution
397 
398 
399 
400 float Jet::constituentEtaPhiSpread() const {
401 
402  Jet::Constituents constituents = this->getJetConstituents();
403 
404 
405  float sum_pt2 = 0.;
406  float sum_pt2deltaR2 = 0.;
407 
408  for( unsigned iConst=0; iConst<constituents.size(); ++iConst ) {
409 
410  LorentzVector thisConstituent = constituents[iConst]->p4();
411 
412  float pt = thisConstituent.Pt();
413  float pt2 = pt*pt;
414  double dR = deltaR (*this, *(constituents[iConst]));
415  float pt2deltaR2 = pt*pt*dR*dR;
416 
417  sum_pt2 += pt2;
418  sum_pt2deltaR2 += pt2deltaR2;
419 
420  } //for constituents
421 
422  float rmsCand_value = (sum_pt2>0.) ? sum_pt2deltaR2/sum_pt2 : 0.;
423 
424  return rmsCand_value;
425 
426 } //constituentEtaPhiSpread
427 
428 
429 
430 std::string Jet::print () const {
431  std::ostringstream out;
432  out << "Jet p/px/py/pz/pt: " << p() << '/' << px () << '/' << py() << '/' << pz() << '/' << pt() << std::endl
433  << " eta/phi: " << eta () << '/' << phi () << std::endl
434  << " # of constituents: " << nConstituents () << std::endl;
435  out << " Constituents:" << std::endl;
436  for (unsigned index = 0; index < numberOfDaughters(); index++) {
437  Constituent constituent = daughterPtr (index); // deref
438  if (constituent.isNonnull()) {
439  out << " #" << index << " p/pt/eta/phi: "
440  << constituent->p() << '/' << constituent->pt() << '/' << constituent->eta() << '/' << constituent->phi()
441  << " productId/index: " << constituent.id() << '/' << constituent.key() << std::endl;
442  }
443  else {
444  out << " #" << index << " constituent is not available in the event" << std::endl;
445  }
446  }
447  return out.str ();
448 }
449 
450 void Jet::scaleEnergy (double fScale) {
451  setP4 (p4() * fScale);
452 }
453 
454 bool Jet::isJet() const {
455  return true;
456 }
int i
Definition: DBlmapReader.cc:9
float etInAnnulus(float fRmin, float fRmax) const
ET in annulus between rmin and rmax around jet direction.
virtual double energy() const =0
energy
static float detectorEta(float fZVertex, float fPhysicsEta)
static function to convert physics eta to detector eta
float constituentEtaPhiSpread() const
virtual void scaleEnergy(double fScale)
scale energy of the jet
CandidatePtr daughterPtr(size_type i) const
reference to daughter at given position
static const double slope[3]
math::XYZTLorentzVector LorentzVector
virtual double phi() const final
momentum azimuthal angle
std::vector< Constituent > Constituents
Definition: Jet.h:23
float etaetaMoment() const
eta-eta second moment, ET weighted
virtual Constituents getJetConstituents() const
list of constituents
std::pair< double, double > Point
Definition: CaloEllipse.h:18
float float float z
tuple result
Definition: mps_fire.py:83
virtual Vector momentum() const =0
spatial momentum vector
virtual size_t numberOfDaughters() const
number of daughters
virtual CandidatePtr sourceCandidatePtr(size_type i) const
int np
Definition: AMPTWrapper.h:33
virtual std::string print() const
Print object.
fixed size vector
Definition: Vector.h:25
T sqrt(T t)
Definition: SSEVec.h:18
math::XYZPoint Point
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
virtual double py() const final
y coordinate of momentum vector
unsigned int index
Definition: LeafCandidate.h:31
virtual double pz() const final
z coordinate of momentum vector
edm::Ptr< Candidate > Constituent
Definition: Jet.h:22
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
EtaPhiMoments etaPhiStatistics() const
eta-phi statistics, ET weighted
static Candidate::LorentzVector detectorP4(const Candidate::Point &vertex, const Candidate &inParticle)
float etaphiMoment() const
eta-phi second moment, ET weighted
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
float maxDistance() const
maximum distance from jet to constituent
bool isJet() const
int nCarrying(float fFraction) const
return # of constituent carrying fraction of energy
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
virtual double p() const final
magnitude of momentum vector
float constituentPtDistribution() const
float phiphiMoment() const
phi-phi second moment, ET weighted
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
Square< F >::type sqr(const F &f)
Definition: Square.h:13
void addDaughter(const CandidatePtr &)
add a daughter via a reference
virtual double px() const final
x coordinate of momentum vector
Jet()
Default constructor.
Definition: Jet.h:36
virtual double et() const final
transverse energy
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
virtual double eta() const final
momentum pseudorapidity
dbl *** dir
Definition: mlp_gen.cc:35
static Candidate::LorentzVector physicsP4(const Candidate::Point &newVertex, const Candidate &inParticle, const Candidate::Point &oldVertex=Candidate::Point(0, 0, 0))
virtual std::vector< const reco::Candidate * > getJetConstituentsQuick() const
quick list of constituents
long double T
static float physicsEta(float fZVertex, float fDetectorEta)
static function to convert detector eta to physics eta
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
*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
virtual double pt() const final
transverse momentum
double p3[4]
Definition: TauolaWrapper.h:91