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():
48 //-----------------------------------------------------------------------
49 AliFlowTrackSimple::AliFlowTrackSimple(Double_t phi, Double_t eta, Double_t pt, Double_t weight, Int_t charge):
62 //-----------------------------------------------------------------------
63 AliFlowTrackSimple::AliFlowTrackSimple( TParticle* p ):
74 TParticlePDG* ppdg = p->GetPDG();
75 fCharge = TMath::Nint(ppdg->Charge()/3.0);
78 //-----------------------------------------------------------------------
79 AliFlowTrackSimple::AliFlowTrackSimple(const AliFlowTrackSimple& aTrack):
84 fTrackWeight(aTrack.fTrackWeight),
85 fCharge(aTrack.fCharge),
86 fFlowBits(aTrack.fFlowBits),
87 fSubEventBits(aTrack.fSubEventBits)
92 //-----------------------------------------------------------------------
93 AliFlowTrackSimple* AliFlowTrackSimple::Clone(const char* /*option*/) const
96 return new AliFlowTrackSimple(*this);
99 //-----------------------------------------------------------------------
100 AliFlowTrackSimple& AliFlowTrackSimple::operator=(const AliFlowTrackSimple& aTrack)
105 fTrackWeight = aTrack.fTrackWeight;
106 fCharge = aTrack.fCharge;
107 fFlowBits = aTrack.fFlowBits;
108 fSubEventBits = aTrack.fSubEventBits;
113 //-----------------------------------------------------------------------
114 AliFlowTrackSimple::~AliFlowTrackSimple()
120 //-----------------------------------------------------------------------
121 void AliFlowTrackSimple::ResolutionPt(Double_t res)
123 //smear the pt by a gaussian with sigma=res
124 fPt += gRandom->Gaus(0.,res);
127 //-----------------------------------------------------------------------
128 void AliFlowTrackSimple::AddV1( Double_t v1,
129 Double_t reactionPlaneAngle,
130 Double_t precisionPhi,
131 Int_t maxNumberOfIterations )
133 //afterburner, adds v1, uses Newton-Raphson iteration
139 for (Int_t i=0; i<maxNumberOfIterations; i++)
141 phiprev=fPhi; //store last value for comparison
142 f = fPhi-phi0+2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle);
143 fp = 1.0+2.0*v1*TMath::Cos(fPhi-reactionPlaneAngle); //first derivative
145 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
149 //-----------------------------------------------------------------------
150 void AliFlowTrackSimple::AddV2( Double_t v2,
151 Double_t reactionPlaneAngle,
152 Double_t precisionPhi,
153 Int_t maxNumberOfIterations )
155 //afterburner, adds v2, uses Newton-Raphson iteration
161 for (Int_t i=0; i<maxNumberOfIterations; i++)
163 phiprev=fPhi; //store last value for comparison
164 f = fPhi-phi0+v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle));
165 fp = 1.0+2.0*v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle)); //first derivative
167 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
171 //-----------------------------------------------------------------------
172 void AliFlowTrackSimple::AddV4( Double_t v4,
173 Double_t reactionPlaneAngle,
174 Double_t precisionPhi,
175 Int_t maxNumberOfIterations )
177 //afterburner, adds v4, uses Newton-Raphson iteration
183 for (Int_t i=0; i<maxNumberOfIterations; i++)
185 phiprev=fPhi; //store last value for comparison
186 f = fPhi-phi0+0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle));
187 fp = 1.0+2.0*v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle)); //first derivative
189 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
193 //______________________________________________________________________________
194 void AliFlowTrackSimple::AddFlow( Double_t v1,
197 Double_t reactionPlaneAngle,
198 Double_t precisionPhi,
199 Int_t maxNumberOfIterations )
201 //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
207 for (Int_t i=0; i<maxNumberOfIterations; i++)
209 phiprev=fPhi; //store last value for comparison
211 +2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle)
212 + v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle))
213 +0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle))
217 +v1*TMath::Cos(fPhi-reactionPlaneAngle)
218 +v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle))
219 +v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle))
220 ); //first derivative
222 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
226 //______________________________________________________________________________
227 void AliFlowTrackSimple::Print( Option_t* /*option*/ ) const
230 printf("Phi: %.3f, Eta: %+.3f, Pt: %.3f, weight: %.3f",fPhi,fEta,fPt,fTrackWeight);
231 if (InRPSelection()) printf(", RP");
232 if (InPOISelection()) printf(", POI");
233 for (Int_t i=0; i<2; i++)
235 if (InSubevent(i)) printf(", subevent %i",i);