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