AliPhysics  fb3b612 (fb3b612)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliFlowTrackSimple.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 /* $Id$ */
17 
18 // AliFlowTrackSimple:
19 // A simple track class to the the AliFlowEventSimple for flow analysis
20 //
21 //
22 // author: N. van der Kolk (kolk@nikhef.nl)
23 // mods: Mikolaj Krzewicki (mikolaj.krzewicki@cern.ch)
24 
25 #include "TObject.h"
26 #include "TParticle.h"
27 #include "TParticlePDG.h"
28 #include "AliFlowTrackSimple.h"
29 #include "TRandom.h"
30 #include "TMath.h"
31 
32 ClassImp(AliFlowTrackSimple)
33 
34 //-----------------------------------------------------------------------
36  TObject(),
37  fEta(0),
38  fPt(0),
39  fPhi(0),
40  fTrackWeight(1.),
41  fCharge(0),
42  fMass(-1),
43  fPOItype(0),
44  fSubEventBits(0),
45  fID(-1),
46  fITStype(0)
47 {
48  //constructor
49 }
50 
51 //-----------------------------------------------------------------------
53  TObject(),
54  fEta(eta),
55  fPt(pt),
56  fPhi(phi),
57  fTrackWeight(weight),
58  fCharge(charge),
59  fMass(mass),
60  fPOItype(0),
61  fSubEventBits(0),
62  fID(-1),
63  fITStype(0)
64 {
65  //constructor
66 }
67 
68 //-----------------------------------------------------------------------
70  TObject(),
71  fEta(p->Eta()),
72  fPt(p->Pt()),
73  fPhi(p->Phi()),
74  fTrackWeight(1.),
75  fCharge(0),
76  fMass(-1),
77  fPOItype(0),
78  fSubEventBits(0),
79  fID(-1),
80  fITStype(0)
81 {
82  //ctor
83  TParticlePDG* ppdg = p->GetPDG();
84  fCharge = TMath::Nint(ppdg->Charge()/3.0);
85  fMass = ppdg->Mass();
86 }
87 
88 //-----------------------------------------------------------------------
89 void AliFlowTrackSimple::Set(TParticle* p)
90 {
91  //set from a TParticle
92  fEta = p->Eta();
93  fPt = p->Pt();
94  fPhi = p->Phi();
95  fTrackWeight = 1.;
96  TParticlePDG* ppdg = p->GetPDG();
97  fCharge = TMath::Nint(ppdg->Charge()/3.0);
98  fMass = ppdg->Mass();
99  fITStype = 0;
100 }
101 
102 //-----------------------------------------------------------------------
104  TObject(aTrack),
105  fEta(aTrack.fEta),
106  fPt(aTrack.fPt),
107  fPhi(aTrack.fPhi),
108  fTrackWeight(aTrack.fTrackWeight),
109  fCharge(aTrack.fCharge),
110  fMass(aTrack.fMass),
111  fPOItype(aTrack.fPOItype),
112  fSubEventBits(aTrack.fSubEventBits),
113  fID(aTrack.fID),
114  fITStype(aTrack.fITStype)
115 {
116  //copy constructor
117 }
118 
119 //-----------------------------------------------------------------------
120 AliFlowTrackSimple* AliFlowTrackSimple::Clone(const char* /*option*/) const
121 {
122  //clone "constructor"
123  return new AliFlowTrackSimple(*this);
124 }
125 
126 //-----------------------------------------------------------------------
128 {
129  //assignment
130  if (&aTrack==this) return *this; //handle self assignmnet
131  fEta = aTrack.fEta;
132  fPt = aTrack.fPt;
133  fPhi = aTrack.fPhi;
134  fTrackWeight = aTrack.fTrackWeight;
135  fCharge = aTrack.fCharge;
136  fMass = aTrack.fMass;
137  fPOItype = aTrack.fPOItype;
138  fSubEventBits = aTrack.fSubEventBits;
139  fID = aTrack.fID;
140  fITStype = aTrack.fITStype;
141 
142  return *this;
143 }
144 
145 //-----------------------------------------------------------------------
147 {
148  //destructor
149 
150 }
151 
152 //-----------------------------------------------------------------------
154 {
155  //smear the pt by a gaussian with sigma=res
156  fPt += gRandom->Gaus(0.,res);
157 }
158 
159 //-----------------------------------------------------------------------
161  Double_t reactionPlaneAngle,
162  Double_t precisionPhi,
163  Int_t maxNumberOfIterations )
164 {
165  //afterburner, adds v1, uses Newton-Raphson iteration
166  Double_t phi0=fPhi;
167  Double_t f=0.;
168  Double_t fp=0.;
169  Double_t phiprev=0.;
170 
171  for (Int_t i=0; i<maxNumberOfIterations; i++)
172  {
173  phiprev=fPhi; //store last value for comparison
174  f = fPhi-phi0+2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle);
175  fp = 1.0+2.0*v1*TMath::Cos(fPhi-reactionPlaneAngle); //first derivative
176  fPhi -= f/fp;
177  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
178  }
179 }
180 
181 //-----------------------------------------------------------------------
183  Double_t reactionPlaneAngle,
184  Double_t precisionPhi,
185  Int_t maxNumberOfIterations )
186 {
187  //afterburner, adds v2, uses Newton-Raphson iteration
188  Double_t phi0=fPhi;
189  Double_t f=0.;
190  Double_t fp=0.;
191  Double_t phiprev=0.;
192 
193  for (Int_t i=0; i<maxNumberOfIterations; i++)
194  {
195  phiprev=fPhi; //store last value for comparison
196  f = fPhi-phi0+v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle));
197  fp = 1.0+2.0*v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle)); //first derivative
198  fPhi -= f/fp;
199  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
200  }
201 }
202 
203 //-----------------------------------------------------------------------
205  Double_t reactionPlaneAngle,
206  Double_t precisionPhi,
207  Int_t maxNumberOfIterations )
208 {
209  //afterburner, adds v3, uses Newton-Raphson iteration
210  Double_t phi0=fPhi;
211  Double_t f=0.;
212  Double_t fp=0.;
213  Double_t phiprev=0.;
214 
215  for (Int_t i=0; i<maxNumberOfIterations; i++)
216  {
217  phiprev=fPhi; //store last value for comparison
218  f = fPhi-phi0+2./3.*v3*TMath::Sin(3.*(fPhi-reactionPlaneAngle));
219  fp = 1.0+2.0*v3*TMath::Cos(3.*(fPhi-reactionPlaneAngle)); //first derivative
220  fPhi -= f/fp;
221  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
222  }
223 }
224 
225 //-----------------------------------------------------------------------
227  Double_t reactionPlaneAngle,
228  Double_t precisionPhi,
229  Int_t maxNumberOfIterations )
230 {
231  //afterburner, adds v4, uses Newton-Raphson iteration
232  Double_t phi0=fPhi;
233  Double_t f=0.;
234  Double_t fp=0.;
235  Double_t phiprev=0.;
236 
237  for (Int_t i=0; i<maxNumberOfIterations; i++)
238  {
239  phiprev=fPhi; //store last value for comparison
240  f = fPhi-phi0+0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle));
241  fp = 1.0+2.0*v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle)); //first derivative
242  fPhi -= f/fp;
243  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
244  }
245 }
246 
247 //-----------------------------------------------------------------------
249  Double_t reactionPlaneAngle,
250  Double_t precisionPhi,
251  Int_t maxNumberOfIterations )
252 {
253  //afterburner, adds v4, uses Newton-Raphson iteration
254  Double_t phi0=fPhi;
255  Double_t f=0.;
256  Double_t fp=0.;
257  Double_t phiprev=0.;
258 
259  for (Int_t i=0; i<maxNumberOfIterations; i++)
260  {
261  phiprev=fPhi; //store last value for comparison
262  f = fPhi-phi0+0.4*v5*TMath::Sin(5.*(fPhi-reactionPlaneAngle));
263  fp = 1.0+2.0*v5*TMath::Cos(5.*(fPhi-reactionPlaneAngle)); //first derivative
264  fPhi -= f/fp;
265  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
266  }
267 }
268 
269 //______________________________________________________________________________
271  Double_t v2,
272  Double_t v3,
273  Double_t v4,
274  Double_t v5,
275  Double_t rp1,
276  Double_t rp2,
277  Double_t rp3,
278  Double_t rp4,
279  Double_t rp5,
280  Double_t precisionPhi,
281  Int_t maxNumberOfIterations )
282 {
283  //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
284  Double_t phi0=fPhi;
285  Double_t f=0.;
286  Double_t fp=0.;
287  Double_t phiprev=0.;
288 
289  for (Int_t i=0; i<maxNumberOfIterations; i++)
290  {
291  phiprev=fPhi; //store last value for comparison
292  f = fPhi-phi0
293  +2.0* v1*TMath::Sin( fPhi-rp1)
294  + v2*TMath::Sin(2.*(fPhi-rp2))
295  +2./3.*v3*TMath::Sin(3.*(fPhi-rp3))
296  +0.5* v4*TMath::Sin(4.*(fPhi-rp4))
297  +0.4* v5*TMath::Sin(5.*(fPhi-rp5))
298  ;
299  fp = 1.0
300  +2.0*(
301  +v1*TMath::Cos( fPhi-rp1)
302  +v2*TMath::Cos(2.*(fPhi-rp2))
303  +v3*TMath::Cos(3.*(fPhi-rp3))
304  +v4*TMath::Cos(4.*(fPhi-rp4))
305  +v5*TMath::Cos(5.*(fPhi-rp5))
306  ); //first derivative
307  fPhi -= f/fp;
308  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
309  }
310 }
311 
312 //______________________________________________________________________________
314  Double_t v2,
315  Double_t v3,
316  Double_t v4,
317  Double_t v5,
318  Double_t rp,
319  Double_t precisionPhi,
320  Int_t maxNumberOfIterations )
321 {
322  //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
323  AddFlow(v1,v2,v3,v4,v5,rp,rp,rp,rp,rp,precisionPhi,maxNumberOfIterations);
324 }
325 
326 //______________________________________________________________________________
328  Double_t v2,
329  Double_t v3,
330  Double_t v4,
331  Double_t reactionPlaneAngle,
332  Double_t precisionPhi,
333  Int_t maxNumberOfIterations )
334 {
335  //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
336  Double_t phi0=fPhi;
337  Double_t f=0.;
338  Double_t fp=0.;
339  Double_t phiprev=0.;
340 
341  for (Int_t i=0; i<maxNumberOfIterations; i++)
342  {
343  phiprev=fPhi; //store last value for comparison
344  f = fPhi-phi0
345  +2.0* v1*TMath::Sin( fPhi-reactionPlaneAngle)
346  + v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle))
347  +2./3.*v3*TMath::Sin(3.*(fPhi-reactionPlaneAngle))
348  +0.5* v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle))
349  ;
350  fp = 1.0
351  +2.0*(
352  +v1*TMath::Cos( fPhi-reactionPlaneAngle)
353  +v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle))
354  +v3*TMath::Cos(3.*(fPhi-reactionPlaneAngle))
355  +v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle))
356  ); //first derivative
357  fPhi -= f/fp;
358  if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
359  }
360 }
361 
362 //______________________________________________________________________________
363 void AliFlowTrackSimple::Print( Option_t* /*option*/ ) const
364 {
365  //print stuff
366  printf("Phi: %.3f, Eta: %+.3f, Pt: %.3f, weight: %.3f",fPhi,fEta,fPt,fTrackWeight);
367  if (InRPSelection()) printf(", RP");
368  if (InPOISelection()) printf(", POI");
369  for (Int_t i=0; i<2; i++)
370  {
371  if (InSubevent(i)) printf(", subevent %i",i);
372  }
373  printf("\n");
374 }
375 
376 //______________________________________________________________________________
378 {
379  //clear track
380  fEta=0.0;
381  fPt=0.0;
382  fPhi=0.0;
383  fTrackWeight=1.0;
384  fCharge=0;
385  fMass=-1;
386  fPOItype.ResetAllBits();
387  fSubEventBits.ResetAllBits();
388  fID=-1;
389  fITStype=0;
390 }
Int_t charge
double Double_t
Definition: External.C:58
virtual void Print(Option_t *option="") const
void AddV1(Double_t v1, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
Bool_t InSubevent(Int_t i) const
void AddV5(Double_t v5, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
Double_t mass
virtual AliFlowTrackSimple * Clone(const char *option="") const
void Set(TParticle *p)
TRandom * gRandom
Bool_t InRPSelection() const
int Int_t
Definition: External.C:63
void AddV4(Double_t v4, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
AliFlowTrackSimple & operator=(const AliFlowTrackSimple &aTrack)
void ResolutionPt(Double_t resolution)
void AddV2(Double_t v2, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
void AddFlow(Double_t v1, Double_t v2, Double_t v3, Double_t v4, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
void AddV3(Double_t v3, Double_t reactionPlaneAngle, Double_t precision, Int_t maxNumberOfIterations=100)
const char Option_t
Definition: External.C:48
Bool_t InPOISelection(Int_t poiType=1) const
virtual void Clear(Option_t *o="")