]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCTrackHits.cxx
0d0aeddd55d43e45486265b394a0c97597f885c1
[u/mrichter/AliRoot.git] / TPC / AliTPCTrackHits.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 /*
17 $Log$
18 */
19 ///////////////////////////////////////////////////////////////////////////////
20 //                                                                           //
21 //  Time Projection Chamber  track hits object                                //
22 //
23 //  Origin: Marian Ivanov , GSI Darmstadt
24 //
25 //  Class for storing simulated AliTPCHits  for given track                  //
26 //             -average compression comparing to classical ClonesArray is    //
27 //              around 5-7 (depending on the required hit precision)       //
28 //                                                                           //
29 //Begin_Html
30 /*
31 <img src="gif/AliTPCTrackHits.gif">
32 */
33 //End_Html
34 //                                                                           //
35 //                                                                          //
36 ///////////////////////////////////////////////////////////////////////////////
37
38 #include "TVector3.h"
39 #include "TClonesArray.h"    
40 #include "AliTPCTrackHits.h"
41 #include "AliTPC.h"
42
43 #include <iostream.h>
44
45
46
47 ClassImp(AliTPCTrackHits) 
48 LClassImp(AliTrackHitsInfo) 
49 LClassImp(AliTrackHitsParam)  
50 LClassImp(AliHitInfo)
51
52 Int_t AliTrackHitsInfo::fgCounter1 =0;
53 Int_t AliTrackHitsInfo::fgCounter2 =0;
54 Int_t AliTrackHitsParam::fgCounter1 =0;
55 Int_t AliTrackHitsParam::fgCounter2 =0;
56 Int_t AliHitInfo::fgCounter1 =0;
57 Int_t AliHitInfo::fgCounter2 =0;
58 Int_t AliTPCTrackHits::fgCounter1 =0;
59 Int_t AliTPCTrackHits::fgCounter2 =0;
60 const Double_t AliTPCTrackHits::fgkPrecision=1e-6;  //precision 
61 const Double_t AliTPCTrackHits::fgkPrecision2=1e-20;  //precision
62
63
64 /************************************************************/
65 //           Interface classes                              // 
66 #include "AliTPCTrackHitsInterfaces.h"
67
68
69
70
71 struct AliTPCCurrentHit {
72   AliTPChit fHit;
73   UInt_t   fInfoIndex;//   - current info pointer 
74   UInt_t   fParamIndex;//  - current param pointer
75   UInt_t   fStackIndex; // - current hit stack index
76   Double_t fR;   //current Radius
77   Bool_t  fStatus; //current status    
78 };   
79
80
81 struct  AliTPCTempHitInfo {
82   enum    { fkStackSize = 10000};
83   AliTPCTempHitInfo();   
84   void     NewParam(Double_t r, Double_t z, Double_t fi, Int_t q);
85   void     SetHit(Double_t r, Double_t z, Double_t fi, Int_t q);
86   Double_t * GetPosition(Int_t index){return &fPositionStack[index*3];}
87   void    UpdateParam(Double_t maxdelta); //recal
88   void   Fit2(Double_t fSumY, Double_t fSumYX, Double_t fSumYX2,
89             Double_t fSumX,  Double_t fSumX2, Double_t fSumX3, 
90             Double_t fSumX4, Int_t n,
91               Double_t &a, Double_t &b, Double_t &c);
92   void  Fit(AliTrackHitsParam * param);
93   Double_t fSumDr;    //
94   Double_t fSumDr2;   //
95   Double_t fSumDr3;   //  
96   Double_t fSumDr4;   //
97   Double_t fSumDFi;  //
98   Double_t fSumDFiDr; //  
99   Double_t fSumDFiDr2;//
100   Double_t fSumDZ;     //
101   Double_t fSumDZDr;  //
102   Double_t fSumDZDr2;  //
103   Double_t fOldR;     //previos r
104   Double_t fPositionStack[3*fkStackSize];  //position stack 
105   UInt_t   fQStack[fkStackSize];           //Q stack
106   UInt_t fStackIndex;   //current stack index 
107   UInt_t fInfoIndex;    //current track info index
108   UInt_t fParamIndex;   //current track parameters index
109   AliTrackHitsInfo  * fInfo; //current track info
110   AliTrackHitsParam * fParam; //current track param
111 };
112
113
114 AliTPCTempHitInfo::AliTPCTempHitInfo()
115 {
116   //
117   //set to default value
118   fSumDr=fSumDr2=fSumDr3=fSumDr4=
119     fSumDFi=fSumDFiDr=fSumDFiDr2=
120     fSumDZ=fSumDZDr=fSumDZDr2=0;  
121   fStackIndex = 0;
122   fInfoIndex  = 0;
123   fParamIndex = 0;
124 }
125
126
127 void AliTPCTempHitInfo::NewParam(Double_t r, Double_t z, Double_t fi, Int_t q)
128 {
129   //
130   //reset stack and sum parameters
131   //store line initial point
132   fSumDr=fSumDr2=fSumDr3=fSumDr4=
133     fSumDFi=fSumDFiDr=fSumDFiDr2=
134     fSumDZ=fSumDZDr=fSumDZDr2=0;  
135   fStackIndex=0;
136   fParam->fR = r;
137   fOldR = r;
138   fParam->fZ = z;
139   fParam->fFi = fi;
140   fParam->fAn = 0.;
141   fParam->fAd = 0.;
142   fParam->fTheta =0.;
143   fParam->fThetaD =0.;
144   SetHit(r,z,fi,q);
145 }
146
147 void AliTPCTempHitInfo::SetHit(Double_t r, Double_t z, Double_t fi, Int_t q)
148 {
149   //
150   //add hit to the stack
151   //recalculate new estimete of line parameters
152   Double_t *f = GetPosition(fStackIndex);  
153   f[0] = r;
154   f[1] = z;
155   f[2] = fi;
156   fQStack[fStackIndex]=q;
157   if (fStackIndex==0) return;
158   Double_t dr  = (r-fParam->fR);
159   if (TMath::Abs(dr)<AliTPCTrackHits::fgkPrecision) dr =AliTPCTrackHits::fgkPrecision;
160   Double_t dfi = fi-fParam->fFi;
161   Double_t dz  = z -fParam->fZ; 
162   Double_t dr2 =dr*dr;
163   Double_t dr3 =dr2*dr;
164   Double_t dr4 =dr3*dr;
165   fSumDr +=dr;
166   fSumDr2+=dr2;
167   fSumDr3+=dr3;
168   fSumDr4+=dr4;
169   fSumDFi +=dfi;
170   fSumDFiDr+=dfi*dr;
171   fSumDFiDr2+=dfi*dr2;
172   fSumDZ +=dz;
173   fSumDZDr+=dz*dr;
174   fSumDZDr2+=dz*dr2;
175   
176   //update fit parameters
177   //
178   Double_t det = fSumDr2*fSumDr4-fSumDr3*fSumDr3;
179   if (TMath::Abs(det)<AliTPCTrackHits::fgkPrecision2) return;
180   if ( ( fStackIndex>1 )  ){
181     fParam->fAn = (fSumDr4*fSumDFiDr-fSumDr3*fSumDFiDr2)/det;
182     fParam->fAd = (fSumDr2*fSumDFiDr2-fSumDr3*fSumDFiDr)/det;
183   }
184   else
185     fParam->fAn = fSumDFiDr/fSumDr2;
186   if ( ( fStackIndex>1 )  ){
187     fParam->fTheta = (fSumDr4*fSumDZDr-fSumDr3*fSumDZDr2)/det;
188     fParam->fThetaD= (fSumDr2*fSumDZDr2-fSumDr3*fSumDZDr)/det;
189   }
190   else
191     fParam->fTheta = fSumDZDr/fSumDr2; 
192 }
193
194
195 void   AliTPCTempHitInfo::UpdateParam(Double_t maxdelta)
196 {
197   //recalc parameters not fixing origin point
198   if (fStackIndex>5){ 
199     Double_t a,b,c;
200     a=b=c=0;
201     Fit2(fSumDFi, fSumDFiDr, fSumDFiDr2, fSumDr,fSumDr2,fSumDr3,fSumDr4,
202          fStackIndex, a,b,c);
203     if (TMath::Abs(a)<maxdelta){
204       fParam->fFi +=a/fParam->fR;    
205       fParam->fAn = b;    
206       fParam->fAd = c;                  
207     }
208     Fit2(fSumDZ, fSumDZDr, fSumDZDr2, fSumDr,fSumDr2,fSumDr3,fSumDr4,
209          fStackIndex, a,b,c) ;   
210     if (TMath::Abs(a)<maxdelta){
211       fParam->fZ +=a;    
212       fParam->fTheta = b;    
213       fParam->fThetaD = c;   
214     }                         
215   }
216       
217 }
218 void   AliTPCTempHitInfo::Fit2(Double_t fSumY, Double_t fSumYX, Double_t fSumYX2,
219             Double_t fSumX,  Double_t fSumX2, Double_t fSumX3, 
220             Double_t fSumX4, Int_t n,
221             Double_t &a, Double_t &b, Double_t &c)
222 {
223   //fit of second order
224   Double_t det = 
225     n* (fSumX2*fSumX4-fSumX3*fSumX3) -
226     fSumX*      (fSumX*fSumX4-fSumX3*fSumX2)+
227     fSumX2*     (fSumX*fSumX3-fSumX2*fSumX2);
228     
229   if (TMath::Abs(det)> AliTPCTrackHits::fgkPrecision) {    
230     a = 
231       (fSumY * (fSumX2*fSumX4-fSumX3*fSumX3)-
232        fSumX *(fSumYX*fSumX4-fSumYX2*fSumX3)+
233        fSumX2*(fSumYX*fSumX3-fSumYX2*fSumX2))/det; 
234     b=
235       (n*(fSumYX*fSumX4-fSumX3*fSumYX2)-
236       fSumY*(fSumX*fSumX4-fSumX3*fSumX2)+
237       fSumX2*(fSumX*fSumYX2-fSumYX*fSumX2))/det;
238     c=
239       (n*(fSumX2*fSumYX2-fSumYX*fSumX3)-
240        fSumX*(fSumX*fSumYX2-fSumYX*fSumX2)+
241        fSumY*(fSumX*fSumX3-fSumX2*fSumX2))/det;  
242   }
243 }
244
245 void   AliTPCTempHitInfo::Fit(AliTrackHitsParam * param)
246 {
247   // fit fixing first and the last point 
248   //result stored in new param
249   Double_t dx2  = (GetPosition(fStackIndex))[0]-fParam->fR;
250   Double_t det = fSumDr4+dx2*fSumDr2-2*dx2*fSumDr3;
251   if ( (TMath::Abs(det)> AliTPCTrackHits::fgkPrecision) &&
252        ((TMath::Abs(dx2)> AliTPCTrackHits::fgkPrecision))){
253     Double_t dfi2 = (GetPosition(fStackIndex))[1]-fParam->fFi;
254     param->fAd = (fSumDFiDr2+dfi2*fSumDr-dx2*fSumDFiDr-dfi2*fSumDr3/dx2)/det;
255     param->fAn  = (dfi2-param->fAd*dx2*dx2)/dx2;
256     
257     Double_t dz2 = (GetPosition(fStackIndex))[1]-fParam->fZ;
258     param->fTheta = (fSumDZDr2+dz2*fSumDr-dx2*fSumDZDr-dz2*fSumDr3/dx2)/det;
259     param->fTheta  = (dz2-param->fAd*dx2*dx2)/dx2;
260   }
261   
262 }
263
264
265
266
267 AliTPCTrackHits::AliTPCTrackHits()
268 {
269   //
270   //default constructor
271   //
272   const Float_t kHitPrecision=0.002; //default precision for hit position in cm
273   const Float_t kStep =0.003;  //30 mum step 
274   const UShort_t kMaxDistance =100;  //maximum distance 100  
275
276   fPrecision=kHitPrecision; //precision in cm
277   fStep = kStep; //step size
278   fMaxDistance = kMaxDistance; //maximum distance
279   fTempInfo =0;
280   fTrackHitsInfo = new AliObjectArray("AliTrackHitsInfo"); 
281   fTrackHitsParam = new AliObjectArray("AliTrackHitsParam");
282   fHitsPosAndQ = new TArrayOfArray_vStack("AliHitInfo");
283   fCurrentHit = new AliTPCCurrentHit;
284   fgCounter1++;
285   fgCounter2++;
286
287
288
289 AliTPCTrackHits::~AliTPCTrackHits()
290 {
291   //
292   //default destructor
293   //
294   if (fTrackHitsInfo) delete fTrackHitsInfo;
295   if (fTrackHitsParam) delete fTrackHitsParam;
296   if (fHitsPosAndQ) delete fHitsPosAndQ;
297   if (fCurrentHit) delete fCurrentHit;
298   if (fTempInfo) delete fTempInfo;
299   fgCounter1--;
300 }
301
302 void AliTPCTrackHits::Clear()
303 {
304   //
305   //clear object 
306   //  fTrackHitsInfo->Clear();
307   //fTrackHitsParam->Clear();
308   
309   fTrackHitsInfo->Resize(0);
310   fTrackHitsParam->Resize(0);
311   fHitsPosAndQ->Clear();
312
313   if (fTempInfo){
314     delete fTempInfo; 
315     fTempInfo =0;
316   } 
317 }
318
319
320 void AliTPCTrackHits::AddHitKartez(Int_t volumeID, Int_t trackID, Double_t x, 
321               Double_t y, Double_t z,Int_t q)
322 {
323   Double_t r = TMath::Sqrt(x*x+y*y);
324   Double_t fi = TMath::ACos(x/r);
325   if (y<0) fi*=-1.;
326     AddHit(volumeID,trackID,r,z,fi,q);
327 }
328
329 void AliTPCTrackHits::AddHit(Int_t volumeID, Int_t trackID, 
330                              Double_t r, Double_t z, Double_t fi, Int_t q)
331 {
332   //
333   Bool_t diff=kFALSE;
334   if (!fTempInfo) { //initialsation of track
335     fTempInfo = new AliTPCTempHitInfo;
336     //
337     if (fTrackHitsInfo->GetCapacity()<10) fTrackHitsInfo->Reserve(10);
338     fTrackHitsInfo->Resize(1);
339     fTempInfo->fInfoIndex =0;
340     if (fTrackHitsParam->GetCapacity()<100) fTrackHitsParam->Reserve(100);    
341     fTrackHitsParam->Resize(1);
342     //
343     fTempInfo->fInfo = 
344       (AliTrackHitsInfo*) (fTrackHitsInfo->At(0));
345     fTempInfo->fInfo->fVolumeID = volumeID;
346     fTempInfo->fInfo->fTrackID = trackID;
347     fTempInfo->fInfo->fHitParamIndex =0;     
348     fTempInfo->fInfoIndex = 0;
349     //
350     fTempInfo->fParam = 
351       (AliTrackHitsParam*) (fTrackHitsParam->At(0));
352     fTempInfo->fParamIndex = 0;
353     fTempInfo->NewParam(r,z,fi,q);
354     return;
355   }
356   
357   Int_t size = fHitsPosAndQ->GetSize();
358   if (size>(Int_t)fTempInfo->fParamIndex) {
359     fTempInfo->fParamIndex++;
360     if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize()) 
361       fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);    
362     fTempInfo->fParam = 
363       (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));  
364     fTempInfo->NewParam(r,z,fi,q);
365     return;
366   }
367   
368
369   // if new volume or new trackID  
370   if ( (volumeID!=fTempInfo->fInfo->fVolumeID) || 
371        (trackID!=fTempInfo->fInfo->fTrackID)){
372     diff=kTRUE;
373
374     FlushHitStack(kTRUE);
375       
376     fTempInfo->fInfoIndex++;
377     if (fTempInfo->fInfoIndex+1>fTrackHitsInfo->GetSize()) 
378       fTrackHitsInfo->Resize(fTempInfo->fInfoIndex+1);
379     fTempInfo->fInfo = 
380       (AliTrackHitsInfo*) (fTrackHitsInfo->At(fTempInfo->fInfoIndex));      
381     fTempInfo->fInfo->fVolumeID = volumeID;
382     fTempInfo->fInfo->fTrackID = trackID;
383     fTempInfo->fInfo->fHitParamIndex =fTempInfo->fParamIndex+1;  
384     // FlushHitStack(kTRUE);
385
386     fTempInfo->fParamIndex++;
387     if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize()) 
388       fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);    
389     fTempInfo->fParam = 
390       (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));  
391     fTempInfo->NewParam(r,z,fi,q);
392     return;
393   }
394      
395   //calculate current fit precission to next point
396   AliTrackHitsParam &param = *(fTempInfo->fParam);
397   Double_t dd=0;
398   Double_t dl=0;
399   Double_t ratio=0;
400   Double_t dr,dz,dfi,ddz,ddfi;
401   Double_t drhit,ddl;
402   dr=dz=dfi=ddz=ddfi=0;
403   drhit = r-fTempInfo->fOldR;
404   { 
405     //Double_t dfi2 = param.fAn+2*param.fAd*(r-param.fR); 
406     Double_t dfi2 = param.fAn;
407     dfi2*=dfi2*fTempInfo->fOldR*fTempInfo->fOldR;
408     //Double_t ddz2 =  param.fTheta+2*param.fThetaD*(r-param.fR);
409     Double_t ddz2 =  param.fTheta;
410     ddz2*=ddz2;
411     ratio = TMath::Sqrt(1.+ dfi2+ ddz2);  
412   }
413   dl = fStep * Short_t(TMath::Nint(drhit*ratio/fStep));
414   ddl = dl - drhit*ratio; 
415   fTempInfo->fOldR += dl/ratio; 
416
417   if (fTempInfo->fStackIndex>2){     
418     dr = r-param.fR;        
419     dz =  z-param.fZ;  
420     dfi = fi-param.fFi;
421     ddz = dr*param.fTheta+dr*dr*param.fThetaD-dz;
422     ddfi= dr*param.fAn+dr*dr*param.fAd-dfi;    
423     dd  = TMath::Sqrt(ddz*ddz+r*r*ddfi*ddfi+ddl*ddl); 
424     //
425   }        
426   //safety factor 1.25
427   if ( ( (dd*1.25>fPrecision) ) ||  
428        (fTempInfo->fStackIndex+4>fTempInfo->fkStackSize) || 
429        (TMath::Abs(dl/fStep)>fMaxDistance)  ) 
430     diff=kTRUE;
431   else{
432     fTempInfo->fStackIndex++;   
433     fTempInfo->SetHit(r,z,fi,q);
434     return;
435   }  
436   //if parameter changed 
437   if (FlushHitStack(kFALSE)){   //if full buffer flushed
438     fTempInfo->fParamIndex++;
439     if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize()) 
440       fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);    
441     fTempInfo->fParam = 
442       (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));  
443     fTempInfo->NewParam(r,z,fi,q);
444   }
445   else{
446     fTempInfo->fStackIndex++;
447     fTempInfo->SetHit(r,z,fi,q);              
448   }
449 }   
450
451 Bool_t AliTPCTrackHits::FlushHitStack(Bool_t force)
452 {
453   //
454   //write fHitsPosAndQ information from the stack to te arrays
455   if (!fTempInfo) return kFALSE; 
456   Int_t size = fHitsPosAndQ->GetSize();
457
458   if ( (size>0)&&(size!=(Int_t)fTempInfo->fParamIndex)) return kFALSE;  
459
460   if (fHitsPosAndQ->Push(fTempInfo->fStackIndex+1)!=fTempInfo->fParamIndex){
461     cout<<"internal error - contact MI\n";
462     return kFALSE;
463   }
464   AliHitInfo * info;
465  
466   AliTrackHitsParam & param = *(fTempInfo->fParam);
467   //recalculate track parameter not fixing first point
468   fTempInfo->UpdateParam(fStep/4.);
469   //fTempInfo->Fit(fTempInfo->fParam);  //- fixing the first and the last point
470
471   Double_t oldr = param.fR; 
472   //cout<<"C3"<<fTempInfo->fStackIndex<<"\n"<<flush;
473   UInt_t i;
474   Double_t dd;
475   for (i=0; i <= fTempInfo->fStackIndex; i++){
476     Double_t * position = fTempInfo->GetPosition(i);
477     Double_t   dr = position[0]-oldr;
478     Double_t   ratio; 
479     { 
480       //Double_t dfi2 = param.fAn+2*param.fAd*(position[0]-param.fR);
481       Double_t dfi2 = param.fAn;
482       dfi2*=dfi2*oldr*oldr;
483       //Double_t ddz2 =  param.fTheta+2*param.fThetaD*(position[0]-param.fR);
484       Double_t ddz2 =  param.fTheta;
485       ddz2*=ddz2;
486       ratio = TMath::Sqrt(1.+ dfi2+ ddz2);  
487     }
488
489     Double_t dl = fStep*(Short_t)TMath::Nint(dr*ratio/fStep);  
490     dr = dl/ratio; 
491     oldr+=dr;
492     //calculate precission
493     AliTrackHitsParam &param = *(fTempInfo->fParam);    
494     //real deltas
495     Double_t dr1=  position[0]-param.fR;
496     Double_t dz =  position[1]-param.fZ;
497
498     Double_t dfi = position[2]-param.fFi;
499     //extrapolated deltas
500     Double_t dr2 = oldr-param.fR; 
501     Double_t ddr = dr2-dr1;
502     Double_t ddz = dr2*param.fTheta+dr2*dr2*param.fThetaD-dz;
503     Double_t ddfi= dr2*param.fAn+dr2*dr2*param.fAd-dfi;    
504     dd = TMath::Sqrt(ddz*ddz+oldr*oldr*ddfi*ddfi+ddr*ddr); 
505
506     if ( (dd>fPrecision) ){ 
507       if (i==0){
508         param.fAn = 0;
509         param.fAd = 0;
510         param.fTheta =0;
511         param.fThetaD =0;
512         Double_t ddz = dr2*param.fTheta+dr2*dr2*param.fThetaD-dz;
513         Double_t ddfi= dr2*param.fAn+dr2*dr2*param.fAd-dfi;    
514         dl = 0;
515         dd = TMath::Sqrt(ddz*ddz+oldr*oldr*ddfi*ddfi+ddr*ddr); 
516       }
517       else
518         break;
519     }
520
521     info = (AliHitInfo*)(fHitsPosAndQ->At(fTempInfo->fParamIndex,i));
522     info->fHitDistance = Short_t(TMath::Nint(dl/fStep));
523     info->fCharge = Short_t(fTempInfo->fQStack[i]);
524     /*
525     cout<<"C2";
526     cout<<" "<<fTempInfo->fStackIndex<<" \t";
527     cout<<" "<<i<<" \t";
528     cout<<position[0]<<"\t";
529     cout<<position[1]<<"\t"; 
530     cout<<position[2]<<"\t";
531     cout<<param.fAn<<"\t";
532     cout<<param.fTheta<<"\t";
533     cout<<dr1<<"\t"<<ddr<<"\t"<<ddz<<"\t"<<ddfi<<"\t"<<dd<<"\n"<<flush;    
534     */
535   }    
536   
537   if (i<=fTempInfo->fStackIndex){ //if previous iteration not succesfull 
538     fHitsPosAndQ->Resize(fTempInfo->fParamIndex,i);
539     //
540     fTempInfo->fParamIndex++;
541     if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize()) 
542       fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);    
543     fTempInfo->fParam = 
544         (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));    
545     Double_t * p = fTempInfo->GetPosition(i);
546     UInt_t index2 = fTempInfo->fStackIndex;
547     fTempInfo->NewParam(p[0],p[1],p[2],fTempInfo->fQStack[i]);        
548     if (i+1<=index2) FlushHitStack2(i+1,index2);
549
550     if (force) return      FlushHitStack(kTRUE);      
551     return kFALSE;
552   }  
553   return kTRUE;
554
555  
556
557 void AliTPCTrackHits::FlushHitStack2(Int_t index1, Int_t index2)
558 {
559   //
560   // second iteration flush stack
561   // call only for hits where first iteration were not succesfully interpolated  
562   Double_t * positionstack = new Double_t[3*(index2-index1+1)];
563   UInt_t   * qstack        = new UInt_t[index2-index1+1];
564   memcpy(positionstack, &fTempInfo->fPositionStack[3*index1],
565          (3*(index2-index1+1))*sizeof(Double_t));
566   memcpy(qstack, &fTempInfo->fQStack[index1],(index2-index1+1)*sizeof(UInt_t));
567   Double_t *p = positionstack;
568   for (Int_t j=0; j<=index2-index1;j++){ 
569     fTempInfo->fStackIndex++;
570     fTempInfo->SetHit(p[3*j+0],p[3*j+1],p[3*j+2],qstack[j]);
571   }  
572   delete []positionstack;
573   delete []qstack;
574 }
575
576
577    
578
579
580
581
582   
583
584 Bool_t AliTPCTrackHits::First()
585 {
586   //
587   //set Current hit for the first hit
588   //
589   AliTrackHitsInfo *info = (AliTrackHitsInfo *)fTrackHitsInfo->At(0);
590   AliTrackHitsParam *param = (AliTrackHitsParam *)fTrackHitsParam->At(0);
591   AliHitInfo * hinfo = (AliHitInfo *)fHitsPosAndQ->At(0,0);
592
593   if (!(info) || !(param) || !(hinfo) ) {
594     fCurrentHit->fStatus = kFALSE;
595     return kFALSE;
596   }
597
598   fCurrentHit->fInfoIndex  = 0;
599   fCurrentHit->fParamIndex = 0;
600   fCurrentHit->fStackIndex = 0;
601
602   fCurrentHit->fHit.fSector = info->fVolumeID;
603   fCurrentHit->fHit.SetTrack(info->fTrackID);
604   fCurrentHit->fHit.SetX(param->fR*TMath::Cos(param->fFi));
605   fCurrentHit->fHit.SetY(param->fR*TMath::Sin(param->fFi));
606   fCurrentHit->fHit.SetZ(param->fZ); 
607   fCurrentHit->fHit.fQ = hinfo->fCharge;
608    
609   fCurrentHit->fR = param->fR;
610   
611   return fCurrentHit->fStatus = kTRUE;
612 }
613
614 Bool_t AliTPCTrackHits::Next()
615 {
616   //
617   //  
618   if (!(fCurrentHit->fStatus)) 
619     return kFALSE;
620
621   fCurrentHit->fStackIndex++;
622   AliHitInfo * hinfo = (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex,
623                                                       fCurrentHit->fStackIndex);
624   AliTrackHitsInfo *info = (AliTrackHitsInfo *)fTrackHitsInfo->At(fCurrentHit->fInfoIndex);
625   AliTrackHitsParam *param =  (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex);
626
627   if (!hinfo) {
628     hinfo = (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex+1, 0);
629     if (!hinfo) 
630       return fCurrentHit->fStatus = kFALSE;
631     if (hinfo){ 
632       fCurrentHit->fParamIndex++;
633       fCurrentHit->fStackIndex = 0;
634       param = (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex);
635       if (!param) 
636         return fCurrentHit->fStatus = kFALSE;     
637       fCurrentHit->fR = param->fR;
638
639       if ((fCurrentHit->fInfoIndex+1<fTrackHitsInfo->GetSize())
640         &&((info+1)->fHitParamIndex<=fCurrentHit->fParamIndex)){
641         fCurrentHit->fInfoIndex++;
642         info = (AliTrackHitsInfo *)fTrackHitsInfo->At(fCurrentHit->fInfoIndex);
643         if (!info) 
644           return fCurrentHit->fStatus = kFALSE;
645         fCurrentHit->fHit.fSector = info->fVolumeID;
646         fCurrentHit->fHit.SetTrack(info->fTrackID);
647       }
648     }  
649   } 
650   Double_t ratio;
651   { 
652     //    Double_t dfi2 = param->fAn+2*param->fAd*(fCurrentHit->fR-param->fR);
653     Double_t dfi2 = param->fAn;
654     dfi2*=dfi2*fCurrentHit->fR*fCurrentHit->fR;
655     //    Double_t ddz2 = param->fTheta+2*param->fThetaD*(fCurrentHit->fR-param->fR);
656     Double_t ddz2 =  param->fTheta;
657     ddz2*=ddz2;
658     ratio = TMath::Sqrt(1.+ dfi2+ ddz2);  
659   }
660
661   fCurrentHit->fHit.fQ = hinfo->fCharge;
662   fCurrentHit->fR += fStep*hinfo->fHitDistance/ratio;
663   Double_t dR = fCurrentHit->fR - param->fR;
664   //Double_t dR =0;
665   Double_t fi = param->fFi + (param->fAn*dR+param->fAd*dR*dR);
666   Double_t z  = param->fZ + (param->fTheta*dR+param->fThetaD*dR*dR);
667   
668   fCurrentHit->fHit.SetX(fCurrentHit->fR*TMath::Cos(fi));
669   fCurrentHit->fHit.SetY(fCurrentHit->fR*TMath::Sin(fi));
670   fCurrentHit->fHit.SetZ(z);  
671   return kTRUE;
672 }
673   
674 AliTPChit * AliTPCTrackHits::GetHit()
675 {
676   //
677    return (fCurrentHit->fStatus)? &fCurrentHit->fHit:0;
678   //return &fCurrentHit->fHit;
679
680 }  
681
682 AliTrackHitsParam * AliTPCTrackHits::GetParam()
683 {
684   //
685   return (fCurrentHit->fStatus)? (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex) :0;
686
687
688 AliHitInfo * AliTPCTrackHits::GetHitInfo()
689 {
690   //
691   return (fCurrentHit->fStatus)? 
692     (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex,fCurrentHit->fStackIndex) :0;
693
694
695
696