]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliFlowTrackSimple.cxx
f8012d8da3e4526527cedd22cbb002a80792ea31
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowTrackSimple.cxx
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 "AliFlowTrackSimple.h"
28 #include "TRandom.h"
29
30 ClassImp(AliFlowTrackSimple)
31
32 //-----------------------------------------------------------------------
33 AliFlowTrackSimple::AliFlowTrackSimple():
34   TObject(),
35   fEta(0),
36   fPt(0),
37   fPhi(0),
38   fTrackWeight(1.),
39   fFlowBits(0),
40   fSubEventBits(0)
41 {
42   //constructor 
43 }
44
45 //-----------------------------------------------------------------------
46 AliFlowTrackSimple::AliFlowTrackSimple(Double_t phi, Double_t eta, Double_t pt, Double_t weight):
47   TObject(),
48   fEta(eta),
49   fPt(pt),
50   fPhi(phi),
51   fTrackWeight(weight),
52   fFlowBits(0),
53   fSubEventBits(0)
54 {
55   //constructor 
56 }
57
58 //-----------------------------------------------------------------------
59 AliFlowTrackSimple::AliFlowTrackSimple(const TParticle* p):
60   TObject(),
61   fEta(p->Eta()),
62   fPt(p->Pt()),
63   fPhi(p->Phi()),
64   fTrackWeight(1.),
65   fFlowBits(0),
66   fSubEventBits(0)
67 {
68   //ctor
69 }
70
71 //-----------------------------------------------------------------------
72 AliFlowTrackSimple::AliFlowTrackSimple(const AliFlowTrackSimple& aTrack):
73   TObject(aTrack),
74   fEta(aTrack.fEta),
75   fPt(aTrack.fPt),
76   fPhi(aTrack.fPhi),
77   fTrackWeight(aTrack.fTrackWeight),
78   fFlowBits(aTrack.fFlowBits),
79   fSubEventBits(aTrack.fSubEventBits)
80 {
81   //copy constructor 
82 }
83
84 //-----------------------------------------------------------------------
85 AliFlowTrackSimple* AliFlowTrackSimple::Clone(const char* /*option*/) const
86 {
87   //clone "constructor"
88   return new AliFlowTrackSimple(*this);
89 }
90
91 //-----------------------------------------------------------------------
92 AliFlowTrackSimple& AliFlowTrackSimple::operator=(const AliFlowTrackSimple& aTrack)
93 {
94   fEta = aTrack.fEta;
95   fPt = aTrack.fPt;
96   fPhi = aTrack.fPhi;
97   fTrackWeight = aTrack.fTrackWeight;
98   fFlowBits = aTrack.fFlowBits;
99   fSubEventBits = aTrack.fSubEventBits;
100
101   return *this;
102 }
103
104 //----------------------------------------------------------------------- 
105 AliFlowTrackSimple::~AliFlowTrackSimple()
106 {
107   //destructor
108   
109 }
110
111 //----------------------------------------------------------------------- 
112 void AliFlowTrackSimple::ResolutionPt(Double_t res)
113 {
114   //smear the pt by a gaussian with sigma=res
115   fPt += gRandom->Gaus(0.,res);
116 }
117
118 //----------------------------------------------------------------------- 
119 void AliFlowTrackSimple::AddV1( Double_t v1,
120                                 Double_t reactionPlaneAngle,
121                                 Double_t precisionPhi,
122                                 Int_t maxNumberOfIterations )
123 {
124   //afterburner, adds v1, uses Newton-Raphson iteration
125   Double_t phi0=fPhi;
126   Double_t f,fp,phiprev;
127
128   for (Int_t i=0; i<maxNumberOfIterations; i++)
129   {
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
133     fPhi -= f/fp;
134     if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
135   }
136 }
137
138 //----------------------------------------------------------------------- 
139 void AliFlowTrackSimple::AddV2( Double_t v2,
140                                 Double_t reactionPlaneAngle,
141                                 Double_t precisionPhi,
142                                 Int_t maxNumberOfIterations )
143 {
144   //afterburner, adds v2, uses Newton-Raphson iteration
145   Double_t phi0=fPhi;
146   Double_t f,fp,phiprev;
147
148   for (Int_t i=0; i<maxNumberOfIterations; i++)
149   {
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
153     fPhi -= f/fp;
154     if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
155   }
156 }
157
158 //----------------------------------------------------------------------- 
159 void AliFlowTrackSimple::AddV4( Double_t v4,
160                                 Double_t reactionPlaneAngle,
161                                 Double_t precisionPhi,
162                                 Int_t maxNumberOfIterations )
163 {
164   //afterburner, adds v4, uses Newton-Raphson iteration
165   Double_t phi0=fPhi;
166   Double_t f,fp,phiprev;
167
168   for (Int_t i=0; i<maxNumberOfIterations; i++)
169   {
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
173     fPhi -= f/fp;
174     if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
175   }
176 }
177
178 //______________________________________________________________________________
179 void AliFlowTrackSimple::AddFlow( Double_t v1,
180                                   Double_t v2,
181                                   Double_t v4,
182                                   Double_t reactionPlaneAngle,
183                                   Double_t precisionPhi,
184                                   Int_t maxNumberOfIterations )
185 {
186   //afterburner, adds v1,v2,v4 uses Newton-Raphson iteration
187   Double_t phi0=fPhi;
188   Double_t f,fp,phiprev;
189
190   for (Int_t i=0; i<maxNumberOfIterations; i++)
191   {
192     phiprev=fPhi; //store last value for comparison
193     f =  fPhi-phi0
194         +2.0*v1*TMath::Sin(fPhi-reactionPlaneAngle)
195         +    v2*TMath::Sin(2.*(fPhi-reactionPlaneAngle))
196         +0.5*v4*TMath::Sin(4.*(fPhi-reactionPlaneAngle))
197         ;
198     fp =  1.0
199          +2.0*(
200            +v1*TMath::Cos(fPhi-reactionPlaneAngle)
201            +v2*TMath::Cos(2.*(fPhi-reactionPlaneAngle))
202            +v4*TMath::Cos(4.*(fPhi-reactionPlaneAngle))
203          ); //first derivative
204     fPhi -= f/fp;
205     if (TMath::AreEqualAbs(phiprev,fPhi,precisionPhi)) break;
206   }
207 }
208
209 //______________________________________________________________________________
210 void AliFlowTrackSimple::Print( Option_t* /*option*/ ) const
211 {
212   //print stuff
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++)
217   {
218     if (InSubevent(i)) printf(", subevent %i",i);
219   }
220   printf("\n");
221 }