1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // AliFlowTrackSimple:
19 // A simple track class to the the AliFlowEventSimple for flow analysis
22 // author: N. van der Kolk (kolk@nikhef.nl)
23 // mods: Mikolaj Krzewicki (mikolaj.krzewicki@cern.ch)
26 #include "TParticle.h"
27 #include "TParticlePDG.h"
28 #include "AliFlowTrackSimple.h"
32 ClassImp(AliFlowTrackSimple)
34 //-----------------------------------------------------------------------
35 AliFlowTrackSimple::AliFlowTrackSimple():
49 //-----------------------------------------------------------------------
50 AliFlowTrackSimple::AliFlowTrackSimple(Double_t phi, Double_t eta, Double_t pt, Double_t weight, Int_t charge):
64 //-----------------------------------------------------------------------
65 AliFlowTrackSimple::AliFlowTrackSimple( TParticle* p ):
77 TParticlePDG* ppdg = p->GetPDG();
78 fCharge = TMath::Nint(ppdg->Charge()/3.0);
81 //-----------------------------------------------------------------------
82 void AliFlowTrackSimple::Set(TParticle* p)
84 //set from a TParticle
89 TParticlePDG* ppdg = p->GetPDG();
90 fCharge = TMath::Nint(ppdg->Charge()/3.0);
93 //-----------------------------------------------------------------------
94 AliFlowTrackSimple::AliFlowTrackSimple(const AliFlowTrackSimple& aTrack):
99 fTrackWeight(aTrack.fTrackWeight),
100 fCharge(aTrack.fCharge),
101 fFlowBits(aTrack.fFlowBits),
102 fSubEventBits(aTrack.fSubEventBits),
108 //-----------------------------------------------------------------------
109 AliFlowTrackSimple* AliFlowTrackSimple::Clone(const char* /*option*/) const
111 //clone "constructor"
112 return new AliFlowTrackSimple(*this);
115 //-----------------------------------------------------------------------
116 AliFlowTrackSimple& AliFlowTrackSimple::operator=(const AliFlowTrackSimple& aTrack)
119 if (&aTrack==this) return *this; //handle self assignmnet
123 fTrackWeight = aTrack.fTrackWeight;
124 fCharge = aTrack.fCharge;
125 fFlowBits = aTrack.fFlowBits;
126 fSubEventBits = aTrack.fSubEventBits;
132 //-----------------------------------------------------------------------
133 AliFlowTrackSimple::~AliFlowTrackSimple()
139 //-----------------------------------------------------------------------
140 void AliFlowTrackSimple::ResolutionPt(Double_t res)
142 //smear the pt by a gaussian with sigma=res
143 fPt += gRandom->Gaus(0.,res);
146 //-----------------------------------------------------------------------
147 void AliFlowTrackSimple::AddV1( Double_t v1,
148 Double_t reactionPlaneAngle,
149 Double_t precisionPhi,
150 Int_t maxNumberOfIterations )
152 //afterburner, adds v1, uses Newton-Raphson iteration
158 for (Int_t i=0; i<maxNumberOfIterations; i++)
160 phiprev=fPhi; //store last value for comparison
161 f = fPhi-phi0+2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle);
162 fp = 1.0+2.0*v1*TMath::Cos(fPhi-reactionPlaneAngle); //first derivative
164 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
168 //-----------------------------------------------------------------------
169 void AliFlowTrackSimple::AddV2( Double_t v2,
170 Double_t reactionPlaneAngle,
171 Double_t precisionPhi,
172 Int_t maxNumberOfIterations )
174 //afterburner, adds v2, uses Newton-Raphson iteration
180 for (Int_t i=0; i<maxNumberOfIterations; i++)
182 phiprev=fPhi; //store last value for comparison
183 f = fPhi-phi0+v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle));
184 fp = 1.0+2.0*v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle)); //first derivative
186 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
190 //-----------------------------------------------------------------------
191 void AliFlowTrackSimple::AddV3( Double_t v3,
192 Double_t reactionPlaneAngle,
193 Double_t precisionPhi,
194 Int_t maxNumberOfIterations )
196 //afterburner, adds v3, uses Newton-Raphson iteration
202 for (Int_t i=0; i<maxNumberOfIterations; i++)
204 phiprev=fPhi; //store last value for comparison
205 f = fPhi-phi0+2./3.*v3*TMath::Sin(3.*(fPhi-reactionPlaneAngle));
206 fp = 1.0+2.0*v3*TMath::Cos(3.*(fPhi-reactionPlaneAngle)); //first derivative
208 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
212 //-----------------------------------------------------------------------
213 void AliFlowTrackSimple::AddV4( Double_t v4,
214 Double_t reactionPlaneAngle,
215 Double_t precisionPhi,
216 Int_t maxNumberOfIterations )
218 //afterburner, adds v4, uses Newton-Raphson iteration
224 for (Int_t i=0; i<maxNumberOfIterations; i++)
226 phiprev=fPhi; //store last value for comparison
227 f = fPhi-phi0+0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle));
228 fp = 1.0+2.0*v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle)); //first derivative
230 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
234 //-----------------------------------------------------------------------
235 void AliFlowTrackSimple::AddV5( Double_t v5,
236 Double_t reactionPlaneAngle,
237 Double_t precisionPhi,
238 Int_t maxNumberOfIterations )
240 //afterburner, adds v4, uses Newton-Raphson iteration
246 for (Int_t i=0; i<maxNumberOfIterations; i++)
248 phiprev=fPhi; //store last value for comparison
249 f = fPhi-phi0+0.4*v5*TMath::Sin(5.*(fPhi-reactionPlaneAngle));
250 fp = 1.0+2.0*v5*TMath::Cos(5.*(fPhi-reactionPlaneAngle)); //first derivative
252 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
256 //______________________________________________________________________________
257 void AliFlowTrackSimple::AddFlow( Double_t v1,
267 Double_t precisionPhi,
268 Int_t maxNumberOfIterations )
270 //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
276 for (Int_t i=0; i<maxNumberOfIterations; i++)
278 phiprev=fPhi; //store last value for comparison
280 +2.0* v1*TMath::Sin( fPhi-rp1)
281 + v2*TMath::Sin(2.*(fPhi-rp2))
282 +2./3.*v3*TMath::Sin(3.*(fPhi-rp3))
283 +0.5* v4*TMath::Sin(4.*(fPhi-rp4))
284 +0.4* v5*TMath::Sin(5.*(fPhi-rp5))
288 +v1*TMath::Cos( fPhi-rp1)
289 +v2*TMath::Cos(2.*(fPhi-rp2))
290 +v3*TMath::Cos(3.*(fPhi-rp3))
291 +v4*TMath::Cos(4.*(fPhi-rp4))
292 +v5*TMath::Cos(5.*(fPhi-rp5))
293 ); //first derivative
295 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
299 //______________________________________________________________________________
300 void AliFlowTrackSimple::AddFlow( Double_t v1,
306 Double_t precisionPhi,
307 Int_t maxNumberOfIterations )
309 //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
310 AddFlow(v1,v2,v3,v4,v5,rp,rp,rp,rp,rp,precisionPhi,maxNumberOfIterations);
313 //______________________________________________________________________________
314 void AliFlowTrackSimple::AddFlow( Double_t v1,
318 Double_t reactionPlaneAngle,
319 Double_t precisionPhi,
320 Int_t maxNumberOfIterations )
322 //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
328 for (Int_t i=0; i<maxNumberOfIterations; i++)
330 phiprev=fPhi; //store last value for comparison
332 +2.0* v1*TMath::Sin( fPhi-reactionPlaneAngle)
333 + v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle))
334 +2./3.*v3*TMath::Sin(3.*(fPhi-reactionPlaneAngle))
335 +0.5* v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle))
339 +v1*TMath::Cos( fPhi-reactionPlaneAngle)
340 +v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle))
341 +v3*TMath::Cos(3.*(fPhi-reactionPlaneAngle))
342 +v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle))
343 ); //first derivative
345 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
349 //______________________________________________________________________________
350 void AliFlowTrackSimple::Print( Option_t* /*option*/ ) const
353 printf("Phi: %.3f, Eta: %+.3f, Pt: %.3f, weight: %.3f",fPhi,fEta,fPt,fTrackWeight);
354 if (InRPSelection()) printf(", RP");
355 if (InPOISelection()) printf(", POI");
356 for (Int_t i=0; i<2; i++)
358 if (InSubevent(i)) printf(", subevent %i",i);
363 //______________________________________________________________________________
364 void AliFlowTrackSimple::Clear(Option_t*)
372 fFlowBits.ResetAllBits();
373 fSubEventBits.ResetAllBits();