Logo ROOT   6.18/05
Reference Guide
TLorentzVector.h
Go to the documentation of this file.
1// @(#)root/physics:$Id$
2// Author: Pasha Murat , Peter Malzacher 12/02/99
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#ifndef ROOT_TLorentzVector
13#define ROOT_TLorentzVector
14
15
16//////////////////////////////////////////////////////////////////////////
17// //
18// TLorentzVector //
19// //
20// Place holder for real lorentz vector class. //
21// //
22//////////////////////////////////////////////////////////////////////////
23
24#include "TMath.h"
25#include "TVector3.h"
26#include "TRotation.h"
27
28
30
31
32class TLorentzVector : public TObject {
33
34private:
35
36 TVector3 fP; // 3 vector component
37 Double_t fE; // time or energy of (x,y,z,t) or (px,py,pz,e)
38
39public:
40
41 typedef Double_t Scalar; // to be able to use it with the ROOT::Math::VectorUtil functions
42
43 enum { kX=0, kY=1, kZ=2, kT=3, kNUM_COORDINATES=4, kSIZE=kNUM_COORDINATES };
44 // Safe indexing of the coordinates when using with matrices, arrays, etc.
45
47
49 // Constructor giving the components x, y, z, t.
50
51 TLorentzVector(const Double_t * carray);
52 TLorentzVector(const Float_t * carray);
53 // Constructor from an array, not checked!
54
55 TLorentzVector(const TVector3 & vector3, Double_t t);
56 // Constructor giving a 3-Vector and a time component.
57
58 TLorentzVector(const TLorentzVector & lorentzvector);
59 // Copy constructor.
60
61 virtual ~TLorentzVector(){};
62 // Destructor
63
64 // inline operator TVector3 () const;
65 // inline operator TVector3 & ();
66 // Conversion (cast) to TVector3.
67
68 inline Double_t X() const;
69 inline Double_t Y() const;
70 inline Double_t Z() const;
71 inline Double_t T() const;
72 // Get position and time.
73
74 inline void SetX(Double_t a);
75 inline void SetY(Double_t a);
76 inline void SetZ(Double_t a);
77 inline void SetT(Double_t a);
78 // Set position and time.
79
80 inline Double_t Px() const;
81 inline Double_t Py() const;
82 inline Double_t Pz() const;
83 inline Double_t P() const;
84 inline Double_t E() const;
85 inline Double_t Energy() const;
86 // Get momentum and energy.
87
88 inline void SetPx(Double_t a);
89 inline void SetPy(Double_t a);
90 inline void SetPz(Double_t a);
91 inline void SetE(Double_t a);
92 // Set momentum and energy.
93
94 inline TVector3 Vect() const ;
95 // Get spatial component.
96
97 inline void SetVect(const TVector3 & vect3);
98 // Set spatial component.
99
100 inline Double_t Theta() const;
101 inline Double_t CosTheta() const;
102 inline Double_t Phi() const; //returns phi from -pi to pi
103 inline Double_t Rho() const;
104 // Get spatial vector components in spherical coordinate system.
105
106 inline void SetTheta(Double_t theta);
107 inline void SetPhi(Double_t phi);
108 inline void SetRho(Double_t rho);
109 // Set spatial vector components in spherical coordinate system.
110
111 inline void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e);
112 inline void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t);
113 inline void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m);
114 inline void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m);
115 inline void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e);
116 // Setters to provide the functionality (but a more meanigful name) of
117 // the previous version eg SetV4... PsetV4...
118
119 inline void GetXYZT(Double_t *carray) const;
120 inline void GetXYZT(Float_t *carray) const;
121 // Getters into an arry
122 // no checking!
123
124 Double_t operator () (int i) const;
125 inline Double_t operator [] (int i) const;
126 // Get components by index.
127
128 Double_t & operator () (int i);
129 inline Double_t & operator [] (int i);
130 // Set components by index.
131
132 inline TLorentzVector & operator = (const TLorentzVector &);
133 // Assignment.
134
135 inline TLorentzVector operator + (const TLorentzVector &) const;
137 // Additions.
138
139 inline TLorentzVector operator - (const TLorentzVector &) const;
141 // Subtractions.
142
143 inline TLorentzVector operator - () const;
144 // Unary minus.
145
146 inline TLorentzVector operator * (Double_t a) const;
148 // Scaling with real numbers.
149
150 inline Bool_t operator == (const TLorentzVector &) const;
151 inline Bool_t operator != (const TLorentzVector &) const;
152 // Comparisons.
153
154 inline Double_t Perp2() const;
155 // Transverse component of the spatial vector squared.
156
157 inline Double_t Pt() const;
158 inline Double_t Perp() const;
159 // Transverse component of the spatial vector (R in cylindrical system).
160
161 inline void SetPerp(Double_t);
162 // Set the transverse component of the spatial vector.
163
164 inline Double_t Perp2(const TVector3 & v) const;
165 // Transverse component of the spatial vector w.r.t. given axis squared.
166
167 inline Double_t Pt(const TVector3 & v) const;
168 inline Double_t Perp(const TVector3 & v) const;
169 // Transverse component of the spatial vector w.r.t. given axis.
170
171 inline Double_t Et2() const;
172 // Transverse energy squared.
173
174 inline Double_t Et() const;
175 // Transverse energy.
176
177 inline Double_t Et2(const TVector3 &) const;
178 // Transverse energy w.r.t. given axis squared.
179
180 inline Double_t Et(const TVector3 &) const;
181 // Transverse energy w.r.t. given axis.
182
183 inline Double_t DeltaPhi(const TLorentzVector &) const;
184 inline Double_t DeltaR(const TLorentzVector &) const;
185 inline Double_t DrEtaPhi(const TLorentzVector &) const;
186 inline TVector2 EtaPhiVector();
187
188 inline Double_t Angle(const TVector3 & v) const;
189 // Angle wrt. another vector.
190
191 inline Double_t Mag2() const;
192 inline Double_t M2() const;
193 // Invariant mass squared.
194
195 inline Double_t Mag() const;
196 inline Double_t M() const;
197 // Invariant mass. If mag2() is negative then -sqrt(-mag2()) is returned.
198
199 inline Double_t Mt2() const;
200 // Transverse mass squared.
201
202 inline Double_t Mt() const;
203 // Transverse mass.
204
205 inline Double_t Beta() const;
206 inline Double_t Gamma() const;
207
208 inline Double_t Dot(const TLorentzVector &) const;
209 inline Double_t operator * (const TLorentzVector &) const;
210 // Scalar product.
211
212 inline void SetVectMag(const TVector3 & spatial, Double_t magnitude);
213 inline void SetVectM(const TVector3 & spatial, Double_t mass);
214 // Copy spatial coordinates, and set energy = sqrt(mass^2 + spatial^2)
215
216 inline Double_t Plus() const;
217 inline Double_t Minus() const;
218 // Returns t +/- z.
219 // Related to the positive/negative light-cone component,
220 // which some define this way and others define as (t +/- z)/sqrt(2)
221
222 inline TVector3 BoostVector() const ;
223 // Returns the spatial components divided by the time component.
224
226 inline void Boost(const TVector3 &);
227 // Lorentz boost.
228
229 Double_t Rapidity() const;
230 // Returns the rapidity, i.e. 0.5*ln((E+pz)/(E-pz))
231
232 inline Double_t Eta() const;
233 inline Double_t PseudoRapidity() const;
234 // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
235
236 inline void RotateX(Double_t angle);
237 // Rotate the spatial component around the x-axis.
238
239 inline void RotateY(Double_t angle);
240 // Rotate the spatial component around the y-axis.
241
242 inline void RotateZ(Double_t angle);
243 // Rotate the spatial component around the z-axis.
244
245 inline void RotateUz(TVector3 & newUzVector);
246 // Rotates the reference frame from Uz to newUz (unit vector).
247
248 inline void Rotate(Double_t, const TVector3 &);
249 // Rotate the spatial component around specified axis.
250
251 inline TLorentzVector & operator *= (const TRotation &);
252 inline TLorentzVector & Transform(const TRotation &);
253 // Transformation with HepRotation.
254
257 // Transformation with HepLorenzRotation.
258
259 virtual void Print(Option_t *option="") const;
260
261 ClassDef(TLorentzVector,4) // A four vector with (-,-,-,+) metric
262};
263
264
265//inline TLorentzVector operator * (const TLorentzVector &, Double_t a);
266// moved to TLorentzVector::operator * (Double_t a)
268// Scaling LorentzVector with a real number
269
270
271inline Double_t TLorentzVector::X() const { return fP.X(); }
272inline Double_t TLorentzVector::Y() const { return fP.Y(); }
273inline Double_t TLorentzVector::Z() const { return fP.Z(); }
274inline Double_t TLorentzVector::T() const { return fE; }
275
279inline void TLorentzVector::SetT(Double_t a) { fE = a; }
280
281inline Double_t TLorentzVector::Px() const { return X(); }
282inline Double_t TLorentzVector::Py() const { return Y(); }
283inline Double_t TLorentzVector::Pz() const { return Z(); }
284inline Double_t TLorentzVector::P() const { return fP.Mag(); }
285inline Double_t TLorentzVector::E() const { return T(); }
286inline Double_t TLorentzVector::Energy() const { return T(); }
287
292
293inline TVector3 TLorentzVector::Vect() const { return fP; }
294
295inline void TLorentzVector::SetVect(const TVector3 &p) { fP = p; }
296
298 return fP.Phi();
299}
300
302 return fP.Theta();
303}
304
306 return fP.CosTheta();
307}
308
309
311 return fP.Mag();
312}
313
315 fP.SetTheta(th);
316}
317
319 fP.SetPhi(phi);
320}
321
323 fP.SetMag(rho);
324}
325
327 fP.SetXYZ(x, y, z);
328 SetT(t);
329}
330
332 SetXYZT(px, py, pz, e);
333}
334
336 if ( m >= 0 )
337 SetXYZT( x, y, z, TMath::Sqrt(x*x+y*y+z*z+m*m) );
338 else
339 SetXYZT( x, y, z, TMath::Sqrt( TMath::Max((x*x+y*y+z*z-m*m), 0. ) ) );
340}
341
343 pt = TMath::Abs(pt);
344 SetXYZM(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,m);
345}
346
348 pt = TMath::Abs(pt);
349 SetXYZT(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,e);
350}
351
352inline void TLorentzVector::GetXYZT(Double_t *carray) const {
353 fP.GetXYZ(carray);
354 carray[3] = fE;
355}
356
357inline void TLorentzVector::GetXYZT(Float_t *carray) const{
358 fP.GetXYZ(carray);
359 carray[3] = fE;
360}
361
362inline Double_t & TLorentzVector::operator [] (int i) { return (*this)(i); }
363inline Double_t TLorentzVector::operator [] (int i) const { return (*this)(i); }
364
366 fP = q.Vect();
367 fE = q.T();
368 return *this;
369}
370
372 return TLorentzVector(fP+q.Vect(), fE+q.T());
373}
374
376 fP += q.Vect();
377 fE += q.T();
378 return *this;
379}
380
382 return TLorentzVector(fP-q.Vect(), fE-q.T());
383}
384
386 fP -= q.Vect();
387 fE -= q.T();
388 return *this;
389}
390
392 return TLorentzVector(-X(), -Y(), -Z(), -T());
393}
394
396 fP *= a;
397 fE *= a;
398 return *this;
399}
400
402 return TLorentzVector(a*X(), a*Y(), a*Z(), a*T());
403}
404
406 return (Vect() == q.Vect() && T() == q.T());
407}
408
410 return (Vect() != q.Vect() || T() != q.T());
411}
412
413inline Double_t TLorentzVector::Perp2() const { return fP.Perp2(); }
414
415inline Double_t TLorentzVector::Perp() const { return fP.Perp(); }
416
417inline Double_t TLorentzVector::Pt() const { return Perp(); }
418
420 fP.SetPerp(r);
421}
422
424 return fP.Perp2(v);
425}
426
428 return fP.Perp(v);
429}
430
431inline Double_t TLorentzVector::Pt(const TVector3 &v) const {
432 return Perp(v);
433}
434
436 Double_t pt2 = fP.Perp2();
437 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+Z()*Z());
438}
439
441 Double_t etet = Et2();
442 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
443}
444
445inline Double_t TLorentzVector::Et2(const TVector3 & v) const {
446 Double_t pt2 = fP.Perp2(v);
447 Double_t pv = fP.Dot(v.Unit());
448 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+pv*pv);
449}
450
451inline Double_t TLorentzVector::Et(const TVector3 & v) const {
452 Double_t etet = Et2(v);
453 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
454}
455
457 return TVector2::Phi_mpi_pi(Phi()-v.Phi());
458}
459
461 return PseudoRapidity();
462}
464 Double_t deta = Eta()-v.Eta();
465 Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
466 return TMath::Sqrt( deta*deta+dphi*dphi );
467}
468
470 return DeltaR(v);
471}
472
474 return TVector2 (Eta(),Phi());
475}
476
477
479 return fP.Angle(v);
480}
481
483 return T()*T() - fP.Mag2();
484}
485
487 Double_t mm = Mag2();
488 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
489}
490
491inline Double_t TLorentzVector::M2() const { return Mag2(); }
492inline Double_t TLorentzVector::M() const { return Mag(); }
493
495 return E()*E() - Z()*Z();
496}
497
499 Double_t mm = Mt2();
500 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
501}
502
504 return fP.Mag() / fE;
505}
506
508 Double_t b = Beta();
509 return 1.0/TMath::Sqrt(1- b*b);
510}
511
512inline void TLorentzVector::SetVectMag(const TVector3 & spatial, Double_t magnitude) {
513 SetXYZM(spatial.X(), spatial.Y(), spatial.Z(), magnitude);
514}
515
516inline void TLorentzVector::SetVectM(const TVector3 & spatial, Double_t mass) {
517 SetVectMag(spatial, mass);
518}
519
521 return T()*q.T() - Z()*q.Z() - Y()*q.Y() - X()*q.X();
522}
523
525 return Dot(q);
526}
527
528//Member functions Plus() and Minus() return the positive and negative
529//light-cone components:
530//
531// Double_t pcone = v.Plus();
532// Double_t mcone = v.Minus();
533//
534//CAVEAT: The values returned are T{+,-}Z. It is known that some authors
535//find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus
536//check what definition is used in the physics you're working in and adapt
537//your code accordingly.
538
540 return T() + Z();
541}
542
544 return T() - Z();
545}
546
548 return TVector3(X()/T(), Y()/T(), Z()/T());
549}
550
551inline void TLorentzVector::Boost(const TVector3 & b) {
552 Boost(b.X(), b.Y(), b.Z());
553}
554
556 return fP.PseudoRapidity();
557}
558
560 fP.RotateX(angle);
561}
562
564 fP.RotateY(angle);
565}
566
568 fP.RotateZ(angle);
569}
570
571inline void TLorentzVector::RotateUz(TVector3 &newUzVector) {
572 fP.RotateUz(newUzVector);
573}
574
576 fP.Rotate(a,v);
577}
578
580 fP *= m;
581 return *this;
582}
583
585 fP.Transform(m);
586 return *this;
587}
588
590 return TLorentzVector(a*p.X(), a*p.Y(), a*p.Z(), a*p.T());
591}
592
594 : fP(), fE(0.0) {}
595
597 : fP(x,y,z), fE(t) {}
598
600 : fP(x0), fE(x0[3]) {}
601
603 : fP(x0), fE(x0[3]) {}
604
606 : fP(p), fE(e) {}
607
609 , fP(p.Vect()), fE(p.T()) {}
610
611
612
614{
615 //dereferencing operator const
616 switch(i) {
617 case kX:
618 return fP.X();
619 case kY:
620 return fP.Y();
621 case kZ:
622 return fP.Z();
623 case kT:
624 return fE;
625 default:
626 Error("operator()()", "bad index (%d) returning 0",i);
627 }
628 return 0.;
629}
630
632{
633 //dereferencing operator
634 switch(i) {
635 case kX:
636 return fP.fX;
637 case kY:
638 return fP.fY;
639 case kZ:
640 return fP.fZ;
641 case kT:
642 return fE;
643 default:
644 Error("operator()()", "bad index (%d) returning &fE",i);
645 }
646 return fE;
647}
648
649#endif
SVector< double, 2 > v
Definition: Dict.h:5
ROOT::R::TRInterface & r
Definition: Object.C:4
#define b(i)
Definition: RSha256.hxx:100
#define e(i)
Definition: RSha256.hxx:103
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
#define ClassDef(name, id)
Definition: Rtypes.h:326
float * q
Definition: THbookFile.cxx:87
TLorentzVector operator*(Double_t a, const TLorentzVector &)
double sinh(double)
double sqrt(double)
The TLorentzRotation class describes Lorentz transformations including Lorentz boosts and rotations (...
Double_t Energy() const
void RotateUz(TVector3 &newUzVector)
Double_t Minus() const
void RotateZ(Double_t angle)
void Rotate(Double_t, const TVector3 &)
Bool_t operator!=(const TLorentzVector &) const
Double_t Rapidity() const
TLorentzVector & operator+=(const TLorentzVector &)
Double_t DeltaR(const TLorentzVector &) const
TLorentzVector operator+(const TLorentzVector &) const
TLorentzVector & operator=(const TLorentzVector &)
void RotateX(Double_t angle)
void SetPy(Double_t a)
TLorentzVector & operator*=(Double_t a)
Bool_t operator==(const TLorentzVector &) const
Double_t Angle(const TVector3 &v) const
void GetXYZT(Double_t *carray) const
Double_t Rho() const
Double_t CosTheta() const
Double_t operator[](int i) const
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Double_t Plus() const
void SetY(Double_t a)
void SetPerp(Double_t)
Double_t Beta() const
Double_t Y() const
TVector3 Vect() const
void SetT(Double_t a)
Double_t Px() const
Double_t Mag() const
Double_t Dot(const TLorentzVector &) const
Double_t M() const
Double_t M2() const
Double_t Pt() const
Double_t Perp2() const
TLorentzVector operator*(Double_t a) const
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
void SetE(Double_t a)
Double_t Pz() const
Double_t X() const
Double_t Mt2() const
Double_t Py() const
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
void SetPhi(Double_t phi)
void RotateY(Double_t angle)
void Boost(Double_t, Double_t, Double_t)
Double_t Eta() const
void SetPz(Double_t a)
TLorentzVector operator-() const
Double_t P() const
TVector2 EtaPhiVector()
void SetRho(Double_t rho)
TLorentzVector & Transform(const TRotation &)
virtual void Print(Option_t *option="") const
Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations.
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
virtual ~TLorentzVector()
Double_t Et2() const
Double_t PseudoRapidity() const
Double_t Theta() const
void SetPx(Double_t a)
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Double_t Mag2() const
Double_t Et() const
Double_t Gamma() const
void SetTheta(Double_t theta)
Double_t operator()(int i) const
Double_t Phi() const
Double_t DeltaPhi(const TLorentzVector &) const
Double_t Mt() const
Double_t Perp() const
void SetZ(Double_t a)
Double_t E() const
void SetVect(const TVector3 &vect3)
void SetVectM(const TVector3 &spatial, Double_t mass)
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
TVector3 BoostVector() const
Double_t DrEtaPhi(const TLorentzVector &) const
TLorentzVector & operator-=(const TLorentzVector &)
Double_t Z() const
Double_t T() const
void SetX(Double_t a)
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
The TRotation class describes a rotation of objects of the TVector3 class.
Definition: TRotation.h:20
TVector2 is a general two vector class, which can be used for the description of different vectors in...
Definition: TVector2.h:18
static Double_t Phi_mpi_pi(Double_t x)
Returns phi angle in the interval [-PI,PI)
Definition: TVector2.cxx:101
TVector3 is a general three vector class, which can be used for the description of different vectors ...
Definition: TVector3.h:22
Double_t Z() const
Definition: TVector3.h:218
void SetY(Double_t)
Definition: TVector3.h:224
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:227
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition: TVector3.cxx:362
Double_t fZ
Definition: TVector3.h:185
Double_t fX
Definition: TVector3.h:185
Double_t Phi() const
Return the azimuth angle. Returns phi from -pi to pi.
Definition: TVector3.cxx:230
Double_t Y() const
Definition: TVector3.h:217
Double_t Dot(const TVector3 &) const
Definition: TVector3.h:331
void RotateZ(Double_t)
Rotate vector around Z.
Definition: TVector3.cxx:279
void SetPerp(Double_t)
Definition: TVector3.h:388
Double_t Mag2() const
Definition: TVector3.h:339
Double_t X() const
Definition: TVector3.h:216
void RotateX(Double_t)
Rotate vector around X.
Definition: TVector3.cxx:257
Double_t Angle(const TVector3 &) const
Return the angle w.r.t. another 3-vector.
Definition: TVector3.cxx:197
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition: TVector3.cxx:320
void Rotate(Double_t, const TVector3 &)
Rotate vector.
Definition: TVector3.cxx:290
void GetXYZ(Double_t *carray) const
Definition: TVector3.h:233
Double_t CosTheta() const
Definition: TVector3.h:371
Double_t Mag() const
Definition: TVector3.h:86
Double_t Theta() const
Return the polar angle.
Definition: TVector3.cxx:238
Double_t fY
Definition: TVector3.h:185
void SetMag(Double_t)
Definition: TVector3.h:376
Double_t Perp() const
Return the transverse component (R in cylindrical coordinate system)
Definition: TVector3.cxx:213
TVector3 & Transform(const TRotation &)
Transform this vector with a TRotation.
Definition: TVector3.cxx:190
void SetX(Double_t)
Definition: TVector3.h:223
Double_t Perp2() const
Definition: TVector3.h:353
void RotateY(Double_t)
Rotate vector around Y.
Definition: TVector3.cxx:268
void SetZ(Double_t)
Definition: TVector3.h:225
void RotateUz(const TVector3 &)
NewUzVector must be normalized !
Definition: TVector3.cxx:299
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition: TVector3.cxx:350
TPaveText * pt
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
double T(double x)
Definition: ChebyshevPol.h:34
static constexpr double mm
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
Double_t Sqrt(Double_t x)
Definition: TMath.h:679
Double_t Cos(Double_t)
Definition: TMath.h:629
Double_t Sin(Double_t)
Definition: TMath.h:625
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * m
Definition: textangle.C:8
auto * a
Definition: textangle.C:12