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