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 "AliFlowTrackSimple.h"
30 ClassImp(AliFlowTrackSimple)
32 //-----------------------------------------------------------------------
33 AliFlowTrackSimple::AliFlowTrackSimple():
45 //-----------------------------------------------------------------------
46 AliFlowTrackSimple::AliFlowTrackSimple(Double_t phi, Double_t eta, Double_t pt, Double_t weight):
58 //-----------------------------------------------------------------------
59 AliFlowTrackSimple::AliFlowTrackSimple(const TParticle* p):
71 //-----------------------------------------------------------------------
72 AliFlowTrackSimple::AliFlowTrackSimple(const AliFlowTrackSimple& aTrack):
77 fTrackWeight(aTrack.fTrackWeight),
78 fFlowBits(aTrack.fFlowBits),
79 fSubEventBits(aTrack.fSubEventBits)
84 //-----------------------------------------------------------------------
85 AliFlowTrackSimple* AliFlowTrackSimple::Clone(const char* /*option*/) const
88 return new AliFlowTrackSimple(*this);
91 //-----------------------------------------------------------------------
92 AliFlowTrackSimple& AliFlowTrackSimple::operator=(const AliFlowTrackSimple& aTrack)
97 fTrackWeight = aTrack.fTrackWeight;
98 fFlowBits = aTrack.fFlowBits;
99 fSubEventBits = aTrack.fSubEventBits;
104 //-----------------------------------------------------------------------
105 AliFlowTrackSimple::~AliFlowTrackSimple()
111 //-----------------------------------------------------------------------
112 void AliFlowTrackSimple::ResolutionPt(Double_t res)
114 //smear the pt by a gaussian with sigma=res
115 fPt += gRandom->Gaus(0.,res);
118 //-----------------------------------------------------------------------
119 void AliFlowTrackSimple::AddV1( Double_t v1,
120 Double_t reactionPlaneAngle,
121 Double_t precisionPhi,
122 Int_t maxNumberOfIterations )
124 //afterburner, adds v1, uses Newton-Raphson iteration
126 Double_t f,fp,phiprev;
128 for (Int_t i=0; i<maxNumberOfIterations; i++)
130 phiprev=fPhi; //store last value for comparison
131 f = fPhi-phi0+2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle);
132 fp = 1.0+2.0*v1*TMath::Cos(fPhi-reactionPlaneAngle); //first derivative
134 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
138 //-----------------------------------------------------------------------
139 void AliFlowTrackSimple::AddV2( Double_t v2,
140 Double_t reactionPlaneAngle,
141 Double_t precisionPhi,
142 Int_t maxNumberOfIterations )
144 //afterburner, adds v2, uses Newton-Raphson iteration
146 Double_t f,fp,phiprev;
148 for (Int_t i=0; i<maxNumberOfIterations; i++)
150 phiprev=fPhi; //store last value for comparison
151 f = fPhi-phi0+v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle));
152 fp = 1.0+2.0*v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle)); //first derivative
154 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
158 //-----------------------------------------------------------------------
159 void AliFlowTrackSimple::AddV4( Double_t v4,
160 Double_t reactionPlaneAngle,
161 Double_t precisionPhi,
162 Int_t maxNumberOfIterations )
164 //afterburner, adds v4, uses Newton-Raphson iteration
166 Double_t f,fp,phiprev;
168 for (Int_t i=0; i<maxNumberOfIterations; i++)
170 phiprev=fPhi; //store last value for comparison
171 f = fPhi-phi0+0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle));
172 fp = 1.0+2.0*v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle)); //first derivative
174 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
178 //______________________________________________________________________________
179 void AliFlowTrackSimple::AddFlow( Double_t v1,
182 Double_t reactionPlaneAngle,
183 Double_t precisionPhi,
184 Int_t maxNumberOfIterations )
186 //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
188 Double_t f,fp,phiprev;
190 for (Int_t i=0; i<maxNumberOfIterations; i++)
192 phiprev=fPhi; //store last value for comparison
194 +2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle)
195 + v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle))
196 +0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle))
200 +v1*TMath::Cos(fPhi-reactionPlaneAngle)
201 +v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle))
202 +v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle))
203 ); //first derivative
205 if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
209 //______________________________________________________________________________
210 void AliFlowTrackSimple::Print( Option_t* /*option*/ ) const
213 printf("Phi: %.3f, Eta: %+.3f, Pt: %.3f, weight: %.3f",fPhi,fEta,fPt,fTrackWeight);
214 if (InRPSelection()) printf(", RP");
215 if (InPOISelection()) printf(", POI");
216 for (Int_t i=0; i<2; i++)
218 if (InSubevent(i)) printf(", subevent %i",i);